ITK  6.0.0
Insight Toolkit
Public Types | Public Member Functions | List of all members
itk::NeighborhoodIterator< TImage, TBoundaryCondition > Class Template Reference

#include <itkNeighborhoodIterator.h>

Detailed Description

template<typename TImage, typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
class itk::NeighborhoodIterator< TImage, TBoundaryCondition >

Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.

This class is a loose extension of the Standard Template Library (STL) bi-directional iterator concept to masks of pixel neighborhoods within itk::Image objects. This NeighborhoodIterator base class defines simple forward and reverse iteration of an N-dimensional neighborhood mask across an image. Elements within the mask can be accessed like elements within an array.

NeighborhoodIterators are designed to encapsulate some of the complexity of working with image neighborhoods, complexity that would otherwise have to be managed at the algorithmic level. Use NeighborhoodIterators to simplify writing algorithms that perform geometrically localized operations on images (for example, convolution and morphological operations).

To motivate the discussion of NeighborhoodIterators and their use in Itk, consider the following code that takes directional derivatives at each point in an image.

operator->SetOrder(1);
operator->SetDirection(0);
operator->CreateDirectional();
iterator(operator->GetRadius(), myImage, myImage->GetRequestedRegion());
iterator.SetToBegin();
while ( ! iterator.IsAtEnd() )
{
std::cout << "Derivative at index " << iterator.GetIndex() << is <<
innerProduct(iterator, operator) << std::endl;
++iterator;
}
A NeighborhoodOperator for taking an n-th order derivative at a pixel.
void SetOrder(const unsigned int order)
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.

Most of the work for the programmer in the code above is in setting up for the iteration. There are three steps. First an inner product function object is created which will be used to effect convolution with the derivative kernel. Setting up the derivative kernel, DerivativeOperator, involves setting the order and direction of the derivative. Finally, we create an iterator over the RequestedRegion of the itk::Image (see Image) using the radius of the derivative kernel as the size.

Itk iterators only loosely follow STL conventions. Notice that instead of asking myImage for myImage.begin() and myImage.end(), iterator.SetToBegin() and iterator.IsAtEnd() are called. Itk iterators are typically more complex objects than traditional, pointer-style STL iterators, and the increased overhead required to conform to the complete STL API is not always justified.

The API for creating and manipulating a NeighborhoodIterator mimics that of the itk::ImageIterators. Like the itk::ImageIterator, a ConstNeighborhoodIterator is defined on a region of interest in an itk::Image. Iteration is constrained within that region of interest.

A NeighborhoodIterator is constructed as a container of pointers (offsets) to a geometric neighborhood of image pixels. As the central pixel position in the mask is moved around the image, the neighboring pixel pointers (offsets) are moved accordingly.

A pixel neighborhood is defined as a central pixel location and an N-dimensional radius extending outward from that location.

Pixels in a neighborhood can be accessed through a NeighborhoodIterator like elements in an array. For example, a 2D neighborhood with radius 2x1 has indices:

0 1 2 3 4
5 6 7 8 9
10 11 12 13 14

Now suppose a NeighborhoodIterator with the above dimensions is constructed and positioned over a neighborhood of values in an Image:

1.2 1.3 1.8 1.4 1.1
1.8 1.1 0.7 1.0 1.0
2.1 1.9 1.7 1.4 2.0

Shown below is some sample pixel access code and the values that it returns.

SizeValueType c = (SizeValueType) (iterator.Size() / 2); // get offset of center pixel
SizeValueType s = iterator.GetStride(1); // y-dimension step size
std::cout << iterator.GetPixel(7) << std::endl;
std::cout << iterator.GetCenterPixel() << std::endl;
std::cout << iterator.GetPixel(c) << std::endl;
std::cout << iterator.GetPixel(c-1) << std::endl;
std::cout << iterator.GetPixel(c-s) << std::endl;
std::cout << iterator.GetPixel(c-s-1) << std::endl;
std::cout << *iterator[c] << std::endl;

Results:

0.7
0.7
0.7
1.1
1.8
1.3
0.7

Use of GetPixel() is preferred over the *iterator[] form, and can be used without loss of efficiency in most cases. Some variations (subclasses) of NeighborhoodIterators may exist which do not support the latter API. Corresponding SetPixel() methods exist to modify pixel values in non-const NeighborhoodIterators.

NeighborhoodIterators are "bidirectional iterators". They move only in two directions through the data set. These directions correspond to the layout of the image data in memory and not to spatial directions of the N-dimensional itk::Image. Iteration always proceeds along the fastest increasing dimension (as defined by the layout of the image data). For itk::Image this is the first dimension specified (i.e. for 3-dimensional (x,y,z) NeighborhoodIterator proceeds along the x-dimension) (For random access iteration through N-dimensional indices, use RandomAccessNeighborhoodIterator).

Each subclass of a ConstNeighborhoodIterator may also define its own mechanism for iteration through an image. In general, the Iterator does not directly keep track of its spatial location in the image, but uses a set of internal loop variables and offsets to trigger wraps at itk::Image region boundaries, and to identify the end of the itk::Image region.

Todo:

Better support for regions with negative indices.

Add Begin() and End() methods?

See also
DerivativeOperator
NeighborhoodInnerProduct
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
ShapedImageNeighborhoodRange
ITK Sphinx Examples:
  • <a href=
  • <a href=
  • <a href=

Definition at line 212 of file itkNeighborhoodIterator.h.

+ Inheritance diagram for itk::NeighborhoodIterator< TImage, TBoundaryCondition >:
+ Collaboration diagram for itk::NeighborhoodIterator< TImage, TBoundaryCondition >:

Public Types

using Self = NeighborhoodIterator
 
using Superclass = ConstNeighborhoodIterator< TImage, TBoundaryCondition >
 
- Public Types inherited from itk::ConstNeighborhoodIterator< TImage, ZeroFluxNeumannBoundaryCondition< TImage > >
using BoundaryConditionType = ZeroFluxNeumannBoundaryCondition< TImage >
 
using ConstIterator = typename AllocatorType::const_iterator
 
using DimensionValueType = unsigned int
 
using ImageBoundaryConditionConstPointerType = const ImageBoundaryCondition< ImageType, OutputImageType > *
 
using ImageBoundaryConditionPointerType = ImageBoundaryCondition< ImageType, OutputImageType > *
 
using ImageType = TImage
 
using IndexType = Index< Self::Dimension >
 
using InternalPixelType = typename TImage::InternalPixelType
 
using Iterator = typename AllocatorType::iterator
 
using NeighborhoodAccessorFunctorType = typename ImageType::NeighborhoodAccessorFunctorType
 
using NeighborhoodType = Neighborhood< PixelType, Self::Dimension >
 
using NeighborIndexType = typename NeighborhoodType::NeighborIndexType
 
using OffsetType = Offset< VDimension >
 
using OutputImageType = typename BoundaryConditionType::OutputImageType
 
using PixelType = typename TImage::PixelType
 
using RadiusType = itk::Size< VDimension >
 
using RegionType = typename TImage::RegionType
 
using Self = ConstNeighborhoodIterator
 
using SizeType = itk::Size< VDimension >
 
using Superclass = Neighborhood< InternalPixelType *, Self::Dimension >
 
- Public Types inherited from itk::Neighborhood< TImage::InternalPixelType *, TImage::ImageDimension >
using AllocatorType = NeighborhoodAllocator< TImage::InternalPixelType * >
 
using ConstIterator = typename AllocatorType::const_iterator
 
using DimensionValueType = unsigned int
 
using Iterator = typename AllocatorType::iterator
 
using NeighborIndexType = SizeValueType
 
using OffsetType = Offset< VDimension >
 
using PixelType = TImage::InternalPixelType *
 
using RadiusType = itk::Size< VDimension >
 
using Self = Neighborhood
 
using SizeType = itk::Size< VDimension >
 
using SizeValueType = typename SizeType::SizeValueType
 
using SliceIteratorType = SliceIterator< TImage::InternalPixelType *, Self >
 

Public Member Functions

InternalPixelTypeGetCenterPointer ()
 
 NeighborhoodIterator ()
 
 NeighborhoodIterator (const NeighborhoodIterator &n)
 
 NeighborhoodIterator (const SizeType &radius, ImageType *ptr, const RegionType &region)
 
void SetCenterPixel (const PixelType &p)
 
void SetNeighborhood (const NeighborhoodType &)
 
void SetNext (const unsigned int axis, const PixelType &v)
 
void SetNext (const unsigned int axis, const unsigned int i, const PixelType &v)
 
void SetPixel (const unsigned int i, const PixelType &v)
 
void SetPixel (const unsigned int i, const PixelType &v, bool &status)
 
void SetPrevious (const unsigned int axis, const unsigned int i, const PixelType &v)
 
Selfoperator= (const Self &orig)
 
void SetPixel (const OffsetType o, const PixelType &v)
 
void SetPrevious (const unsigned int axis, const PixelType &v)
 
- Public Member Functions inherited from itk::ConstNeighborhoodIterator< TImage, ZeroFluxNeumannBoundaryCondition< TImage > >
OffsetType ComputeInternalIndex (const NeighborIndexType n) const
 
 ConstNeighborhoodIterator ()=default
 
 ConstNeighborhoodIterator (const ConstNeighborhoodIterator &)
 
IndexType GetBeginIndex () const
 
IndexType GetBound () const
 
IndexValueType GetBound (NeighborIndexType n) const
 
ImageBoundaryConditionPointerType GetBoundaryCondition () const
 
RegionType GetBoundingBoxAsImageRegion () const
 
PixelType GetCenterPixel () const
 
const InternalPixelTypeGetCenterPointer () const
 
IndexType GetFastIndexPlusOffset (const OffsetType &o) const
 
const ImageTypeGetImagePointer () const
 
IndexType GetIndex () const
 
IndexType GetIndex (const OffsetType &o) const
 
IndexType GetIndex (NeighborIndexType i) const
 
bool GetNeedToUseBoundaryCondition () const
 
NeighborhoodType GetNeighborhood () const
 
PixelType GetNext (const unsigned int axis) const
 
PixelType GetNext (const unsigned int axis, NeighborIndexType i) const
 
PixelType GetPixel (const NeighborIndexType i) const
 
PixelType GetPixel (const OffsetType &o) const
 
PixelType GetPixel (const OffsetType &o, bool &IsInBounds) const
 
PixelType GetPixel (NeighborIndexType i, bool &IsInBounds) const
 
PixelType GetPrevious (const unsigned int axis) const
 
PixelType GetPrevious (const unsigned int axis, NeighborIndexType i) const
 
RegionType GetRegion () const
 
OffsetType GetWrapOffset () const
 
OffsetValueType GetWrapOffset (NeighborIndexType n) const
 
void GoToBegin ()
 
void GoToEnd ()
 
bool InBounds () const
 
bool IndexInBounds (const NeighborIndexType n) const
 
bool IndexInBounds (const NeighborIndexType n, OffsetType &internalIndex, OffsetType &offset) const
 
void Initialize (const SizeType &radius, const ImageType *ptr, const RegionType &region)
 
bool IsAtBegin () const
 
 ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION (Self)
 
void NeedToUseBoundaryConditionOff ()
 
void NeedToUseBoundaryConditionOn ()
 
Selfoperator++ ()
 
Selfoperator+= (const OffsetType &)
 
OffsetType operator- (const Self &b) const
 
Selfoperator-- ()
 
Selfoperator-= (const OffsetType &)
 
bool operator< (const Self &it) const
 
bool operator<= (const Self &it) const
 
Selfoperator= (const Self &orig)
 
bool operator== (const Self &it) const
 
bool operator> (const Self &it) const
 
bool operator>= (const Self &it) const
 
void OverrideBoundaryCondition (const ImageBoundaryConditionPointerType i)
 
void PrintSelf (std::ostream &, Indent) const override
 
void ResetBoundaryCondition ()
 
void SetBoundaryCondition (const ZeroFluxNeumannBoundaryCondition< TImage > &c)
 
void SetNeedToUseBoundaryCondition (bool b)
 
void SetRegion (const RegionType &region)
 
 ~ConstNeighborhoodIterator () override=default
 
 ConstNeighborhoodIterator (const SizeType &radius, const ImageType *ptr, const RegionType &region)
 
bool IsAtEnd () const
 
void SetLocation (const IndexType &position)
 
- Public Member Functions inherited from itk::Neighborhood< TImage::InternalPixelType *, TImage::ImageDimension >
NeighborIndexType GetCenterNeighborhoodIndex () const
 
TImage::InternalPixelType * GetCenterValue () const
 
virtual const char * GetNameOfClass () const
 
virtual NeighborIndexType GetNeighborhoodIndex (const OffsetType &) const
 
OffsetType GetOffset (NeighborIndexType i) const
 
const SizeType GetRadius () const
 
SizeValueType GetRadius (DimensionValueType n) const
 
SizeType GetSize () const
 
SizeValueType GetSize (DimensionValueType n) const
 
std::slice GetSlice (unsigned int) const
 
OffsetValueType GetStride (DimensionValueType axis) const
 
 ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION (Self)
 
 Neighborhood ()=default
 
 Neighborhood (const Self &)=default
 
 Neighborhood (Self &&)=default
 
Selfoperator= (const Self &)=default
 
Selfoperator= (Self &&)=default
 
bool operator== (const Self &other) const
 
void Print (std::ostream &os) const
 
void SetRadius (const SizeType &)
 
void SetRadius (const SizeValueType)
 
NeighborIndexType Size () const
 
virtual ~Neighborhood ()=default
 
Iterator End ()
 
ConstIterator End () const
 
Iterator Begin ()
 
ConstIterator Begin () const
 
TImage::InternalPixelType * & operator[] (NeighborIndexType i)
 
const TImage::InternalPixelType * & operator[] (NeighborIndexType i) const
 
TImage::InternalPixelType * & GetElement (NeighborIndexType i)
 
TImage::InternalPixelType * & operator[] (const OffsetType &o)
 
const TImage::InternalPixelType * & operator[] (const OffsetType &o) const
 
void SetRadius (const SizeValueType *rad)
 
AllocatorTypeGetBufferReference ()
 
const AllocatorTypeGetBufferReference () const
 

Additional Inherited Members

- Static Public Attributes inherited from itk::ConstNeighborhoodIterator< TImage, ZeroFluxNeumannBoundaryCondition< TImage > >
static constexpr DimensionValueType Dimension
 
- Static Public Attributes inherited from itk::Neighborhood< TImage::InternalPixelType *, TImage::ImageDimension >
static constexpr unsigned int NeighborhoodDimension
 
- Protected Member Functions inherited from itk::ConstNeighborhoodIterator< TImage, ZeroFluxNeumannBoundaryCondition< TImage > >
void SetBeginIndex (const IndexType &start)
 
void SetBound (const SizeType &)
 
void SetEndIndex ()
 
void SetLoop (const IndexType &p)
 
void SetPixelPointers (const IndexType &)
 
- Protected Member Functions inherited from itk::Neighborhood< TImage::InternalPixelType *, TImage::ImageDimension >
virtual void Allocate (NeighborIndexType i)
 
virtual void ComputeNeighborhoodOffsetTable ()
 
virtual void ComputeNeighborhoodStrideTable ()
 
virtual void PrintSelf (std::ostream &, Indent) const
 
void SetSize ()
 
- Protected Attributes inherited from itk::ConstNeighborhoodIterator< TImage, ZeroFluxNeumannBoundaryCondition< TImage > >
const InternalPixelTypem_Begin
 
IndexType m_BeginIndex
 
IndexType m_Bound
 
ImageBoundaryConditionPointerType m_BoundaryCondition
 
ImageType::ConstWeakPointer m_ConstImage
 
const InternalPixelTypem_End
 
IndexType m_EndIndex
 
bool m_InBounds [Dimension]
 
IndexType m_InnerBoundsHigh
 
IndexType m_InnerBoundsLow
 
ZeroFluxNeumannBoundaryCondition< TImage > m_InternalBoundaryCondition
 
bool m_IsInBounds
 
bool m_IsInBoundsValid
 
IndexType m_Loop
 
bool m_NeedToUseBoundaryCondition
 
NeighborhoodAccessorFunctorType m_NeighborhoodAccessorFunctor
 
RegionType m_Region
 
OffsetType m_WrapOffset
 

Member Typedef Documentation

◆ Self

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
using itk::NeighborhoodIterator< TImage, TBoundaryCondition >::Self = NeighborhoodIterator

Standard class type aliases.

Definition at line 216 of file itkNeighborhoodIterator.h.

◆ Superclass

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
using itk::NeighborhoodIterator< TImage, TBoundaryCondition >::Superclass = ConstNeighborhoodIterator<TImage, TBoundaryCondition>

Definition at line 217 of file itkNeighborhoodIterator.h.

Constructor & Destructor Documentation

◆ NeighborhoodIterator() [1/3]

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
itk::NeighborhoodIterator< TImage, TBoundaryCondition >::NeighborhoodIterator ( )
inline

Default constructor.

Definition at line 234 of file itkNeighborhoodIterator.h.

◆ NeighborhoodIterator() [2/3]

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
itk::NeighborhoodIterator< TImage, TBoundaryCondition >::NeighborhoodIterator ( const NeighborhoodIterator< TImage, TBoundaryCondition > &  n)
inline

Copy constructor

Definition at line 239 of file itkNeighborhoodIterator.h.

◆ NeighborhoodIterator() [3/3]

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
itk::NeighborhoodIterator< TImage, TBoundaryCondition >::NeighborhoodIterator ( const SizeType radius,
ImageType ptr,
const RegionType region 
)
inline

Constructor which establishes the region size, neighborhood, and image over which to walk.

Definition at line 254 of file itkNeighborhoodIterator.h.

Member Function Documentation

◆ GetCenterPointer()

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
InternalPixelType * itk::NeighborhoodIterator< TImage, TBoundaryCondition >::GetCenterPointer ( )
inline

Returns the central memory pointer of the neighborhood.

Definition at line 260 of file itkNeighborhoodIterator.h.

◆ operator=()

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
Self & itk::NeighborhoodIterator< TImage, TBoundaryCondition >::operator= ( const Self orig)
inline

Assignment operator

Definition at line 245 of file itkNeighborhoodIterator.h.

◆ SetCenterPixel()

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetCenterPixel ( const PixelType p)
inline

Returns the central pixel of the neighborhood.

Definition at line 267 of file itkNeighborhoodIterator.h.

◆ SetNeighborhood()

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetNeighborhood ( const NeighborhoodType )

Virtual function that replaces the pixel values in the image neighborhood that are pointed to by this NeighborhoodIterator with the pixel values contained in a Neighborhood.

◆ SetNext() [1/2]

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetNext ( const unsigned int  axis,
const PixelType v 
)
inline

Sets the pixel value located one pixel distant from the neighborhood center in the specified positive axis direction. No bounds checking is done on the size of the neighborhood.

Definition at line 311 of file itkNeighborhoodIterator.h.

◆ SetNext() [2/2]

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetNext ( const unsigned int  axis,
const unsigned int  i,
const PixelType v 
)
inline

Sets the pixel value located i pixels distant from the neighborhood center in the positive specified "axis" direction. No bounds checking is done on the size of the neighborhood.

Definition at line 302 of file itkNeighborhoodIterator.h.

◆ SetPixel() [1/3]

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetPixel ( const OffsetType  o,
const PixelType v 
)
inline

Set the pixel at offset o from the neighborhood center

Definition at line 291 of file itkNeighborhoodIterator.h.

◆ SetPixel() [2/3]

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetPixel ( const unsigned int  i,
const PixelType v 
)

Set the pixel at the ith location.

◆ SetPixel() [3/3]

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetPixel ( const unsigned int  i,
const PixelType v,
bool &  status 
)

Special SetPixel method which quietly ignores out-of-bounds attempts. Sets status TRUE if pixel has been set, FALSE otherwise.

◆ SetPrevious() [1/2]

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetPrevious ( const unsigned int  axis,
const PixelType v 
)
inline

Sets the pixel value located one pixel distant from the neighborhood center in the specified negative axis direction. No bounds checking is done on the size of the neighborhood.

Definition at line 329 of file itkNeighborhoodIterator.h.

◆ SetPrevious() [2/2]

template<typename TImage , typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetPrevious ( const unsigned int  axis,
const unsigned int  i,
const PixelType v 
)
inline

Sets the pixel value located i pixels distant from the neighborhood center in the negative specified "axis" direction. No bounds checking is done on the size of the neighborhood.

Definition at line 320 of file itkNeighborhoodIterator.h.


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