RTK  2.5.0
Reconstruction Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage > Class Template Reference

#include <rtkReconstructionConjugateGradientOperator.h>

+ Inheritance diagram for rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >:
+ Collaboration diagram for rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >:

Public Types

using AddFilterType = itk::AddImageFilter< TOutputImage >
 
using BackProjectionFilterPointer = typename BackProjectionFilterType::Pointer
 
using BackProjectionFilterType = rtk::BackProjectionImageFilter< TOutputImage, TOutputImage >
 
using ConstantSourceType = rtk::ConstantImageSource< TOutputImage >
 
using ForwardProjectionFilterPointer = typename ForwardProjectionFilterType::Pointer
 
using ForwardProjectionFilterType = rtk::ForwardProjectionImageFilter< TOutputImage, TOutputImage >
 
using GradientImageType = itk::Image< itk::CovariantVector< typename TOutputImage::PixelType, TOutputImage::ImageDimension >, TOutputImage::ImageDimension >
 
using MatrixVectorMultiplyFilterType = rtk::BlockDiagonalMatrixVectorMultiplyImageFilter< TOutputImage, TWeightsImage >
 
using MultiplyFilterType = itk::MultiplyImageFilter< TOutputImage, TSingleComponentImage >
 
typedef std::conditional< std::is_same< TSingleComponentImage, TOutputImage >::value, PlainMultiplyFilterType, MatrixVectorMultiplyFilterType >::type MultiplyWithWeightsFilterType
 
using OutputImagePointer = typename TOutputImage::Pointer
 
using PlainMultiplyFilterType = itk::MultiplyImageFilter< TOutputImage, TOutputImage, TOutputImage >
 
using Pointer = itk::SmartPointer< Self >
 
using Self = ReconstructionConjugateGradientOperator
 
using Superclass = ConjugateGradientOperator< TOutputImage >
 
- Public Types inherited from rtk::ConjugateGradientOperator< TOutputImage >
using Pointer = itk::SmartPointer< Self >
 
using Self = ConjugateGradientOperator
 
using Superclass = itk::ImageToImageFilter< TOutputImage, TOutputImage >
 

Public Member Functions

virtual ::itk::LightObject::Pointer CreateAnother () const
 
void SetBackProjectionFilter (const BackProjectionFilterPointer _arg)
 
void SetForwardProjectionFilter (const ForwardProjectionFilterPointer _arg)
 
virtual void SetGeometry (const ThreeDCircularProjectionGeometry *_arg)
 
void SetInputVolume (const TOutputImage *vol)
 
void SetInputProjectionStack (const TOutputImage *projs)
 
void SetInputWeights (const TWeightsImage *weights)
 
virtual const char * GetNameOfClass () const
 
void SetSupportMask (const TSingleComponentImage *SupportMask)
 
TSingleComponentImage::ConstPointer GetSupportMask ()
 
void SetLocalRegularizationWeights (const TSingleComponentImage *localRegularizationWeights)
 
TSingleComponentImage::ConstPointer GetLocalRegularizationWeights ()
 
virtual void SetGamma (float _arg)
 
virtual float GetGamma ()
 
virtual void SetTikhonov (float _arg)
 
virtual float GetTikhonov ()
 
- Public Member Functions inherited from rtk::ConjugateGradientOperator< TOutputImage >
virtual ::itk::LightObject::Pointer CreateAnother () const
 
virtual void SetX (const TOutputImage *OutputImage)
 

Static Public Member Functions

static Pointer New ()
 
- Static Public Member Functions inherited from rtk::ConjugateGradientOperator< TOutputImage >
static Pointer New ()
 

Protected Member Functions

template<typename ImageType >
std::enable_if< std::is_same< TSingleComponentImage, ImageType >::value, ImageType >::type::Pointer ConnectGradientRegularization ()
 
template<typename ImageType >
std::enable_if<!std::is_same< TSingleComponentImage, ImageType >::value, ImageType >::type::Pointer ConnectGradientRegularization ()
 
void GenerateData () override
 
 ReconstructionConjugateGradientOperator ()
 
void VerifyInputInformation () const override
 
void VerifyPreconditions () ITKv5_CONST override
 
 ~ReconstructionConjugateGradientOperator () override=default
 
void GenerateInputRequestedRegion () override
 
void GenerateOutputInformation () override
 
TOutputImage::ConstPointer GetInputVolume ()
 
TOutputImage::ConstPointer GetInputProjectionStack ()
 
TWeightsImage::ConstPointer GetInputWeights ()
 
- Protected Member Functions inherited from rtk::ConjugateGradientOperator< TOutputImage >
 ConjugateGradientOperator ()
 
 ~ConjugateGradientOperator () override=default
 

Protected Attributes

AddFilterType::Pointer m_AddLaplacianFilter
 
AddFilterType::Pointer m_AddTikhonovFilter
 
BackProjectionFilterPointer m_BackProjectionFilter
 
ConstantSourceType::Pointer m_ConstantProjectionsSource
 
ConstantSourceType::Pointer m_ConstantVolumeSource
 
TOutputImage::Pointer m_FloatingInputPointer
 
TOutputImage::Pointer m_FloatingOutputPointer
 
ForwardProjectionFilterPointer m_ForwardProjectionFilter
 
float m_Gamma { 0 }
 
rtk::ThreeDCircularProjectionGeometry::ConstPointer m_Geometry { nullptr }
 
itk::ImageToImageFilter< TOutputImage, TOutputImage >::Pointer m_LaplacianFilter
 
MultiplyFilterType::Pointer m_MultiplyInputVolumeFilter
 
MultiplyFilterType::Pointer m_MultiplyLaplacianFilter
 
MultiplyFilterType::Pointer m_MultiplyOutputVolumeFilter
 
MultiplyFilterType::Pointer m_MultiplyTikhonovFilter
 
MultiplyFilterType::Pointer m_MultiplyTikhonovWeightsFilter
 
MultiplyWithWeightsFilterType::Pointer m_MultiplyWithWeightsFilter
 
float m_Tikhonov { 0 }
 

Detailed Description

template<typename TOutputImage, typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
class rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >

Implements the operator A used in conjugate gradient reconstruction.

This filter implements the operator A used in the conjugate gradient reconstruction method, which attempts to find the f that minimizes || sqrt(D) (Rf -p) ||_2^2 + gamma || grad f ||_2^2 + Tikhonov || f ||_2^2, with R the forward projection operator, p the measured projections, and D the displaced detector weighting operator.

With gamma=0, this it is similar to the ART and SART methods. The difference lies in the algorithm employed to minimize this cost function. ART uses the Kaczmarz method (projects and back projects one ray at a time), SART the block-Kaczmarz method (projects and back projects one projection at a time), and ConjugateGradient a conjugate gradient method (projects and back projects all projections together).

This filter takes in input f and outputs R_t D R f + gamma Laplacian f + Tikhonov f

Test:
rtkconjugategradienttest.cxx
Author
Cyril Mory

Definition at line 118 of file rtkReconstructionConjugateGradientOperator.h.

Member Typedef Documentation

◆ AddFilterType

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::AddFilterType = itk::AddImageFilter<TOutputImage>

◆ BackProjectionFilterPointer

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::BackProjectionFilterPointer = typename BackProjectionFilterType::Pointer

◆ BackProjectionFilterType

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::BackProjectionFilterType = rtk::BackProjectionImageFilter<TOutputImage, TOutputImage>

◆ ConstantSourceType

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::ConstantSourceType = rtk::ConstantImageSource<TOutputImage>

◆ ForwardProjectionFilterPointer

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::ForwardProjectionFilterPointer = typename ForwardProjectionFilterType::Pointer

◆ ForwardProjectionFilterType

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::ForwardProjectionFilterType = rtk::ForwardProjectionImageFilter<TOutputImage, TOutputImage>

◆ GradientImageType

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GradientImageType = itk::Image<itk::CovariantVector<typename TOutputImage::PixelType, TOutputImage::ImageDimension>, TOutputImage::ImageDimension>

◆ MatrixVectorMultiplyFilterType

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::MatrixVectorMultiplyFilterType = rtk::BlockDiagonalMatrixVectorMultiplyImageFilter<TOutputImage, TWeightsImage>

◆ MultiplyFilterType

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::MultiplyFilterType = itk::MultiplyImageFilter<TOutputImage, TSingleComponentImage>

◆ MultiplyWithWeightsFilterType

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
typedef std::conditional<std::is_same<TSingleComponentImage, TOutputImage>::value, PlainMultiplyFilterType, MatrixVectorMultiplyFilterType>::type rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::MultiplyWithWeightsFilterType

◆ OutputImagePointer

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::OutputImagePointer = typename TOutputImage::Pointer

◆ PlainMultiplyFilterType

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::PlainMultiplyFilterType = itk::MultiplyImageFilter<TOutputImage, TOutputImage, TOutputImage>

◆ Pointer

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::Pointer = itk::SmartPointer<Self>

◆ Self

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::Self = ReconstructionConjugateGradientOperator

Standard class type alias.

Definition at line 124 of file rtkReconstructionConjugateGradientOperator.h.

◆ Superclass

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::Superclass = ConjugateGradientOperator<TOutputImage>

Constructor & Destructor Documentation

◆ ReconstructionConjugateGradientOperator()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::ReconstructionConjugateGradientOperator ( )
protected

◆ ~ReconstructionConjugateGradientOperator()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::~ReconstructionConjugateGradientOperator ( )
overrideprotecteddefault

Member Function Documentation

◆ ConnectGradientRegularization() [1/2]

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
template<typename ImageType >
std::enable_if<std::is_same<TSingleComponentImage, ImageType>::value, ImageType>::type::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::ConnectGradientRegularization ( )
protected

◆ ConnectGradientRegularization() [2/2]

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
template<typename ImageType >
std::enable_if<!std::is_same<TSingleComponentImage, ImageType>::value, ImageType>::type::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::ConnectGradientRegularization ( )
protected

◆ CreateAnother()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
virtual::itk::LightObject::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::CreateAnother ( ) const
virtual

Reimplemented from itk::Object.

◆ GenerateData()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GenerateData ( )
overrideprotectedvirtual

Does the real work.

Reimplemented from itk::ImageSource< TOutputImage >.

◆ GenerateInputRequestedRegion()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GenerateInputRequestedRegion ( )
overrideprotectedvirtual

The volume and the projections must have different requested regions

Reimplemented from itk::ProcessObject.

◆ GenerateOutputInformation()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GenerateOutputInformation ( )
overrideprotectedvirtual

The volume and the projections must have different requested regions

Reimplemented from itk::ProcessObject.

◆ GetGamma()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
virtual float rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GetGamma ( )
virtual

Perform laplacian-based and/or Tikhonov regularization during reconstruction (gamma is the strength of laplacian the regularization)

◆ GetInputProjectionStack()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
TOutputImage::ConstPointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GetInputProjectionStack ( )
protected

Getters for the inputs

◆ GetInputVolume()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
TOutputImage::ConstPointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GetInputVolume ( )
protected

Getters for the inputs

◆ GetInputWeights()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
TWeightsImage::ConstPointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GetInputWeights ( )
protected

Getters for the inputs

◆ GetLocalRegularizationWeights()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
TSingleComponentImage::ConstPointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GetLocalRegularizationWeights ( )

Set local regularization weights. The map should have the same information (size, spacing, origin etc.) as the reconstructed volume. The same map is used in Laplacian and Tikhonov regularization.

◆ GetNameOfClass()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
virtual const char* rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GetNameOfClass ( ) const
virtual

Run-time type information (and related methods).

Reimplemented from rtk::ConjugateGradientOperator< TOutputImage >.

◆ GetSupportMask()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
TSingleComponentImage::ConstPointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GetSupportMask ( )

Set the support mask, if any, for support constraint in reconstruction

◆ GetTikhonov()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
virtual float rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GetTikhonov ( )
virtual

Perform laplacian-based and/or Tikhonov regularization during reconstruction (gamma is the strength of laplacian the regularization)

◆ New()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
static Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::New ( )
static

Method for creation through the object factory.

◆ SetBackProjectionFilter()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetBackProjectionFilter ( const BackProjectionFilterPointer  _arg)

Set the backprojection filter

◆ SetForwardProjectionFilter()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetForwardProjectionFilter ( const ForwardProjectionFilterPointer  _arg)

Set the forward projection filter

◆ SetGamma()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
virtual void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetGamma ( float  _arg)
virtual

Perform laplacian-based and/or Tikhonov regularization during reconstruction (gamma is the strength of laplacian the regularization)

◆ SetGeometry()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
virtual void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetGeometry ( const ThreeDCircularProjectionGeometry _arg)
virtual

Set the geometry of both m_BackProjectionFilter and m_ForwardProjectionFilter

◆ SetInputProjectionStack()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetInputProjectionStack ( const TOutputImage *  projs)

Setters for the inputs

◆ SetInputVolume()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetInputVolume ( const TOutputImage *  vol)

Setters for the inputs

◆ SetInputWeights()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetInputWeights ( const TWeightsImage *  weights)

Setters for the inputs

◆ SetLocalRegularizationWeights()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetLocalRegularizationWeights ( const TSingleComponentImage *  localRegularizationWeights)

Set local regularization weights. The map should have the same information (size, spacing, origin etc.) as the reconstructed volume. The same map is used in Laplacian and Tikhonov regularization.

◆ SetSupportMask()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetSupportMask ( const TSingleComponentImage *  SupportMask)

Set the support mask, if any, for support constraint in reconstruction

◆ SetTikhonov()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
virtual void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetTikhonov ( float  _arg)
virtual

Perform laplacian-based and/or Tikhonov regularization during reconstruction (gamma is the strength of laplacian the regularization)

◆ VerifyInputInformation()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::VerifyInputInformation ( ) const
inlineoverrideprotected

When the inputs have the same type, ITK checks whether they occupy the same physical space or not. Obviously they dont, so we have to remove this check

Definition at line 260 of file rtkReconstructionConjugateGradientOperator.h.

◆ VerifyPreconditions()

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::VerifyPreconditions ( )
overrideprotectedvirtual

Checks that inputs are correctly set.

Reimplemented from itk::ProcessObject.

Member Data Documentation

◆ m_AddLaplacianFilter

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
AddFilterType::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_AddLaplacianFilter
protected

◆ m_AddTikhonovFilter

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
AddFilterType::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_AddTikhonovFilter
protected

◆ m_BackProjectionFilter

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
BackProjectionFilterPointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_BackProjectionFilter
protected

Member pointers to the filters used internally (for convenience)

Definition at line 234 of file rtkReconstructionConjugateGradientOperator.h.

◆ m_ConstantProjectionsSource

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
ConstantSourceType::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_ConstantProjectionsSource
protected

◆ m_ConstantVolumeSource

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
ConstantSourceType::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_ConstantVolumeSource
protected

◆ m_FloatingInputPointer

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
TOutputImage::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_FloatingInputPointer
protected

Pointers to intermediate images, used to simplify complex branching

Definition at line 255 of file rtkReconstructionConjugateGradientOperator.h.

◆ m_FloatingOutputPointer

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
TOutputImage::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_FloatingOutputPointer
protected

◆ m_ForwardProjectionFilter

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
ForwardProjectionFilterPointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_ForwardProjectionFilter
protected

◆ m_Gamma

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
float rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_Gamma { 0 }
protected

◆ m_Geometry

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
rtk::ThreeDCircularProjectionGeometry::ConstPointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_Geometry { nullptr }
protected

Member attributes

Definition at line 250 of file rtkReconstructionConjugateGradientOperator.h.

◆ m_LaplacianFilter

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
itk::ImageToImageFilter<TOutputImage, TOutputImage>::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_LaplacianFilter
protected

◆ m_MultiplyInputVolumeFilter

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
MultiplyFilterType::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_MultiplyInputVolumeFilter
protected

◆ m_MultiplyLaplacianFilter

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
MultiplyFilterType::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_MultiplyLaplacianFilter
protected

◆ m_MultiplyOutputVolumeFilter

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
MultiplyFilterType::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_MultiplyOutputVolumeFilter
protected

◆ m_MultiplyTikhonovFilter

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
MultiplyFilterType::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_MultiplyTikhonovFilter
protected

◆ m_MultiplyTikhonovWeightsFilter

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
MultiplyFilterType::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_MultiplyTikhonovWeightsFilter
protected

◆ m_MultiplyWithWeightsFilter

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
MultiplyWithWeightsFilterType::Pointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_MultiplyWithWeightsFilter
protected

◆ m_Tikhonov

template<typename TOutputImage , typename TSingleComponentImage = TOutputImage, typename TWeightsImage = TOutputImage>
float rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::m_Tikhonov { 0 }
protected

The documentation for this class was generated from the following file: