RTK  2.6.0
Reconstruction Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType > Class Template Reference

#include <rtkFourDToProjectionStackImageFilter.h>

+ Inheritance diagram for rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >:
+ Collaboration diagram for rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >:

Public Types

using ConstantProjectionStackSourceType = rtk::ConstantImageSource< ProjectionStackType >
 
using ConstantVolumeSourceType = rtk::ConstantImageSource< VolumeType >
 
using ForwardProjectionFilterType = rtk::ForwardProjectionImageFilter< ProjectionStackType, ProjectionStackType >
 
using GeometryType = rtk::ThreeDCircularProjectionGeometry
 
using InterpolatorFilterType = rtk::InterpolatorWithKnownWeightsImageFilter< VolumeType, VolumeSeriesType >
 
using PasteFilterType = itk::PasteImageFilter< ProjectionStackType, ProjectionStackType >
 
using Pointer = itk::SmartPointer< Self >
 
using Self = FourDToProjectionStackImageFilter
 
using Superclass = itk::ImageToImageFilter< ProjectionStackType, ProjectionStackType >
 
using VolumeType = ProjectionStackType
 

Public Member Functions

virtual ::itk::LightObject::Pointer CreateAnother () const
 
void InitializeConstantVolumeSource ()
 
void SetForwardProjectionFilter (const typename ForwardProjectionFilterType::Pointer _arg)
 
virtual void SetGeometry (GeometryType::Pointer _arg)
 
void SetInputProjectionStack (const ProjectionStackType *Projection)
 
void SetInputVolumeSeries (const VolumeSeriesType *VolumeSeries)
 
virtual void SetSignal (const std::vector< double > signal)
 
void SetWeights (const itk::Array2D< float > _arg)
 
virtual const char * GetNameOfClass () const
 

Static Public Member Functions

static Pointer New ()
 

Protected Member Functions

 FourDToProjectionStackImageFilter ()
 
void GenerateData () override
 
void GenerateInputRequestedRegion () override
 
void GenerateOutputInformation () override
 
ProjectionStackType::Pointer GetInputProjectionStack ()
 
VolumeSeriesType::ConstPointer GetInputVolumeSeries ()
 
 ~FourDToProjectionStackImageFilter () override=default
 

Protected Attributes

ConstantProjectionStackSourceType::Pointer m_ConstantProjectionStackSource
 
ConstantVolumeSourceType::Pointer m_ConstantVolumeSource
 
ForwardProjectionFilterType::Pointer m_ForwardProjectionFilter
 
GeometryType::Pointer m_Geometry
 
InterpolatorFilterType::Pointer m_InterpolationFilter
 
PasteFilterType::Pointer m_PasteFilter
 
ConstantProjectionStackSourceType::OutputImageRegionType m_PasteRegion
 
std::vector< double > m_Signal
 
itk::Array2D< float > m_Weights
 

Detailed Description

template<typename ProjectionStackType, typename VolumeSeriesType>
class rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >

Implements part of the 4D reconstruction by conjugate gradient.

See the reference paper: "Cardiac C-arm computed tomography using a 3D + time ROI reconstruction method with spatial and temporal regularization" by Mory et al.

4D conjugate gradient reconstruction consists in iteratively minimizing the following cost function:

Sum_over_theta || R_theta S_theta f - p_theta ||_2^2

with

Computing the gradient of this cost function yields:

S_theta^T R_theta^T R_theta S_theta f - S_theta^T R_theta^T p_theta

where A^T means the adjoint of operator A.

FourDToProjectionStackImageFilter implements R_theta S_theta.

Test:
rtkfourdconjugategradienttest.cxx
Author
Cyril Mory

Definition at line 100 of file rtkFourDToProjectionStackImageFilter.h.

Member Typedef Documentation

◆ ConstantProjectionStackSourceType

template<typename ProjectionStackType , typename VolumeSeriesType >
using rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::ConstantProjectionStackSourceType = rtk::ConstantImageSource<ProjectionStackType>

Definition at line 140 of file rtkFourDToProjectionStackImageFilter.h.

◆ ConstantVolumeSourceType

template<typename ProjectionStackType , typename VolumeSeriesType >
using rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::ConstantVolumeSourceType = rtk::ConstantImageSource<VolumeType>

Definition at line 139 of file rtkFourDToProjectionStackImageFilter.h.

◆ ForwardProjectionFilterType

template<typename ProjectionStackType , typename VolumeSeriesType >
using rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::ForwardProjectionFilterType = rtk::ForwardProjectionImageFilter<ProjectionStackType, ProjectionStackType>

Typedefs for the sub filters

Definition at line 136 of file rtkFourDToProjectionStackImageFilter.h.

◆ GeometryType

template<typename ProjectionStackType , typename VolumeSeriesType >
using rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::GeometryType = rtk::ThreeDCircularProjectionGeometry

Definition at line 141 of file rtkFourDToProjectionStackImageFilter.h.

◆ InterpolatorFilterType

template<typename ProjectionStackType , typename VolumeSeriesType >
using rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::InterpolatorFilterType = rtk::InterpolatorWithKnownWeightsImageFilter<VolumeType, VolumeSeriesType>

Definition at line 138 of file rtkFourDToProjectionStackImageFilter.h.

◆ PasteFilterType

template<typename ProjectionStackType , typename VolumeSeriesType >
using rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::PasteFilterType = itk::PasteImageFilter<ProjectionStackType, ProjectionStackType>

Definition at line 137 of file rtkFourDToProjectionStackImageFilter.h.

◆ Pointer

template<typename ProjectionStackType , typename VolumeSeriesType >
using rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::Pointer = itk::SmartPointer<Self>

Definition at line 109 of file rtkFourDToProjectionStackImageFilter.h.

◆ Self

template<typename ProjectionStackType , typename VolumeSeriesType >
using rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::Self = FourDToProjectionStackImageFilter

Standard class type alias.

Definition at line 107 of file rtkFourDToProjectionStackImageFilter.h.

◆ Superclass

template<typename ProjectionStackType , typename VolumeSeriesType >
using rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::Superclass = itk::ImageToImageFilter<ProjectionStackType, ProjectionStackType>

Definition at line 108 of file rtkFourDToProjectionStackImageFilter.h.

◆ VolumeType

template<typename ProjectionStackType , typename VolumeSeriesType >
using rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::VolumeType = ProjectionStackType

Convenient type alias

Definition at line 112 of file rtkFourDToProjectionStackImageFilter.h.

Constructor & Destructor Documentation

◆ FourDToProjectionStackImageFilter()

template<typename ProjectionStackType , typename VolumeSeriesType >
rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::FourDToProjectionStackImageFilter ( )
protected

◆ ~FourDToProjectionStackImageFilter()

template<typename ProjectionStackType , typename VolumeSeriesType >
rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::~FourDToProjectionStackImageFilter ( )
overrideprotecteddefault

Member Function Documentation

◆ CreateAnother()

template<typename ProjectionStackType , typename VolumeSeriesType >
virtual::itk::LightObject::Pointer rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::CreateAnother ( ) const
virtual

◆ GenerateData()

template<typename ProjectionStackType , typename VolumeSeriesType >
void rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::GenerateData ( )
overrideprotectedvirtual

◆ GenerateInputRequestedRegion()

template<typename ProjectionStackType , typename VolumeSeriesType >
void rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::GenerateInputRequestedRegion ( )
overrideprotectedvirtual

◆ GenerateOutputInformation()

template<typename ProjectionStackType , typename VolumeSeriesType >
void rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::GenerateOutputInformation ( )
overrideprotectedvirtual

◆ GetInputProjectionStack()

template<typename ProjectionStackType , typename VolumeSeriesType >
ProjectionStackType::Pointer rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::GetInputProjectionStack ( )
protected

◆ GetInputVolumeSeries()

template<typename ProjectionStackType , typename VolumeSeriesType >
VolumeSeriesType::ConstPointer rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::GetInputVolumeSeries ( )
protected

◆ GetNameOfClass()

template<typename ProjectionStackType , typename VolumeSeriesType >
virtual const char* rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::GetNameOfClass ( ) const
virtual

Run-time type information (and related methods).

Reimplemented from itk::ImageSource< TOutputImage >.

Reimplemented in rtk::WarpFourDToProjectionStackImageFilter< VolumeSeriesType, ProjectionStackType >.

◆ InitializeConstantVolumeSource()

template<typename ProjectionStackType , typename VolumeSeriesType >
void rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::InitializeConstantVolumeSource ( )

Initializes the empty volume source, set it and update it

◆ New()

template<typename ProjectionStackType , typename VolumeSeriesType >
static Pointer rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::New ( )
static

Method for creation through the object factory.

◆ SetForwardProjectionFilter()

template<typename ProjectionStackType , typename VolumeSeriesType >
void rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::SetForwardProjectionFilter ( const typename ForwardProjectionFilterType::Pointer  _arg)

Set the ForwardProjection filter

◆ SetGeometry()

template<typename ProjectionStackType , typename VolumeSeriesType >
virtual void rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::SetGeometry ( GeometryType::Pointer  _arg)
virtual

Pass the geometry to SingleProjectionToFourDFilter

◆ SetInputProjectionStack()

template<typename ProjectionStackType , typename VolumeSeriesType >
void rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::SetInputProjectionStack ( const ProjectionStackType *  Projection)

The image that will be backprojected, then added, with coefficients, to each 3D volume of the 4D image. It is 3D because the ForwardProjection filters need it, but the third dimension, which is the number of projections, is 1

◆ SetInputVolumeSeries()

template<typename ProjectionStackType , typename VolumeSeriesType >
void rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::SetInputVolumeSeries ( const VolumeSeriesType *  VolumeSeries)

The 4D image to be updated.

◆ SetSignal()

template<typename ProjectionStackType , typename VolumeSeriesType >
virtual void rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::SetSignal ( const std::vector< double >  signal)
virtual

Store the phase signal in a member variable

Reimplemented in rtk::WarpFourDToProjectionStackImageFilter< VolumeSeriesType, ProjectionStackType >.

◆ SetWeights()

template<typename ProjectionStackType , typename VolumeSeriesType >
void rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::SetWeights ( const itk::Array2D< float >  _arg)

Pass the interpolation weights to SingleProjectionToFourDFilter

Member Data Documentation

◆ m_ConstantProjectionStackSource

template<typename ProjectionStackType , typename VolumeSeriesType >
ConstantProjectionStackSourceType::Pointer rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::m_ConstantProjectionStackSource
protected

Definition at line 186 of file rtkFourDToProjectionStackImageFilter.h.

◆ m_ConstantVolumeSource

template<typename ProjectionStackType , typename VolumeSeriesType >
ConstantVolumeSourceType::Pointer rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::m_ConstantVolumeSource
protected

Definition at line 185 of file rtkFourDToProjectionStackImageFilter.h.

◆ m_ForwardProjectionFilter

template<typename ProjectionStackType , typename VolumeSeriesType >
ForwardProjectionFilterType::Pointer rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::m_ForwardProjectionFilter
protected

Definition at line 187 of file rtkFourDToProjectionStackImageFilter.h.

◆ m_Geometry

template<typename ProjectionStackType , typename VolumeSeriesType >
GeometryType::Pointer rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::m_Geometry
protected

Definition at line 191 of file rtkFourDToProjectionStackImageFilter.h.

◆ m_InterpolationFilter

template<typename ProjectionStackType , typename VolumeSeriesType >
InterpolatorFilterType::Pointer rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::m_InterpolationFilter
protected

Definition at line 184 of file rtkFourDToProjectionStackImageFilter.h.

◆ m_PasteFilter

template<typename ProjectionStackType , typename VolumeSeriesType >
PasteFilterType::Pointer rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::m_PasteFilter
protected

Member pointers to the filters used internally (for convenience)

Definition at line 183 of file rtkFourDToProjectionStackImageFilter.h.

◆ m_PasteRegion

template<typename ProjectionStackType , typename VolumeSeriesType >
ConstantProjectionStackSourceType::OutputImageRegionType rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::m_PasteRegion
protected

Definition at line 192 of file rtkFourDToProjectionStackImageFilter.h.

◆ m_Signal

template<typename ProjectionStackType , typename VolumeSeriesType >
std::vector<double> rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::m_Signal
protected

Definition at line 193 of file rtkFourDToProjectionStackImageFilter.h.

◆ m_Weights

template<typename ProjectionStackType , typename VolumeSeriesType >
itk::Array2D<float> rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType >::m_Weights
protected

Other member variables

Definition at line 190 of file rtkFourDToProjectionStackImageFilter.h.


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