ITK  5.2.0
Insight Toolkit
SphinxExamples/src/Filtering/FFT/ComputeImageSpectralDensity/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.
*
*=========================================================================*/
int
main(int argc, char * argv[])
{
if (argc != 3)
{
std::cerr << "Usage: " << std::endl;
std::cerr << argv[0];
std::cerr << " <InputFileName> <OutputFileName>";
std::cerr << std::endl;
return EXIT_FAILURE;
}
const char * inputFileName = argv[1];
const char * outputFileName = argv[2];
constexpr unsigned int Dimension = 3;
using PixelType = float;
using RealImageType = itk::Image<PixelType, Dimension>;
ReaderType::Pointer reader = ReaderType::New();
reader->SetFileName(inputFileName);
// 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);
using ForwardFFTFilterType = itk::ForwardFFTImageFilter<RealImageType>;
using ComplexImageType = ForwardFFTFilterType::OutputImageType;
ForwardFFTFilterType::Pointer forwardFFTFilter = ForwardFFTFilterType::New();
forwardFFTFilter->SetInput(padFilter->GetOutput());
ComplexToModulusFilterType::Pointer complexToModulusFilter = ComplexToModulusFilterType::New();
complexToModulusFilter->SetInput(forwardFFTFilter->GetOutput());
// Window and shift the output for visualization.
using OutputPixelType = unsigned short;
using OutputImageType = itk::Image<OutputPixelType, Dimension>;
WindowingFilterType::Pointer windowingFilter = WindowingFilterType::New();
windowingFilter->SetInput(complexToModulusFilter->GetOutput());
windowingFilter->SetWindowMinimum(0);
windowingFilter->SetWindowMaximum(20000);
FFTShiftFilterType::Pointer fftShiftFilter = FFTShiftFilterType::New();
fftShiftFilter->SetInput(windowingFilter->GetOutput());
WriterType::Pointer writer = WriterType::New();
writer->SetFileName(outputFileName);
writer->SetInput(fftShiftFilter->GetOutput());
try
{
writer->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
itkIntensityWindowingImageFilter.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
itkWrapPadImageFilter.h
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itkForwardFFTImageFilter.h
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:88
itkFFTShiftImageFilter.h
itkComplexToModulusImageFilter.h
itk::IntensityWindowingImageFilter
Applies a linear transformation to the intensity levels of the input Image that are inside a user-def...
Definition: itkIntensityWindowingImageFilter.h:149
itkImageFileWriter.h
itk::FFTShiftImageFilter
Shift the zero-frequency components of a Fourier transform to the center of the image.
Definition: itkFFTShiftImageFilter.h:49
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