RTK  2.6.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  itkOverrideGetNameOfClassMacro(ConjugateGradientConeBeamReconstructionFilter);
126 
128  void
129  SetInputVolume(const TOutputImage * vol);
130  void
131  SetInputProjectionStack(const TOutputImage * projs);
132  void
133  SetInputWeights(const TWeightsImage * weights);
134  void
135  SetLocalRegularizationWeights(const TSingleComponentImage * weights);
137 
144  using CGOperatorFilterType =
148  using OutputImagePointer = typename TOutputImage::Pointer;
150 
151  using ForwardProjectionType = typename Superclass::ForwardProjectionType;
152  using BackProjectionType = typename Superclass::BackProjectionType;
153 
154  // If TOutputImage is an itk::Image of floats or double, so are the weights, and a simple Multiply filter is required
155  // If TOutputImage is an itk::Image of itk::Vector<float (or double)>, a BlockDiagonalMatrixVectorMultiply filter
156  // is needed. Thus the meta-programming construct
159  typedef typename std::conditional<std::is_same<TSingleComponentImage, TOutputImage>::value,
163 #ifdef RTK_USE_CUDA
164  typedef typename std::conditional<!std::is_same<TOutputImage, CPUOutputImageType>::value &&
165  std::is_same<TSingleComponentImage, TOutputImage>::value,
166  CudaDisplacedDetectorImageFilter,
168  typedef typename std::conditional<!std::is_same<TOutputImage, CPUOutputImageType>::value &&
169  std::is_same<TSingleComponentImage, TOutputImage>::value,
170  CudaConstantVolumeSource,
172  typedef typename std::conditional<!std::is_same<TOutputImage, CPUOutputImageType>::value &&
173  std::is_same<TSingleComponentImage, TOutputImage>::value,
174  CudaConstantVolumeSource,
176 #else
180 #endif
181 
183  void
184  SetSupportMask(const TSingleComponentImage * SupportMask);
185  typename TSingleComponentImage::ConstPointer
186  GetSupportMask();
188 
190  itkSetConstObjectMacro(Geometry, ThreeDCircularProjectionGeometry);
191 
192  itkSetMacro(NumberOfIterations, int);
193  itkGetMacro(NumberOfIterations, int);
194 
196  itkSetMacro(DisableDisplacedDetectorFilter, bool);
197  itkGetMacro(DisableDisplacedDetectorFilter, bool);
199 
202  itkSetMacro(Tikhonov, float);
203  itkGetMacro(Tikhonov, float);
204  itkSetMacro(Gamma, float);
205  itkGetMacro(Gamma, float);
207 
209  itkGetMacro(CudaConjugateGradient, bool);
210  itkSetMacro(CudaConjugateGradient, bool);
212 
213 protected:
215  ~ConjugateGradientConeBeamReconstructionFilter() override = default;
216 
218  void
219  VerifyPreconditions() const override;
220 
222  void
223  GenerateData() override;
224 
236  typename MultiplyWithWeightsFilterType::Pointer m_MultiplyWithWeightsFilter;
237 
241  void
242  VerifyInputInformation() const override
243  {}
244 
247  void
248  GenerateInputRequestedRegion() override;
249  void
250  GenerateOutputInformation() override;
252 
254  typename TOutputImage::ConstPointer
255  GetInputVolume();
256  typename TOutputImage::ConstPointer
257  GetInputProjectionStack();
258  typename TWeightsImage::ConstPointer
259  GetInputWeights();
260  typename TSingleComponentImage::ConstPointer
261  GetLocalRegularizationWeights();
263 
264  template <typename ImageType,
265  typename IterativeConeBeamReconstructionFilter<TOutputImage>::template EnableCudaScalarAndVectorType<
266  ImageType> * = nullptr>
267  ConjugateGradientFilterPointer
268  InstantiateCudaConjugateGradientImageFilter();
269 
270  template <typename ImageType,
271  typename IterativeConeBeamReconstructionFilter<TOutputImage>::template DisableCudaScalarAndVectorType<
272  ImageType> * = nullptr>
273  ConjugateGradientFilterPointer
274  InstantiateCudaConjugateGradientImageFilter();
275 
277  void
278  ReportProgress(itk::Object *, const itk::EventObject &);
279 
280 private:
282 
284  float m_Gamma;
285  float m_Tikhonov;
288 
289  // Iteration reporting
291 };
292 } // namespace rtk
293 
294 
295 #ifndef ITK_MANUAL_INSTANTIATION
296 # include "rtkConjugateGradientConeBeamReconstructionFilter.hxx"
297 #endif
298 
299 #endif
typename OutputImageType::Pointer OutputImagePointer
Implements the operator A used in conjugate gradient reconstruction.
BackProjectionImageFilter< TOutputImage, TOutputImage >::Pointer m_BackProjectionFilter
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
BackProjectionImageFilter< TOutputImage, TOutputImage >::Pointer m_BackProjectionFilterForB
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.
ForwardProjectionImageFilter< TOutputImage, TOutputImage >::Pointer m_ForwardProjectionFilter