RTK  2.6.0
Reconstruction Toolkit
rtkFourDConjugateGradientConeBeamReconstructionFilter.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 rtkFourDConjugateGradientConeBeamReconstructionFilter_h
20 #define rtkFourDConjugateGradientConeBeamReconstructionFilter_h
21 
29 
30 #include <itkExtractImageFilter.h>
31 #include <itkSubtractImageFilter.h>
32 #include <itkMultiplyImageFilter.h>
33 #include <itkIterationReporter.h>
34 #ifdef RTK_USE_CUDA
36 #endif
37 
38 namespace rtk
39 {
40 
95 template <typename VolumeSeriesType, typename ProjectionStackType>
97  : public rtk::IterativeConeBeamReconstructionFilter<VolumeSeriesType, ProjectionStackType>
98 {
99 public:
100  ITK_DISALLOW_COPY_AND_MOVE(FourDConjugateGradientConeBeamReconstructionFilter);
101 
107 
109  using InputImageType = VolumeSeriesType;
110  using OutputImageType = VolumeSeriesType;
111  using VolumeType = ProjectionStackType;
112 
120 
121  using ForwardProjectionType = typename Superclass::ForwardProjectionType;
122  using BackProjectionType = typename Superclass::BackProjectionType;
123 
125  using CPUVolumeSeriesType =
127 #ifdef RTK_USE_CUDA
128  typedef typename std::conditional<std::is_same<VolumeSeriesType, CPUVolumeSeriesType>::value,
130  CudaConjugateGradientImageFilter<VolumeSeriesType>>::type
132 #else
134 #endif
135 
137  itkNewMacro(Self);
138 
140  itkOverrideGetNameOfClassMacro(FourDConjugateGradientConeBeamReconstructionFilter);
141 
143  itkGetConstObjectMacro(Geometry, ThreeDCircularProjectionGeometry);
144  itkSetConstObjectMacro(Geometry, ThreeDCircularProjectionGeometry);
146 
148  itkGetMacro(NumberOfIterations, unsigned int);
149  itkSetMacro(NumberOfIterations, unsigned int);
151 
153  itkGetMacro(CudaConjugateGradient, bool);
154  itkSetMacro(CudaConjugateGradient, bool);
156 
158  void
159  SetInputVolumeSeries(const VolumeSeriesType * VolumeSeries);
160  typename VolumeSeriesType::ConstPointer
161  GetInputVolumeSeries();
163 
165  void
166  SetInputProjectionStack(const ProjectionStackType * Projections);
167  typename ProjectionStackType::ConstPointer
168  GetInputProjectionStack();
170 
172  void
173  SetWeights(const itk::Array2D<float> _arg);
174 
176  virtual void
177  SetSignal(const std::vector<double> signal);
178 
180  itkSetMacro(DisableDisplacedDetectorFilter, bool);
181  itkGetMacro(DisableDisplacedDetectorFilter, bool);
183 
184 protected:
187 
189  void
190  VerifyPreconditions() const override;
191 
192  void
193  GenerateOutputInformation() override;
194 
195  void
196  GenerateInputRequestedRegion() override;
197 
198  void
199  GenerateData() override;
200 
203  void
204  VerifyInputInformation() const override
205  {}
206 
215 
217  std::vector<double> m_Signal;
219 
220  // Iteration reporting
222 
223 private:
226 
228  unsigned int m_NumberOfIterations;
229 
231  void
232  ReportProgress(itk::Object *, const itk::EventObject &);
233 
234 }; // end of class
235 
236 } // end namespace rtk
237 
238 #ifndef ITK_MANUAL_INSTANTIATION
239 # include "rtkFourDConjugateGradientConeBeamReconstructionFilter.hxx"
240 #endif
241 
242 #endif
Base class for forward projection, i.e. accumulation along x-ray lines.
Weigting for displaced detectors.
typename itk::Image< typename VolumeSeriesType::PixelType, VolumeSeriesType::ImageDimension > CPUVolumeSeriesType
Projection geometry for a source and a 2-D flat panel.
#define itkSetMacro(name, type)
Implements part of the 4D reconstruction by conjugate gradient.
TOutputImage OutputImageType
Mother class for cone beam reconstruction filters which need runtime selection of their forward and b...
Solves AX = B by conjugate gradient.
Implements part of the 4D reconstruction by conjugate gradient.