Visualize an Evolving Dense 2D Level Set as Elevation Map#

Synopsis#

Visualize an evolving dense level-set function 2D rendered as an elevation map.

Results#

Input image (cells)

Input image#

Evolving level-sets rendered as an elevation map

Evolving level-sets#

Evolving level-sets animation

Code#

C++#

#include "itkBinaryImageToLevelSetImageAdaptor.h"
#include "itkImageFileReader.h"
#include "itkLevelSetIterationUpdateCommand.h"
#include "itkLevelSetContainer.h"
#include "itkLevelSetEquationChanAndVeseInternalTerm.h"
#include "itkLevelSetEquationChanAndVeseExternalTerm.h"
#include "itkLevelSetEquationContainer.h"
#include "itkLevelSetEquationTermContainer.h"
#include "itkLevelSetEvolution.h"
#include "itkLevelSetEvolutionNumberOfIterationsStoppingCriterion.h"
#include "itkLevelSetDenseImage.h"
#include "itkVTKVisualize2DLevelSetAsElevationMap.h"
#include "itkSinRegularizedHeavisideStepFunction.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 LevelSetType = itk::LevelSetDenseImage<LevelSetImageType>;

  using LevelSetOutputType = LevelSetType::OutputType;
  using LevelSetRealType = LevelSetType::OutputRealType;

  // Generate a binary mask that will be used as initialization for the level
  // set evolution.
  using BinaryImageType = itk::Image<LevelSetOutputType, Dimension>;
  auto binary = BinaryImageType::New();
  binary->SetRegions(input->GetLargestPossibleRegion());
  binary->CopyInformation(input);
  binary->Allocate();
  binary->FillBuffer(itk::NumericTraits<LevelSetOutputType>::Zero);

  BinaryImageType::RegionType region;
  BinaryImageType::IndexType  index;
  BinaryImageType::SizeType   size;

  index.Fill(5);
  size.Fill(120);

  region.SetIndex(index);
  region.SetSize(size);

  using InputIteratorType = itk::ImageRegionIteratorWithIndex<BinaryImageType>;
  InputIteratorType iIt(binary, region);
  iIt.GoToBegin();
  while (!iIt.IsAtEnd())
  {
    iIt.Set(itk::NumericTraits<LevelSetOutputType>::One);
    ++iIt;
  }

  // convert a binary mask to a level-set function
  using BinaryImageToLevelSetType = itk::BinaryImageToLevelSetImageAdaptor<BinaryImageType, LevelSetType>;

  auto adaptor = BinaryImageToLevelSetType::New();
  adaptor->SetInputImage(binary);
  adaptor->Initialize();
  LevelSetType::Pointer levelSet = adaptor->GetLevelSet();

  // The Heaviside function
  using HeavisideFunctionType = itk::SinRegularizedHeavisideStepFunction<LevelSetRealType, LevelSetRealType>;
  auto heaviside = HeavisideFunctionType::New();
  heaviside->SetEpsilon(1.5);

  // Create the level set container
  using LevelSetContainerType = itk::LevelSetContainer<itk::IdentifierType, LevelSetType>;
  auto levelSetContainer = LevelSetContainerType::New();
  levelSetContainer->SetHeaviside(heaviside);
  levelSetContainer->AddLevelSet(0, levelSet);

  // Create the terms.
  //
  // // Chan and Vese internal term
  using ChanAndVeseInternalTermType =
    itk::LevelSetEquationChanAndVeseInternalTerm<InputImageType, LevelSetContainerType>;
  auto cvInternalTerm = ChanAndVeseInternalTermType::New();
  cvInternalTerm->SetInput(input);
  cvInternalTerm->SetCoefficient(0.5);

  // // Chan and Vese external term
  using ChanAndVeseExternalTermType =
    itk::LevelSetEquationChanAndVeseExternalTerm<InputImageType, LevelSetContainerType>;
  auto cvExternalTerm = ChanAndVeseExternalTermType::New();
  cvExternalTerm->SetInput(input);

  // Create term container (equation rhs)
  using TermContainerType = itk::LevelSetEquationTermContainer<InputImageType, LevelSetContainerType>;
  auto termContainer = TermContainerType::New();
  termContainer->SetLevelSetContainer(levelSetContainer);
  termContainer->SetInput(input);
  termContainer->AddTerm(0, cvInternalTerm);
  termContainer->AddTerm(1, cvExternalTerm);

  // Create equation container
  using EquationContainerType = itk::LevelSetEquationContainer<TermContainerType>;
  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
  using LevelSetEvolutionType = itk::LevelSetEvolution<EquationContainerType, LevelSetType>;
  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;
}

Classes demonstrated#