RTK  2.6.0
Reconstruction Toolkit
rtkMechlemOneStepSpectralReconstructionFilter.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 rtkMechlemOneStepSpectralReconstructionFilter_h
20 #define rtkMechlemOneStepSpectralReconstructionFilter_h
21 
26 #include "rtkConstantImageSource.h"
31 
32 #include <itkExtractImageFilter.h>
33 #include <itkAddImageFilter.h>
34 #include <itkMultiplyImageFilter.h>
35 
36 #ifdef RTK_USE_CUDA
38 #endif
39 
40 namespace rtk
41 {
152 template <typename TOutputImage, typename TPhotonCounts, typename TSpectrum>
154  : public rtk::IterativeConeBeamReconstructionFilter<TOutputImage, TOutputImage>
155 {
156 public:
157  ITK_DISALLOW_COPY_AND_MOVE(MechlemOneStepSpectralReconstructionFilter);
158 
163 
165  itkNewMacro(Self);
166 
168 #ifdef itkOverrideGetNameOfClassMacro
169  itkOverrideGetNameOfClassMacro(MechlemOneStepSpectralReconstructionFilter);
170 #else
172 #endif
173 
174 
176  static constexpr unsigned int nBins = TPhotonCounts::PixelType::Dimension;
177  static constexpr unsigned int nMaterials = TOutputImage::PixelType::Dimension;
178  using dataType = typename TOutputImage::PixelType::ValueType;
179 
182 #ifdef RTK_USE_CUDA
183  typedef
184  typename std::conditional<std::is_same<TOutputImage, CPUOutputImageType>::value,
185  itk::Image<itk::Vector<dataType, nMaterials * nMaterials>, TOutputImage::ImageDimension>,
186  itk::CudaImage<itk::Vector<dataType, nMaterials * nMaterials>,
187  TOutputImage::ImageDimension>>::type HessiansImageType;
188  typedef
189  typename std::conditional<std::is_same<TOutputImage, CPUOutputImageType>::value,
191  itk::CudaImage<dataType, TOutputImage::ImageDimension>>::type SingleComponentImageType;
192 #else
193  using HessiansImageType =
194  typename itk::Image<itk::Vector<dataType, nMaterials * nMaterials>, TOutputImage::ImageDimension>;
196 #endif
197 
198 #if !defined(ITK_WRAPPING_PARSER)
199 # ifdef RTK_USE_CUDA
200  typedef typename std::conditional<std::is_same<TOutputImage, CPUOutputImageType>::value,
202  CudaWeidingerForwardModelImageFilter<TOutputImage, TPhotonCounts, TSpectrum>>::type
204  typedef
205  typename std::conditional<std::is_same<TOutputImage, CPUOutputImageType>::value,
207  CudaForwardProjectionImageFilter<SingleComponentImageType, SingleComponentImageType>>::
209  typedef typename std::conditional<std::is_same<TOutputImage, CPUOutputImageType>::value,
211  CudaBackProjectionImageFilter<HessiansImageType>>::type
213 # else
218 # endif
219  using GradientsImageType = TOutputImage;
220 #endif
221 
222  using ForwardProjectionType = typename Superclass::ForwardProjectionType;
223  using BackProjectionType = typename Superclass::BackProjectionType;
224 
225 #if !defined(ITK_WRAPPING_PARSER)
226 
228  using AddFilterType = itk::AddImageFilter<GradientsImageType>;
246 #endif
247 
249  itkSetConstObjectMacro(Geometry, ThreeDCircularProjectionGeometry);
250 
251  itkSetMacro(NumberOfIterations, int);
252  itkGetMacro(NumberOfIterations, int);
253 
255  itkSetMacro(NumberOfSubsets, int);
256  itkGetMacro(NumberOfSubsets, int);
258 
262  itkSetMacro(ResetNesterovEvery, int);
263  itkGetMacro(ResetNesterovEvery, int);
265 
267  void
268  SetInputMaterialVolumes(const TOutputImage * materialVolumes);
269  void
270  SetInputPhotonCounts(const TPhotonCounts * photonCounts);
271  void
272  SetInputSpectrum(const TSpectrum * spectrum);
273  void
274  SetSupportMask(const SingleComponentImageType * support);
275  void
276  SetSpatialRegularizationWeights(const SingleComponentImageType * regweights);
277  void
278  SetProjectionWeights(const SingleComponentImageType * weiprojections);
280 
282  itkSetMacro(RegularizationWeights, typename TOutputImage::PixelType);
283  itkGetMacro(RegularizationWeights, typename TOutputImage::PixelType);
285 
287  itkSetMacro(RegularizationRadius, typename TOutputImage::RegionType::SizeType);
288  itkGetMacro(RegularizationRadius, typename TOutputImage::RegionType::SizeType);
290 
293  using BinnedDetectorResponseType = vnl_matrix<dataType>;
294  using MaterialAttenuationsType = vnl_matrix<dataType>;
295  virtual void
296  SetBinnedDetectorResponse(const BinnedDetectorResponseType & detResp);
297  virtual void
298  SetMaterialAttenuations(const MaterialAttenuationsType & matAtt);
300 
301 protected:
303  ~MechlemOneStepSpectralReconstructionFilter() override = default;
304 
306  void
307  VerifyPreconditions() const override;
308 
310  void
311  GenerateData() override;
312 
313 #if !defined(ITK_WRAPPING_PARSER)
314 
316  typename AddFilterType::Pointer m_AddGradients;
337 #endif
338 
342  void
343  VerifyInputInformation() const override
344  {}
345 
348  void
349  GenerateInputRequestedRegion() override;
350  void
351  GenerateOutputInformation() override;
353 
355  typename TOutputImage::ConstPointer
356  GetInputMaterialVolumes();
357  typename TPhotonCounts::ConstPointer
358  GetInputPhotonCounts();
359  typename TSpectrum::ConstPointer
360  GetInputSpectrum();
361  typename SingleComponentImageType::ConstPointer
362  GetSupportMask();
363  typename SingleComponentImageType::ConstPointer
364  GetSpatialRegularizationWeights();
365  typename SingleComponentImageType::ConstPointer
366  GetProjectionWeights();
368 
369 #if !defined(ITK_WRAPPING_PARSER)
370 
372  typename SingleComponentForwardProjectionFilterType::Pointer
373  InstantiateSingleComponentForwardProjectionFilter(int fwtype);
374  typename HessiansBackProjectionFilterType::Pointer
375  InstantiateHessiansBackProjectionFilter(int bptype);
376 #endif
377 
378 
380 
387 
388  typename TOutputImage::PixelType m_RegularizationWeights;
389  typename TOutputImage::RegionType::SizeType m_RegularizationRadius;
390 };
391 } // namespace rtk
392 
393 
394 #ifndef ITK_MANUAL_INSTANTIATION
395 # include "rtkMechlemOneStepSpectralReconstructionFilter.hxx"
396 #endif
397 
398 #endif
typename itk::Image< dataType, TOutputImage::ImageDimension > SingleComponentImageType
Base class for forward projection, i.e. accumulation along x-ray lines.
Applies Nesterov&#39;s momentum technique.
Generate an n-dimensional image with constant pixel values.
ReorderProjectionsFilterProjectionsWeightsType::Pointer m_ReorderProjectionsWeightsFilter
typename itk::Image< itk::Vector< dataType, nMaterials *nMaterials >, TOutputImage::ImageDimension > HessiansImageType
ReorderProjectionsFilterPhotonCountsType::Pointer m_ReorderPhotonCountsFilter
Projection geometry for a source and a 2-D flat panel.
#define itkSetMacro(name, type)
Computes update from gradient and Hessian in Newton&#39;s method.
Performs intermediate computations in Weidinger2016.
typename itk::Image< typename TOutputImage::PixelType, TOutputImage::ImageDimension > CPUOutputImageType
For each vector-valued pixel, adds a vector to the diagonal of a matrix.
For one-step inversion of spectral CT data by the method Mechlem2017, computes regularization term&#39;s ...
Mother class for cone beam reconstruction filters which need runtime selection of their forward and b...
Implements the one-step spectral CT inversion method described by Mechlem et al.
Sorts or shuffle projections and geometry inputs.
SingleComponentForwardProjectionFilterType::Pointer m_SingleComponentForwardProjectionFilter