ITK  6.0.0
Insight Toolkit
itkFFTWInverse1DFFTImageFilter.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright NumFOCUS
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License");
6 * you may not use this file except in compliance with the License.
7 * You may obtain a copy of the License at
8 *
9 * https://www.apache.org/licenses/LICENSE-2.0.txt
10 *
11 * Unless required by applicable law or agreed to in writing, software
12 * distributed under the License is distributed on an "AS IS" BASIS,
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 * See the License for the specific language governing permissions and
15 * limitations under the License.
16 *
17 *=========================================================================*/
18#ifndef itkFFTWInverse1DFFTImageFilter_h
19#define itkFFTWInverse1DFFTImageFilter_h
20
24
26
27#include <vector>
28
29namespace itk
30{
38template <typename TInputImage,
39 typename TOutputImage =
41class ITK_TEMPLATE_EXPORT FFTWInverse1DFFTImageFilter : public Inverse1DFFTImageFilter<TInputImage, TOutputImage>
42{
43public:
44 ITK_DISALLOW_COPY_AND_MOVE(FFTWInverse1DFFTImageFilter);
45
50
52 using InputImageType = typename Superclass::InputImageType;
53 using OutputImageType = typename Superclass::OutputImageType;
55
64 using PlanArrayType = typename std::vector<typename FFTW1DProxyType::PlanType>;
65 using PlanBufferPointerType = typename std::vector<typename FFTW1DProxyType::ComplexType *>;
66
68 itkNewMacro(Self);
69
71 itkOverrideGetNameOfClassMacro(FFTWInverse1DFFTImageFilter);
72
73
74protected:
77
78 void
80 void
81 ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadID) override;
82
86 GetImageRegionSplitter() const override;
87
88private:
90
92 void
94
95 bool m_PlanComputed{ false };
96 PlanArrayType m_PlanArray{};
97 unsigned int m_LastImageSize{ 0 };
98 PlanBufferPointerType m_InputBufferArray{};
99 PlanBufferPointerType m_OutputBufferArray{};
100};
101
102
103// Describe whether input/output are real- or complex-valued
104// for factory registration
105template <>
107{
108 template <typename TUnderlying>
109 using InputPixelType = std::complex<TUnderlying>;
110 template <typename TUnderlying>
111 using OutputPixelType = TUnderlying;
112 using FilterDimensions = std::integer_sequence<unsigned int, 4, 3, 2, 1>;
113};
114
115} // namespace itk
116
117#ifndef ITK_MANUAL_INSTANTIATION
118# include "itkFFTWInverse1DFFTImageFilter.hxx"
119#endif
120
121#endif // itkFFTWInverse1DFFTImageFilter_h
only do FFT along one dimension using FFTW as a backend.
typename std::vector< typename FFTW1DProxyType::ComplexType * > PlanBufferPointerType
typename fftw::ComplexToComplexProxy< typename TOutputImage::PixelType > FFTW1DProxyType
const ImageRegionSplitterBase * GetImageRegionSplitter() const override
typename Superclass::InputImageType InputImageType
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadID) override
void BeforeThreadedGenerateData() override
typename Superclass::OutputImageType OutputImageType
typename std::vector< typename FFTW1DProxyType::PlanType > PlanArrayType
typename OutputImageType::RegionType OutputImageRegionType
Divide an image region into several pieces.
TPixel ValueType
Definition: itkImage.h:111
Perform the Fast Fourier Transform, in the reverse direction, with real output, but only along one di...
Light weight base class for most itk classes.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType
Definition: itkIntTypes.h:102
std::integer_sequence< unsigned int, 4, 3, 2, 1 > FilterDimensions
Helper defining pixel traits for templated FFT image filters.