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

#include <itkNeighborhoodIterator.h>

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

Public Types

using ConstIterator = typename Superclass::ConstIterator
 
using ImageBoundaryConditionPointerType = typename Superclass::ImageBoundaryConditionPointerType
 
using ImageType = typename Superclass::ImageType
 
using IndexType = typename Superclass::IndexType
 
using InternalPixelType = typename Superclass::InternalPixelType
 
using Iterator = typename Superclass::Iterator
 
using NeighborhoodType = typename Superclass::NeighborhoodType
 
using OffsetType = typename Superclass::OffsetType
 
using PixelType = typename Superclass::PixelType
 
using RadiusType = typename Superclass::RadiusType
 
using RegionType = typename Superclass::RegionType
 
using Self = NeighborhoodIterator
 
using SizeType = typename Superclass::SizeType
 
using Superclass = ConstNeighborhoodIterator< TImage, TBoundaryCondition >
 
- Public Types inherited from itk::ConstNeighborhoodIterator< TImage, TBoundaryCondition >
using BoundaryConditionType = TBoundaryCondition
 
using ConstIterator = typename Superclass::ConstIterator
 
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 Superclass::Iterator
 
using NeighborhoodAccessorFunctorType = typename ImageType::NeighborhoodAccessorFunctorType
 
using NeighborhoodType = Neighborhood< PixelType, Self::Dimension >
 
using NeighborIndexType = typename NeighborhoodType::NeighborIndexType
 
using OffsetType = typename Superclass::OffsetType
 
using OutputImageType = typename BoundaryConditionType::OutputImageType
 
using PixelType = typename TImage::PixelType
 
using RadiusType = typename Superclass::RadiusType
 
using RegionType = typename TImage::RegionType
 
using Self = ConstNeighborhoodIterator
 
using SizeType = typename Superclass::SizeType
 
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

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

Additional Inherited Members

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

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;
}

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:
Examples
Examples/Iterators/NeighborhoodIterators6.cxx, SphinxExamples/src/Core/Common/IterateRegionWithNeighborhood/Code.cxx, and SphinxExamples/src/Core/Common/NeighborhoodIteratorOnVectorImage/Code.cxx.

Definition at line 212 of file itkNeighborhoodIterator.h.

Member Typedef Documentation

◆ ConstIterator

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

Definition at line 230 of file itkNeighborhoodIterator.h.

◆ ImageBoundaryConditionPointerType

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

Definition at line 231 of file itkNeighborhoodIterator.h.

◆ ImageType

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

Definition at line 223 of file itkNeighborhoodIterator.h.

◆ IndexType

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

Definition at line 225 of file itkNeighborhoodIterator.h.

◆ InternalPixelType

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

Extract type alias from superclass.

Definition at line 220 of file itkNeighborhoodIterator.h.

◆ Iterator

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

Definition at line 229 of file itkNeighborhoodIterator.h.

◆ NeighborhoodType

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

Definition at line 228 of file itkNeighborhoodIterator.h.

◆ OffsetType

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

Definition at line 226 of file itkNeighborhoodIterator.h.

◆ PixelType

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

Definition at line 221 of file itkNeighborhoodIterator.h.

◆ RadiusType

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

Definition at line 227 of file itkNeighborhoodIterator.h.

◆ RegionType

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

Definition at line 224 of file itkNeighborhoodIterator.h.

◆ 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.

◆ SizeType

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

Definition at line 222 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 264 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.

◆ PrintSelf()

template<typename TImage, typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::PrintSelf ( std::ostream &  ,
Indent   
) const
overridevirtual

◆ SetCenterPixel()

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

Returns the central pixel of the neighborhood.

Definition at line 271 of file itkNeighborhoodIterator.h.

◆ SetNeighborhood()

template<typename TImage, typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
ITK_ITERATOR_VIRTUAL 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>>
ITK_ITERATOR_VIRTUAL void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetNext ( const unsigned  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 315 of file itkNeighborhoodIterator.h.

◆ SetNext() [2/2]

template<typename TImage, typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
ITK_ITERATOR_VIRTUAL void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetNext ( const unsigned  axis,
const unsigned  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 306 of file itkNeighborhoodIterator.h.

◆ SetPixel() [1/3]

template<typename TImage, typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
ITK_ITERATOR_VIRTUAL 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 295 of file itkNeighborhoodIterator.h.

◆ SetPixel() [2/3]

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

Set the pixel at the ith location.

◆ SetPixel() [3/3]

template<typename TImage, typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
ITK_ITERATOR_VIRTUAL void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetPixel ( const unsigned  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>>
ITK_ITERATOR_VIRTUAL void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetPrevious ( const unsigned  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 333 of file itkNeighborhoodIterator.h.

◆ SetPrevious() [2/2]

template<typename TImage, typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
ITK_ITERATOR_VIRTUAL void itk::NeighborhoodIterator< TImage, TBoundaryCondition >::SetPrevious ( const unsigned  axis,
const unsigned  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 324 of file itkNeighborhoodIterator.h.


The documentation for this class was generated from the following file:
itk::Neighborhood< TImage::InternalPixelType *, TImage::ImageDimension >::SizeValueType
typename SizeType::SizeValueType SizeValueType
Definition: itkNeighborhood.h:80
itk::ConstNeighborhoodIterator::GetIndex
ITK_ITERATOR_VIRTUAL IndexType GetIndex() const ITK_ITERATOR_FINAL
Definition: itkConstNeighborhoodIterator.h:177
itk::NeighborhoodIterator
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Definition: itkNeighborhoodIterator.h:212
itk::DerivativeOperator
A NeighborhoodOperator for taking an n-th order derivative at a pixel.
Definition: itkDerivativeOperator.h:69
itk::DerivativeOperator::SetOrder
void SetOrder(const unsigned int &order)
Definition: itkDerivativeOperator.h:82
itk::NeighborhoodInnerProduct< ImageType >