ITK  6.0.0
Insight Toolkit
Public Types | Public Member Functions | Private Attributes | List of all members
itk::ImageLinearConstIteratorWithIndex< TImage > Class Template Reference

#include <itkImageLinearConstIteratorWithIndex.h>

Detailed Description

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

A multi-dimensional image iterator that visits image pixels within a region in a "scan-line" order.

ImageLinearConstIteratorWithIndex is templated over image type and is constrained to walk within a specified image region. It is designed for line-by-line processing of images. This iterator walks a linear path along a selected image direction that is parallel to one of the coordinate axes of the image. The iterator conceptually breaks the image into a set of parallel lines that span the selected image dimension.

ImageLinearConstIteratorWithIndex 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 preselected direction constraining the movement to within a region of image. The user can verify when the iterator reaches the boundary of the region along this direction, by calling the IsAtEndOfLine() method. Then it is possible to pass to the next line starting at the first pixel in the row that is part of the region by calling the NextLine() method.

This is the typical use of this iterator in a loop:

ImageLinearConstIteratorWithIndex<ImageType> it( image, image->GetRequestedRegion() );
it.SetDirection(2);
it.GoToBegin();
while( !it.IsAtEnd() )
{
while( !it.IsAtEndOfLine() )
{
value = it.Get(); // it.Set() doesn't exist in the Const Iterator
++it;
}
it.NextLine();
}
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 https://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

Definition at line 101 of file itkImageLinearConstIteratorWithIndex.h.

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

Public Types

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

Public Member Functions

unsigned int GetDirection ()
 
void GoToBeginOfLine ()
 
void GoToEndOfLine ()
 
void GoToReverseBeginOfLine ()
 
 ImageLinearConstIteratorWithIndex ()
 
 ImageLinearConstIteratorWithIndex (const ImageConstIteratorWithIndex< TImage > &it)
 
 ImageLinearConstIteratorWithIndex (const ImageType *ptr, const RegionType &region)
 
bool IsAtEndOfLine () const
 
bool IsAtReverseEndOfLine () const
 
void NextLine ()
 
void PreviousLine ()
 
void SetDirection (unsigned int direction)
 
Selfoperator++ ()
 
Selfoperator-- ()
 
- Public Member Functions inherited from itk::ImageConstIteratorWithIndex< TImage >
PixelType Get () const
 
const IndexTypeGetIndex () const
 
const RegionTypeGetRegion () const
 
void GoToBegin ()
 
void GoToReverseBegin ()
 
 ImageConstIteratorWithIndex ()=default
 
 ImageConstIteratorWithIndex (const Self &it)
 
 ImageConstIteratorWithIndex (const TImage *ptr, const RegionType &region)
 
bool IsAtEnd () const
 
bool IsAtReverseEnd () const
 
 ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION (Self)
 
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 () const
 
virtual ~ImageConstIteratorWithIndex ()=default
 
void SetIndex (const IndexType &ind)
 

Private Attributes

unsigned int m_Direction { 0 }
 
OffsetValueType m_Jump { 0 }
 

Additional Inherited Members

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

Member Typedef Documentation

◆ ImageType

template<typename TImage >
using itk::ImageLinearConstIteratorWithIndex< TImage >::ImageType = TImage

Image type alias support While this was already typedef'ed in the superclass, it needs to be redone here for this subclass to compile properly with gcc. Note that we have to rescope Index back to itk::Index so that it is not confused with ImageIterator::Index.

Definition at line 124 of file itkImageLinearConstIteratorWithIndex.h.

◆ IndexType

template<typename TImage >
using itk::ImageLinearConstIteratorWithIndex< TImage >::IndexType = typename TImage::IndexType

Index type alias support While this was already typedef'ed in the superclass, it needs to be redone here for this subclass to compile properly with gcc. Note that we have to rescope Index back to itk::Index so that it is not confused with ImageIterator::Index.

Definition at line 112 of file itkImageLinearConstIteratorWithIndex.h.

◆ PixelContainer

template<typename TImage >
using itk::ImageLinearConstIteratorWithIndex< TImage >::PixelContainer = typename TImage::PixelContainer

PixelContainer type alias support Used to refer to the container for the pixel data. While this was already typedef'ed in the superclass, it needs to be redone here for this subclass to compile properly with gcc.

Definition at line 129 of file itkImageLinearConstIteratorWithIndex.h.

◆ PixelContainerPointer

template<typename TImage >
using itk::ImageLinearConstIteratorWithIndex< TImage >::PixelContainerPointer = typename PixelContainer::Pointer

Definition at line 130 of file itkImageLinearConstIteratorWithIndex.h.

◆ RegionType

template<typename TImage >
using itk::ImageLinearConstIteratorWithIndex< TImage >::RegionType = typename TImage::RegionType

Region type alias support While this was already typedef'ed in the superclass, it needs to be redone here for this subclass to compile properly with gcc. Note that we have to rescope Region back to itk::ImageRegion so that it is not confused with ImageIterator::Index.

Definition at line 118 of file itkImageLinearConstIteratorWithIndex.h.

◆ Self

Standard class type aliases.

Definition at line 105 of file itkImageLinearConstIteratorWithIndex.h.

◆ Superclass

template<typename TImage >
using itk::ImageLinearConstIteratorWithIndex< TImage >::Superclass = ImageConstIteratorWithIndex<TImage>

Definition at line 106 of file itkImageLinearConstIteratorWithIndex.h.

Constructor & Destructor Documentation

◆ ImageLinearConstIteratorWithIndex() [1/3]

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

Default constructor. Needed since we provide a cast constructor.

Definition at line 133 of file itkImageLinearConstIteratorWithIndex.h.

◆ ImageLinearConstIteratorWithIndex() [2/3]

template<typename TImage >
itk::ImageLinearConstIteratorWithIndex< TImage >::ImageLinearConstIteratorWithIndex ( const ImageType ptr,
const RegionType region 
)

Constructor establishes an iterator to walk a particular image and a particular region of that image. Initializes the iterator at the begin of the region.

◆ ImageLinearConstIteratorWithIndex() [3/3]

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

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

Definition at line 148 of file itkImageLinearConstIteratorWithIndex.h.

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

Member Function Documentation

◆ GetDirection()

template<typename TImage >
unsigned int itk::ImageLinearConstIteratorWithIndex< TImage >::GetDirection ( )
inline

get the direction of movement

Definition at line 208 of file itkImageLinearConstIteratorWithIndex.h.

◆ GoToBeginOfLine()

template<typename TImage >
void itk::ImageLinearConstIteratorWithIndex< TImage >::GoToBeginOfLine ( )

Go to the beginning pixel of the current line.

See also
GoToReverseBeginOfLine
operator++
operator--
NextLine
IsAtEndOfLine

◆ GoToEndOfLine()

template<typename TImage >
void itk::ImageLinearConstIteratorWithIndex< TImage >::GoToEndOfLine ( )

Go to the past end pixel of the current line.

See also
GoToBeginOfLine
operator++
operator--
NextLine
IsAtEndOfLine

◆ GoToReverseBeginOfLine()

template<typename TImage >
void itk::ImageLinearConstIteratorWithIndex< TImage >::GoToReverseBeginOfLine ( )

Go to the beginning pixel of the current line.

See also
GoToBeginOfLine
operator++
operator--
NextLine
IsAtEndOfLine

◆ IsAtEndOfLine()

template<typename TImage >
bool itk::ImageLinearConstIteratorWithIndex< TImage >::IsAtEndOfLine ( ) const
inline

Test if the index is at the end of line

Definition at line 180 of file itkImageLinearConstIteratorWithIndex.h.

◆ IsAtReverseEndOfLine()

template<typename TImage >
bool itk::ImageLinearConstIteratorWithIndex< TImage >::IsAtReverseEndOfLine ( ) const
inline

Test if the index is at the begin of line

Definition at line 187 of file itkImageLinearConstIteratorWithIndex.h.

◆ NextLine()

template<typename TImage >
void itk::ImageLinearConstIteratorWithIndex< TImage >::NextLine
inline

Go to the next line.

See also
operator++
operator--
IsAtEndOfLine
PreviousLine
End

Definition at line 245 of file itkImageLinearConstIteratorWithIndex.h.

◆ operator++()

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

Increment (prefix) the selected dimension. No bounds checking is performed.

See also
GetIndex
operator--

Definition at line 216 of file itkImageLinearConstIteratorWithIndex.h.

◆ operator--()

template<typename TImage >
Self & itk::ImageLinearConstIteratorWithIndex< TImage >::operator-- ( )
inline

Decrement (prefix) the selected dimension. No bounds checking is performed.

See also
GetIndex
operator++

Definition at line 227 of file itkImageLinearConstIteratorWithIndex.h.

◆ PreviousLine()

template<typename TImage >
void itk::ImageLinearConstIteratorWithIndex< TImage >::PreviousLine
inline

Go to the previous line.

See also
operator++
operator--
IsAtEndOfLine
NextLine
End

Definition at line 281 of file itkImageLinearConstIteratorWithIndex.h.

◆ SetDirection()

template<typename TImage >
void itk::ImageLinearConstIteratorWithIndex< TImage >::SetDirection ( unsigned int  direction)
inline

Set the direction of movement

Definition at line 194 of file itkImageLinearConstIteratorWithIndex.h.

Member Data Documentation

◆ m_Direction

template<typename TImage >
unsigned int itk::ImageLinearConstIteratorWithIndex< TImage >::m_Direction { 0 }
private

Definition at line 237 of file itkImageLinearConstIteratorWithIndex.h.

◆ m_Jump

template<typename TImage >
OffsetValueType itk::ImageLinearConstIteratorWithIndex< TImage >::m_Jump { 0 }
private

Definition at line 236 of file itkImageLinearConstIteratorWithIndex.h.


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