ITK  5.2.0
Insight Toolkit
SphinxExamples/src/Core/Common/IterateOverSpecificRegion/Code.cxx
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* http://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#include "itkImage.h"
namespace
{
using ImageType = itk::Image<int, 2>;
}
static void CreateImage(ImageType::Pointer);
int
main(int, char *[])
{
ImageType::SizeType exclusionRegionSize;
exclusionRegionSize.Fill(2);
ImageType::IndexType exclusionRegionIndex;
exclusionRegionIndex.Fill(0);
ImageType::RegionType exclusionRegion(exclusionRegionIndex, exclusionRegionSize);
ImageType::Pointer image = ImageType::New();
CreateImage(image);
itk::ImageRegionExclusionConstIteratorWithIndex<ImageType> imageIterator(image, image->GetLargestPossibleRegion());
imageIterator.SetExclusionRegion(exclusionRegion);
unsigned int numberVisited = 0;
while (!imageIterator.IsAtEnd())
{
std::cout << imageIterator.Get() << std::endl;
++imageIterator;
++numberVisited;
}
std::cout << "Visited " << numberVisited << std::endl;
return EXIT_SUCCESS;
}
void
CreateImage(ImageType::Pointer image)
{
ImageType::SizeType regionSize;
regionSize.Fill(3);
ImageType::IndexType regionIndex;
regionIndex.Fill(0);
ImageType::RegionType region(regionIndex, regionSize);
image->SetRegions(region);
image->Allocate();
image->FillBuffer(0);
}
itkImageFileReader.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::ImageRegionExclusionConstIteratorWithIndex
A multi-dimensional image iterator that walks an image region, excluding a second region that may par...
Definition: itkImageRegionExclusionConstIteratorWithIndex.h:131
itkImage.h
itk::Index::Fill
void Fill(IndexValueType value)
Definition: itkIndex.h:270
itk::Size::Fill
void Fill(SizeValueType value)
Definition: itkSize.h:211
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ImageRegionExclusionConstIteratorWithIndex::SetExclusionRegion
void SetExclusionRegion(const RegionType &region)
itkImageRegionExclusionConstIteratorWithIndex.h
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86