ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkDirectFourierReconstructionImageToImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software 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  * http://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 
21 #include "itkImageToImageFilter.h"
22 #include "itkImage.h"
23 
26 
29 
31 
32 #include <cmath>
33 
34 namespace itk
35 {
50 template< typename TInputImage, typename TOutputImage=TInputImage >
52  public ImageToImageFilter< TInputImage, TOutputImage >
53 {
54 public:
55 
58 
59  typedef TInputImage InputImageType;
60  typedef typename InputImageType::PixelType InputPixelType;
61  typedef TOutputImage OutputImageType;
62  typedef typename OutputImageType::PixelType OutputPixelType;
63 
66 
69 
72 
73  itkNewMacro(Self);
75 
77  typedef typename InputImageType::RegionType RegionType;
78 
80  typedef typename InputImageType::IndexType IndexType;
81 
83  typedef typename InputImageType::SizeType SizeType;
84 
86  typedef typename InputImageType::PointType PointType;
87 
89  typedef typename InputImageType::SpacingType SpacingType;
90 
92  typedef typename InputImageType::ConstPointer ConstInputImagePointer;
93 
95  typedef typename InputImageType::Pointer InputImagePointer;
96 
98  typedef typename OutputImageType::Pointer OutputImagePointer;
99 
100  itkSetMacro(ZeroPadding, unsigned short int);
101  itkGetConstMacro(ZeroPadding, unsigned short int);
102 
103  itkSetMacro(OverSampling, unsigned short int);
104  itkGetConstMacro(OverSampling, unsigned short int);
105 
106  itkSetMacro(Cutoff, double);
107  itkGetConstMacro(Cutoff, double);
108 
109  itkSetMacro(AlphaRange, double);
110  itkGetConstMacro(AlphaRange, double);
111 
112  itkSetMacro(AlphaDirection, unsigned short int);
113  itkGetConstMacro(AlphaDirection, unsigned short int);
114 
115  itkSetMacro(ZDirection, unsigned short int);
116  itkGetConstMacro(ZDirection, unsigned short int);
117 
118  itkSetMacro(RDirection, unsigned short int);
119  itkGetConstMacro(RDirection, unsigned short int);
120 
121  itkSetMacro(RadialSplineOrder, unsigned short int);
122  itkGetConstMacro(RadialSplineOrder, unsigned short int);
123 
124 protected:
127 
130 
132  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
133 
135  void GenerateOutputInformation() ITK_OVERRIDE;
136 
138  void GenerateInputRequestedRegion() ITK_OVERRIDE;
139 
141  void GenerateData() ITK_OVERRIDE;
142 
143 private:
146 
148  typedef Image< double, 1 > LineImageType;
149  typedef VnlForwardFFTImageFilter< LineImageType > FFTLineFilterType;
150 
152  typedef FFTLineFilterType::OutputImageType FFTLineType;
153 
155  typedef FFTLineFilterType::InputImageType ProjectionLineType;
156 
159 
161  typedef ComplexBSplineInterpolateImageFunction< FFTLineType, double, double > FFTLineInterpolatorType;
162 
164  typedef Image< std::complex<double>, 2> IFFTImageType;
166 
168  typedef IFFTSliceFilterType::InputImageType FFTSliceType;
169 
171  typedef IFFTSliceFilterType::OutputImageType OutputSliceType;
172 
175 
178 
204  DirectFourierReconstructionImageToImageFilter(const Self &) ITK_DELETE_FUNCTION;
205  void operator=(const Self &) ITK_DELETE_FUNCTION;
206 };
207 } // namespace itk
208 
209 #ifndef ITK_MANUAL_INSTANTIATION
210 #include "itkDirectFourierReconstructionImageToImageFilter.hxx"
211 #endif
212 
213 #endif /* itkDirectFourierReconstructionImageToImageFilter_h */
void PrintSelf(std::ostream &os, Indent indent) const override
Base class for all process objects that output image data.
Complex wrapper around BSplineInterpolateImageFunction.
Direct fourier reconstruction filter of a tomographic volume.
VNL based forward Fast Fourier Transform.
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
VNL-based reverse Fast Fourier Transform.
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
TOutputImage OutputImageType
Multi-dimensional image iterator which only walks a region.
Templated n-dimensional image class.
Definition: itkImage.h:75