ITK  6.0.0
Insight Toolkit
SphinxExamples/src/Filtering/ImageGradient/ApplyGradientRecursiveGaussian/Code.cxx
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* https://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#include "itkImage.h"
int
main(int argc, char * argv[])
{
if (argc != 5)
{
std::cerr << "Usage: " << std::endl;
std::cerr << argv[0];
std::cerr << " <InputFileName> <OutputFileNameX> <OutputFileNameY> <OutputFileNameMagnitude>";
std::cerr << std::endl;
return EXIT_FAILURE;
}
const char * inputFileName = argv[1];
const char * outputFileNameX = argv[2];
const char * outputFileNameY = argv[3];
const char * outputFileNameMagnitude = argv[4];
const char * filenames[2];
filenames[0] = outputFileNameX;
filenames[1] = outputFileNameY;
constexpr unsigned int Dimension = 2;
// Input and output are png files, use unsigned char
using PixelType = unsigned char;
// Double type for GradientRecursiveGaussianImageFilter
using DoublePixelType = double;
using DoubleImageType = itk::Image<DoublePixelType, Dimension>;
// The output of GradientRecursiveGaussianImageFilter
// are images of the gradient along X and Y, so the type of
// the output is a covariant vector of dimension 2 (X, Y)
const auto input = itk::ReadImage<ImageType>(inputFileName);
auto filter = FilterType::New();
filter->SetInput(input);
// Allows to select the X or Y output images
auto indexSelectionFilter = IndexSelectionType::New();
indexSelectionFilter->SetInput(filter->GetOutput());
// Rescale for png output
auto rescaler = RescalerType::New();
rescaler->SetOutputMinimum(itk::NumericTraits<PixelType>::min());
rescaler->SetOutputMaximum(itk::NumericTraits<PixelType>::max());
rescaler->SetInput(indexSelectionFilter->GetOutput());
// Write the X and Y images
for (int i = 0; i < 2; ++i)
{
indexSelectionFilter->SetIndex(i);
try
{
itk::WriteImage(rescaler->GetOutput(), filenames[i]);
}
catch (const itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
}
// Compute the magnitude of the vector and output the image
auto magnitudeFilter = MagnitudeType::New();
magnitudeFilter->SetInput(filter->GetOutput());
// Rescale for png output
rescaler->SetInput(magnitudeFilter->GetOutput());
try
{
itk::WriteImage(rescaler->GetOutput(), outputFileNameMagnitude);
}
catch (const itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
itk::VectorMagnitudeImageFilter
Take an image of vectors as input and produce an image with the magnitude of those vectors.
Definition: itkVectorMagnitudeImageFilter.h:68
itkVectorIndexSelectionCastImageFilter.h
itkImageFileReader.h
itkImage.h
itk::VectorIndexSelectionCastImageFilter
Extracts the selected index of the vector that is the input pixel type.
Definition: itkVectorIndexSelectionCastImageFilter.h:88
itkGradientRecursiveGaussianImageFilter.h
itk::GradientRecursiveGaussianImageFilter
Computes the gradient of an image by convolution with the first derivative of a Gaussian.
Definition: itkGradientRecursiveGaussianImageFilter.h:59
itkRescaleIntensityImageFilter.h
itkVectorMagnitudeImageFilter.h
itkImageFileWriter.h
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:60
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::RescaleIntensityImageFilter
Applies a linear transformation to the intensity levels of the input Image.
Definition: itkRescaleIntensityImageFilter.h:133
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itk::WriteImage
ITK_TEMPLATE_EXPORT void WriteImage(TImagePointer &&image, const std::string &filename, bool compress=false)
Definition: itkImageFileWriter.h:256