ITK  6.0.0
Insight Toolkit
itkFFTWForward1DFFTImageFilter.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 itkFFTWForward1DFFTImageFilter_h
19#define itkFFTWForward1DFFTImageFilter_h
20
24
26
27#include <vector>
28
29
30namespace itk
31{
32
39template <typename TInputImage,
40 typename TOutputImage = Image<std::complex<typename TInputImage::PixelType>, TInputImage::ImageDimension>>
41class ITK_TEMPLATE_EXPORT FFTWForward1DFFTImageFilter : public Forward1DFFTImageFilter<TInputImage, TOutputImage>
42{
43public:
44 ITK_DISALLOW_COPY_AND_MOVE(FFTWForward1DFFTImageFilter);
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(FFTWForward1DFFTImageFilter);
72
73
74protected:
77
78 void
80 void
82
86 GetImageRegionSplitter() const override;
87
88
89private:
91
93 void
95
96 bool m_PlanComputed{ false };
97 PlanArrayType m_PlanArray{};
98 unsigned int m_LastImageSize{ 0 };
99 PlanBufferPointerType m_InputBufferArray{};
100 PlanBufferPointerType m_OutputBufferArray{};
101};
102
103
104// Describe whether input/output are real- or complex-valued
105// for factory registration
106template <>
108{
109 template <typename TUnderlying>
110 using InputPixelType = TUnderlying;
111 template <typename TUnderlying>
112 using OutputPixelType = std::complex<TUnderlying>;
113 using FilterDimensions = std::integer_sequence<unsigned int, 4, 3, 2, 1>;
114};
115
116} // namespace itk
117
118#ifndef ITK_MANUAL_INSTANTIATION
119# include "itkFFTWForward1DFFTImageFilter.hxx"
120#endif
121
122#endif // itkFFTWForward1DFFTImageFilter_h
only do FFT along one dimension using FFTW as a backend.
void ThreadedGenerateData(const OutputImageRegionType &, ThreadIdType threadID) override
typename Superclass::OutputImageType OutputImageType
void BeforeThreadedGenerateData() override
const ImageRegionSplitterBase * GetImageRegionSplitter() const override
typename std::vector< typename FFTW1DProxyType::ComplexType * > PlanBufferPointerType
typename std::vector< typename FFTW1DProxyType::PlanType > PlanArrayType
typename Superclass::InputImageType InputImageType
typename OutputImageType::RegionType OutputImageRegionType
typename fftw::ComplexToComplexProxy< typename TInputImage::PixelType > FFTW1DProxyType
Perform the Fast Fourier Transform, in the forward direction, with real inputs, but only along one di...
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
std::integer_sequence< unsigned int, 4, 3, 2, 1 > FilterDimensions
Helper defining pixel traits for templated FFT image filters.