Iterate Region in Image With Access to Current Index With Write Access#
Synopsis#
Iterate over a region of an image with efficient access to the current index (with write access).
Results#
Yinyang.png#
Yinyang.png In VTK Window With Index Access#
Output:
Indices:
Index: [0, 0] value: ?
Index: [1, 0] value: ?
Index: [2, 0] value: ?
Index: [3, 0] value: ?
Index: [4, 0] value: ?
Index: [0, 1] value: ?
Index: [1, 1] value: ?
Index: [2, 1] value: ?
Index: [3, 1] value: ?
Index: [4, 1] value: ?
Index: [0, 2] value: ?
Index: [1, 2] value: ?
Index: [2, 2] value: ?
Index: [3, 2] value: ?
Index: [4, 2] value: ?
Index: [0, 3] value: ?
Index: [1, 3] value: ?
Index: [2, 3] value: ?
Index: [3, 3] value: ?
Index: [4, 3] value: ?
Code#
C++#
#include "itkImage.h"
#include "itkImageFileReader.h"
#include "itkImageRegionIteratorWithIndex.h"
#ifdef ENABLE_QUICKVIEW
# include "QuickView.h"
#endif
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::ImageRegionIteratorWithIndex<ImageType> imageIterator(image, region);
while (!imageIterator.IsAtEnd())
{
std::cout << "Index: " << imageIterator.GetIndex() << " value: " << (int)imageIterator.Get() << std::endl;
// Set the current pixel to white
imageIterator.Set(255);
++imageIterator;
}
#ifdef ENABLE_QUICKVIEW
// Visualize
QuickView viewer;
viewer.AddImage<ImageType>(image);
viewer.Visualize();
#endif
return EXIT_SUCCESS;
}
Classes demonstrated#
-
template<typename TImage>
class ImageRegionIteratorWithIndex : public itk::ImageRegionConstIteratorWithIndex<TImage> A multi-dimensional iterator templated over image type that walks pixels within a region and is specialized to keep track of its image index location.
This class is a specialization of ImageRegionConstIteratorWithIndex that adds write-access (the Set() method). Please see ImageRegionConstIteratorWithIndex for more information.
See also
ImageConstIterator
See also
ConditionalConstIterator
See also
ConstNeighborhoodIterator
See also
ConstShapedNeighborhoodIterator
See also
ConstSliceIterator
See also
CorrespondenceDataStructureIterator
See also
FloodFilledFunctionConditionalConstIterator
See also
FloodFilledImageFunctionConditionalConstIterator
See also
FloodFilledImageFunctionConditionalIterator
See also
FloodFilledSpatialFunctionConditionalConstIterator
See also
FloodFilledSpatialFunctionConditionalIterator
See also
ImageConstIterator
See also
ImageConstIteratorWithIndex
See also
ImageIterator
See also
ImageIteratorWithIndex
See also
ImageLinearConstIteratorWithIndex
See also
ImageLinearIteratorWithIndex
See also
ImageRandomConstIteratorWithIndex
See also
ImageRandomIteratorWithIndex
See also
ImageRegionConstIterator
See also
ImageRegionConstIteratorWithIndex
See also
ImageRegionExclusionConstIteratorWithIndex
See also
ImageRegionExclusionIteratorWithIndex
See also
ImageRegionIterator
See also
ImageRegionIteratorWithIndex
See also
ImageRegionReverseConstIterator
See also
ImageRegionReverseIterator
See also
ImageReverseConstIterator
See also
ImageReverseIterator
See also
ImageSliceConstIteratorWithIndex
See also
ImageSliceIteratorWithIndex
See also
NeighborhoodIterator
See also
PathConstIterator
See also
PathIterator
See also
ShapedNeighborhoodIterator
See also
SliceIterator
See also
ImageConstIteratorWithIndex
See also
ImageRegionRange
See also
ImageRegionIndexRange
- 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.
- ITK Sphinx Examples:
