ITK  5.4.0
Insight Toolkit
SphinxExamples/src/Filtering/FFT/FilterImageInFourierDomain/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.
*
*=========================================================================*/
int
main(int argc, char * argv[])
{
if (argc < 3)
{
std::cerr << "Usage: " << std::endl;
std::cerr << argv[0];
std::cerr << " <InputFileName> <OutputFileName>";
std::cerr << " [Sigma]";
std::cerr << std::endl;
return EXIT_FAILURE;
}
const char * inputFileName = argv[1];
const char * outputFileName = argv[2];
double sigmaValue = 0.5;
if (argc > 3)
{
sigmaValue = std::stod(argv[3]);
}
constexpr unsigned int Dimension = 3;
using PixelType = float;
using RealImageType = itk::Image<PixelType, Dimension>;
const auto input = itk::ReadImage<RealImageType>(inputFileName);
// Some FFT filter implementations, like VNL's, need the image size to be a
// multiple of small prime numbers.
auto padFilter = PadFilterType::New();
padFilter->SetInput(input);
// 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;
auto forwardFFTFilter = ForwardFFTFilterType::New();
forwardFFTFilter->SetInput(padFilter->GetOutput());
try
{
forwardFFTFilter->UpdateOutputInformation();
}
catch (const itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
// A Gaussian is used here to create a low-pass filter.
using GaussianSourceType = itk::GaussianImageSource<RealImageType>;
auto gaussianSource = GaussianSourceType::New();
gaussianSource->SetNormalized(true);
ComplexImageType::ConstPointer transformedInput = forwardFFTFilter->GetOutput();
const ComplexImageType::RegionType inputRegion(transformedInput->GetLargestPossibleRegion());
const ComplexImageType::SizeType inputSize = inputRegion.GetSize();
const ComplexImageType::SpacingType inputSpacing = transformedInput->GetSpacing();
const ComplexImageType::PointType inputOrigin = transformedInput->GetOrigin();
const ComplexImageType::DirectionType inputDirection = transformedInput->GetDirection();
gaussianSource->SetSize(inputSize);
gaussianSource->SetSpacing(inputSpacing);
gaussianSource->SetOrigin(inputOrigin);
gaussianSource->SetDirection(inputDirection);
GaussianSourceType::ArrayType sigma;
sigma.Fill(sigmaValue);
for (unsigned int ii = 0; ii < Dimension; ++ii)
{
const double halfLength = inputSize[ii] * inputSpacing[ii] / 2.0;
sigma[ii] *= halfLength;
mean[ii] = inputOrigin[ii] + halfLength;
}
mean = inputDirection * mean;
gaussianSource->SetSigma(sigma);
gaussianSource->SetMean(mean);
auto fftShiftFilter = FFTShiftFilterType::New();
fftShiftFilter->SetInput(gaussianSource->GetOutput());
auto multiplyFilter = MultiplyFilterType::New();
multiplyFilter->SetInput1(forwardFFTFilter->GetOutput());
multiplyFilter->SetInput2(fftShiftFilter->GetOutput());
auto inverseFFTFilter = InverseFilterType::New();
inverseFFTFilter->SetInput(multiplyFilter->GetOutput());
try
{
itk::WriteImage(inverseFFTFilter->GetOutput(), outputFileName);
}
catch (const itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itkInverseFFTImageFilter.h
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itkImageFileReader.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::ForwardFFTImageFilter
Base class for forward Fast Fourier Transform.
Definition: itkForwardFFTImageFilter.h:65
itk::WrapPadImageFilter
Increase the image size by padding with replicants of the input image value.
Definition: itkWrapPadImageFilter.h:54
itkWrapPadImageFilter.h
itkForwardFFTImageFilter.h
itkFFTShiftImageFilter.h
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkMultiplyImageFilter.h
itk::InverseFFTImageFilter
Base class for inverse Fast Fourier Transform.
Definition: itkInverseFFTImageFilter.h:51
itk::GaussianImageSource
Generate an n-dimensional image of a Gaussian.
Definition: itkGaussianImageSource.h:45
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:88
New
static Pointer New()
itkGaussianImageSource.h
itk::MultiplyImageFilter
Pixel-wise multiplication of two images.
Definition: itkMultiplyImageFilter.h:43
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:254
itk::FixedArray::Fill
void Fill(const ValueType &)
itk::Size::GetSize
const SizeValueType * GetSize() const
Definition: itkSize.h:171