ITK  4.9.0
Insight Segmentation and Registration Toolkit
Public Types | Public Member Functions | List of all members
itk::ImageRegionConstIteratorWithIndex< TImage > Class Template Reference

#include <itkImageRegionConstIteratorWithIndex.h>

+ Inheritance diagram for itk::ImageRegionConstIteratorWithIndex< TImage >:
+ Collaboration diagram for itk::ImageRegionConstIteratorWithIndex< TImage >:

Detailed Description

template<typename TImage>
class itk::ImageRegionConstIteratorWithIndex< TImage >

A multi-dimensional iterator templated over image type that walks an image region and is specialized to keep track of its index location.

The "WithIndex" family of iteators was designed for algorithms that use both the values and locations of image pixels in calculations. Unlike ImageRegionIterator, which calculates an index only when requested, ImageRegionIteratorWithIndex maintains its index location as a member variable that is updated during increment and decrement operations. Iteration speed is penalized, but index queries become more efficient.

ImageRegionConstIteratorWithIndex is a multi-dimensional iterator, requiring more information be specified before the iterator can be used than conventional iterators. Whereas the std::vector::iterator from the STL only needs to be passed a pointer to establish the iterator, the multi-dimensional image iterator needs a pointer, the size of the buffer, the size of the region, the start index of the buffer, and the start index of the region. To gain access to this information, ImageRegionConstIteratorWithIndex holds a reference to the image over which it is traversing.

ImageRegionConstIteratorWithIndex assumes a particular layout of the image data. The is arranged in a 1D array as if it were [][][][slice][row][col] with Index[0] = col, Index[1] = row, Index[2] = slice, etc.

operator++ provides a simple syntax for walking around a region of a multidimensional image. operator++ iterates across a row, constraining the movement to within a region of image. When the iterator reaches the boundary of the region along a row, the iterator automatically wraps to the next row, starting at the first pixel in the row that is part of the region. This allows for simple processing loops of the form:

*
* IteratorType it( image, image->GetRequestedRegion() );
*
* it.Begin();
*
* while( ! it.IsAtEnd() )
* {
* it.Set( 100.0 + it.Get() );
* ++it;
* }
*
*

It also can be used for walking in the reverse direction like

*
* IteratorType it( image, image->GetRequestedRegion() );
*
* it.End();
*
* while( !it.IsAtBegin() )
* {
* it.Set( 100.0 );
* --it;
* }
*
*
MORE INFORMATION

For a complete description of the ITK Image Iterators and their API, please see the Iterators chapter in the ITK Software Guide. The ITK Software Guide is available in print and as a free .pdf download from http://www.itk.org.

See Also
ImageConstIterator
ConditionalConstIterator
ConstNeighborhoodIterator
ConstShapedNeighborhoodIterator
ConstSliceIterator
CorrespondenceDataStructureIterator
FloodFilledFunctionConditionalConstIterator
FloodFilledImageFunctionConditionalConstIterator
FloodFilledImageFunctionConditionalIterator
FloodFilledSpatialFunctionConditionalConstIterator
FloodFilledSpatialFunctionConditionalIterator
ImageConstIterator
ImageConstIteratorWithIndex
ImageIterator
ImageIteratorWithIndex
ImageLinearConstIteratorWithIndex
ImageLinearIteratorWithIndex
ImageRandomConstIteratorWithIndex
ImageRandomIteratorWithIndex
ImageRegionConstIterator
ImageRegionConstIteratorWithIndex
ImageRegionExclusionConstIteratorWithIndex
ImageRegionExclusionIteratorWithIndex
ImageRegionIterator
ImageRegionIteratorWithIndex
ImageRegionReverseConstIterator
ImageRegionReverseIterator
ImageReverseConstIterator
ImageReverseIterator
ImageSliceConstIteratorWithIndex
ImageSliceIteratorWithIndex
NeighborhoodIterator
PathConstIterator
PathIterator
ShapedNeighborhoodIterator
SliceIterator
ImageConstIteratorWithIndex
Wiki Examples:
Examples:
Examples/RegistrationITKv3/ModelToImageRegistration1.cxx, Examples/RegistrationITKv4/ModelToImageRegistration1.cxx, and WikiExamples/Iterators/ImageRegionConstIteratorWithIndex.cxx.

Definition at line 127 of file itkImageRegionConstIteratorWithIndex.h.

Public Types

typedef Superclass::AccessorType AccessorType
 
typedef Superclass::ImageType ImageType
 
typedef Superclass::IndexType IndexType
 
typedef
Superclass::InternalPixelType 
InternalPixelType
 
typedef Superclass::OffsetType OffsetType
 
typedef Superclass::PixelContainer PixelContainer
 
typedef
Superclass::PixelContainerPointer 
PixelContainerPointer
 
typedef Superclass::PixelType PixelType
 
typedef Superclass::RegionType RegionType
 
typedef
ImageRegionConstIteratorWithIndex 
Self
 
typedef Superclass::SizeType SizeType
 
typedef
ImageConstIteratorWithIndex
< TImage > 
Superclass
 
- Public Types inherited from itk::ImageConstIteratorWithIndex< TImage >
typedef TImage::AccessorFunctorType AccessorFunctorType
 
typedef TImage::AccessorType AccessorType
 
typedef TImage ImageType
 
typedef TImage::IndexType IndexType
 
typedef IndexType::IndexValueType IndexValueType
 
typedef TImage::InternalPixelType InternalPixelType
 
typedef TImage::OffsetType OffsetType
 
typedef OffsetType::OffsetValueType OffsetValueType
 
typedef TImage::PixelContainer PixelContainer
 
typedef PixelContainer::Pointer PixelContainerPointer
 
typedef TImage::PixelType PixelType
 
typedef TImage::RegionType RegionType
 
typedef ImageConstIteratorWithIndex Self
 
typedef TImage::SizeType SizeType
 
typedef SizeType::SizeValueType SizeValueType
 

Public Member Functions

 ImageRegionConstIteratorWithIndex ()
 
 ImageRegionConstIteratorWithIndex (const TImage *ptr, const RegionType &region)
 
 ImageRegionConstIteratorWithIndex (const ImageConstIteratorWithIndex< TImage > &it)
 
Selfoperator++ ()
 
Selfoperator-- ()
 
- Public Member Functions inherited from itk::ImageConstIteratorWithIndex< TImage >
PixelType Get (void) const
 
const IndexTypeGetIndex () const
 
const RegionTypeGetRegion () const
 
void GoToBegin ()
 
void GoToReverseBegin ()
 
 ImageConstIteratorWithIndex ()
 
 ImageConstIteratorWithIndex (const Self &it)
 
 ImageConstIteratorWithIndex (const TImage *ptr, const RegionType &region)
 
bool IsAtEnd (void) const
 
bool IsAtReverseEnd (void) const
 
 itkLegacyMacro (Self Begin(void) const)
 
 itkLegacyMacro (Self End(void) const)
 
bool operator!= (const Self &it) const
 
bool operator< (const Self &it) const
 
bool operator<= (const Self &it) const
 
Selfoperator= (const Self &it)
 
bool operator== (const Self &it) const
 
bool operator> (const Self &it) const
 
bool operator>= (const Self &it) const
 
bool Remaining ()
 
const PixelTypeValue (void) const
 
virtual ~ImageConstIteratorWithIndex ()
 
void SetIndex (const IndexType &ind)
 

Additional Inherited Members

- Static Public Member Functions inherited from itk::ImageConstIteratorWithIndex< TImage >
static unsigned int GetImageDimension ()
 
- Static Public Attributes inherited from itk::ImageConstIteratorWithIndex< TImage >
static const unsigned int ImageDimension = TImage::ImageDimension
 
- Protected Attributes inherited from itk::ImageConstIteratorWithIndex< TImage >
const InternalPixelTypem_Begin
 
IndexType m_BeginIndex
 
const InternalPixelTypem_End
 
IndexType m_EndIndex
 
TImage::ConstWeakPointer m_Image
 
OffsetValueType m_OffsetTable [ImageDimension+1]
 
AccessorType m_PixelAccessor
 
AccessorFunctorType m_PixelAccessorFunctor
 
const InternalPixelTypem_Position
 
IndexType m_PositionIndex
 
RegionType m_Region
 
bool m_Remaining
 

Member Typedef Documentation

Definition at line 149 of file itkImageRegionConstIteratorWithIndex.h.

template<typename TImage>
typedef Superclass::ImageType itk::ImageRegionConstIteratorWithIndex< TImage >::ImageType

Definition at line 144 of file itkImageRegionConstIteratorWithIndex.h.

template<typename TImage>
typedef Superclass::IndexType itk::ImageRegionConstIteratorWithIndex< TImage >::IndexType

Index typedef support. While these were already typdef'ed in the superclass they need to be redone here for this subclass to compile properly with gcc.Types inherited from the Superclass

Definition at line 140 of file itkImageRegionConstIteratorWithIndex.h.

Definition at line 147 of file itkImageRegionConstIteratorWithIndex.h.

template<typename TImage>
typedef Superclass::OffsetType itk::ImageRegionConstIteratorWithIndex< TImage >::OffsetType

Definition at line 142 of file itkImageRegionConstIteratorWithIndex.h.

Definition at line 145 of file itkImageRegionConstIteratorWithIndex.h.

Definition at line 146 of file itkImageRegionConstIteratorWithIndex.h.

template<typename TImage>
typedef Superclass::PixelType itk::ImageRegionConstIteratorWithIndex< TImage >::PixelType

Definition at line 148 of file itkImageRegionConstIteratorWithIndex.h.

template<typename TImage>
typedef Superclass::RegionType itk::ImageRegionConstIteratorWithIndex< TImage >::RegionType

Definition at line 143 of file itkImageRegionConstIteratorWithIndex.h.

Standard class typedefs.

Definition at line 131 of file itkImageRegionConstIteratorWithIndex.h.

template<typename TImage>
typedef Superclass::SizeType itk::ImageRegionConstIteratorWithIndex< TImage >::SizeType

Definition at line 141 of file itkImageRegionConstIteratorWithIndex.h.

template<typename TImage>
typedef ImageConstIteratorWithIndex< TImage > itk::ImageRegionConstIteratorWithIndex< TImage >::Superclass

Definition at line 132 of file itkImageRegionConstIteratorWithIndex.h.

Constructor & Destructor Documentation

template<typename TImage>
itk::ImageRegionConstIteratorWithIndex< TImage >::ImageRegionConstIteratorWithIndex ( )
inline

Default constructor. Needed since we provide a cast constructor.

Definition at line 152 of file itkImageRegionConstIteratorWithIndex.h.

template<typename TImage>
itk::ImageRegionConstIteratorWithIndex< TImage >::ImageRegionConstIteratorWithIndex ( const TImage *  ptr,
const RegionType region 
)
inline

Constructor establishes an iterator to walk a particular image and a particular region of that image.

Definition at line 156 of file itkImageRegionConstIteratorWithIndex.h.

template<typename TImage>
itk::ImageRegionConstIteratorWithIndex< TImage >::ImageRegionConstIteratorWithIndex ( const ImageConstIteratorWithIndex< TImage > &  it)
inline

Constructor that can be used to cast from an ImageIterator to an ImageRegionConstIteratorWithIndex. Many routines return an ImageIterator but for a particular task, you may want an ImageRegionConstIteratorWithIndex. Rather than provide overloaded APIs that return different types of Iterators, itk returns ImageIterators and uses constructors to cast from an ImageIterator to a ImageRegionConstIteratorWithIndex.

Definition at line 166 of file itkImageRegionConstIteratorWithIndex.h.

References itk::ImageConstIteratorWithIndex< TImage >::operator=().

Member Function Documentation

template<typename TImage>
Self& itk::ImageRegionConstIteratorWithIndex< TImage >::operator++ ( )

Increment (prefix) the fastest moving dimension of the iterator's index. This operator will constrain the iterator within the region (i.e. the iterator will automatically wrap from the end of the row of the region to the beginning of the next row of the region) up until the iterator tries to moves past the last pixel of the region. Here, the iterator will be set to be one pixel past the end of the region.

See Also
operator--
template<typename TImage>
Self& itk::ImageRegionConstIteratorWithIndex< TImage >::operator-- ( )

Decrement (prefix) the fastest moving dimension of the iterator's index. This operator will constrain the iterator within the region (i.e. the iterator will automatically wrap from the beginning of the row of the region to the end of the previous row of the region) up until the iterator tries to moves past the first pixel of the region. Here, the iterator will be set to be one pixel past the beginning of the region.

See Also
operator++

The documentation for this class was generated from the following file: