ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkFFTWInverseFFTImageFilter.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 itkFFTWInverseFFTImageFilter_h
19#define itkFFTWInverseFFTImageFilter_h
20
22#include "itkFFTWCommon.h"
23// Include fftw3.h directly: the proxy header skips it when first included before ITK_USE_FFTW* is defined.
24#if defined(ITK_USE_FFTWF) || defined(ITK_USE_FFTWD)
25# include "fftw3.h"
26#endif
27
29
30namespace itk
31{
53template <typename TInputImage,
55class ITK_TEMPLATE_EXPORT FFTWInverseFFTImageFilter : public InverseFFTImageFilter<TInputImage, TOutputImage>
56{
57public:
58 ITK_DISALLOW_COPY_AND_MOVE(FFTWInverseFFTImageFilter);
59
61 using InputImageType = TInputImage;
62 using InputPixelType = typename InputImageType::PixelType;
63 using InputSizeType = typename InputImageType::SizeType;
64 using OutputImageType = TOutputImage;
65 using OutputPixelType = typename OutputImageType::PixelType;
66 using OutputSizeType = typename OutputImageType::SizeType;
67
72
79
80 using OutputImageRegionType = typename OutputImageType::RegionType;
81
83 itkNewMacro(Self);
84
86 itkOverrideGetNameOfClassMacro(FFTWInverseFFTImageFilter);
87
89 static constexpr unsigned int ImageDimension = InputImageType::ImageDimension;
90
101 virtual void
102 SetPlanRigor(const int & value)
103 {
104#ifndef ITK_USE_CUFFTW
105 // Use that method to check the value.
107#endif
108 if (m_PlanRigor != value)
109 {
110 m_PlanRigor = value;
111 this->Modified();
112 }
113 }
114 itkGetConstReferenceMacro(PlanRigor, int);
115 void
116 SetPlanRigor(const std::string & name)
117 {
118#ifndef ITK_USE_CUFFTW
120#endif
121 }
122
125
126protected:
128 ~FFTWInverseFFTImageFilter() override = default;
129
130 void
132
133 void
134 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
135
136 void
137 PrintSelf(std::ostream & os, Indent indent) const override;
138
139private:
141};
142
143
144// Describe whether input/output are real- or complex-valued
145// for factory registration
146template <>
148{
149 template <typename TUnderlying>
150 using InputPixelType = std::complex<TUnderlying>;
151 template <typename TUnderlying>
152 using OutputPixelType = TUnderlying;
153 using FilterDimensions = std::integer_sequence<unsigned int, 4, 3, 2, 1>;
154};
155
156// Disable the precision whose FFTW backend is absent (avoids instantiating an undefined fftw proxy).
157#if !defined(ITK_USE_FFTWF)
158template <>
160{};
161#endif
162#if !defined(ITK_USE_FFTWD)
163template <>
164struct FFTImageFilterEnableDouble<FFTWInverseFFTImageFilter> : std::false_type
165{};
166#endif
167
168} // namespace itk
169
170#ifndef ITK_MANUAL_INSTANTIATION
171# include "itkFFTWInverseFFTImageFilter.hxx"
172#endif
173
174#endif // itkFFTWInverseFFTImageFilter_h
static int GetPlanRigorValue(const std::string &name)
static std::string GetPlanRigorName(const int value)
FFTW-based inverse Fast Fourier Transform.
~FFTWInverseFFTImageFilter() override=default
void PrintSelf(std::ostream &os, Indent indent) const override
typename OutputImageType::RegionType OutputImageRegionType
static constexpr unsigned int ImageDimension
typename fftw::Proxy< OutputPixelType > FFTWProxyType
SizeValueType GetSizeGreatestPrimeFactor() const override
typename OutputImageType::SizeType OutputSizeType
void BeforeThreadedGenerateData() override
InverseFFTImageFilter< InputImageType, OutputImageType > Superclass
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
typename InputImageType::SizeType InputSizeType
virtual void SetPlanRigor(const int &value)
typename InputImageType::PixelType InputPixelType
void SetPlanRigor(const std::string &name)
typename OutputImageType::PixelType OutputPixelType
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:51
virtual void Modified() const
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86
std::integer_sequence< unsigned int, 4, 3, 2, 1 > FilterDimensions
Helper defining pixel traits for templated FFT image filters.