Segment With Geodesic Active Contour Level Set¶
Synopsis¶
Segments structures in images based on a user supplied edge potential map.
This example is the same as the one in the ITK software guide. It was re-organized and does not contain any comments about the executed filters. Please refer to the guide for more details.
Results¶
![Input image](../../../../_images/BrainProtonDensitySlice8.png)
Input image¶
![Evolution of the segmentation depending on the number of iterations (10 to 500).](../../../../_images/Evolution.png)
Used parameters:
Left ventricle: 81 114 5.0 1.0 -0.5 3.0 2.0
Right ventricle: 99 114 5.0 1.0 -0.5 3.0 2.0
White matter: 56 92 5.0 1.0 -0.3 2.0 10.0
Gray matter: 40 90 5.0 .5 -0.3 2.0 10.0
Code¶
Python¶
#!/usr/bin/env python
import itk
import argparse
parser = argparse.ArgumentParser(
description="Segment With Geodesic Active Contour Level Set."
)
parser.add_argument("input_image")
parser.add_argument("output_image")
parser.add_argument("seed_x", type=int)
parser.add_argument("seed_y", type=int)
parser.add_argument("initial_distance", type=float)
parser.add_argument("sigma", type=float)
parser.add_argument("sigmoid_alpha", type=float)
parser.add_argument("sigmoid_beta", type=float)
parser.add_argument("propagation_scaling", type=float)
parser.add_argument("number_of_iterations", type=int)
args = parser.parse_args()
seedValue = -args.initial_distance
Dimension = 2
InputPixelType = itk.F
OutputPixelType = itk.UC
InputImageType = itk.Image[InputPixelType, Dimension]
OutputImageType = itk.Image[OutputPixelType, Dimension]
ReaderType = itk.ImageFileReader[InputImageType]
WriterType = itk.ImageFileWriter[OutputImageType]
reader = ReaderType.New()
reader.SetFileName(args.input_image)
SmoothingFilterType = itk.CurvatureAnisotropicDiffusionImageFilter[
InputImageType, InputImageType
]
smoothing = SmoothingFilterType.New()
smoothing.SetTimeStep(0.125)
smoothing.SetNumberOfIterations(5)
smoothing.SetConductanceParameter(9.0)
smoothing.SetInput(reader.GetOutput())
GradientFilterType = itk.GradientMagnitudeRecursiveGaussianImageFilter[
InputImageType, InputImageType
]
gradientMagnitude = GradientFilterType.New()
gradientMagnitude.SetSigma(args.sigma)
gradientMagnitude.SetInput(smoothing.GetOutput())
SigmoidFilterType = itk.SigmoidImageFilter[InputImageType, InputImageType]
sigmoid = SigmoidFilterType.New()
sigmoid.SetOutputMinimum(0.0)
sigmoid.SetOutputMaximum(1.0)
sigmoid.SetAlpha(args.sigmoid_alpha)
sigmoid.SetBeta(args.sigmoid_beta)
sigmoid.SetInput(gradientMagnitude.GetOutput())
FastMarchingFilterType = itk.FastMarchingImageFilter[InputImageType, InputImageType]
fastMarching = FastMarchingFilterType.New()
GeoActiveContourFilterType = itk.GeodesicActiveContourLevelSetImageFilter[
InputImageType, InputImageType, InputPixelType
]
geodesicActiveContour = GeoActiveContourFilterType.New()
geodesicActiveContour.SetPropagationScaling(args.propagation_scaling)
geodesicActiveContour.SetCurvatureScaling(1.0)
geodesicActiveContour.SetAdvectionScaling(1.0)
geodesicActiveContour.SetMaximumRMSError(0.02)
geodesicActiveContour.SetNumberOfIterations(args.number_of_iterations)
geodesicActiveContour.SetInput(fastMarching.GetOutput())
geodesicActiveContour.SetFeatureImage(sigmoid.GetOutput())
ThresholdingFilterType = itk.BinaryThresholdImageFilter[InputImageType, OutputImageType]
thresholder = ThresholdingFilterType.New()
thresholder.SetLowerThreshold(-1000.0)
thresholder.SetUpperThreshold(0.0)
thresholder.SetOutsideValue(itk.NumericTraits[OutputPixelType].min())
thresholder.SetInsideValue(itk.NumericTraits[OutputPixelType].max())
thresholder.SetInput(geodesicActiveContour.GetOutput())
seedPosition = itk.Index[Dimension]()
seedPosition[0] = args.seed_x
seedPosition[1] = args.seed_y
node = itk.LevelSetNode[InputPixelType, Dimension]()
node.SetValue(seedValue)
node.SetIndex(seedPosition)
seeds = itk.VectorContainer[itk.UI, itk.LevelSetNode[InputPixelType, Dimension]].New()
seeds.Initialize()
seeds.InsertElement(0, node)
fastMarching.SetTrialPoints(seeds)
fastMarching.SetSpeedConstant(1.0)
CastFilterType = itk.RescaleIntensityImageFilter[InputImageType, OutputImageType]
caster1 = CastFilterType.New()
caster2 = CastFilterType.New()
caster3 = CastFilterType.New()
caster4 = CastFilterType.New()
writer1 = WriterType.New()
writer2 = WriterType.New()
writer3 = WriterType.New()
writer4 = WriterType.New()
caster1.SetInput(smoothing.GetOutput())
writer1.SetInput(caster1.GetOutput())
writer1.SetFileName("GeodesicActiveContourImageFilterOutput1.png")
caster1.SetOutputMinimum(itk.NumericTraits[OutputPixelType].min())
caster1.SetOutputMaximum(itk.NumericTraits[OutputPixelType].max())
writer1.Update()
caster2.SetInput(gradientMagnitude.GetOutput())
writer2.SetInput(caster2.GetOutput())
writer2.SetFileName("GeodesicActiveContourImageFilterOutput2.png")
caster2.SetOutputMinimum(itk.NumericTraits[OutputPixelType].min())
caster2.SetOutputMaximum(itk.NumericTraits[OutputPixelType].max())
writer2.Update()
caster3.SetInput(sigmoid.GetOutput())
writer3.SetInput(caster3.GetOutput())
writer3.SetFileName("GeodesicActiveContourImageFilterOutput3.png")
caster3.SetOutputMinimum(itk.NumericTraits[OutputPixelType].min())
caster3.SetOutputMaximum(itk.NumericTraits[OutputPixelType].max())
writer3.Update()
caster4.SetInput(fastMarching.GetOutput())
writer4.SetInput(caster4.GetOutput())
writer4.SetFileName("GeodesicActiveContourImageFilterOutput4.png")
caster4.SetOutputMinimum(itk.NumericTraits[OutputPixelType].min())
caster4.SetOutputMaximum(itk.NumericTraits[OutputPixelType].max())
fastMarching.SetOutputSize(reader.GetOutput().GetBufferedRegion().GetSize())
writer = WriterType.New()
writer.SetFileName(args.output_image)
writer.SetInput(thresholder.GetOutput())
writer.Update()
print(
"Max. no. iterations: " + str(geodesicActiveContour.GetNumberOfIterations()) + "\n"
)
print("Max. RMS error: " + str(geodesicActiveContour.GetMaximumRMSError()) + "\n")
print(
"No. elpased iterations: "
+ str(geodesicActiveContour.GetElapsedIterations())
+ "\n"
)
print("RMS change: " + str(geodesicActiveContour.GetRMSChange()) + "\n")
writer4.Update()
InternalWriterType = itk.ImageFileWriter[InputImageType]
mapWriter = InternalWriterType.New()
mapWriter.SetInput(fastMarching.GetOutput())
mapWriter.SetFileName("GeodesicActiveContourImageFilterOutput4.mha")
mapWriter.Update()
speedWriter = InternalWriterType.New()
speedWriter.SetInput(sigmoid.GetOutput())
speedWriter.SetFileName("GeodesicActiveContourImageFilterOutput3.mha")
speedWriter.Update()
gradientWriter = InternalWriterType.New()
gradientWriter.SetInput(gradientMagnitude.GetOutput())
gradientWriter.SetFileName("GeodesicActiveContourImageFilterOutput2.mha")
gradientWriter.Update()
C++¶
#include "itkImageFileReader.h"
#include "itkImageFileWriter.h"
#include "itkGeodesicActiveContourLevelSetImageFilter.h"
#include "itkCurvatureAnisotropicDiffusionImageFilter.h"
#include "itkGradientMagnitudeRecursiveGaussianImageFilter.h"
#include "itkSigmoidImageFilter.h"
#include "itkFastMarchingImageFilter.h"
#include "itkRescaleIntensityImageFilter.h"
#include "itkBinaryThresholdImageFilter.h"
int
main(int argc, char * argv[])
{
if (argc != 11)
{
std::cerr << "Usage: " << argv[0];
std::cerr << " <InputFileName> <OutputFileName>";
std::cerr << " <seedX> <seedY> <InitialDistance>";
std::cerr << " <Sigma> <SigmoidAlpha> <SigmoidBeta>";
std::cerr << " <PropagationScaling> <NumberOfIterations>" << std::endl;
return EXIT_FAILURE;
}
const char * inputFileName = argv[1];
const char * outputFileName = argv[2];
const int seedPosX = std::stoi(argv[3]);
const int seedPosY = std::stoi(argv[4]);
const double initialDistance = std::stod(argv[5]);
const double sigma = std::stod(argv[6]);
const double alpha = std::stod(argv[7]);
const double beta = std::stod(argv[8]);
const double propagationScaling = std::stod(argv[9]);
const double numberOfIterations = std::stoi(argv[10]);
const double seedValue = -initialDistance;
constexpr unsigned int Dimension = 2;
using InputPixelType = float;
using InputImageType = itk::Image<InputPixelType, Dimension>;
using OutputPixelType = unsigned char;
using OutputImageType = itk::Image<OutputPixelType, Dimension>;
using ReaderType = itk::ImageFileReader<InputImageType>;
using WriterType = itk::ImageFileWriter<OutputImageType>;
ReaderType::Pointer reader = ReaderType::New();
reader->SetFileName(inputFileName);
using SmoothingFilterType = itk::CurvatureAnisotropicDiffusionImageFilter<InputImageType, InputImageType>;
SmoothingFilterType::Pointer smoothing = SmoothingFilterType::New();
smoothing->SetTimeStep(0.125);
smoothing->SetNumberOfIterations(5);
smoothing->SetConductanceParameter(9.0);
smoothing->SetInput(reader->GetOutput());
using GradientFilterType = itk::GradientMagnitudeRecursiveGaussianImageFilter<InputImageType, InputImageType>;
GradientFilterType::Pointer gradientMagnitude = GradientFilterType::New();
gradientMagnitude->SetSigma(sigma);
gradientMagnitude->SetInput(smoothing->GetOutput());
using SigmoidFilterType = itk::SigmoidImageFilter<InputImageType, InputImageType>;
SigmoidFilterType::Pointer sigmoid = SigmoidFilterType::New();
sigmoid->SetOutputMinimum(0.0);
sigmoid->SetOutputMaximum(1.0);
sigmoid->SetAlpha(alpha);
sigmoid->SetBeta(beta);
sigmoid->SetInput(gradientMagnitude->GetOutput());
using FastMarchingFilterType = itk::FastMarchingImageFilter<InputImageType, InputImageType>;
FastMarchingFilterType::Pointer fastMarching = FastMarchingFilterType::New();
using GeodesicActiveContourFilterType = itk::GeodesicActiveContourLevelSetImageFilter<InputImageType, InputImageType>;
GeodesicActiveContourFilterType::Pointer geodesicActiveContour = GeodesicActiveContourFilterType::New();
geodesicActiveContour->SetPropagationScaling(propagationScaling);
geodesicActiveContour->SetCurvatureScaling(1.0);
geodesicActiveContour->SetAdvectionScaling(1.0);
geodesicActiveContour->SetMaximumRMSError(0.02);
geodesicActiveContour->SetNumberOfIterations(numberOfIterations);
geodesicActiveContour->SetInput(fastMarching->GetOutput());
geodesicActiveContour->SetFeatureImage(sigmoid->GetOutput());
using ThresholdingFilterType = itk::BinaryThresholdImageFilter<InputImageType, OutputImageType>;
ThresholdingFilterType::Pointer thresholder = ThresholdingFilterType::New();
thresholder->SetLowerThreshold(-1000.0);
thresholder->SetUpperThreshold(0.0);
thresholder->SetOutsideValue(itk::NumericTraits<OutputPixelType>::min());
thresholder->SetInsideValue(itk::NumericTraits<OutputPixelType>::max());
thresholder->SetInput(geodesicActiveContour->GetOutput());
using NodeContainer = FastMarchingFilterType::NodeContainer;
using NodeType = FastMarchingFilterType::NodeType;
InputImageType::IndexType seedPosition;
seedPosition[0] = seedPosX;
seedPosition[1] = seedPosY;
NodeContainer::Pointer seeds = NodeContainer::New();
NodeType node;
node.SetValue(seedValue);
node.SetIndex(seedPosition);
seeds->Initialize();
seeds->InsertElement(0, node);
fastMarching->SetTrialPoints(seeds);
fastMarching->SetSpeedConstant(1.0);
using CastFilterType = itk::RescaleIntensityImageFilter<InputImageType, OutputImageType>;
CastFilterType::Pointer caster1 = CastFilterType::New();
CastFilterType::Pointer caster2 = CastFilterType::New();
CastFilterType::Pointer caster3 = CastFilterType::New();
CastFilterType::Pointer caster4 = CastFilterType::New();
WriterType::Pointer writer1 = WriterType::New();
WriterType::Pointer writer2 = WriterType::New();
WriterType::Pointer writer3 = WriterType::New();
WriterType::Pointer writer4 = WriterType::New();
caster1->SetInput(smoothing->GetOutput());
writer1->SetInput(caster1->GetOutput());
writer1->SetFileName("GeodesicActiveContourImageFilterOutput1.png");
caster1->SetOutputMinimum(itk::NumericTraits<OutputPixelType>::min());
caster1->SetOutputMaximum(itk::NumericTraits<OutputPixelType>::max());
writer1->Update();
caster2->SetInput(gradientMagnitude->GetOutput());
writer2->SetInput(caster2->GetOutput());
writer2->SetFileName("GeodesicActiveContourImageFilterOutput2.png");
caster2->SetOutputMinimum(itk::NumericTraits<OutputPixelType>::min());
caster2->SetOutputMaximum(itk::NumericTraits<OutputPixelType>::max());
writer2->Update();
caster3->SetInput(sigmoid->GetOutput());
writer3->SetInput(caster3->GetOutput());
writer3->SetFileName("GeodesicActiveContourImageFilterOutput3.png");
caster3->SetOutputMinimum(itk::NumericTraits<OutputPixelType>::min());
caster3->SetOutputMaximum(itk::NumericTraits<OutputPixelType>::max());
writer3->Update();
caster4->SetInput(fastMarching->GetOutput());
writer4->SetInput(caster4->GetOutput());
writer4->SetFileName("GeodesicActiveContourImageFilterOutput4.png");
caster4->SetOutputMinimum(itk::NumericTraits<OutputPixelType>::min());
caster4->SetOutputMaximum(itk::NumericTraits<OutputPixelType>::max());
InputImageType::Pointer inImage = reader->GetOutput();
fastMarching->SetOutputDirection(inImage->GetDirection());
fastMarching->SetOutputOrigin(inImage->GetOrigin());
fastMarching->SetOutputRegion(inImage->GetBufferedRegion());
fastMarching->SetOutputSpacing(inImage->GetSpacing());
WriterType::Pointer writer = WriterType::New();
writer->SetFileName(outputFileName);
writer->SetInput(thresholder->GetOutput());
try
{
writer->Update();
}
catch (itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
std::cout << std::endl;
std::cout << "Max. no. iterations: " << geodesicActiveContour->GetNumberOfIterations() << std::endl;
std::cout << "Max. RMS error: " << geodesicActiveContour->GetMaximumRMSError() << std::endl;
std::cout << std::endl;
std::cout << "No. elpased iterations: " << geodesicActiveContour->GetElapsedIterations() << std::endl;
std::cout << "RMS change: " << geodesicActiveContour->GetRMSChange() << std::endl;
try
{
writer4->Update();
}
catch (itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
using InternalWriterType = itk::ImageFileWriter<InputImageType>;
InternalWriterType::Pointer mapWriter = InternalWriterType::New();
mapWriter->SetInput(fastMarching->GetOutput());
mapWriter->SetFileName("GeodesicActiveContourImageFilterOutput4.mha");
try
{
mapWriter->Update();
}
catch (itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
InternalWriterType::Pointer speedWriter = InternalWriterType::New();
speedWriter->SetInput(sigmoid->GetOutput());
speedWriter->SetFileName("GeodesicActiveContourImageFilterOutput3.mha");
try
{
speedWriter->Update();
}
catch (itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
InternalWriterType::Pointer gradientWriter = InternalWriterType::New();
gradientWriter->SetInput(gradientMagnitude->GetOutput());
gradientWriter->SetFileName("GeodesicActiveContourImageFilterOutput2.mha");
try
{
gradientWriter->Update();
}
catch (itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
Classes demonstrated¶
-
template<typename
TInputImage
, typenameTFeatureImage
, typenameTOutputPixelType
= float>
classGeodesicActiveContourLevelSetImageFilter
: public itk::SegmentationLevelSetImageFilter<TInputImage, TFeatureImage, TOutputPixelType> Segments structures in images based on a user supplied edge potential map.
- IMPORTANT
The SegmentationLevelSetImageFilter class and the GeodesicActiveContourLevelSetFunction class contain additional information necessary to gain full understanding of how to use this filter.
- OVERVIEW
This class is a level set method segmentation filter. An initial contour is propagated outwards (or inwards) until it ‘’sticks’’ to the shape boundaries. This is done by using a level set speed function based on a user supplied edge potential map.
- INPUTS
This filter requires two inputs. The first input is a initial level set. The initial level set is a real image which contains the initial contour/surface as the zero level set. For example, a signed distance function from the initial contour/surface is typically used. Unlike the simpler ShapeDetectionLevelSetImageFilter the initial contour does not have to lie wholly within the shape to be segmented. The initial contour is allow to overlap the shape boundary. The extra advection term in the update equation behaves like a doublet and attracts the contour to the boundary. This approach for segmentation follows that of Caselles et al (1997).
The second input is the feature image. For this filter, this is the edge potential map. General characteristics of an edge potential map is that it has values close to zero in regions near the edges and values close to one inside the shape itself. Typically, the edge potential map is compute from the image gradient, for example:
where
is image intensity and
is the derivative of Gaussian operator.
This implementation allows the user to set the weights between the propagation, advection and curvature term using methods
SetPropagationScaling(), SetAdvectionScaling(), SetCurvatureScaling(). In general, the larger the CurvatureScaling, the smoother the resulting contour. To follow the implementation in Caselles et al paper, set the PropagationScaling to(the inflation or ballon force) and AdvectionScaling and CurvatureScaling both to 1.0.
See SegmentationLevelSetImageFilter and SparseFieldLevelSetImageFilter for more information on Inputs.
- PARAMETERS
The PropagationScaling parameter can be used to switch from propagation outwards (POSITIVE scaling parameter) versus propagating inwards (NEGATIVE scaling parameter).
- OUTPUTS
The filter outputs a single, scalar, real-valued image. Negative values in the output image represent the inside of the segmented region and positive values in the image represent the outside of the segmented region. The zero crossings of the image correspond to the position of the propagating front.
See SparseFieldLevelSetImageFilter and SegmentationLevelSetImageFilter for more information.
- REFERENCES
”Geodesic Active Contours”, V. Caselles, R. Kimmel and G. Sapiro. International Journal on Computer Vision, Vol 22, No. 1, pp 61-97, 1997
- See
SegmentationLevelSetImageFilter
- See
GeodesicActiveContourLevelSetFunction
- See
SparseFieldLevelSetImageFilter