RTK  2.6.0
Reconstruction Toolkit
rtkFFTRampImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright RTK Consortium
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 
19 #ifndef rtkFFTRampImageFilter_h
20 #define rtkFFTRampImageFilter_h
21 
22 #include <itkConceptChecking.h>
23 #include "rtkConfiguration.h"
25 #include "rtkMacro.h"
26 
27 // The Set macro is redefined to clear the current FFT kernel when a parameter
28 // is modified.
29 // clang-format off
30 #ifndef ITK_GCC_PRAGMA_PUSH
31 #define ITK_GCC_PRAGMA_PUSH CLANG_PRAGMA_PUSH
32 #endif
33 #ifndef ITK_GCC_SUPPRESS_Wfloat_equal
34 #define ITK_GCC_SUPPRESS_Wfloat_equal CLANG_SUPPRESS_Wfloat_equal
35 #endif
36 #ifndef ITK_GCC_PRAGMA_POP
37 #define ITK_GCC_PRAGMA_POP CLANG_PRAGMA_POP
38 #endif
39 #undef itkSetMacro
40 #define itkSetMacro(name, type) \
41  virtual void Set##name(type _arg) \
42  { \
43  itkDebugMacro("setting " #name " to " << _arg); \
44  ITK_GCC_PRAGMA_PUSH \
45  ITK_GCC_SUPPRESS_Wfloat_equal \
46  if (this->m_##name != _arg) \
47  { \
48  this->m_##name = std::move(_arg); \
49  this->Modified(); \
50  this->m_KernelFFT = nullptr; \
51  } \
52  ITK_GCC_PRAGMA_POP \
53  } \
54  ITK_MACROEND_NOOP_STATEMENT
55 // clang-format on
56 
57 namespace rtk
58 {
59 
73 template <class TInputImage, class TOutputImage = TInputImage, class TFFTPrecision = double>
74 class ITK_TEMPLATE_EXPORT FFTRampImageFilter
75  : public rtk::FFTProjectionsConvolutionImageFilter<TInputImage, TOutputImage, TFFTPrecision>
76 {
77 public:
78  ITK_DISALLOW_COPY_AND_MOVE(FFTRampImageFilter);
79 
85 
87  using InputImageType = TInputImage;
88  using OutputImageType = TOutputImage;
89  using FFTPrecisionType = TFFTPrecision;
90  using IndexType = typename InputImageType::IndexType;
91  using SizeType = typename InputImageType::SizeType;
92 
93  using FFTInputImageType = typename Superclass::FFTInputImageType;
94  using FFTInputImagePointer = typename FFTInputImageType::Pointer;
95  using FFTOutputImageType = typename Superclass::FFTOutputImageType;
96  using FFTOutputImagePointer = typename FFTOutputImageType::Pointer;
97 
99  itkNewMacro(Self);
100 
102  itkOverrideGetNameOfClassMacro(FFTRampImageFilter);
103 
105  itkGetConstMacro(HannCutFrequency, double);
106  itkSetMacro(HannCutFrequency, double);
108 
110  itkGetConstMacro(CosineCutFrequency, double);
111  itkSetMacro(CosineCutFrequency, double);
113 
115  itkGetConstMacro(HammingFrequency, double);
116  itkSetMacro(HammingFrequency, double);
118 
120  itkGetConstMacro(HannCutFrequencyY, double);
121  itkSetMacro(HannCutFrequencyY, double);
123 
131  itkGetConstMacro(RamLakCutFrequency, double);
132  itkSetMacro(RamLakCutFrequency, double);
134 
142  itkGetConstMacro(SheppLoganCutFrequency, double);
143  itkSetMacro(SheppLoganCutFrequency, double);
145 
146 protected:
148  ~FFTRampImageFilter() override = default;
149 
150  void
151  GenerateInputRequestedRegion() override;
152 
155  void
156  UpdateFFTProjectionsConvolutionKernel(const SizeType s) override;
157 
158 private:
163  double m_HannCutFrequency{ 0. };
164  double m_CosineCutFrequency{ 0. };
165  double m_HammingFrequency{ 0. };
166  double m_HannCutFrequencyY{ 0. };
167 
170  double m_RamLakCutFrequency{ 0. };
171  double m_SheppLoganCutFrequency{ 0. };
172 
174 }; // end of class
175 
176 } // end namespace rtk
177 
178 #ifndef ITK_MANUAL_INSTANTIATION
179 # include "rtkFFTRampImageFilter.hxx"
180 #endif
181 
182 // Rollback to the original definition of the Set macro
183 // clang-format off
184 #undef itkSetMacro
185 #define itkSetMacro(name, type) \
186  virtual void Set##name(type _arg) \
187  { \
188  itkDebugMacro("setting " #name " to " << _arg); \
189  ITK_GCC_PRAGMA_PUSH \
190  ITK_GCC_SUPPRESS_Wfloat_equal \
191  if (this->m_##name != _arg) \
192  { \
193  this->m_##name = std::move(_arg); \
194  this->Modified(); \
195  } \
196  ITK_GCC_PRAGMA_POP \
197  } \
198  ITK_MACROEND_NOOP_STATEMENT
199 // clang-format on
200 #endif
Base class for 1D or 2D FFT based convolution of projections.
#define itkSetMacro(name, type)
TOutputImage OutputImageType
typename itk::Image< TFFTPrecision, TInputImage::ImageDimension > FFTInputImageType
Implements the ramp image filter of the filtered backprojection algorithm.
typename itk::Image< std::complex< TFFTPrecision >, TInputImage::ImageDimension > FFTOutputImageType