ITK  5.4.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 
22 #include "itkFFTWCommonExtended.h"
24 
26 
27 #include <vector>
28 
29 
30 namespace itk
31 {
32 
39 template <typename TInputImage, typename TOutputImage = TInputImage>
40 class ITK_TEMPLATE_EXPORT FFTWComplexToComplex1DFFTImageFilter
41  : public ComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>
42 {
43 public:
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 
72 protected:
75 
76  void
77  BeforeThreadedGenerateData() override;
78  void
79  ThreadedGenerateData(const OutputImageRegionType &, ThreadIdType threadID) override;
80 
84  GetImageRegionSplitter() const override;
85 
86 private:
87  FFTWComplexToComplex1DFFTImageFilter(const Self &); // purposely not implemented
88  void
89  operator=(const Self &); // purposely not implemented
90 
91  ImageRegionSplitterDirection::Pointer m_ImageRegionSplitter{};
92 
94  void
95  DestroyPlans();
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
107 template <>
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
itk::FFTImageFilterTraits< FFTWComplexToComplex1DFFTImageFilter >::OutputPixelType
std::complex< TUnderlying > OutputPixelType
Definition: itkFFTWComplexToComplex1DFFTImageFilter.h:113
itk::FFTImageFilterTraits< FFTWComplexToComplex1DFFTImageFilter >::InputPixelType
std::complex< TUnderlying > InputPixelType
Definition: itkFFTWComplexToComplex1DFFTImageFilter.h:111
itk::fftw::ComplexToComplexProxy
Definition: itkFFTWCommonExtended.h:45
itkImageRegionSplitterDirection.h
itk::FFTImageFilterTraits
Helper defining pixel traits for templated FFT image filters.
Definition: itkFFTImageFilterFactory.h:42
itk::SmartPointer< Self >
itk::FFTWComplexToComplex1DFFTImageFilter::PlanBufferPointerType
typename std::vector< typename FFTW1DProxyType::ComplexType * > PlanBufferPointerType
Definition: itkFFTWComplexToComplex1DFFTImageFilter.h:63
itkComplexToComplex1DFFTImageFilter.h
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::ImageRegionSplitterBase
Divide an image region into several pieces.
Definition: itkImageRegionSplitterBase.h:58
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ComplexToComplex1DFFTImageFilter
Perform the Fast Fourier Transform, complex input to complex output, but only along one dimension.
Definition: itkComplexToComplex1DFFTImageFilter.h:43
itkFFTWCommonExtended.h
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itk::FFTWComplexToComplex1DFFTImageFilter
only do FFT along one dimension using FFTW as a backend.
Definition: itkFFTWComplexToComplex1DFFTImageFilter.h:40
itk::FFTWComplexToComplex1DFFTImageFilter::PlanArrayType
typename std::vector< typename FFTW1DProxyType::PlanType > PlanArrayType
Definition: itkFFTWComplexToComplex1DFFTImageFilter.h:62
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::FFTWComplexToComplex1DFFTImageFilter::FFTW1DProxyType
typename fftw::ComplexToComplexProxy< typename TInputImage::PixelType::value_type > FFTW1DProxyType
Definition: itkFFTWComplexToComplex1DFFTImageFilter.h:61
itk::FFTImageFilterTraits< FFTWComplexToComplex1DFFTImageFilter >::FilterDimensions
std::integer_sequence< unsigned int, 4, 3, 2, 1 > FilterDimensions
Definition: itkFFTWComplexToComplex1DFFTImageFilter.h:114
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itkFFTImageFilterFactory.h
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90