Random Selection of Pixels From Region#

Synopsis#

Randomly select pixels from a region of an image.

Results#

Output:

[4, 2]
[1, 1]
[1, 3]
[4, 0]
[3, 0]
[4, 3]
[1, 0]
[1, 1]
[4, 2]
[4, 2]
[1, 2]
[2, 1]
[1, 2]
[2, 3]
[2, 3]
[4, 1]
[4, 2]
[3, 1]
[2, 2]
[2, 2]
[1, 1]
[2, 2]
[1, 1]
[1, 1]
[1, 0]
[3, 1]
[2, 1]
[4, 0]
[2, 2]
[4, 1]
[0, 2]
[2, 1]
[0, 2]
[1, 0]
[1, 0]
[2, 3]
[1, 1]
[1, 3]
[0, 1]
[3, 0]
[4, 2]
[4, 0]
[3, 1]
[0, 1]
[0, 3]
[0, 0]
[1, 2]
[3, 2]
[2, 0]
[2, 2]
[1, 1]
[3, 2]
[1, 1]
[2, 1]
[2, 1]
[4, 2]
[1, 3]
[4, 1]
[0, 3]
[1, 2]
[4, 3]
[0, 0]
[2, 1]
[1, 0]
[4, 3]
[4, 3]
[2, 3]
[2, 3]
[0, 3]
[1, 2]
[1, 3]
[2, 0]
[1, 0]
[2, 3]
[3, 3]
[3, 2]
[0, 0]
[3, 3]
[1, 2]
[1, 0]
[3, 0]
[3, 3]
[4, 0]
[1, 2]
[0, 2]
[1, 3]
[3, 3]
[1, 1]
[0, 3]
[0, 0]
[0, 2]
[0, 2]
[2, 2]
[2, 2]
[4, 1]
[2, 2]
[1, 1]
[1, 0]
[4, 1]
[1, 0]
[3, 1]
[4, 3]
[0, 0]
[3, 2]
[1, 2]
[1, 2]
[2, 3]
[4, 2]
[2, 2]
[3, 0]
[0, 3]
[2, 2]
[1, 1]
[2, 3]
[4, 0]
[0, 3]
[4, 2]
[3, 2]
[3, 3]
[1, 2]
[4, 3]
[0, 0]
[1, 3]
[0, 0]
[4, 2]
[0, 2]
[0, 2]
[0, 3]
[4, 2]
[4, 2]
[4, 2]
[3, 2]
[4, 0]
[0, 0]
[0, 1]
[3, 0]
[0, 0]
[0, 0]
[0, 0]
[0, 0]
[4, 2]
[1, 2]
[0, 3]
[3, 3]
[4, 2]
[2, 2]
[2, 2]
[3, 3]
[1, 2]
[2, 1]
[0, 2]
[1, 2]
[4, 0]
[4, 3]
[3, 3]
[1, 2]
[0, 1]
[3, 2]
[0, 2]
[4, 1]
[1, 2]
[3, 0]
[4, 3]
[1, 2]
[2, 2]
[2, 3]
[0, 0]
[2, 3]
[4, 1]
[4, 2]
[4, 2]
[1, 2]
[3, 3]
[1, 3]
[0, 2]
[3, 2]
[3, 0]
[2, 3]
[4, 2]
[2, 3]
[4, 2]
[4, 0]
[2, 0]
[1, 0]
[1, 2]
[1, 0]
[3, 2]
[0, 0]
[3, 1]
[4, 3]
[0, 3]
[0, 2]
[2, 0]
[4, 2]
[3, 2]
[1, 0]
[1, 0]
[4, 0]
[1, 2]
[0, 2]

Code#

C++#

#include "itkImage.h"
#include "itkImageRandomConstIteratorWithIndex.h"

int
main()
{
  using ImageType = itk::Image<unsigned char, 2>;
  auto image = ImageType::New();

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

  image->SetRegions(region);
  image->Allocate();
  image->FillBuffer(0);

  itk::ImageRandomConstIteratorWithIndex<ImageType> imageIterator(image, image->GetLargestPossibleRegion());
  imageIterator.SetNumberOfSamples(200);
  imageIterator.GoToBegin();

  while (!imageIterator.IsAtEnd())
  {
    std::cout << imageIterator.GetIndex() << std::endl;

    ++imageIterator;
  }

  return EXIT_SUCCESS;
}

Classes demonstrated#

template<typename TImage>
class ImageRandomConstIteratorWithIndex : public itk::ImageConstIteratorWithIndex<TImage>

A multi-dimensional image iterator that visits a random set of pixels within an image region.

ImageRandomConstIteratorWithIndex is a multi-dimensional iterator class that is templated over image type. ImageRandomConstIteratorWithIndex is constrained to walk only within the specified region. It samples random pixel positions at each increment or decrement.

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

The operator++ method provides a simple syntax for walking around a region of a multidimensional image. operator++ performs a jump to a random position within the specified image region. This is designed to facilitate the extraction of random samples from the image.

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

ImageRandomConstIteratorWithIndex<ImageType> it( image, image->GetRequestedRegion() );

it.SetNumberOfSamples(200);
it.GoToBegin();
while( !it.IsAtEnd() )
{
  it.Get();
  ++it;  // here it jumps to another random position inside the region
 }

or

ImageRandomConstIteratorWithIndex<ImageType> it( image, image->GetRequestedRegion() );

it.SetNumberOfSamples(200);
it.GoToEnd();
while( !it.IsAtBegin() )
{
  it.Get();
  --it;  // here it jumps to another random position inside the region
 }

Warning

Incrementing the iterator (++it) followed by a decrement (it) or vice versa does not in general return the iterator to the same position.

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

ITK Sphinx Examples:

Subclassed by itk::ImageRandomIteratorWithIndex< TImage >

See itk::ImageRandomConstIteratorWithIndex for additional documentation.