RTK  2.7.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 #include <itkCastImageFilter.h>
37 
38 #ifdef RTK_USE_CUDA
40 #endif
41 
42 namespace rtk
43 {
154 template <typename TOutputImage, typename TMeasuredProjections, typename TIncidentSpectrum>
156  : public rtk::IterativeConeBeamReconstructionFilter<TOutputImage, TOutputImage>
157 {
158 public:
159  ITK_DISALLOW_COPY_AND_MOVE(MechlemOneStepSpectralReconstructionFilter);
160 
165 
167  itkNewMacro(Self);
168 
170  itkOverrideGetNameOfClassMacro(MechlemOneStepSpectralReconstructionFilter);
171 
173  static constexpr unsigned int nBins = TMeasuredProjections::PixelType::Dimension;
174  static constexpr unsigned int nMaterials = TOutputImage::PixelType::Dimension;
175  using dataType = typename TOutputImage::PixelType::ValueType;
176 
179 
182 #ifdef RTK_USE_CUDA
183  using HessiansImageType = typename std::conditional_t<
184  std::is_same_v<TOutputImage, CPUOutputImageType>,
185  itk::Image<itk::Vector<dataType, nMaterials * nMaterials>, TOutputImage::ImageDimension>,
186  itk::CudaImage<itk::Vector<dataType, nMaterials * nMaterials>, TOutputImage::ImageDimension>>;
187  using SingleComponentImageType = typename std::conditional_t<std::is_same_v<TOutputImage, CPUOutputImageType>,
189  itk::CudaImage<dataType, TOutputImage::ImageDimension>>;
190 #else
191  using HessiansImageType =
192  typename itk::Image<itk::Vector<dataType, nMaterials * nMaterials>, TOutputImage::ImageDimension>;
194 #endif
195 
196 #if !defined(ITK_WRAPPING_PARSER)
197 # ifdef RTK_USE_CUDA
198  using WeidingerForwardModelType = typename std::conditional_t<
199  std::is_same_v<TOutputImage, CPUOutputImageType>,
203  typename std::conditional_t<std::is_same_v<TOutputImage, CPUOutputImageType>,
207  typename std::conditional_t<std::is_same_v<TOutputImage, CPUOutputImageType>,
210 # else
216 # endif
217  using GradientsImageType = TOutputImage;
218 #endif
219 
220  using ForwardProjectionType = typename Superclass::ForwardProjectionType;
221  using BackProjectionType = typename Superclass::BackProjectionType;
222 
223 #if !defined(ITK_WRAPPING_PARSER)
224 
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  SetInputMaterialVolumes(const VectorImageType * materialVolumes);
271  void
272  SetInputMeasuredProjections(const TMeasuredProjections * measuredProjections);
273  void
274  SetInputMeasuredProjections(const VectorImageType * measuredProjections);
275  void
276  SetInputIncidentSpectrum(const TIncidentSpectrum * incidentSpectrum);
277 #ifndef ITK_FUTURE_LEGACY_REMOVE
278  void
279  SetInputPhotonCounts(const TMeasuredProjections * measuredProjections);
280  void
281  SetInputSpectrum(const TIncidentSpectrum * incidentSpectrum);
282 #endif
283  void
284  SetSupportMask(const SingleComponentImageType * support);
285  void
286  SetSpatialRegularizationWeights(const SingleComponentImageType * regweights);
287  void
288  SetProjectionWeights(const SingleComponentImageType * weiprojections);
290 
292  itkSetMacro(RegularizationWeights, typename TOutputImage::PixelType);
293  itkGetMacro(RegularizationWeights, typename TOutputImage::PixelType);
295 
297  itkSetMacro(RegularizationRadius, typename TOutputImage::RegionType::SizeType);
298  itkGetMacro(RegularizationRadius, typename TOutputImage::RegionType::SizeType);
300 
303  using BinnedDetectorResponseType = vnl_matrix<dataType>;
304  using MaterialAttenuationsType = vnl_matrix<dataType>;
305  virtual void
306  SetBinnedDetectorResponse(const BinnedDetectorResponseType & detResp);
307  virtual void
308  SetMaterialAttenuations(const MaterialAttenuationsType & matAtt);
310 
311 protected:
313  ~MechlemOneStepSpectralReconstructionFilter() override = default;
314 
316  void
317  VerifyPreconditions() const override;
318 
320  void
321  GenerateData() override;
322 
323 #if !defined(ITK_WRAPPING_PARSER)
324 
335  typename WeidingerForwardModelType::Pointer m_WeidingerForward;
349 #endif
350 
354  void
355  VerifyInputInformation() const override
356  {}
357 
360  void
361  GenerateInputRequestedRegion() override;
362  void
363  GenerateOutputInformation() override;
365 
367  typename TOutputImage::ConstPointer
368  GetInputMaterialVolumes();
369  typename TMeasuredProjections::ConstPointer
370  GetInputMeasuredProjections();
371  typename TIncidentSpectrum::ConstPointer
372  GetInputIncidentSpectrum();
373 #ifndef ITK_FUTURE_LEGACY_REMOVE
374  typename TMeasuredProjections::ConstPointer
375  GetInputPhotonCounts();
376  typename TIncidentSpectrum::ConstPointer
377  GetInputSpectrum();
378 #endif
379  typename SingleComponentImageType::ConstPointer
380  GetSupportMask();
381  typename SingleComponentImageType::ConstPointer
382  GetSpatialRegularizationWeights();
383  typename SingleComponentImageType::ConstPointer
384  GetProjectionWeights();
386 
387 #if !defined(ITK_WRAPPING_PARSER)
388 
390  typename SingleComponentForwardProjectionFilterType::Pointer
391  InstantiateSingleComponentForwardProjectionFilter(int fwtype);
392  typename HessiansBackProjectionFilterType::Pointer
393  InstantiateHessiansBackProjectionFilter(int bptype);
394 #endif
395 
396 
398 
405 
406  typename TOutputImage::PixelType m_RegularizationWeights;
407  typename TOutputImage::RegionType::SizeType m_RegularizationRadius;
408 };
409 } // namespace rtk
410 
411 
412 #ifndef ITK_MANUAL_INSTANTIATION
413 # include "rtkMechlemOneStepSpectralReconstructionFilter.hxx"
414 #endif
415 
416 #endif
ReorderProjectionsWeightsFilterType::Pointer m_ReorderProjectionsWeightsFilter
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.
typename itk::VectorImage< dataType, TOutputImage::ImageDimension > VectorImageType
Projection geometry for a source and a 2-D flat panel.
typename std::conditional_t< std::is_same_v< TOutputImage, CPUOutputImageType >, BackProjectionImageFilter< HessiansImageType, HessiansImageType >, CudaBackProjectionImageFilter< HessiansImageType > > CudaHessiansBackProjectionImageFilterType
#define itkSetMacro(name, type)
Computes update from gradient and Hessian in Newton&#39;s method.
ExtractMeasuredProjectionsFilterType::Pointer m_ExtractMeasuredProjectionsFilter
Performs intermediate computations in Weidinger2016.
ReorderMeasuredProjectionsFilterType::Pointer m_ReorderMeasuredProjectionsFilter
typename std::conditional_t< std::is_same_v< TOutputImage, CPUOutputImageType >, itk::Image< itk::Vector< dataType, nMaterials *nMaterials >, TOutputImage::ImageDimension >, itk::CudaImage< itk::Vector< dataType, nMaterials *nMaterials >, TOutputImage::ImageDimension > > HessiansImageType
Cuda version of the backprojection.
typename std::conditional_t< std::is_same_v< TOutputImage, CPUOutputImageType >, JosephForwardProjectionImageFilter< SingleComponentImageType, SingleComponentImageType >, CudaForwardProjectionImageFilter< SingleComponentImageType, SingleComponentImageType > > CudaSingleComponentForwardProjectionImageFilterType
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.
typename itk::Image< typename TOutputImage::PixelType, TOutputImage::ImageDimension > CPUOutputImageType
typename std::conditional_t< std::is_same_v< TOutputImage, CPUOutputImageType >, WeidingerForwardModelImageFilter< TOutputImage, TMeasuredProjections, TIncidentSpectrum >, CudaWeidingerForwardModelImageFilter< TOutputImage, TMeasuredProjections, TIncidentSpectrum > > WeidingerForwardModelType
SingleComponentForwardProjectionFilterType::Pointer m_SingleComponentForwardProjectionFilter
Sorts or shuffle projections and geometry inputs.
typename std::conditional_t< std::is_same_v< TOutputImage, CPUOutputImageType >, itk::Image< dataType, TOutputImage::ImageDimension >, itk::CudaImage< dataType, TOutputImage::ImageDimension > > SingleComponentImageType