RTK  2.5.0
Reconstruction Toolkit
rtkConjugateGradientConeBeamReconstructionFilter.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 rtkConjugateGradientConeBeamReconstructionFilter_h
20 #define rtkConjugateGradientConeBeamReconstructionFilter_h
21 
22 #include <itkMultiplyImageFilter.h>
24 #include <itkProcessObject.h>
25 #include <itkObject.h>
26 #include <itkCommand.h>
27 #include <itkIterationReporter.h>
28 
34 #include "rtkConstantImageSource.h"
37 
38 #ifdef RTK_USE_CUDA
42 #endif
43 
44 namespace rtk
45 {
109 template <typename TOutputImage, typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
111  : public IterativeConeBeamReconstructionFilter<TOutputImage>
112 {
113 public:
114  ITK_DISALLOW_COPY_AND_MOVE(ConjugateGradientConeBeamReconstructionFilter);
115 
120 
122  itkNewMacro(Self);
123 
125 #ifdef itkOverrideGetNameOfClassMacro
126  itkOverrideGetNameOfClassMacro(ConjugateGradientConeBeamReconstructionFilter);
127 #else
129 #endif
130 
131 
133  void
134  SetInputVolume(const TOutputImage * vol);
135  void
136  SetInputProjectionStack(const TOutputImage * projs);
137  void
138  SetInputWeights(const TWeightsImage * weights);
139  void
140  SetLocalRegularizationWeights(const TSingleComponentImage * weights);
142 
149  using CGOperatorFilterType =
153  using OutputImagePointer = typename TOutputImage::Pointer;
155 
156  using ForwardProjectionType = typename Superclass::ForwardProjectionType;
157  using BackProjectionType = typename Superclass::BackProjectionType;
158 
159  // If TOutputImage is an itk::Image of floats or double, so are the weights, and a simple Multiply filter is required
160  // If TOutputImage is an itk::Image of itk::Vector<float (or double)>, a BlockDiagonalMatrixVectorMultiply filter
161  // is needed. Thus the meta-programming construct
164  typedef typename std::conditional<std::is_same<TSingleComponentImage, TOutputImage>::value,
168 #ifdef RTK_USE_CUDA
169  typedef typename std::conditional<!std::is_same<TOutputImage, CPUOutputImageType>::value &&
170  std::is_same<TSingleComponentImage, TOutputImage>::value,
171  CudaDisplacedDetectorImageFilter,
173  typedef typename std::conditional<!std::is_same<TOutputImage, CPUOutputImageType>::value &&
174  std::is_same<TSingleComponentImage, TOutputImage>::value,
175  CudaConstantVolumeSource,
177  typedef typename std::conditional<!std::is_same<TOutputImage, CPUOutputImageType>::value &&
178  std::is_same<TSingleComponentImage, TOutputImage>::value,
179  CudaConstantVolumeSource,
181 #else
185 #endif
186 
188  void
189  SetSupportMask(const TSingleComponentImage * SupportMask);
190  typename TSingleComponentImage::ConstPointer
191  GetSupportMask();
193 
195  itkSetConstObjectMacro(Geometry, ThreeDCircularProjectionGeometry);
196 
197  itkSetMacro(NumberOfIterations, int);
198  itkGetMacro(NumberOfIterations, int);
199 
201  itkSetMacro(DisableDisplacedDetectorFilter, bool);
202  itkGetMacro(DisableDisplacedDetectorFilter, bool);
204 
207  itkSetMacro(Tikhonov, float);
208  itkGetMacro(Tikhonov, float);
209  itkSetMacro(Gamma, float);
210  itkGetMacro(Gamma, float);
212 
214  itkGetMacro(CudaConjugateGradient, bool);
215  itkSetMacro(CudaConjugateGradient, bool);
217 
218 protected:
220  ~ConjugateGradientConeBeamReconstructionFilter() override = default;
221 
223  void
224  VerifyPreconditions() ITKv5_CONST override;
225 
227  void
228  GenerateData() override;
229 
231  typename MultiplyFilterType::Pointer m_MultiplyProjectionsFilter;
232  typename MultiplyFilterType::Pointer m_MultiplyVolumeFilter;
233  typename MultiplyFilterType::Pointer m_MultiplyOutputFilter;
234  ConjugateGradientFilterPointer m_ConjugateGradientFilter;
235  typename CGOperatorFilterType::Pointer m_CGOperator;
236  typename ForwardProjectionImageFilter<TOutputImage, TOutputImage>::Pointer m_ForwardProjectionFilter;
237  typename BackProjectionImageFilter<TOutputImage, TOutputImage>::Pointer m_BackProjectionFilter;
238  typename BackProjectionImageFilter<TOutputImage, TOutputImage>::Pointer m_BackProjectionFilterForB;
239  typename DisplacedDetectorFilterType::Pointer m_DisplacedDetectorFilter;
240  typename ConstantImageSourceType::Pointer m_ConstantVolumeSource;
241  typename MultiplyWithWeightsFilterType::Pointer m_MultiplyWithWeightsFilter;
242 
246  void
247  VerifyInputInformation() const override
248  {}
249 
252  void
253  GenerateInputRequestedRegion() override;
254  void
255  GenerateOutputInformation() override;
257 
259  typename TOutputImage::ConstPointer
260  GetInputVolume();
261  typename TOutputImage::ConstPointer
262  GetInputProjectionStack();
263  typename TWeightsImage::ConstPointer
264  GetInputWeights();
265  typename TSingleComponentImage::ConstPointer
266  GetLocalRegularizationWeights();
268 
269  template <typename ImageType,
270  typename IterativeConeBeamReconstructionFilter<TOutputImage>::template EnableCudaScalarAndVectorType<
271  ImageType> * = nullptr>
272  ConjugateGradientFilterPointer
273  InstantiateCudaConjugateGradientImageFilter();
274 
275  template <typename ImageType,
276  typename IterativeConeBeamReconstructionFilter<TOutputImage>::template DisableCudaScalarAndVectorType<
277  ImageType> * = nullptr>
278  ConjugateGradientFilterPointer
279  InstantiateCudaConjugateGradientImageFilter();
280 
282  void
283  ReportProgress(itk::Object *, const itk::EventObject &);
284 
285 private:
287 
289  float m_Gamma;
290  float m_Tikhonov;
293 
294  // Iteration reporting
296 };
297 } // namespace rtk
298 
299 
300 #ifndef ITK_MANUAL_INSTANTIATION
301 # include "rtkConjugateGradientConeBeamReconstructionFilter.hxx"
302 #endif
303 
304 #endif
typename OutputImageType::Pointer OutputImagePointer
Implements the operator A used in conjugate gradient reconstruction.
typename itk::Image< typename TOutputImage::PixelType, TOutputImage::ImageDimension > CPUOutputImageType
Generate an n-dimensional image with constant pixel values.
std::conditional< std::is_same< TSingleComponentImage, TOutputImage >::value, PlainMultiplyFilterType, MatrixVectorMultiplyFilterType >::type MultiplyWithWeightsFilterType
Weigting for displaced detectors.
Projection geometry for a source and a 2-D flat panel.
#define itkSetMacro(name, type)
Mother class for cone beam reconstruction filters which need runtime selection of their forward and b...
Solves AX = B by conjugate gradient.