ITK  6.0.0
Insight Toolkit
itkFFTWComplexToComplex1DFFTImageFilter.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 itkFFTWComplexToComplex1DFFTImageFilter_h
19#define itkFFTWComplexToComplex1DFFTImageFilter_h
20
24
26
27#include <vector>
28
29
30namespace itk
31{
32
39template <typename TInputImage, typename TOutputImage = TInputImage>
40class ITK_TEMPLATE_EXPORT FFTWComplexToComplex1DFFTImageFilter
41 : public ComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>
42{
43public:
48
50 using InputImageType = typename Superclass::InputImageType;
51 using OutputImageType = typename Superclass::OutputImageType;
53
62 using PlanArrayType = typename std::vector<typename FFTW1DProxyType::PlanType>;
63 using PlanBufferPointerType = typename std::vector<typename FFTW1DProxyType::ComplexType *>;
64
66 itkNewMacro(Self);
67
69 itkOverrideGetNameOfClassMacro(FFTWComplexToComplex1DFFTImageFilter);
70
71
72protected:
75
76 void
78 void
80
84 GetImageRegionSplitter() const override;
85
86private:
87 FFTWComplexToComplex1DFFTImageFilter(const Self &); // purposely not implemented
88 void
89 operator=(const Self &); // purposely not implemented
90
92
94 void
96
97 bool m_PlanComputed{ false };
98 PlanArrayType m_PlanArray{};
99 unsigned int m_LastImageSize{ 0 };
100 PlanBufferPointerType m_InputBufferArray{};
101 PlanBufferPointerType m_OutputBufferArray{};
102};
103
104
105// Describe whether input/output are real- or complex-valued
106// for factory registration
107template <>
109{
110 template <typename TUnderlying>
111 using InputPixelType = std::complex<TUnderlying>;
112 template <typename TUnderlying>
113 using OutputPixelType = std::complex<TUnderlying>;
114 using FilterDimensions = std::integer_sequence<unsigned int, 4, 3, 2, 1>;
115};
116
117} // namespace itk
118
119#ifndef ITK_MANUAL_INSTANTIATION
120# include "itkFFTWComplexToComplex1DFFTImageFilter.hxx"
121#endif
122
123#endif // itkFFTWComplexToComplex1DFFTImageFilter_h
Perform the Fast Fourier Transform, complex input to complex output, but only along one dimension.
only do FFT along one dimension using FFTW as a backend.
void ThreadedGenerateData(const OutputImageRegionType &, ThreadIdType threadID) override
typename std::vector< typename FFTW1DProxyType::PlanType > PlanArrayType
typename fftw::ComplexToComplexProxy< typename TInputImage::PixelType::value_type > FFTW1DProxyType
const ImageRegionSplitterBase * GetImageRegionSplitter() const override
typename std::vector< typename FFTW1DProxyType::ComplexType * > PlanBufferPointerType
Divide an image region into several pieces.
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
Helper defining pixel traits for templated FFT image filters.