![]() |
RTK
2.6.0
Reconstruction Toolkit
|
#include <rtkReconstructionConjugateGradientOperator.h>
Public Member Functions | |
virtual ::itk::LightObject::Pointer | CreateAnother () const |
const char * | GetNameOfClass () const override |
void | SetBackProjectionFilter (BackProjectionFilterType *_arg) |
void | SetForwardProjectionFilter (ForwardProjectionFilterType *_arg) |
virtual void | SetGeometry (const ThreeDCircularProjectionGeometry *_arg) |
void | SetInputVolume (const TOutputImage *vol) |
void | SetInputProjectionStack (const TOutputImage *projs) |
void | SetInputWeights (const TWeightsImage *weights) |
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 () |
![]() | |
virtual ::itk::LightObject::Pointer | CreateAnother () const |
const char * | GetNameOfClass () const override |
virtual void | SetX (const TOutputImage *OutputImage) |
Static Public Member Functions | |
static Pointer | New () |
![]() | |
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 () const override |
~ReconstructionConjugateGradientOperator () override=default | |
void | GenerateInputRequestedRegion () override |
void | GenerateOutputInformation () override |
TOutputImage::ConstPointer | GetInputVolume () |
TOutputImage::ConstPointer | GetInputProjectionStack () |
TWeightsImage::ConstPointer | GetInputWeights () |
![]() | |
ConjugateGradientOperator () | |
~ConjugateGradientOperator () override=default | |
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
Definition at line 118 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::AddFilterType = itk::AddImageFilter<TOutputImage> |
Definition at line 160 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::BackProjectionFilterPointer = typename BackProjectionFilterType::Pointer |
Definition at line 153 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::BackProjectionFilterType = rtk::BackProjectionImageFilter<TOutputImage, TOutputImage> |
Definition at line 152 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::ConstantSourceType = rtk::ConstantImageSource<TOutputImage> |
Definition at line 158 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::ForwardProjectionFilterPointer = typename ForwardProjectionFilterType::Pointer |
Definition at line 156 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::ForwardProjectionFilterType = rtk::ForwardProjectionImageFilter<TOutputImage, TOutputImage> |
Definition at line 155 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GradientImageType = itk::Image<itk::CovariantVector<typename TOutputImage::PixelType, TOutputImage::ImageDimension>, TOutputImage::ImageDimension> |
Definition at line 134 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::MatrixVectorMultiplyFilterType = rtk::BlockDiagonalMatrixVectorMultiplyImageFilter<TOutputImage, TWeightsImage> |
Definition at line 165 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::MultiplyFilterType = itk::MultiplyImageFilter<TOutputImage, TSingleComponentImage> |
Definition at line 159 of file rtkReconstructionConjugateGradientOperator.h.
typedef std::conditional<std::is_same<TSingleComponentImage, TOutputImage>::value, PlainMultiplyFilterType, MatrixVectorMultiplyFilterType>::type rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::MultiplyWithWeightsFilterType |
Definition at line 169 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::OutputImagePointer = typename TOutputImage::Pointer |
Definition at line 171 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::PlainMultiplyFilterType = itk::MultiplyImageFilter<TOutputImage, TOutputImage, TOutputImage> |
Definition at line 166 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::Pointer = itk::SmartPointer<Self> |
Definition at line 126 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::Self = ReconstructionConjugateGradientOperator |
Standard class type alias.
Definition at line 124 of file rtkReconstructionConjugateGradientOperator.h.
using rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::Superclass = ConjugateGradientOperator<TOutputImage> |
Definition at line 125 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
|
overrideprotecteddefault |
|
protected |
|
protected |
|
virtual |
Reimplemented from itk::Object.
|
overrideprotectedvirtual |
Does the real work.
Reimplemented from itk::ImageSource< typename TOutputImage >.
|
overrideprotectedvirtual |
The volume and the projections must have different requested regions
Reimplemented from itk::ProcessObject.
|
overrideprotectedvirtual |
The volume and the projections must have different requested regions
Reimplemented from itk::ProcessObject.
|
virtual |
Perform laplacian-based and/or Tikhonov regularization during reconstruction (gamma is the strength of laplacian the regularization)
|
protected |
Getters for the inputs
|
protected |
Getters for the inputs
|
protected |
Getters for the inputs
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.
|
overridevirtual |
Run-time type information (and related methods).
Reimplemented from itk::ImageSource< typename TOutputImage >.
TSingleComponentImage::ConstPointer rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::GetSupportMask | ( | ) |
Set the support mask, if any, for support constraint in reconstruction
|
virtual |
Perform laplacian-based and/or Tikhonov regularization during reconstruction (gamma is the strength of laplacian the regularization)
|
static |
Method for creation through the object factory.
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetBackProjectionFilter | ( | BackProjectionFilterType * | _arg | ) |
Set the backprojection filter
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetForwardProjectionFilter | ( | ForwardProjectionFilterType * | _arg | ) |
Set the forward projection filter
|
virtual |
Perform laplacian-based and/or Tikhonov regularization during reconstruction (gamma is the strength of laplacian the regularization)
|
virtual |
Set the geometry of both m_BackProjectionFilter and m_ForwardProjectionFilter
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetInputProjectionStack | ( | const TOutputImage * | projs | ) |
Setters for the inputs
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetInputVolume | ( | const TOutputImage * | vol | ) |
Setters for the inputs
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetInputWeights | ( | const TWeightsImage * | weights | ) |
Setters for the inputs
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.
void rtk::ReconstructionConjugateGradientOperator< TOutputImage, TSingleComponentImage, TWeightsImage >::SetSupportMask | ( | const TSingleComponentImage * | SupportMask | ) |
Set the support mask, if any, for support constraint in reconstruction
|
virtual |
Perform laplacian-based and/or Tikhonov regularization during reconstruction (gamma is the strength of laplacian the regularization)
|
inlineoverrideprotectedvirtual |
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
Reimplemented from itk::ProcessObject.
Definition at line 255 of file rtkReconstructionConjugateGradientOperator.h.
|
overrideprotectedvirtual |
Checks that inputs are correctly set.
Reimplemented from itk::ProcessObject.
|
protected |
Definition at line 239 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 240 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Member pointers to the filters used internally (for convenience)
Definition at line 229 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 232 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 233 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Pointers to intermediate images, used to simplify complex branching
Definition at line 250 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 250 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 230 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 246 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Member attributes
Definition at line 245 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 241 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 235 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 236 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 234 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 237 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 238 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 242 of file rtkReconstructionConjugateGradientOperator.h.
|
protected |
Definition at line 247 of file rtkReconstructionConjugateGradientOperator.h.