ITK  5.2.0
Insight Toolkit
SphinxExamples/src/Filtering/FFT/ComputeForwardFFT/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
*
* http://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> <Real filename> <Imaginary filename> <Modulus filename>";
std::cerr << std::endl;
return EXIT_FAILURE;
}
const char * inputFileName = argv[1];
const char * realFileName = argv[2];
const char * imaginaryFileName = argv[3];
const char * modulusFileName = argv[4];
constexpr unsigned int Dimension = 3;
using FloatPixelType = float;
using FloatImageType = itk::Image<FloatPixelType, Dimension>;
ReaderType::Pointer reader = ReaderType::New();
reader->SetFileName(inputFileName);
using UnsignedCharPixelType = unsigned char;
using UnsignedCharImageType = itk::Image<UnsignedCharPixelType, Dimension>;
// Some FFT filter implementations, like VNL's, need the image size to be a
// multiple of small prime numbers.
PadFilterType::Pointer padFilter = PadFilterType::New();
padFilter->SetInput(reader->GetOutput());
// Input size is [48, 62, 42]. Pad to [48, 64, 48].
padding[0] = 0;
padding[1] = 2;
padding[2] = 6;
padFilter->SetPadUpperBound(padding);
FFTType::Pointer fftFilter = FFTType::New();
fftFilter->SetInput(padFilter->GetOutput());
using FFTOutputImageType = FFTType::OutputImageType;
// Extract the real part
RealFilterType::Pointer realFilter = RealFilterType::New();
realFilter->SetInput(fftFilter->GetOutput());
RescaleFilterType::Pointer realRescaleFilter = RescaleFilterType::New();
realRescaleFilter->SetInput(realFilter->GetOutput());
realRescaleFilter->SetOutputMinimum(itk::NumericTraits<UnsignedCharPixelType>::min());
realRescaleFilter->SetOutputMaximum(itk::NumericTraits<UnsignedCharPixelType>::max());
WriterType::Pointer realWriter = WriterType::New();
realWriter->SetFileName(realFileName);
realWriter->SetInput(realRescaleFilter->GetOutput());
try
{
realWriter->Update();
}
catch (itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
// Extract the imaginary part
ImaginaryFilterType::Pointer imaginaryFilter = ImaginaryFilterType::New();
imaginaryFilter->SetInput(fftFilter->GetOutput());
RescaleFilterType::Pointer imaginaryRescaleFilter = RescaleFilterType::New();
imaginaryRescaleFilter->SetInput(imaginaryFilter->GetOutput());
imaginaryRescaleFilter->SetOutputMinimum(itk::NumericTraits<UnsignedCharPixelType>::min());
imaginaryRescaleFilter->SetOutputMaximum(itk::NumericTraits<UnsignedCharPixelType>::max());
WriterType::Pointer complexWriter = WriterType::New();
complexWriter->SetFileName(imaginaryFileName);
complexWriter->SetInput(imaginaryRescaleFilter->GetOutput());
try
{
complexWriter->Update();
}
catch (itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
// Compute the magnitude
ModulusFilterType::Pointer modulusFilter = ModulusFilterType::New();
modulusFilter->SetInput(fftFilter->GetOutput());
RescaleFilterType::Pointer magnitudeRescaleFilter = RescaleFilterType::New();
magnitudeRescaleFilter->SetInput(modulusFilter->GetOutput());
magnitudeRescaleFilter->SetOutputMinimum(itk::NumericTraits<UnsignedCharPixelType>::min());
magnitudeRescaleFilter->SetOutputMaximum(itk::NumericTraits<UnsignedCharPixelType>::max());
WriterType::Pointer magnitudeWriter = WriterType::New();
magnitudeWriter->SetFileName(modulusFileName);
magnitudeWriter->SetInput(magnitudeRescaleFilter->GetOutput());
try
{
magnitudeWriter->Update();
}
catch (itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
itkImageFileReader.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itk::ForwardFFTImageFilter
Base class for forward Fast Fourier Transform.
Definition: itkForwardFFTImageFilter.h:62
itk::WrapPadImageFilter
Increase the image size by padding with replicants of the input image value.
Definition: itkWrapPadImageFilter.h:54
itkComplexToImaginaryImageFilter.h
itkWrapPadImageFilter.h
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itkForwardFFTImageFilter.h
itk::ComplexToImaginaryImageFilter
Computes pixel-wise the imaginary part of a complex image.
Definition: itkComplexToImaginaryImageFilter.h:63
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:88
itkComplexToModulusImageFilter.h
itkComplexToRealImageFilter.h
itkRescaleIntensityImageFilter.h
itkImageFileWriter.h
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:58
itk::RescaleIntensityImageFilter
Applies a linear transformation to the intensity levels of the input Image.
Definition: itkRescaleIntensityImageFilter.h:154
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itk::ComplexToModulusImageFilter
Computes pixel-wise the Modulus of a complex image.
Definition: itkComplexToModulusImageFilter.h:61
itk::ComplexToRealImageFilter
Computes pixel-wise the real(x) part of a complex image.
Definition: itkComplexToRealImageFilter.h:61