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 
22 #include "itkFFTWCommonExtended.h"
24 
26 
27 #include <vector>
28 
29 namespace itk
30 {
38 template <typename TInputImage,
39  typename TOutputImage =
41 class ITK_TEMPLATE_EXPORT FFTWInverse1DFFTImageFilter : public Inverse1DFFTImageFilter<TInputImage, TOutputImage>
42 {
43 public:
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 
74 protected:
76  ~FFTWInverse1DFFTImageFilter() override;
77 
78  void
79  BeforeThreadedGenerateData() override;
80  void
81  ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadID) override;
82 
86  GetImageRegionSplitter() const override;
87 
88 private:
89  ImageRegionSplitterDirection::Pointer m_ImageRegionSplitter{};
90 
92  void
93  DestroyPlans();
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
105 template <>
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
itk::FFTWInverse1DFFTImageFilter
only do FFT along one dimension using FFTW as a backend.
Definition: itkFFTWInverse1DFFTImageFilter.h:41
itk::fftw::ComplexToComplexProxy
Definition: itkFFTWCommonExtended.h:45
itk::FFTImageFilterTraits< FFTWInverse1DFFTImageFilter >::InputPixelType
std::complex< TUnderlying > InputPixelType
Definition: itkFFTWInverse1DFFTImageFilter.h:109
itkImageRegionSplitterDirection.h
itkInverse1DFFTImageFilter.h
itk::FFTImageFilterTraits
Helper defining pixel traits for templated FFT image filters.
Definition: itkFFTImageFilterFactory.h:42
itk::SmartPointer< Self >
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:102
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::Image::ValueType
TPixel ValueType
Definition: itkImage.h:111
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
itkFFTWCommonExtended.h
itk::FFTWInverse1DFFTImageFilter::FFTW1DProxyType
typename fftw::ComplexToComplexProxy< typename TOutputImage::PixelType > FFTW1DProxyType
Definition: itkFFTWInverse1DFFTImageFilter.h:63
itk::FFTWInverse1DFFTImageFilter::PlanArrayType
typename std::vector< typename FFTW1DProxyType::PlanType > PlanArrayType
Definition: itkFFTWInverse1DFFTImageFilter.h:64
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itk::FFTImageFilterTraits< FFTWInverse1DFFTImageFilter >::OutputPixelType
TUnderlying OutputPixelType
Definition: itkFFTWInverse1DFFTImageFilter.h:111
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itk::Inverse1DFFTImageFilter
Perform the Fast Fourier Transform, in the reverse direction, with real output, but only along one di...
Definition: itkInverse1DFFTImageFilter.h:38
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::FFTImageFilterTraits< FFTWInverse1DFFTImageFilter >::FilterDimensions
std::integer_sequence< unsigned int, 4, 3, 2, 1 > FilterDimensions
Definition: itkFFTWInverse1DFFTImageFilter.h:112
itk::FFTWInverse1DFFTImageFilter::PlanBufferPointerType
typename std::vector< typename FFTW1DProxyType::ComplexType * > PlanBufferPointerType
Definition: itkFFTWInverse1DFFTImageFilter.h:65
itkFFTImageFilterFactory.h
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90