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  using HessiansImageType =
183  typename TOutputImage::template RebindImageType<itk::Vector<dataType, nMaterials * nMaterials>,
184  TOutputImage::ImageDimension>;
186  typename TOutputImage::template RebindImageType<dataType, TOutputImage::ImageDimension>;
187 
188 #if !defined(ITK_WRAPPING_PARSER)
189 # ifdef RTK_USE_CUDA
190  using WeidingerForwardModelType = typename std::conditional_t<
191  std::is_same_v<TOutputImage, CPUOutputImageType>,
195  typename std::conditional_t<std::is_same_v<TOutputImage, CPUOutputImageType>,
199  typename std::conditional_t<std::is_same_v<TOutputImage, CPUOutputImageType>,
202 # else
208 # endif
209  using GradientsImageType = TOutputImage;
210 #endif
211 
212  using ForwardProjectionType = typename Superclass::ForwardProjectionType;
213  using BackProjectionType = typename Superclass::BackProjectionType;
214 
215 #if !defined(ITK_WRAPPING_PARSER)
216 
238 #endif
239 
241  itkSetConstObjectMacro(Geometry, ThreeDCircularProjectionGeometry);
242 
243  itkSetMacro(NumberOfIterations, int);
244  itkGetMacro(NumberOfIterations, int);
245 
247  itkSetMacro(NumberOfSubsets, int);
248  itkGetMacro(NumberOfSubsets, int);
250 
254  itkSetMacro(ResetNesterovEvery, int);
255  itkGetMacro(ResetNesterovEvery, int);
257 
259  void
260  SetInputMaterialVolumes(const TOutputImage * materialVolumes);
261  void
262  SetInputMaterialVolumes(const VectorImageType * materialVolumes);
263  void
264  SetInputMeasuredProjections(const TMeasuredProjections * measuredProjections);
265  void
266  SetInputMeasuredProjections(const VectorImageType * measuredProjections);
267  void
268  SetInputIncidentSpectrum(const TIncidentSpectrum * incidentSpectrum);
269 #ifndef ITK_FUTURE_LEGACY_REMOVE
270  void
271  SetInputPhotonCounts(const TMeasuredProjections * measuredProjections);
272  void
273  SetInputSpectrum(const TIncidentSpectrum * incidentSpectrum);
274 #endif
275  void
276  SetSupportMask(const SingleComponentImageType * support);
277  void
278  SetSpatialRegularizationWeights(const SingleComponentImageType * regweights);
279  void
280  SetProjectionWeights(const SingleComponentImageType * weiprojections);
282 
284  itkSetMacro(RegularizationWeights, typename TOutputImage::PixelType);
285  itkGetMacro(RegularizationWeights, typename TOutputImage::PixelType);
287 
289  itkSetMacro(RegularizationRadius, typename TOutputImage::RegionType::SizeType);
290  itkGetMacro(RegularizationRadius, typename TOutputImage::RegionType::SizeType);
292 
295  using BinnedDetectorResponseType = vnl_matrix<dataType>;
296  using MaterialAttenuationsType = vnl_matrix<dataType>;
297  virtual void
298  SetBinnedDetectorResponse(const BinnedDetectorResponseType & detResp);
299  virtual void
300  SetMaterialAttenuations(const MaterialAttenuationsType & matAtt);
302 
303 protected:
305  ~MechlemOneStepSpectralReconstructionFilter() override = default;
306 
308  void
309  VerifyPreconditions() const override;
310 
312  void
313  GenerateData() override;
314 
315 #if !defined(ITK_WRAPPING_PARSER)
316 
327  typename WeidingerForwardModelType::Pointer m_WeidingerForward;
341 #endif
342 
346  void
347  VerifyInputInformation() const override
348  {}
349 
352  void
353  GenerateInputRequestedRegion() override;
354  void
355  GenerateOutputInformation() override;
357 
359  typename TOutputImage::ConstPointer
360  GetInputMaterialVolumes();
361  typename TMeasuredProjections::ConstPointer
362  GetInputMeasuredProjections();
363  typename TIncidentSpectrum::ConstPointer
364  GetInputIncidentSpectrum();
365 #ifndef ITK_FUTURE_LEGACY_REMOVE
366  typename TMeasuredProjections::ConstPointer
367  GetInputPhotonCounts();
368  typename TIncidentSpectrum::ConstPointer
369  GetInputSpectrum();
370 #endif
371  typename SingleComponentImageType::ConstPointer
372  GetSupportMask();
373  typename SingleComponentImageType::ConstPointer
374  GetSpatialRegularizationWeights();
375  typename SingleComponentImageType::ConstPointer
376  GetProjectionWeights();
378 
379 #if !defined(ITK_WRAPPING_PARSER)
380 
382  typename SingleComponentForwardProjectionFilterType::Pointer
383  InstantiateSingleComponentForwardProjectionFilter(int fwtype);
384  typename HessiansBackProjectionFilterType::Pointer
385  InstantiateHessiansBackProjectionFilter(int bptype);
386 #endif
387 
388 
390 
397 
398  typename TOutputImage::PixelType m_RegularizationWeights;
399  typename TOutputImage::RegionType::SizeType m_RegularizationRadius;
400 };
401 } // namespace rtk
402 
403 
404 #ifndef ITK_MANUAL_INSTANTIATION
405 # include "rtkMechlemOneStepSpectralReconstructionFilter.hxx"
406 #endif
407 
408 #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 TOutputImage::template RebindImageType< dataType, TOutputImage::ImageDimension > SingleComponentImageType
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
typename TOutputImage::template RebindImageType< itk::Vector< dataType, nMaterials *nMaterials >, TOutputImage::ImageDimension > HessiansImageType
Sorts or shuffle projections and geometry inputs.