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 #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  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<
201  std::is_same<TOutputImage, CPUOutputImageType>::value,
203  CudaWeidingerForwardModelImageFilter<TOutputImage, TMeasuredProjections, TIncidentSpectrum>>::type
205  typedef
206  typename std::conditional<std::is_same<TOutputImage, CPUOutputImageType>::value,
208  CudaForwardProjectionImageFilter<SingleComponentImageType, SingleComponentImageType>>::
210  typedef typename std::conditional<std::is_same<TOutputImage, CPUOutputImageType>::value,
212  CudaBackProjectionImageFilter<HessiansImageType>>::type
214 # else
220 # endif
221  using GradientsImageType = TOutputImage;
222 #endif
223 
224  using ForwardProjectionType = typename Superclass::ForwardProjectionType;
225  using BackProjectionType = typename Superclass::BackProjectionType;
226 
227 #if !defined(ITK_WRAPPING_PARSER)
228 
232  using AddFilterType = itk::AddImageFilter<GradientsImageType>;
250 #endif
251 
253  itkSetConstObjectMacro(Geometry, ThreeDCircularProjectionGeometry);
254 
255  itkSetMacro(NumberOfIterations, int);
256  itkGetMacro(NumberOfIterations, int);
257 
259  itkSetMacro(NumberOfSubsets, int);
260  itkGetMacro(NumberOfSubsets, int);
262 
266  itkSetMacro(ResetNesterovEvery, int);
267  itkGetMacro(ResetNesterovEvery, int);
269 
271  void
272  SetInputMaterialVolumes(const TOutputImage * materialVolumes);
273  void
274  SetInputMaterialVolumes(const VectorImageType * materialVolumes);
275  void
276  SetInputMeasuredProjections(const TMeasuredProjections * measuredProjections);
277  void
278  SetInputMeasuredProjections(const VectorImageType * measuredProjections);
279  void
280  SetInputIncidentSpectrum(const TIncidentSpectrum * incidentSpectrum);
281 #ifndef ITK_FUTURE_LEGACY_REMOVE
282  void
283  SetInputPhotonCounts(const TMeasuredProjections * measuredProjections);
284  void
285  SetInputSpectrum(const TIncidentSpectrum * incidentSpectrum);
286 #endif
287  void
288  SetSupportMask(const SingleComponentImageType * support);
289  void
290  SetSpatialRegularizationWeights(const SingleComponentImageType * regweights);
291  void
292  SetProjectionWeights(const SingleComponentImageType * weiprojections);
294 
296  itkSetMacro(RegularizationWeights, typename TOutputImage::PixelType);
297  itkGetMacro(RegularizationWeights, typename TOutputImage::PixelType);
299 
301  itkSetMacro(RegularizationRadius, typename TOutputImage::RegionType::SizeType);
302  itkGetMacro(RegularizationRadius, typename TOutputImage::RegionType::SizeType);
304 
307  using BinnedDetectorResponseType = vnl_matrix<dataType>;
308  using MaterialAttenuationsType = vnl_matrix<dataType>;
309  virtual void
310  SetBinnedDetectorResponse(const BinnedDetectorResponseType & detResp);
311  virtual void
312  SetMaterialAttenuations(const MaterialAttenuationsType & matAtt);
314 
315 protected:
317  ~MechlemOneStepSpectralReconstructionFilter() override = default;
318 
320  void
321  VerifyPreconditions() const override;
322 
324  void
325  GenerateData() override;
326 
327 #if !defined(ITK_WRAPPING_PARSER)
328 
332  typename AddFilterType::Pointer m_AddGradients;
353 #endif
354 
358  void
359  VerifyInputInformation() const override
360  {}
361 
364  void
365  GenerateInputRequestedRegion() override;
366  void
367  GenerateOutputInformation() override;
369 
371  typename TOutputImage::ConstPointer
372  GetInputMaterialVolumes();
373  typename TMeasuredProjections::ConstPointer
374  GetInputMeasuredProjections();
375  typename TIncidentSpectrum::ConstPointer
376  GetInputIncidentSpectrum();
377 #ifndef ITK_FUTURE_LEGACY_REMOVE
378  typename TMeasuredProjections::ConstPointer
379  GetInputPhotonCounts();
380  typename TIncidentSpectrum::ConstPointer
381  GetInputSpectrum();
382 #endif
383  typename SingleComponentImageType::ConstPointer
384  GetSupportMask();
385  typename SingleComponentImageType::ConstPointer
386  GetSpatialRegularizationWeights();
387  typename SingleComponentImageType::ConstPointer
388  GetProjectionWeights();
390 
391 #if !defined(ITK_WRAPPING_PARSER)
392 
394  typename SingleComponentForwardProjectionFilterType::Pointer
395  InstantiateSingleComponentForwardProjectionFilter(int fwtype);
396  typename HessiansBackProjectionFilterType::Pointer
397  InstantiateHessiansBackProjectionFilter(int bptype);
398 #endif
399 
400 
402 
409 
410  typename TOutputImage::PixelType m_RegularizationWeights;
411  typename TOutputImage::RegionType::SizeType m_RegularizationRadius;
412 };
413 } // namespace rtk
414 
415 
416 #ifndef ITK_MANUAL_INSTANTIATION
417 # include "rtkMechlemOneStepSpectralReconstructionFilter.hxx"
418 #endif
419 
420 #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.
#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
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...
typename itk::Image< dataType, TOutputImage::ImageDimension > SingleComponentImageType
Implements the one-step spectral CT inversion method described by Mechlem et al.
typename itk::Image< typename TOutputImage::PixelType, TOutputImage::ImageDimension > CPUOutputImageType
SingleComponentForwardProjectionFilterType::Pointer m_SingleComponentForwardProjectionFilter
typename itk::Image< itk::Vector< dataType, nMaterials *nMaterials >, TOutputImage::ImageDimension > HessiansImageType
Sorts or shuffle projections and geometry inputs.