ITK  5.4.0
Insight Toolkit
SphinxExamples/src/Core/Common/MakeOutOfBoundsPixelsReturnConstValue/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
*
* https://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"
using ImageType = itk::Image<unsigned char, 2>;
static void
CreateImage(ImageType::Pointer image);
int
main()
{
auto image = ImageType::New();
CreateImage(image);
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;
using BoundaryConditionType = itk::ConstantBoundaryCondition<ImageType>;
while (!iterator.IsAtEnd())
{
for (unsigned int i = 0; i < 9; ++i)
{
ImageType::IndexType index = iterator.GetIndex(i);
std::cout << "Index: " << index << " Pixel: " << i << " = " << (int)iterator.GetPixel(i) << std::endl;
}
++iterator;
}
// Visualize
auto 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;
}
void
CreateImage(ImageType::Pointer image)
{
// Create an image with 2 connected components
start[0] = 0;
start[1] = 0;
unsigned int NumRows = 5;
unsigned int NumCols = 5;
size[0] = NumRows;
size[1] = NumCols;
region.SetSize(size);
region.SetIndex(start);
image->SetRegions(region);
image->Allocate();
itk::ImageRegionIterator<ImageType> imageIterator(image, region);
// Set all pixels to white
while (!imageIterator.IsAtEnd())
{
imageIterator.Set(255);
++imageIterator;
}
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itkConstNeighborhoodIterator.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itkImageRegionIterator.h
itkImageToVTKImageFilter.h
itkConstantBoundaryCondition.h
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ConstantBoundaryCondition
This boundary condition returns a constant value for out-of-bounds image pixels.
Definition: itkConstantBoundaryCondition.h:68
itk::Index::GetIndex
const IndexValueType * GetIndex() const
Definition: itkIndex.h:232
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ImageToVTKImageFilter
Converts an ITK Image into a VTK image and plugs a ITK data pipeline to a VTK data pipeline.
Definition: itkImageToVTKImageFilter.h:47
itk::Size::SetSize
void SetSize(const SizeValueType val[VDimension])
Definition: itkSize.h:181
itk::ConstNeighborhoodIterator
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Definition: itkConstNeighborhoodIterator.h:51
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itk::ImageRegion::SetSize
void SetSize(const SizeType &size)
Definition: itkImageRegion.h:202