KWStyle - itkImageSliceConstIteratorWithIndex.h
 
Matrix View
Description

1 /*=========================================================================
2
3   Program:   Insight Segmentation & Registration Toolkit
4   Module:    $RCSfile: itkImageSliceConstIteratorWithIndex.h.html,v $
5   Language:  C++
6   Date:      $Date: 2006/01/17 19:15:40 $
7   Version:   $Revision: 1.4 $
8
9   Copyright (c) Insight Software Consortium. All rights reserved.
10   See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
11
12      This software is distributed WITHOUT ANY WARRANTY; without even 
13      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 
14      PURPOSE.  See the above copyright notices for more information.
15
16 =========================================================================*/
17 #ifndef __itkImageSliceConstIteratorWithIndex_h
18 #define __itkImageSliceConstIteratorWithIndex_h
19
20 #include "itkImageConstIteratorWithIndex.h"
21
22 namespace itk
23 {
24
25 /** \class ImageSliceConstIteratorWithIndex
26  * \brief Multi-dimensional image iterator which only walks a region.
27  *
28  * \brief A multi-dimensional image iterator that extends the
29  * ImageLinearConstIteratorWithIndex from iteration along lines in an image to
30  * iteration along both lines and planes (slices) within an image.  A slice is
31  * defined as a 2D plane spanned by two vectors pointing along orthogonal
32  * coordinate axes. The slice orientation of the iterator is defined by
33  * specifying its two spanning axes using the methods:
34  * \code
35  * SetFirstDirection(n)
36  * SetSecondDirection(n)
37  * \endcode
38  * where n is the number of the axis.
39  *
40  * Use the following methods to move the iterator between slices:
41  * \code
42  * NextSlice()
43  * PreviousSlice()
44  * \endcode
45  *
46  * To test the position of the iterator with respect to the end or beginning of
47  * the slice use the following methods:
48  * \code
49  * IsAtReverseEndOfSlice()
50  * IsAtEndOfSlice()
51  * \endcode
52  *
53  * The following code, for example, illustrates the typical use of this
54  * iterator.  For more information please see the Software Guide.
55  *
56  * \code
57  *  
58 LEN  * ImageSliceConstIteratorWithIndex<ImageType> it( image, image->GetRequestedRegion() );
59  * 
60  * it.SetFirstDirection(2);
61  * it.SetSecondDirection(0);
62  *
63  * it.GoToBegin();
64  * while( !it.IsAtEnd() )
65  * {
66 IND **   while( !it.IsAtEndOfSlice() )
67 IND **   {
68 IND **     while( !it.IsAtEndOfLine() )
69 IND **     {
70 IND **        value = it.Get();  // it.Set() doesn't exist in the Const Iterator
71 IND **        ++it;
72 IND **     }
73 IND **     it.NextLine();
74 IND **   }
75 IND **   it.NextSlice();
76 IND **  } 
77  *
78  *  \endcode
79  *
80  * \example  Common/itkImageSliceIteratorTest.cxx
81 IND ***
82  * \par MORE INFORMATION
83  * For a complete description of the ITK Image Iterators and their API, please
84  * see the Iterators chapter in the ITK Software Guide.  The ITK Software Guide
85  * is available in print and as a free .pdf download from http://www.itk.org.
86  *
87  * \ingroup ImageIterators
88  *
89  * \sa ImageConstIterator \sa ConditionalConstIterator
90  * \sa ConstNeighborhoodIterator \sa ConstShapedNeighborhoodIterator
91  * \sa ConstSliceIterator  \sa CorrespondenceDataStructureIterator 
92  * \sa FloodFilledFunctionConditionalConstIterator 
93  * \sa FloodFilledImageFunctionConditionalConstIterator 
94  * \sa FloodFilledImageFunctionConditionalIterator 
95  * \sa FloodFilledSpatialFunctionConditionalConstIterator 
96  * \sa FloodFilledSpatialFunctionConditionalIterator 
97  * \sa ImageConstIterator \sa ImageConstIteratorWithIndex 
98  * \sa ImageIterator \sa ImageIteratorWithIndex
99  * \sa ImageLinearConstIteratorWithIndex  \sa ImageLinearIteratorWithIndex 
100  * \sa ImageRandomConstIteratorWithIndex  \sa ImageRandomIteratorWithIndex 
101  * \sa ImageRegionConstIterator \sa ImageRegionConstIteratorWithIndex 
102  * \sa ImageRegionExclusionConstIteratorWithIndex 
103  * \sa ImageRegionExclusionIteratorWithIndex 
104  * \sa ImageRegionIterator  \sa ImageRegionIteratorWithIndex 
105  * \sa ImageRegionReverseConstIterator  \sa ImageRegionReverseIterator 
106  * \sa ImageReverseConstIterator  \sa ImageReverseIterator 
107  * \sa ImageSliceConstIteratorWithIndex  \sa ImageSliceIteratorWithIndex 
108  * \sa NeighborhoodIterator \sa PathConstIterator  \sa PathIterator 
109  * \sa ShapedNeighborhoodIterator  \sa SliceIterator 
110  * \sa ImageConstIteratorWithIndex */
111 template<typename TImage>
112 LEN class ITK_EXPORT ImageSliceConstIteratorWithIndex : public ImageConstIteratorWithIndex<TImage>
113 {
114 public:
115   /** Standard class typedefs. */
116   typedef ImageSliceConstIteratorWithIndex Self;
117 TDA   typedef ImageConstIteratorWithIndex<TImage>  Superclass;
118   
119   /** Index typedef support. While this was already typdef'ed in the superclass
120    * it needs to be redone here for this subclass to compile properly with gcc.
121    * Note that we have to rescope Index back to itk::Index to that is it not
122    * confused with ImageIterator::Index. */
123   typedef typename TImage::IndexType IndexType;
124
125   /** Image typedef support. While this was already typdef'ed in the superclass
126    * it needs to be redone here for this subclass to compile properly with gcc.
127    * Note that we have to rescope Index back to itk::Index to that is it not
128    * confused with ImageIterator::Index. */
129   typedef TImage ImageType;
130
131   /** Region typedef support. */
132   typedef typename TImage::RegionType   RegionType;
133
134   /** PixelContainer typedef support. Used to refer to the container for
135    * the pixel data. While this was already typdef'ed in the superclass
136 LEN    * it needs to be redone here for this subclass to compile properly with gcc. */
137   typedef typename TImage::PixelContainer PixelContainer;
138 TDA   typedef typename PixelContainer::Pointer PixelContainerPointer;
139   
140   /** Default constructor. Needed since we provide a cast constructor. */
141   ImageSliceConstIteratorWithIndex() : ImageConstIteratorWithIndex<TImage>() {}
142   
143   /** Constructor establishes an iterator to walk a particular image and a
144    * particular region of that image. */
145   ImageSliceConstIteratorWithIndex( const ImageType *ptr,
146                       const RegionType & region)
147 IND ****: ImageConstIteratorWithIndex<TImage>(ptr, region) 
148     {
149 IND ******m_Direction_A = 0;
150 IND ******m_Direction_B = 1;
151     }
152
153   /** Constructor that can be used to cast from an ImageIterator to an
154 LEN    * ImageSliceConstIteratorWithIndex. Many routines return an ImageIterator but for a
155 LEN    * particular task, you may want an ImageSliceConstIteratorWithIndex.  Rather than
156    * provide overloaded APIs that return different types of Iterators, itk
157    * returns ImageIterators and uses constructors to cast from an
158    * ImageIterator to a ImageSliceConstIteratorWithIndex. */
159 LEN   ImageSliceConstIteratorWithIndex( const ImageConstIteratorWithIndex<TImage> &it)
160     { this->ImageConstIteratorWithIndex<TImage>::operator=(it); }
161
162   /** Go to the next line
163    * \sa operator++ \sa EndOfLine \sa End \sa NextSlice */
164   void NextLine(void);
165   
166   /** Go to the next slice
167    * \sa operator++ \sa EndOfLine \sa End */
168   void NextSlice(void);
169
170   /** Go to the next line
171    * \sa operator-- \sa BeginOfLine \sa BeginOfSlice \sa Begin */
172   void PreviousLine(void);
173   
174   /** Go to the next slice
175    * \sa operator-- \sa BeginOfLine \sa BeginOfSlice \sa Begin */
176   void PreviousSlice(void);
177
178   /** Test if the index is at the end of line */
179   bool IsAtEndOfLine(void);
180
181   /** Test if the index is at the end of the slice */
182   bool IsAtEndOfSlice(void);
183
184   /** Test if the index is at the begin of line */
185   bool IsAtReverseEndOfLine(void);
186
187   /** Test if the index is at the begin of the slice */
188   bool IsAtReverseEndOfSlice(void);
189
190   /** Set the fastest direction of movement */
191   void SetFirstDirection(unsigned int direction);
192
193   /** Set the second fastest direction of movement */
194   void SetSecondDirection(unsigned int direction);
195
196   /** Increment (prefix) the selected dimension.
197    * No bounds checking is performed. 
198    * \sa operator-- \sa GetIndex */
199   inline Self & operator++();
200
201   /** Decrement (prefix) the selected dimension.
202    * No bounds checking is performed. 
203    * \sa operator++ \sa GetIndex */
204   inline Self & operator--();
205
206 private:
207   unsigned long  m_PixelJump;
208   unsigned long  m_LineJump;
209   unsigned int   m_Direction_A;
210   unsigned int   m_Direction_B;
211 };
212
213 // end namespace itk
214
215 #ifndef ITK_MANUAL_INSTANTIATION
216 #include "itkImageSliceConstIteratorWithIndex.txx"
217 #endif
218
219 #endif 
220

Generated by KWStyle 1.0b on Tuesday January,17 at 02:14:15PM
© Kitware Inc.