1 /*=========================================================================
2  *
3  *  Copyright Insight Software Consortium
4  *
5  *  Licensed under the Apache License, Version 2.0 (the "License");
6  *  you may not use this file except in compliance with the License.
7  *  You may obtain a copy of the License at
8  *
9  *         http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  *  Unless required by applicable law or agreed to in writing, software
12  *  distributed under the License is distributed on an "AS IS" BASIS,
13  *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  *  See the License for the specific language governing permissions and
15  *  limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkImageLinearConstIteratorWithIndex_h
19 #define itkImageLinearConstIteratorWithIndex_h
20 
21 #include "itkImageConstIteratorWithIndex.h"
22 
23 namespace itk
24 {
25 /** \class ImageLinearConstIteratorWithIndex
26  * \brief A multi-dimensional image iterator that visits image pixels within a
27  * region in a "scan-line" order.
28  *
29  * ImageLinearConstIteratorWithIndex is templated over image type and is
30  * constrained to walk within a specified image region. It is designed for
31  * line-by-line processing of images.  This iterator walks a linear path along
32  * a selected image direction that is parallel to one of the coordinate axes
33  * of the image.  The iterator conceptually breaks the image into a set of
34  * parallel lines that span the selected image dimension.
35  *
36  * ImageLinearConstIteratorWithIndex assumes a particular layout of the image
37  * data. The is arranged in a 1D array as if it were [][][][slice][row][col]
38  * with Index[0] = col, Index[1] = row, Index[2] = slice, etc.
39  *
40  * operator++ provides a simple syntax for walking around a region of a
41  * multidimensional image. operator++ iterates across a preselected direction
42  * constraining the movement to within a region of image. The user can verify
43  * when the iterator reaches the boundary of the region along this direction,
44  * by calling the IsAtEndOfLine() method. Then it is possible to pass to the
45  * next line starting at the first pixel in the row that is part of the region
46  * by calling the NextLine() method.
47  *
48  * This is the typical use of this iterator in a loop:
49  *
50    \code
51 
52    ImageLinearConstIteratorWithIndex<ImageType> it( image, image->GetRequestedRegion() );
53 
54    it.SetDirection(2);
55    it.GoToBegin();
56    while( !it.IsAtEnd() )
57    {
58      while( !it.IsAtEndOfLine() )
59      {
60         value = it.Get();  // it.Set() doesn't exist in the Const Iterator
61         ++it;
62      }
63      it.NextLine();
64     }
65 
66     \endcode
67  *
68  * \par MORE INFORMATION
69  * For a complete description of the ITK Image Iterators and their API, please
70  * see the Iterators chapter in the ITK Software Guide.  The ITK Software Guide
71  * is available in print and as a free .pdf download from https://www.itk.org.
72  *
73  * \ingroup ImageIterators
74  *
75  * \sa ImageConstIterator \sa ConditionalConstIterator
76  * \sa ConstNeighborhoodIterator \sa ConstShapedNeighborhoodIterator
77  * \sa ConstSliceIterator  \sa CorrespondenceDataStructureIterator
78  * \sa FloodFilledFunctionConditionalConstIterator
79  * \sa FloodFilledImageFunctionConditionalConstIterator
80  * \sa FloodFilledImageFunctionConditionalIterator
81  * \sa FloodFilledSpatialFunctionConditionalConstIterator
82  * \sa FloodFilledSpatialFunctionConditionalIterator
83  * \sa ImageConstIterator \sa ImageConstIteratorWithIndex
84  * \sa ImageIterator \sa ImageIteratorWithIndex
85  * \sa ImageLinearConstIteratorWithIndex  \sa ImageLinearIteratorWithIndex
86  * \sa ImageRandomConstIteratorWithIndex  \sa ImageRandomIteratorWithIndex
87  * \sa ImageRegionConstIterator \sa ImageRegionConstIteratorWithIndex
88  * \sa ImageRegionExclusionConstIteratorWithIndex
89  * \sa ImageRegionExclusionIteratorWithIndex
90  * \sa ImageRegionIterator  \sa ImageRegionIteratorWithIndex
91  * \sa ImageRegionReverseConstIterator  \sa ImageRegionReverseIterator
92  * \sa ImageReverseConstIterator  \sa ImageReverseIterator
93  * \sa ImageSliceConstIteratorWithIndex  \sa ImageSliceIteratorWithIndex
94  * \sa NeighborhoodIterator \sa PathConstIterator  \sa PathIterator
95  * \sa ShapedNeighborhoodIterator  \sa SliceIterator
96  * \sa ImageConstIteratorWithIndex
97  *
98  * \ingroup ITKCommon
99  */
100 template< typename TImage >
101 class ITK_TEMPLATE_EXPORT ImageLinearConstIteratorWithIndex:public ImageConstIteratorWithIndex< TImage >
102 {
103 public:
104   /** Standard class type aliases. */
105   using Self = ImageLinearConstIteratorWithIndex;
106   using Superclass = ImageConstIteratorWithIndex< TImage >;
107 
108   /** Index type alias support While this was already typdef'ed in the superclass,
109    * it needs to be redone here for this subclass to compile properly with gcc.
110    * Note that we have to rescope Index back to itk::Index so that it is not
111    * confused with ImageIterator::Index. */
112   using IndexType = typename TImage::IndexType;
113 
114   /** Region type alias support While this was already typdef'ed in the superclass,
115    * it needs to be redone here for this subclass to compile properly with gcc.
116    * Note that we have to rescope Region back to itk::ImageRegion so that it
117    * is not confused with ImageIterator::Index. */
118   using RegionType = typename TImage::RegionType;
119 
120   /** Image type alias support While this was already typdef'ed in the superclass,
121    * it needs to be redone here for this subclass to compile properly with gcc.
122    * Note that we have to rescope Index back to itk::Index so that it is not
123    * confused with ImageIterator::Index. */
124   using ImageType = TImage;
125 
126   /** PixelContainer type alias support Used to refer to the container for
127    * the pixel data. While this was already typdef'ed in the superclass,
128    * it needs to be redone here for this subclass to compile properly with gcc. */
129   using PixelContainer = typename TImage::PixelContainer;
130   using PixelContainerPointer = typename PixelContainer::Pointer;
131 
132   /** Default constructor. Needed since we provide a cast constructor. */
ImageLinearConstIteratorWithIndex()133   ImageLinearConstIteratorWithIndex() :
134     ImageConstIteratorWithIndex< TImage >()
135 
136   {}
137 
138   /** Constructor establishes an iterator to walk a particular image and a
139    * particular region of that image. */
140   ImageLinearConstIteratorWithIndex(const ImageType *ptr, const RegionType & region);
141 
142   /** Constructor that can be used to cast from an ImageIterator to an
143    * ImageLinearConstIteratorWithIndex. Many routines return an ImageIterator but for a
144    * particular task, you may want an ImageLinearConstIteratorWithIndex.  Rather than
145    * provide overloaded APIs that return different types of Iterators, itk
146    * returns ImageIterators and uses constructors to cast from an
147    * ImageIterator to a ImageLinearConstIteratorWithIndex. */
ImageLinearConstIteratorWithIndex(const ImageConstIteratorWithIndex<TImage> & it)148   ImageLinearConstIteratorWithIndex(const ImageConstIteratorWithIndex< TImage > & it) : m_Direction(0)
149   { this->ImageConstIteratorWithIndex< TImage >::operator=(it); }
150 
151   /** Go to the next line.
152    * \sa operator++  \sa operator-- \sa IsAtEndOfLine \sa PreviousLine \sa End */
153   inline void NextLine();
154 
155   /** Go to the previous line.
156    * \sa operator++ \sa operator-- \sa IsAtEndOfLine \sa NextLine \sa End */
157   inline void PreviousLine();
158 
159   /** Go to the beginning pixel of the current line.
160    * \sa GoToReverseBeginOfLine \sa operator++ \sa operator-- \sa NextLine \sa IsAtEndOfLine */
161   void GoToBeginOfLine();
162 
163   /** Go to the beginning pixel of the current line.
164    * \sa GoToBeginOfLine \sa operator++ \sa operator-- \sa NextLine \sa IsAtEndOfLine */
165   void GoToReverseBeginOfLine();
166 
167   /** Go to the past end pixel of the current line.
168    * \sa GoToBeginOfLine \sa operator++ \sa operator-- \sa NextLine \sa IsAtEndOfLine */
169   void GoToEndOfLine();
170 
171   /** Test if the index is at the end of line */
IsAtEndOfLine()172   inline bool IsAtEndOfLine() const
173   {
174     return this->m_PositionIndex[m_Direction] >= this->m_EndIndex[m_Direction];
175   }
176 
177   /** Test if the index is at the begin of line */
IsAtReverseEndOfLine()178   inline bool IsAtReverseEndOfLine() const
179   {
180     return this->m_PositionIndex[m_Direction] < this->m_BeginIndex[m_Direction];
181   }
182 
183   /** Set the direction of movement */
SetDirection(unsigned int direction)184   inline void SetDirection(unsigned int direction)
185   {
186     if ( direction >= TImage::ImageDimension )
187       {
188       itkGenericExceptionMacro(
189         << "In image of dimension " << TImage::ImageDimension << " Direction " << direction << " sas selected");
190       }
191     m_Direction = direction;
192     m_Jump = this->m_OffsetTable[m_Direction];
193   }
194 
195   /** get the direction of movement */
GetDirection()196   unsigned int GetDirection()
197   {
198     return m_Direction;
199   }
200 
201   /** Increment (prefix) the selected dimension.
202    * No bounds checking is performed. \sa GetIndex \sa operator-- */
203   inline Self & operator++()
204   {
205     this->m_PositionIndex[m_Direction]++;
206     this->m_Position += m_Jump;
207     return *this;
208   }
209 
210   /** Decrement (prefix) the selected dimension.
211    * No bounds checking is performed.  \sa GetIndex \sa operator++ */
212   inline Self & operator--()
213   {
214     this->m_PositionIndex[m_Direction]--;
215     this->m_Position -= m_Jump;
216     return *this;
217   }
218 
219 private:
220   OffsetValueType m_Jump{0};
221   unsigned int    m_Direction{0};
222 };
223 
224 //----------------------------------------------------------------------
225 //  Go to next line
226 //----------------------------------------------------------------------
227 template< typename TImage >
228 inline
229 void
230 ImageLinearConstIteratorWithIndex< TImage >
NextLine()231 ::NextLine()
232 {
233   this->m_Position -= this->m_OffsetTable[m_Direction]
234                       * ( this->m_PositionIndex[m_Direction] - this->m_BeginIndex[m_Direction] );
235 
236   this->m_PositionIndex[m_Direction] = this->m_BeginIndex[m_Direction];
237 
238   for ( unsigned int n = 0; n < TImage::ImageDimension; n++ )
239   {
240     this->m_Remaining = false;
241 
242     if ( n == m_Direction )
243     {
244       continue;
245     }
246 
247     this->m_PositionIndex[n]++;
248     if ( this->m_PositionIndex[n] <  this->m_EndIndex[n] )
249     {
250       this->m_Position += this->m_OffsetTable[n];
251       this->m_Remaining = true;
252       break;
253     }
254     else
255     {
256       this->m_Position -= this->m_OffsetTable[n] * ( this->m_Region.GetSize()[n] - 1 );
257       this->m_PositionIndex[n] = this->m_BeginIndex[n];
258     }
259   }
260 }
261 
262 //----------------------------------------------------------------------
263 //  Pass to the last pixel on the previous line
264 //----------------------------------------------------------------------
265 template< typename TImage >
266 inline
267 void
268 ImageLinearConstIteratorWithIndex< TImage >
PreviousLine()269 ::PreviousLine()
270 {
271   this->m_Position += this->m_OffsetTable[m_Direction]
272                       * ( this->m_EndIndex[m_Direction] - 1 - this->m_PositionIndex[m_Direction] );
273 
274   this->m_PositionIndex[m_Direction] = this->m_EndIndex[m_Direction] - 1;
275 
276   for ( unsigned int n = 0; n < TImage::ImageDimension; n++ )
277   {
278     this->m_Remaining = false;
279 
280     if ( n == m_Direction )
281     {
282       continue;
283     }
284 
285     this->m_PositionIndex[n]--;
286     if ( this->m_PositionIndex[n] >=  this->m_BeginIndex[n] )
287     {
288       this->m_Position -= this->m_OffsetTable[n];
289       this->m_Remaining = true;
290       break;
291     }
292     else
293     {
294       this->m_Position += this->m_OffsetTable[n] * ( this->m_Region.GetSize()[n] - 1 );
295       this->m_PositionIndex[n] = this->m_EndIndex[n] - 1;
296     }
297   }
298 }
299 } // end namespace itk
300 
301 #ifndef ITK_MANUAL_INSTANTIATION
302 #include "itkImageLinearConstIteratorWithIndex.hxx"
303 #endif
304 
305 #endif
306