Visualize a Static Sparse Shi 2D Level-Set Layers#
Synopsis#
Visualize a static sparse Shi level-set function 2D’s layers. From the input image, first an otsu thresholding technique is used to get a binary mask, which is then converted to a sparse level-set function.
Note that Shi’s representation is composed of 4 layers (values = {-3, -1, +1, +3})
Results#
Code#
C++#
#include "itkVTKVisualize2DSparseLevelSetLayers.h"
#include "itkBinaryImageToLevelSetImageAdaptor.h"
#include "itkImageFileReader.h"
#include "itkShiSparseLevelSetImage.h"
#include "itkOtsuMultipleThresholdsImageFilter.h"
#include "itkRescaleIntensityImageFilter.h"
#include "vtkRenderWindowInteractor.h"
int
main(int argc, char * argv[])
{
if (argc != 3)
{
std::cerr << "Missing Arguments" << std::endl;
std::cerr << argv[0] << std::endl;
std::cerr << "<Input Image> <Interactive (0 or 1)>" << 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]);
using LevelSetType = itk::ShiSparseLevelSetImage<Dimension>;
// Generate a binary mask that will be used as initialization for the level
// set evolution.
using OtsuFilterType = itk::OtsuMultipleThresholdsImageFilter<InputImageType, InputImageType>;
auto otsu = OtsuFilterType::New();
otsu->SetInput(input);
otsu->SetNumberOfHistogramBins(256);
otsu->SetNumberOfThresholds(1);
using RescaleType = itk::RescaleIntensityImageFilter<InputImageType, InputImageType>;
auto rescaler = RescaleType::New();
rescaler->SetInput(otsu->GetOutput());
rescaler->SetOutputMinimum(0);
rescaler->SetOutputMaximum(1);
// convert a binary mask to a level-set function
using BinaryImageToLevelSetType = itk::BinaryImageToLevelSetImageAdaptor<InputImageType, LevelSetType>;
auto adaptor = BinaryImageToLevelSetType::New();
adaptor->SetInputImage(rescaler->GetOutput());
adaptor->Initialize();
LevelSetType::Pointer levelSet = adaptor->GetModifiableLevelSet();
// Create the visualizer
using VisualizationType = itk::VTKVisualize2DSparseLevelSetLayers<InputImageType, LevelSetType>;
auto visualizer = VisualizationType::New();
visualizer->SetInputImage(input);
visualizer->SetLevelSet(levelSet);
visualizer->SetScreenCapture(true);
vtkSmartPointer<vtkRenderWindowInteractor> renderWindowInteractor = vtkSmartPointer<vtkRenderWindowInteractor>::New();
renderWindowInteractor->SetRenderWindow(visualizer->GetRenderWindow());
try
{
visualizer->Update();
}
catch (const itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
bool interactive = (std::stoi(argv[2]) != 0);
if (interactive)
{
renderWindowInteractor->Start();
}
return EXIT_SUCCESS;
}