[Insight-users] watershed segmentation without RGB?

Miller, James V (Research) millerjv at crd.ge.com
Thu, 29 Apr 2004 09:40:33 -0400


Watershed can be run on grayscale images.  The example
(WatershedSegmentation1.cxx) is set up to process RGB pixels. But you can
easily run watershed on a grayscale image.

You'll want to create a similar program to WatershedSegmentation1.cxx that
reads a scalar image, produces a gradient magnitude image
(GradientMagnitudeImageFilter or RecursiveGradientMagnitudeImageFilter), and
runs watershed on the gradient magnitude image.

So you just need to change the pixel type, delete the word "Vector" from any
class type in example, and change the CastType to a CastImageFilter that
converts a scalar image into an image of floating point pixels (alternative,
you can remove the cast filter altogether and have the reader return an
image of floating point pixels.  Here is a quick stab at it (haven't tried
to compile it....)

#ifdef _MSC_VER
#pragma warning ( disable : 4786 )
#endif

#include <iostream>

#include "itkGradientAnisotropicDiffusionImageFilter.h"
#include "itkGradientMagnitudeImageFilter.h"
#include "itkWatershedImageFilter.h"

#include "itkImageFileReader.h"
#include "itkImageFileWriter.h"
#include "itkUnaryFunctorImageFilter.h"
#include "itkScalarToRGBPixelFunctor.h"

int main( int argc, char *argv[] )
{
  if (argc < 8 )
    {
    std::cerr << "Missing Parameters " << std::endl;
    std::cerr << "Usage: " << argv[0];
    std::cerr << " inputImage outputImage conductanceTerm
diffusionIterations lowerThreshold outputScaleLevel gradientMode " <<
std::endl;
    return 1;
    }
  
  typedef itk::RGBPixel<unsigned char>   RGBPixelType;
  typedef itk::Image<RGBPixelType, 2>    RGBImageType;
  typedef itk::Image<unsigned long, 2>   LabeledImageType;
  typedef itk::Image<float, 2>           ScalarImageType;

  typedef itk::ImageFileReader<ScalarImageType> FileReaderType;
  typedef itk::GradientAnisotropicDiffusionImageFilter<ScalarImageType,
    ScalarImageType>  DiffusionFilterType;
  typedef itk::GradientMagnitudeImageFilter<ScalarImageType>
    ScalarMagnitudeFilterType; 
  typedef itk::WatershedImageFilter<ScalarImageType> WatershedFilterType;
  // Software Guide : EndCodeSnippet

  typedef itk::ImageFileWriter<RGBImageType> FileWriterType;

  FileReaderType::Pointer reader = FileReaderType::New();
  reader->SetFileName(argv[1]);
  
  DiffusionFilterType::Pointer diffusion = DiffusionFilterType::New();
  diffusion->SetNumberOfIterations( atoi(argv[4]) );
  diffusion->SetConductanceParameter( atof(argv[3]) );
  diffusion->SetTimeStep(0.125);


  GradientMagnitudeFilterType::Pointer gradient =
GradientMagnitudeFilterType::New();
  gradient->SetUsePrincipleComponents(atoi(argv[7]));


  WatershedFilterType::Pointer watershed = WatershedFilterType::New();
  watershed->SetLevel( atof(argv[6]) );
  watershed->SetThreshold( atof(argv[5]) );


  typedef itk::Functor::ScalarToRGBPixelFunctor<unsigned long>
    ColorMapFunctorType;
  typedef itk::UnaryFunctorImageFilter<LabeledImageType,
    RGBImageType, ColorMapFunctorType> ColorMapFilterType;
  ColorMapFilterType::Pointer colormapper = ColorMapFilterType::New();
  
  FileWriterType::Pointer writer = FileWriterType::New();
  writer->SetFileName(argv[2]);

  diffusion->SetInput(reader->GetOutput());
  gradient->SetInput(diffusion->GetOutput());
  watershed->SetInput(gradient->GetOutput());
  colormapper->SetInput(watershed->GetOutput());
  writer->SetInput(colormapper->GetOutput());

  try 
    {
    writer->Update();
    }
  catch (itk::ExceptionObject &e)
    {
    std::cerr << e << std::endl;
    }
    
  return 0;
}





-----Original Message-----
From: Laurent Mundeleer [mailto:lmundele at ulb.ac.be]
Sent: Thursday, April 29, 2004 6:10 AM
To: insight-users at itk.org
Subject: [Insight-users] watershed segmentation without RGB?


Hi all,

I'd like to apply the watershed filter to CT images, is it possible? In 
the documentation, it 's specified that RGB values are required...
How can I do?

Thanks in advance

Regards,

Laurent

-- 
********************************************
Laurent Mundeleer
Universite Libre de Bruxelles (ULB)
Service des Systemes Logiques et Numeriques (SLN) CP165/57
50, Av. F.Roosevelt
1050 Bruxelles
Belgium
tel : ++32.2.650.22.97
fax : ++32.2.650.22.98
e-mail : lmundele at ulb.ac.be
********************************************

_______________________________________________
Insight-users mailing list
Insight-users at itk.org
http://www.itk.org/mailman/listinfo/insight-users