Perona Malik Anisotropic Diffusion on Grayscale Image¶
Synopsis¶
Perona Malik Anisotropic Diffusion for scalar valued images.
Results¶
![Anisotropic diffusion comparison.](../../../../_images/ComputePeronaMalikAnisotropicDiffusionComparison.png)
Before anisotropic diffusion (left) and after anisotropic diffusion (right).¶
Code¶
Python¶
#!/usr/bin/env python
import itk
import argparse
parser = argparse.ArgumentParser(
description="Compute Perona Malik Anisotropic Diffusion."
)
parser.add_argument("input_image")
parser.add_argument("output_image")
parser.add_argument("number_of_iterations", type=int)
parser.add_argument("conductance", type=float)
args = parser.parse_args()
Dimension = 2
InputPixelType = itk.UC
InputImageType = itk.Image[InputPixelType, Dimension]
OutputPixelType = itk.F
OutputImageType = itk.Image[OutputPixelType, Dimension]
ReaderType = itk.ImageFileReader[InputImageType]
reader = ReaderType.New()
reader.SetFileName(args.input_image)
CastFilterType = itk.CastImageFilter[InputImageType, OutputImageType]
castfilter = CastFilterType.New()
castfilter.SetInput(reader)
FilterType = itk.GradientAnisotropicDiffusionImageFilter[
OutputImageType, OutputImageType
]
gradientfilter = FilterType.New()
gradientfilter.SetInput(castfilter.GetOutput())
gradientfilter.SetNumberOfIterations(args.number_of_iterations)
gradientfilter.SetTimeStep(0.125)
gradientfilter.SetConductanceParameter(args.conductance)
WriterType = itk.ImageFileWriter[OutputImageType]
writer = WriterType.New()
writer.SetFileName(args.output_image)
writer.SetInput(gradientfilter.GetOutput())
writer.Update()
C++¶
#include "itkImage.h"
#include "itkImageFileReader.h"
#include "itkImageFileWriter.h"
#include "itkGradientAnisotropicDiffusionImageFilter.h"
int
main(int argc, char * argv[])
{
if (argc != 5)
{
std::cerr << "Usage: " << std::endl;
std::cerr << argv[0];
std::cerr << " <InputFileName>";
std::cerr << " <OutputFileName>";
std::cerr << " <NumberOfIterations> ";
std::cerr << " <Conductance>" << std::endl;
return EXIT_FAILURE;
}
constexpr unsigned int Dimension = 2;
using InputPixelType = unsigned char;
using InputImageType = itk::Image<InputPixelType, Dimension>;
using ReaderType = itk::ImageFileReader<InputImageType>;
ReaderType::Pointer reader = ReaderType::New();
reader->SetFileName(argv[1]);
using OutputPixelType = float;
using OutputImageType = itk::Image<OutputPixelType, Dimension>;
using FilterType = itk::GradientAnisotropicDiffusionImageFilter<InputImageType, OutputImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput(reader->GetOutput());
filter->SetNumberOfIterations(std::stoi(argv[3]));
filter->SetTimeStep(0.125);
filter->SetConductanceParameter(std::stod(argv[4]));
using WriterType = itk::ImageFileWriter<OutputImageType>;
WriterType::Pointer writer = WriterType::New();
writer->SetFileName(argv[2]);
writer->SetInput(filter->GetOutput());
try
{
writer->Update();
}
catch (itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
Classes demonstrated¶
-
template<typename
TInputImage
, typenameTOutputImage
>
classGradientAnisotropicDiffusionImageFilter
: public itk::AnisotropicDiffusionImageFilter<TInputImage, TOutputImage> This filter performs anisotropic diffusion on a scalar itk::Image using the classic Perona-Malik, gradient magnitude based equation.
For detailed information on anisotropic diffusion, see itkAnisotropicDiffusionFunction and itkGradientNDAnisotropicDiffusionFunction.
- Inputs and Outputs
The input to this filter should be a scalar itk::Image of any dimensionality. The output image will be a diffused copy of the input.
- Parameters
Please see the description of parameters given in itkAnisotropicDiffusionImageFilter.
- See
AnisotropicDiffusionImageFilter
- See
AnisotropicDiffusionFunction
- See
GradientAnisotropicDiffusionFunction