RTK  2.6.0
Reconstruction Toolkit
rtkFourDReconstructionConjugateGradientOperator.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 #ifndef rtkFourDReconstructionConjugateGradientOperator_h
19 #define rtkFourDReconstructionConjugateGradientOperator_h
20 
22 
23 #include <itkArray2D.h>
24 #include <itkMultiplyImageFilter.h>
25 
26 #include "rtkConstantImageSource.h"
33 
34 #ifdef RTK_USE_CUDA
36 # include "rtkCudaSplatImageFilter.h"
40 #endif
41 
42 namespace rtk
43 {
129 template <typename VolumeSeriesType, typename ProjectionStackType>
131  : public ConjugateGradientOperator<VolumeSeriesType>
132 {
133 public:
134  ITK_DISALLOW_COPY_AND_MOVE(FourDReconstructionConjugateGradientOperator);
135 
140 
142  using VolumeType = ProjectionStackType;
143 
145  itkNewMacro(Self);
146 
148  itkOverrideGetNameOfClassMacro(FourDReconstructionConjugateGradientOperator);
149 
151  void
152  SetInputVolumeSeries(const VolumeSeriesType * VolumeSeries);
153  typename VolumeSeriesType::ConstPointer
154  GetInputVolumeSeries();
156 
158  void
159  SetInputProjectionStack(const ProjectionStackType * Projections);
160  typename ProjectionStackType::ConstPointer
161  GetInputProjectionStack();
163 
171 
173  using CPUProjectionStackType =
175 #ifdef RTK_USE_CUDA
176  typedef typename std::conditional<std::is_same<ProjectionStackType, CPUProjectionStackType>::value,
178  CudaDisplacedDetectorImageFilter>::type DisplacedDetectorFilterType;
179  typedef typename std::conditional<std::is_same<ProjectionStackType, CPUProjectionStackType>::value,
181  CudaInterpolateImageFilter>::type CudaInterpolateImageFilterType;
182  typedef typename std::conditional<std::is_same<ProjectionStackType, CPUProjectionStackType>::value,
184  CudaSplatImageFilter>::type CudaSplatImageFilterType;
185  typedef typename std::conditional<std::is_same<ProjectionStackType, CPUProjectionStackType>::value,
187  CudaConstantVolumeSource>::type CudaConstantVolumeSourceType;
188  typedef typename std::conditional<std::is_same<ProjectionStackType, CPUProjectionStackType>::value,
190  CudaConstantVolumeSeriesSource>::type CudaConstantVolumeSeriesSourceType;
191 #else
197 #endif
198 
200  void
201  SetBackProjectionFilter(BackProjectionFilterType * _arg);
202 
204  void
205  SetForwardProjectionFilter(ForwardProjectionFilterType * _arg);
206 
208  itkSetConstObjectMacro(Geometry, ThreeDCircularProjectionGeometry);
209 
211  itkSetMacro(UseCudaInterpolation, bool);
212  itkGetMacro(UseCudaInterpolation, bool);
213  itkSetMacro(UseCudaSplat, bool);
214  itkGetMacro(UseCudaSplat, bool);
215  itkSetMacro(UseCudaSources, bool);
216  itkGetMacro(UseCudaSources, bool);
218 
220  itkGetMacro(Weights, itk::Array2D<float>);
223 
225  virtual void
226  SetSignal(const std::vector<double> signal);
227 
229  itkSetMacro(DisableDisplacedDetectorFilter, bool);
230  itkGetMacro(DisableDisplacedDetectorFilter, bool);
232 
233 protected:
235  ~FourDReconstructionConjugateGradientOperator() override = default;
236 
238  void
239  VerifyPreconditions() const override;
240 
242  void
243  GenerateOutputInformation() override;
244 
246  void
247  GenerateInputRequestedRegion() override;
248 
250  void
251  GenerateData() override;
252 
254  void
255  InitializeConstantSources();
256 
267 
273  std::vector<double> m_Signal;
275 };
276 } // namespace rtk
277 
278 
279 #ifndef ITK_MANUAL_INSTANTIATION
280 # include "rtkFourDReconstructionConjugateGradientOperator.hxx"
281 #endif
282 
283 #endif
Base class for forward projection, i.e. accumulation along x-ray lines.
typename itk::Image< typename ProjectionStackType::PixelType, ProjectionStackType::ImageDimension > CPUProjectionStackType
Generate an n-dimensional image with constant pixel values.
Weigting for displaced detectors.
Splats (linearly) a 3D volume into a 3D+t sequence of volumes.
Projection geometry for a source and a 2-D flat panel.
#define itkSetMacro(name, type)
Implements part of the 4D reconstruction by conjugate gradient.
Interpolates (linearly) in a 3D+t sequence of volumes to get a 3D volume.