ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkDirectFourierReconstructionImageToImageFilter.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 itkDirectFourierReconstructionImageToImageFilter_h
19#define itkDirectFourierReconstructionImageToImageFilter_h
20
22#include "itkImage.h"
23
26
29
31
32#include <cmath>
33
34namespace itk
35{
50template <typename TInputImage, typename TOutputImage = TInputImage>
52 : public ImageToImageFilter<TInputImage, TOutputImage>
53{
54public:
55 ITK_DISALLOW_COPY_AND_MOVE(DirectFourierReconstructionImageToImageFilter);
56
59
60 using InputImageType = TInputImage;
61 using InputPixelType = typename InputImageType::PixelType;
62 using OutputImageType = TOutputImage;
63 using OutputPixelType = typename OutputImageType::PixelType;
64
67
72
73 itkNewMacro(Self);
74 itkOverrideGetNameOfClassMacro(DirectFourierReconstructionImageToImageFilter);
75
77 using RegionType = typename InputImageType::RegionType;
79 using IndexType = typename InputImageType::IndexType;
81 using SizeType = typename InputImageType::SizeType;
83 using PointType = typename InputImageType::PointType;
85 using SpacingType = typename InputImageType::SpacingType;
86
88 using ConstInputImagePointer = typename InputImageType::ConstPointer;
90 using InputImagePointer = typename InputImageType::Pointer;
92 using OutputImagePointer = typename OutputImageType::Pointer;
93
94 itkSetMacro(ZeroPadding, unsigned short);
95 itkGetConstMacro(ZeroPadding, unsigned short);
96
97 itkSetMacro(OverSampling, unsigned short);
98 itkGetConstMacro(OverSampling, unsigned short);
99
100 itkSetMacro(Cutoff, double);
101 itkGetConstMacro(Cutoff, double);
102
103 itkSetMacro(AlphaRange, double);
104 itkGetConstMacro(AlphaRange, double);
105
106 itkSetMacro(AlphaDirection, unsigned short);
107 itkGetConstMacro(AlphaDirection, unsigned short);
108
109 itkSetMacro(ZDirection, unsigned short);
110 itkGetConstMacro(ZDirection, unsigned short);
111
112 itkSetMacro(RDirection, unsigned short);
113 itkGetConstMacro(RDirection, unsigned short);
114
115 itkSetMacro(RadialSplineOrder, unsigned short);
116 itkGetConstMacro(RadialSplineOrder, unsigned short);
117
118protected:
123
125 void
126 PrintSelf(std::ostream & os, Indent indent) const override;
127
129 void
131
133 void
135
137 void
138 GenerateData() override;
139
140private:
143
155
167
168 unsigned short m_ZeroPadding{};
169 unsigned short m_OverSampling{};
170
171 double m_Cutoff{};
173 double m_AlphaRange{};
174
175 unsigned short m_ZDirection{};
176 unsigned short m_AlphaDirection{};
178 unsigned short m_RDirection{};
180 unsigned short m_RadialSplineOrder{};
182
183 double m_PI{};
184
187};
188} // namespace itk
189
190#ifndef ITK_MANUAL_INSTANTIATION
191# include "itkDirectFourierReconstructionImageToImageFilter.hxx"
192#endif
193
194#endif /* itkDirectFourierReconstructionImageToImageFilter_h */
Complex wrapper around BSplineInterpolateImageFunction.
void PrintSelf(std::ostream &os, Indent indent) const override
ComplexBSplineInterpolateImageFunction< FFTLineType, double, double > FFTLineInterpolatorType
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
Multi-dimensional image iterator which only walks a region.
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
Implements transparent reference counting.
VNL based forward Fast Fourier Transform.
Image< std::complex< typename LineImageType::PixelType >, LineImageType::ImageDimension > OutputImageType
VNL-based reverse Fast Fourier Transform.
Image< typename IFFTImageType::PixelType::value_type, IFFTImageType::ImageDimension > OutputImageType
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....