Watershed Image Filter

Synopsis

This example illustrates how to segment an image using the watershed method.

Results

Input image

Input image

Segmented image

Segmented image

Code

C++

#include "itkImageFileReader.h"
#include "itkImageFileWriter.h"
#include "itkScalarToRGBPixelFunctor.h"
#include "itkUnaryFunctorImageFilter.h"
#include "itkVectorGradientAnisotropicDiffusionImageFilter.h"
#include "itkWatershedImageFilter.h"
#include "itkRescaleIntensityImageFilter.h"
#include "itkScalarToRGBColormapImageFilter.h"
#include "itkGradientMagnitudeImageFilter.h"

// Run with:
// ./SegmentWithWatershedImageFilter inputImageFile outputImageFile threshold level
// e.g.
// ./SegmentWithWatershedImageFilter BrainProtonDensitySlice.png OutBrainWatershed.png 0.005 .5
// (A rule of thumb is to set the Threshold to be about 1 / 100 of the Level.)

int
main(int argc, char * argv[])
{
  if (argc < 5)
  {
    std::cerr << "Missing parameters." << std::endl;
    std::cerr << "Usage: " << argv[0] << " inputImageFile outputImageFile threshold level" << std::endl;
    return EXIT_FAILURE;
  }

  constexpr unsigned int Dimension = 2;

  using InputImageType = itk::Image<unsigned char, Dimension>;
  using FloatImageType = itk::Image<float, Dimension>;
  using RGBPixelType = itk::RGBPixel<unsigned char>;
  using RGBImageType = itk::Image<RGBPixelType, Dimension>;
  using LabeledImageType = itk::Image<itk::IdentifierType, Dimension>;

  using FileReaderType = itk::ImageFileReader<InputImageType>;
  FileReaderType::Pointer reader = FileReaderType::New();
  reader->SetFileName(argv[1]);

  using GradientMagnitudeImageFilterType = itk::GradientMagnitudeImageFilter<InputImageType, FloatImageType>;
  GradientMagnitudeImageFilterType::Pointer gradientMagnitudeImageFilter = GradientMagnitudeImageFilterType::New();

  gradientMagnitudeImageFilter->SetInput(reader->GetOutput());
  gradientMagnitudeImageFilter->Update();

  using WatershedFilterType = itk::WatershedImageFilter<FloatImageType>;
  WatershedFilterType::Pointer watershed = WatershedFilterType::New();

  float threshold = std::stod(argv[3]);
  float level = std::stod(argv[4]);

  watershed->SetThreshold(threshold);
  watershed->SetLevel(level);

  watershed->SetInput(gradientMagnitudeImageFilter->GetOutput());
  watershed->Update();

  using RGBFilterType = itk::ScalarToRGBColormapImageFilter<LabeledImageType, RGBImageType>;
  RGBFilterType::Pointer colormapImageFilter = RGBFilterType::New();
  colormapImageFilter->SetColormap(itk::RGBColormapFilterEnumType::Jet);
  colormapImageFilter->SetInput(watershed->GetOutput());
  colormapImageFilter->Update();

  using FileWriterType = itk::ImageFileWriter<RGBImageType>;
  FileWriterType::Pointer writer = FileWriterType::New();
  writer->SetFileName(argv[2]);
  writer->SetInput(colormapImageFilter->GetOutput());
  writer->Update();

  return EXIT_SUCCESS;
}

Python

#!/usr/bin/env python
# ./WatershedImageFilter.py <InputFileName> <OutputFileName> <Threshold> <Level>
# e.g.
# ./WatershedImageFilter.py BrainProtonDensitySlice.png OutBrainWatershed.png 0.005 .5
# (A rule of thumb is to set the Threshold to be about 1 / 100 of the Level.)
#
#  threshold: absolute minimum height value used during processing.
#        Raising this threshold percentage effectively decreases the number of local minima in the input, 
#        resulting in an initial segmentation with fewer regions. 
#        The assumption is that the shallow regions that thresholding removes are of of less interest.
#  level: parameter controls the depth of metaphorical flooding of the image. 
#        That is, it sets the maximum saliency value of interest in the result. 
#        Raising and lowering the Level influences the number of segments 
#        in the basic segmentation that are merged to produce the final output. 
#        A level of 1.0 is analogous to flooding the image up to a 
#        depth that is 100 percent of the maximum value in the image. 
#        A level of 0.0 produces the basic segmentation, which will typically be very oversegmented. 
#        Level values of interest are typically low (i.e. less than about 0.40 or 40%),
#        since higher values quickly start to undersegment the image.

import sys
import itk

if len(sys.argv) != 5:
    print('Usage: ' + sys.argv[0] +
          ' <InputFileName> <OutputFileName> <Threshold> <Level>')
    sys.exit(1)

inputFileName = sys.argv[1]
outputFileName = sys.argv[2]

Dimension = 2

FloatPixelType = itk.ctype('float')
FloatImageType = itk.Image[FloatPixelType, Dimension]

reader = itk.ImageFileReader[FloatImageType].New()
reader.SetFileName(inputFileName)

gradientMagnitude = \
    itk.GradientMagnitudeImageFilter.New(Input=reader.GetOutput())

watershed = \
    itk.WatershedImageFilter.New(Input=gradientMagnitude.GetOutput())

threshold = float(sys.argv[3])
level = float(sys.argv[4])
watershed.SetThreshold(threshold)
watershed.SetLevel(level)

LabeledImageType = type(watershed.GetOutput())

PixelType = itk.ctype('unsigned char')
RGBPixelType = itk.RGBPixel[PixelType]
RGBImageType = itk.Image[RGBPixelType, Dimension]

ScalarToRGBColormapFilterType = \
    itk.ScalarToRGBColormapImageFilter[LabeledImageType, RGBImageType]
colormapImageFilter = ScalarToRGBColormapFilterType.New()
colormapImageFilter.SetColormap(ScalarToRGBColormapFilterType.Jet)
colormapImageFilter.SetInput(watershed.GetOutput())
colormapImageFilter.Update()

WriterType = itk.ImageFileWriter[RGBImageType]
writer = WriterType.New()
writer.SetFileName(outputFileName)
writer.SetInput(colormapImageFilter.GetOutput())
writer.Update()

Classes demonstrated