RTK  2.5.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 #undef itkSetMacro
31 #define itkSetMacro(name, type) \
32  virtual void Set##name(type _arg) \
33  { \
34  itkDebugMacro("setting " #name " to " << _arg); \
35  CLANG_PRAGMA_PUSH \
36  CLANG_SUPPRESS_Wfloat_equal \
37  if (this->m_##name != _arg) \
38  { \
39  this->m_##name = std::move(_arg); \
40  this->Modified(); \
41  this->m_KernelFFT = nullptr; \
42  } \
43  CLANG_PRAGMA_POP \
44  } \
45  ITK_MACROEND_NOOP_STATEMENT
46 // clang-format on
47 
48 namespace rtk
49 {
50 
64 template <class TInputImage, class TOutputImage = TInputImage, class TFFTPrecision = double>
65 class ITK_TEMPLATE_EXPORT FFTRampImageFilter
66  : public rtk::FFTProjectionsConvolutionImageFilter<TInputImage, TOutputImage, TFFTPrecision>
67 {
68 public:
69  ITK_DISALLOW_COPY_AND_MOVE(FFTRampImageFilter);
70 
76 
78  using InputImageType = TInputImage;
79  using OutputImageType = TOutputImage;
80  using FFTPrecisionType = TFFTPrecision;
81  using IndexType = typename InputImageType::IndexType;
82  using SizeType = typename InputImageType::SizeType;
83 
84  using FFTInputImageType = typename Superclass::FFTInputImageType;
85  using FFTInputImagePointer = typename FFTInputImageType::Pointer;
86  using FFTOutputImageType = typename Superclass::FFTOutputImageType;
87  using FFTOutputImagePointer = typename FFTOutputImageType::Pointer;
88 
90  itkNewMacro(Self);
91 
93 #ifdef itkOverrideGetNameOfClassMacro
94  itkOverrideGetNameOfClassMacro(FFTRampImageFilter);
95 #else
97 #endif
98 
99 
101  itkGetConstMacro(HannCutFrequency, double);
102  itkSetMacro(HannCutFrequency, double);
104 
106  itkGetConstMacro(CosineCutFrequency, double);
107  itkSetMacro(CosineCutFrequency, double);
109 
111  itkGetConstMacro(HammingFrequency, double);
112  itkSetMacro(HammingFrequency, double);
114 
116  itkGetConstMacro(HannCutFrequencyY, double);
117  itkSetMacro(HannCutFrequencyY, double);
119 
127  itkGetConstMacro(RamLakCutFrequency, double);
128  itkSetMacro(RamLakCutFrequency, double);
130 
138  itkGetConstMacro(SheppLoganCutFrequency, double);
139  itkSetMacro(SheppLoganCutFrequency, double);
141 
142 protected:
144  ~FFTRampImageFilter() override = default;
145 
146  void
147  GenerateInputRequestedRegion() override;
148 
151  void
152  UpdateFFTProjectionsConvolutionKernel(const SizeType s) override;
153 
154 private:
159  double m_HannCutFrequency{ 0. };
160  double m_CosineCutFrequency{ 0. };
161  double m_HammingFrequency{ 0. };
162  double m_HannCutFrequencyY{ 0. };
163 
166  double m_RamLakCutFrequency{ 0. };
167  double m_SheppLoganCutFrequency{ 0. };
168 
170 }; // end of class
171 
172 } // end namespace rtk
173 
174 #ifndef ITK_MANUAL_INSTANTIATION
175 # include "rtkFFTRampImageFilter.hxx"
176 #endif
177 
178 // Rollback to the original definition of the Set macro
179 // clang-format off
180 #undef itkSetMacro
181 #define itkSetMacro(name, type) \
182  virtual void Set##name(type _arg) \
183  { \
184  itkDebugMacro("setting " #name " to " << _arg); \
185  CLANG_PRAGMA_PUSH \
186  CLANG_SUPPRESS_Wfloat_equal \
187  if (this->m_##name != _arg) \
188  { \
189  this->m_##name = std::move(_arg); \
190  this->Modified(); \
191  } \
192  CLANG_PRAGMA_POP \
193  } \
194  ITK_MACROEND_NOOP_STATEMENT
195 // clang-format on
196 #endif
Base class for 1D or 2D FFT based convolution of projections.
TInputImage InputImageType
#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