ITK  6.0.0
Insight Toolkit
SphinxExamples/src/Bridge/VtkGlue/VisualizeEvolvingDense2DLevelSetAsElevationMap/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 "itkLevelSetIterationUpdateCommand.h"
#include "itkVTKVisualize2DLevelSetAsElevationMap.h"
int
main(int argc, char * argv[])
{
if (argc != 3)
{
std::cerr << "Missing Arguments" << std::endl;
std::cerr << argv[0] << std::endl;
std::cerr << "1- Input Image" << std::endl;
std::cerr << "2- Number of Iterations" << std::endl;
return EXIT_FAILURE;
}
// Image Dimension
constexpr unsigned int Dimension = 2;
using InputPixelType = unsigned char;
using InputImageType = itk::Image<InputPixelType, Dimension>;
InputImageType::Pointer input = itk::ReadImage<InputImageType>(argv[1]);
int numberOfIterations = std::stoi(argv[2]);
using LevelSetPixelType = float;
using LevelSetImageType = itk::Image<LevelSetPixelType, Dimension>;
using LevelSetOutputType = LevelSetType::OutputType;
using LevelSetRealType = LevelSetType::OutputRealType;
// Generate a binary mask that will be used as initialization for the level
// set evolution.
auto binary = BinaryImageType::New();
binary->SetRegions(input->GetLargestPossibleRegion());
binary->CopyInformation(input);
binary->Allocate();
index.Fill(5);
size.Fill(120);
region.SetIndex(index);
region.SetSize(size);
InputIteratorType iIt(binary, region);
iIt.GoToBegin();
while (!iIt.IsAtEnd())
{
++iIt;
}
// convert a binary mask to a level-set function
adaptor->SetInputImage(binary);
adaptor->Initialize();
LevelSetType::Pointer levelSet = adaptor->GetLevelSet();
// The Heaviside function
auto heaviside = HeavisideFunctionType::New();
heaviside->SetEpsilon(1.5);
// Create the level set container
auto levelSetContainer = LevelSetContainerType::New();
levelSetContainer->SetHeaviside(heaviside);
levelSetContainer->AddLevelSet(0, levelSet);
// Create the terms.
//
// // Chan and Vese internal term
using ChanAndVeseInternalTermType =
auto cvInternalTerm = ChanAndVeseInternalTermType::New();
cvInternalTerm->SetInput(input);
cvInternalTerm->SetCoefficient(0.5);
// // Chan and Vese external term
using ChanAndVeseExternalTermType =
auto cvExternalTerm = ChanAndVeseExternalTermType::New();
cvExternalTerm->SetInput(input);
// Create term container (equation rhs)
auto termContainer = TermContainerType::New();
termContainer->SetLevelSetContainer(levelSetContainer);
termContainer->SetInput(input);
termContainer->AddTerm(0, cvInternalTerm);
termContainer->AddTerm(1, cvExternalTerm);
// Create equation container
auto equationContainer = EquationContainerType::New();
equationContainer->SetLevelSetContainer(levelSetContainer);
equationContainer->AddEquation(0, termContainer);
// Create stopping criteria
using StoppingCriterionType = itk::LevelSetEvolutionNumberOfIterationsStoppingCriterion<LevelSetContainerType>;
auto criterion = StoppingCriterionType::New();
criterion->SetNumberOfIterations(numberOfIterations);
// Create the visualizer
using VisualizationType = itk::VTKVisualize2DLevelSetAsElevationMap<InputImageType, LevelSetType>;
auto visualizer = VisualizationType::New();
visualizer->SetInputImage(input);
visualizer->SetLevelSet(levelSet);
visualizer->SetScreenCapture(true);
// Create evolution class
auto evolution = LevelSetEvolutionType::New();
evolution->SetEquationContainer(equationContainer);
evolution->SetStoppingCriterion(criterion);
evolution->SetLevelSetContainer(levelSetContainer);
using IterationUpdateCommandType = itk::LevelSetIterationUpdateCommand<LevelSetEvolutionType, VisualizationType>;
auto iterationUpdateCommand = IterationUpdateCommandType::New();
iterationUpdateCommand->SetFilterToUpdate(visualizer);
iterationUpdateCommand->SetUpdatePeriod(5);
evolution->AddObserver(itk::IterationEvent(), iterationUpdateCommand);
evolution->Update();
return EXIT_SUCCESS;
}
itk::BinaryImageToLevelSetImageAdaptor
Converts one binary image to the appropriate level-set type provided by the template argument TLevelS...
Definition: itkBinaryImageToLevelSetImageAdaptor.h:51
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itkBinaryImageToLevelSetImageAdaptor.h
itkLevelSetEquationChanAndVeseExternalTerm.h
itk::LevelSetDenseImage
Base class for the "dense" representation of a level-set function on one image.
Definition: itkLevelSetDenseImage.h:41
itkImageFileReader.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkLevelSetEquationContainer.h
itkLevelSetEquationTermContainer.h
itk::Size::Fill
void Fill(SizeValueType value)
Definition: itkSize.h:211
itkSinRegularizedHeavisideStepFunction.h
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::SinRegularizedHeavisideStepFunction
Sin-based implementation of the Regularized (smoothed) Heaviside functions.
Definition: itkSinRegularizedHeavisideStepFunction.h:50
itk::LevelSetEquationChanAndVeseExternalTerm
Class to represent the external energy Chan And Vese term.
Definition: itkLevelSetEquationChanAndVeseExternalTerm.h:50
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkLevelSetEquationChanAndVeseInternalTerm.h
itk::ImageRegionIteratorWithIndex
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
Definition: itkImageRegionIteratorWithIndex.h:73
itk::LevelSetContainer
Container of Level-Sets.
Definition: itkLevelSetContainer.h:39
itkLevelSetEvolutionNumberOfIterationsStoppingCriterion.h
itkLevelSetDenseImage.h
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:60
itkLevelSetEvolution.h
itk::LevelSetEquationChanAndVeseInternalTerm
Class to represent the internal energy Chan And Vese term.
Definition: itkLevelSetEquationChanAndVeseInternalTerm.h:48
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::LevelSetEvolution
Class for iterating and evolving the level-set function.
Definition: itkLevelSetEvolution.h:48
New
static Pointer New()
itk::LevelSetEquationContainer
Class for holding a set of level set equations (PDEs).
Definition: itkLevelSetEquationContainer.h:58
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itk::LevelSetEquationTermContainer
Class for container holding the terms of a given level set update equation.
Definition: itkLevelSetEquationTermContainer.h:42
itkLevelSetContainer.h