ITK  6.0.0
Insight Toolkit
Examples/Filtering/FFTDirectInverse2.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 "itkConfigure.h"
//
// This example is based on the on that was contributed by Stephan in the
// users list
//
// https://public.kitware.com/pipermail/insight-users/2005-June/013482.html
//
//
// Software Guide : BeginLatex
//
// This example illustrates how to compute the direct Fourier transform
// followed by the inverse Fourier transform using the FFTW library.
//
// Software Guide : EndLatex
#include "itkImage.h"
#if !defined(ITK_USE_FFTWF)
// #error "This example only works when single precision FFTW is used
// Changing WorkPixeltype to double and changing this conditional to
// ITK_USE_FFTWD will also work.
#endif
int
main(int argc, char * argv[])
{
if (argc != 3)
{
std::cerr << "Usage: " << argv[0] << " input output" << std::endl;
return EXIT_FAILURE;
}
// Software Guide : BeginLatex
//
// First we set up the types of the input and output images.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
constexpr unsigned int Dimension = 2;
// using OutputPixelType = unsigned char;
using OutputPixelType = unsigned short;
using WorkPixelType = float;
using InputImageType = itk::Image<WorkPixelType, Dimension>;
using WorkImageType = itk::Image<WorkPixelType, Dimension>;
using OutputImageType = itk::Image<OutputPixelType, Dimension>;
// Software Guide : EndCodeSnippet
// File handling
auto inputreader = ReaderType::New();
auto writer = WriterType::New();
inputreader->SetFileName(argv[1]);
writer->SetFileName(argv[2]);
// Read the image and get its size
inputreader->Update();
// Forward FFT filter
auto fftinput = FFTFilterType::New();
fftinput->SetInput(inputreader->GetOutput());
fftinput->Update();
// This is the output type from the FFT filters
using ComplexImageType = FFTFilterType::OutputImageType;
// Do the inverse transform = forward transform + flip all axes
auto fftoutput = invFFTFilterType::New();
fftoutput->SetInput(
fftinput->GetOutput()); // try to recover the input image
fftoutput->Update();
// Rescale the output to suit the output image type
using RescaleFilterType =
auto intensityrescaler = RescaleFilterType::New();
std::cout << fftoutput->GetOutput()->GetLargestPossibleRegion().GetSize()
<< std::endl;
intensityrescaler->SetInput(fftoutput->GetOutput());
intensityrescaler->SetOutputMinimum(0);
intensityrescaler->SetOutputMaximum(65535);
// Write the output
writer->SetInput(intensityrescaler->GetOutput());
writer->Update();
// DEBUG: std::cout << "inputreader
// "<<inputreader->GetOutput()->GetLargestPossibleRegion().GetSize() <<
// std::endl; DEBUG: std::cout << "fftinput "
// <<fftinput->GetOutput()->GetLargestPossibleRegion().GetSize() <<
// std::endl; DEBUG: std::cout << "fftoutput "
// <<fftoutput->GetOutput()->GetLargestPossibleRegion().GetSize() <<
// std::endl; DEBUG: std::cout << "intensityrescaler "
// <<intensityrescaler->GetOutput()->GetLargestPossibleRegion().GetSize() <<
// std::endl;
return EXIT_SUCCESS;
}
itkFFTWForwardFFTImageFilter.h
itkFlipImageFilter.h
itkImageFileReader.h
itkImage.h
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:90
itk::FFTWForwardFFTImageFilter
FFTW-based forward Fast Fourier Transform.
Definition: itkFFTWForwardFFTImageFilter.h:57
itkRescaleIntensityImageFilter.h
itkImageFileWriter.h
itkFFTWInverseFFTImageFilter.h
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()
itkResampleImageFilter.h
itk::FFTWInverseFFTImageFilter
FFTW-based inverse Fast Fourier Transform.
Definition: itkFFTWInverseFFTImageFilter.h:51
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44