- 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 http://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
#if defined(_MSC_VER)
#pragma warning ( disable : 4786 )
#endif
#include <iostream>
#include "itkImage.h"
#include "itkImageSliceIteratorWithIndex.h"
#include "itkImageSliceConstIteratorWithIndex.h"
int itkImageSliceIteratorTest(int, char* [] )
{
std::cout << "Creating an image of indices" << std::endl;
const unsigned int ImageDimension = 3;
typedef itk::Index< ImageDimension > PixelType;
typedef itk::Image< PixelType, ImageDimension > ImageType;
ImageType::Pointer myImage = ImageType::New();
ImageType::SizeType size;
size[0] = 100;
size[1] = 100;
size[2] = 100;
ImageType::IndexType start;
start.Fill(0);
ImageType::RegionType region;
region.SetIndex( start );
region.SetSize( size );
myImage->SetLargestPossibleRegion( region );
myImage->SetBufferedRegion( region );
myImage->SetRequestedRegion( region );
myImage->Allocate();
typedef itk::ImageSliceIteratorWithIndex< ImageType > IteratorType;
typedef itk::ImageSliceConstIteratorWithIndex< ImageType > ConstIteratorType;
IteratorType it( myImage, region );
it.GoToBegin();
it.SetFirstDirection( 0 );
it.SetSecondDirection( 1 );
ImageType::IndexType index;
while( !it.IsAtEnd() )
{
while( !it.IsAtEndOfSlice() )
{
while( !it.IsAtEndOfLine() )
{
index = it.GetIndex();
it.Set( index );
++it;
}
it.NextLine();
}
it.NextSlice();
}
{
std::cout << "Verifying for iterator...";
IteratorType ot( myImage, region );
ot.GoToBegin();
ot.SetFirstDirection( 0 );
ot.SetSecondDirection( 1 );
while( !ot.IsAtEnd() )
{
while( !ot.IsAtEndOfSlice() )
{
while( !ot.IsAtEndOfLine() )
{
index = ot.GetIndex();
if( ot.Get() != index )
{
std::cerr << "Values don't correspond to what was stored "
<< std::endl;
std::cerr << "Test failed at index ";
std::cerr << index << std::endl;
std::cerr << "It should be ";
std::cerr << ot.Get() << std::endl;
return EXIT_FAILURE;
}
++ot;
}
ot.NextLine();
}
ot.NextSlice();
}
std::cout << " Done !" << std::endl;
std::cout << "Verifying for const iterator...";
ConstIteratorType cot( myImage, region );
cot.GoToBegin();
cot.SetFirstDirection( 0 );
cot.SetSecondDirection( 1 );
while( !cot.IsAtEnd() )
{
while( !cot.IsAtEndOfSlice() )
{
while( !cot.IsAtEndOfLine() )
{
index = cot.GetIndex();
if( cot.Get() != index )
{
std::cerr << "Values don't correspond to what was stored "
<< std::endl;
std::cerr << "Test failed at index ";
std::cerr << index << std::endl;
std::cerr << "It should be ";
std::cerr << cot.Get() << std::endl;
return EXIT_FAILURE;
}
++cot;
}
cot.NextLine();
}
cot.NextSlice();
}
std::cout << " Done !" << std::endl;
std::cout << "Verifying iterator in reverse direction... ";
IteratorType ior( myImage, region );
ior.GoToReverseBegin();
ior.SetFirstDirection( 0 );
ior.SetSecondDirection( 1 );
while( !ior.IsAtReverseEnd() )
{
while( !ior.IsAtReverseEndOfSlice() )
{
while( !ior.IsAtReverseEndOfLine() )
{
index = ior.GetIndex();
if( ior.Get() != index )
{
std::cerr << "Values don't correspond to what was stored "
<< std::endl;
std::cerr << "Test failed at index ";
std::cerr << index << " value is " << ior.Get() << std::endl;
return EXIT_FAILURE;
}
--ior;
}
ior.PreviousLine();
}
ior.PreviousSlice();
}
std::cout << " Done ! " << std::endl;
std::cout << "Verifying const iterator in reverse direction... ";
ConstIteratorType cor( myImage, region );
cor.GoToReverseBegin();
cor.SetFirstDirection( 0 );
cor.SetSecondDirection( 1 );
while( !cor.IsAtReverseEnd() )
{
while( !cor.IsAtReverseEndOfSlice() )
{
while( !cor.IsAtReverseEndOfLine() )
{
index = cor.GetIndex();
if( cor.Get() != index )
{
std::cerr << "Values don't correspond to what was stored "
<< std::endl;
std::cerr << "Test failed at index ";
std::cerr << index << " value is " << cor.Get() << std::endl;
return EXIT_FAILURE;
}
--cor;
}
cor.PreviousLine();
}
cor.PreviousSlice();
}
std::cout << " Done ! " << std::endl;
}
{
std::cout << "Test in a region < LargestPossibleRegion" << std::endl;
std::cout << "Verifying for iterator...";
size[0] = 50;
size[1] = 50;
size[2] = 50;
start[0] = 25;
start[1] = 25;
start[2] = 25;
region.SetIndex( start );
region.SetSize( size );
IteratorType ot( myImage, region );
ot.GoToBegin();
ot.SetFirstDirection( 0 );
ot.SetSecondDirection( 1 );
while( !ot.IsAtEnd() )
{
while( !ot.IsAtEndOfSlice() )
{
while( !ot.IsAtEndOfLine() )
{
index = ot.GetIndex();
if( ot.Get() != index )
{
std::cerr << "Values don't correspond to what was stored "
<< std::endl;
std::cerr << "Test failed at index ";
std::cerr << index << std::endl;
std::cerr << "It should be ";
std::cerr << ot.Get() << std::endl;
return EXIT_FAILURE;
}
++ot;
}
ot.NextLine();
}
ot.NextSlice();
}
std::cout << " Done !" << std::endl;
std::cout << "Verifying for const iterator...";
ConstIteratorType cot( myImage, region );
cot.GoToBegin();
cot.SetFirstDirection( 0 );
cot.SetSecondDirection( 1 );
while( !cot.IsAtEnd() )
{
while( !cot.IsAtEndOfSlice() )
{
while( !cot.IsAtEndOfLine() )
{
index = cot.GetIndex();
if( cot.Get() != index )
{
std::cerr << "Values don't correspond to what was stored "
<< std::endl;
std::cerr << "Test failed at index ";
std::cerr << index << std::endl;
std::cerr << "It should be ";
std::cerr << cot.Get() << std::endl;
return EXIT_FAILURE;
}
++cot;
}
cot.NextLine();
}
cot.NextSlice();
}
std::cout << " Done !" << std::endl;
std::cout << "Verifying iterator in reverse direction... ";
IteratorType ior( myImage, region );
ior.GoToReverseBegin();
ior.SetFirstDirection( 0 );
ior.SetSecondDirection( 1 );
while( !ior.IsAtReverseEnd() )
{
while( !ior.IsAtReverseEndOfSlice() )
{
while( !ior.IsAtReverseEndOfLine() )
{
index = ior.GetIndex();
if( ior.Get() != index )
{
std::cerr << "Values don't correspond to what was stored "
<< std::endl;
std::cerr << "Test failed at index ";
std::cerr << index << " value is " << ior.Get() << std::endl;
return EXIT_FAILURE;
}
--ior;
}
ior.PreviousLine();
}
ior.PreviousSlice();
}
std::cout << " Done ! " << std::endl;
std::cout << "Verifying const iterator in reverse direction... ";
ConstIteratorType cor( myImage, region );
cor.GoToReverseBegin();
cor.SetFirstDirection( 0 );
cor.SetSecondDirection( 1 );
while( !cor.IsAtReverseEnd() )
{
while( !cor.IsAtReverseEndOfSlice() )
{
while( !cor.IsAtReverseEndOfLine() )
{
index = cor.GetIndex();
if( cor.Get() != index )
{
std::cerr << "Values don't correspond to what was stored "
<< std::endl;
std::cerr << "Test failed at index ";
std::cerr << index << " value is " << cor.Get() << std::endl;
return EXIT_FAILURE;
}
--cor;
}
cor.PreviousLine();
}
cor.PreviousSlice();
}
std::cout << " Done ! " << std::endl;
}
std::cout << "Test passed" << std::endl;
return EXIT_SUCCESS;
}