ITK  5.2.0
Insight Toolkit
SphinxExamples/src/Core/Common/IterateRegionWithNeighborhood/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"
#include "vtkVersion.h"
#include "vtkImageViewer.h"
#include "vtkImageMapper3D.h"
#include "vtkRenderWindowInteractor.h"
#include "vtkSmartPointer.h"
#include "vtkImageActor.h"
#include "vtkInteractorStyleImage.h"
#include "vtkRenderer.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>;
using ReaderType = itk::ImageFileReader<ImageType>;
ReaderType::Pointer reader = ReaderType::New();
reader->SetFileName(argv[1]);
reader->Update();
ImageType::Pointer image = reader->GetOutput();
ImageType::SizeType regionSize;
regionSize[0] = 50;
regionSize[1] = 1;
ImageType::IndexType regionIndex;
regionIndex[0] = 0;
regionIndex[1] = 0;
region.SetSize(regionSize);
region.SetIndex(regionIndex);
radius[0] = 1;
radius[1] = 1;
itk::NeighborhoodIterator<ImageType> iterator(radius, image, region);
while (!iterator.IsAtEnd())
{
// Set the current pixel to white
iterator.SetCenterPixel(255);
for (unsigned int i = 0; i < 9; i++)
{
ImageType::IndexType index = iterator.GetIndex(i);
std::cout << index[0] << " " << index[1] << std::endl;
bool IsInBounds;
iterator.GetPixel(i, IsInBounds);
if (IsInBounds)
{
iterator.SetPixel(i, 255);
}
}
++iterator;
}
// Visualize
ConnectorType::Pointer connector = ConnectorType::New();
connector->SetInput(image);
vtkSmartPointer<vtkImageActor> actor = vtkSmartPointer<vtkImageActor>::New();
#if VTK_MAJOR_VERSION <= 5
actor->SetInput(connector->GetOutput());
#else
connector->Update();
actor->GetMapper()->SetInputData(connector->GetOutput());
#endif
vtkSmartPointer<vtkRenderWindow> renderWindow = vtkSmartPointer<vtkRenderWindow>::New();
vtkSmartPointer<vtkRenderWindowInteractor> interactor = vtkSmartPointer<vtkRenderWindowInteractor>::New();
interactor->SetRenderWindow(renderWindow);
vtkSmartPointer<vtkRenderer> renderer = vtkSmartPointer<vtkRenderer>::New();
renderWindow->AddRenderer(renderer);
renderer->AddActor(actor);
renderer->ResetCamera();
renderWindow->Render();
vtkSmartPointer<vtkInteractorStyleImage> style = vtkSmartPointer<vtkInteractorStyleImage>::New();
interactor->SetInteractorStyle(style);
interactor->Start();
return EXIT_SUCCESS;
}
itk::ConstNeighborhoodIterator::IsAtEnd
ITK_ITERATOR_VIRTUAL bool IsAtEnd() const ITK_ITERATOR_FINAL
Definition: itkConstNeighborhoodIterator.h:354
itkNeighborhoodIterator.h
itk::ConstNeighborhoodIterator::GetPixel
ITK_ITERATOR_VIRTUAL PixelType GetPixel(const NeighborIndexType i) const ITK_ITERATOR_FINAL
Definition: itkConstNeighborhoodIterator.h:195
itkImageFileReader.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itkImageToVTKImageFilter.h
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
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::ConstNeighborhoodIterator::GetIndex
ITK_ITERATOR_VIRTUAL IndexType GetIndex() const ITK_ITERATOR_FINAL
Definition: itkConstNeighborhoodIterator.h:177
itk::ImageToVTKImageFilter
Converts an ITK image into a VTK image and plugs a itk data pipeline to a VTK datapipeline.
Definition: itkImageToVTKImageFilter.h:47
itk::NeighborhoodIterator::SetPixel
ITK_ITERATOR_VIRTUAL void SetPixel(const unsigned i, const PixelType &v, bool &status) ITK_ITERATOR_FINAL
itk::NeighborhoodIterator
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Definition: itkNeighborhoodIterator.h:212
itk::NeighborhoodIterator::SetCenterPixel
ITK_ITERATOR_VIRTUAL void SetCenterPixel(const PixelType &p) ITK_ITERATOR_FINAL
Definition: itkNeighborhoodIterator.h:271
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::ImageRegion::SetSize
void SetSize(const SizeType &size)
Definition: itkImageRegion.h:174