Iterate Region in Image Without Write Access#

Synopsis#

Iterate over a region of an image (without write access).

Results#

../../../../_images/Yinyang7.png

Yinyang.png#

Output:

255
255
255
255
255
255
255
255
255
255
255
255
255
255
255
255
255
255
255
255

Code#

C++#

#include "itkImage.h"
#include "itkImageFileReader.h"
#include "itkImageRegionConstIterator.h"

int
main(int argc, char * argv[])
{
  if (argc < 2)
  {
    std::cerr << "Required: filename" << std::endl;
    return EXIT_FAILURE;
  }

  using ImageType = itk::Image<unsigned char, 2>;

  ImageType::Pointer image = itk::ReadImage<ImageType>(argv[1]);

  ImageType::SizeType regionSize;
  regionSize[0] = 5;
  regionSize[1] = 4;

  ImageType::IndexType regionIndex;
  regionIndex[0] = 0;
  regionIndex[1] = 0;

  ImageType::RegionType region;
  region.SetSize(regionSize);
  region.SetIndex(regionIndex);

  itk::ImageRegionConstIterator<ImageType> imageIterator(image, region);
  // itk::ImageRegionConstIterator<ImageType> imageIterator(image,image->GetLargestPossibleRegion());

  while (!imageIterator.IsAtEnd())
  {
    // Get the value of the current pixel
    unsigned char val = imageIterator.Get();
    std::cout << (int)val << std::endl;

    ++imageIterator;
  }

  return EXIT_SUCCESS;
}

Classes demonstrated#

template<typename TImage>
class ImageRegionConstIterator : public itk::ImageConstIterator<TImage>

A multi-dimensional iterator templated over image type that walks a region of pixels.

The itk::ImageRegionConstIterator is optimized for iteration speed and is the first choice for iterative, pixel-wise operations on an image. ImageRegionIterator is the least specialized of the ITK image iterator classes. ImageRegionConstIterator is templated over the image type, and is constrained to walk only within the specified region and along a line parallel to one of the coordinate axes, “wrapping” to the next line as it reaches the boundary of the image. To walk the entire image, specify the region to be image->GetRequestedRegion().

ImageRegionConstIterator provides read-only access to image data. It is the base class for the read/write access ImageRegionIterator.

ImageRegionConstIterator 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, ImageRegionConstIterator holds a reference to the image over which it is traversing.

ImageRegionConstIterator 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:

it = it.Begin();
for (; !it.IsAtEnd(); ++it)
   {
   *it += 100.0;
   }

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

ImageConstIterator

See

ConditionalConstIterator

See

ConstNeighborhoodIterator

See

ConstShapedNeighborhoodIterator

See

ConstSliceIterator

See

CorrespondenceDataStructureIterator

See

FloodFilledFunctionConditionalConstIterator

See

FloodFilledImageFunctionConditionalConstIterator

See

FloodFilledImageFunctionConditionalIterator

See

FloodFilledSpatialFunctionConditionalConstIterator

See

FloodFilledSpatialFunctionConditionalIterator

See

ImageConstIterator

See

ImageConstIteratorWithIndex

See

ImageIterator

See

ImageIteratorWithIndex

See

ImageLinearConstIteratorWithIndex

See

ImageLinearIteratorWithIndex

See

ImageRandomConstIteratorWithIndex

See

ImageRandomIteratorWithIndex

See

ImageRegionConstIterator

See

ImageRegionConstIteratorWithIndex

See

ImageRegionExclusionConstIteratorWithIndex

See

ImageRegionExclusionIteratorWithIndex

See

ImageRegionIterator

See

ImageRegionIteratorWithIndex

See

ImageRegionReverseConstIterator

See

ImageRegionReverseIterator

See

ImageReverseConstIterator

See

ImageReverseIterator

See

ImageSliceConstIteratorWithIndex

See

ImageSliceIteratorWithIndex

See

NeighborhoodIterator

See

PathConstIterator

See

PathIterator

See

ShapedNeighborhoodIterator

See

SliceIterator

See

ImageConstIteratorWithIndex

See

ImageRegionRange

See

ImageRegionIndexRange

ITK Sphinx Examples:

Subclassed by itk::ImageRegionIterator< TImage >

See itk::ImageRegionConstIterator for additional documentation.