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

#include <rtkProjectionStackToFourDImageFilter.h>

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

Public Types

using BackProjectionFilterType = rtk::BackProjectionImageFilter< VolumeType, VolumeType >
 
using ConstantVolumeSeriesSourceType = rtk::ConstantImageSource< VolumeSeriesType >
 
using ConstantVolumeSourceType = rtk::ConstantImageSource< VolumeType >
 
using CPUVolumeSeriesType = typename itk::Image< typename VolumeSeriesType::PixelType, VolumeSeriesType::ImageDimension >
 
using CudaConstantVolumeSeriesSourceType = ConstantVolumeSeriesSourceType
 
using CudaConstantVolumeSourceType = ConstantVolumeSourceType
 
using CudaSplatImageFilterType = SplatFilterType
 
using ExtractFilterType = itk::ExtractImageFilter< ProjectionStackType, ProjectionStackType >
 
using GeometryType = rtk::ThreeDCircularProjectionGeometry
 
using Pointer = itk::SmartPointer< Self >
 
using Self = ProjectionStackToFourDImageFilter
 
using SplatFilterType = rtk::SplatWithKnownWeightsImageFilter< VolumeSeriesType, VolumeType >
 
using Superclass = itk::ImageToImageFilter< VolumeSeriesType, VolumeSeriesType >
 
using VolumeType = ProjectionStackType
 

Public Member Functions

virtual ::itk::LightObject::Pointer CreateAnother () const
 
virtual const char * GetNameOfClass () const
 
virtual bool GetUseCudaSources ()
 
virtual bool GetUseCudaSplat ()
 
virtual itk::Array2D< float > GetWeights ()
 
void SetBackProjectionFilter (const typename BackProjectionFilterType::Pointer _arg)
 
virtual void SetGeometry (const GeometryType *_arg)
 
virtual void SetSignal (const std::vector< double > signal)
 
virtual void SetUseCudaSources (bool _arg)
 
virtual void SetUseCudaSplat (bool _arg)
 
virtual void SetWeights (itk::Array2D< float > _arg)
 
void SetInputVolumeSeries (const VolumeSeriesType *VolumeSeries)
 
VolumeSeriesType::ConstPointer GetInputVolumeSeries ()
 
void SetInputProjectionStack (const ProjectionStackType *Projections)
 
ProjectionStackType::ConstPointer GetInputProjectionStack ()
 

Static Public Member Functions

static Pointer New ()
 

Protected Member Functions

void GenerateData () override
 
void GenerateInputRequestedRegion () override
 
void GenerateOutputInformation () override
 
void InitializeConstantSource ()
 
 ProjectionStackToFourDImageFilter ()
 
 ~ProjectionStackToFourDImageFilter () override=default
 

Protected Attributes

BackProjectionFilterType::Pointer m_BackProjectionFilter
 
ConstantVolumeSeriesSourceType::Pointer m_ConstantVolumeSeriesSource
 
ConstantVolumeSourceType::Pointer m_ConstantVolumeSource
 
ExtractFilterType::Pointer m_ExtractFilter
 
GeometryType::ConstPointer m_Geometry
 
std::vector< double > m_Signal
 
SplatFilterType::Pointer m_SplatFilter
 
bool m_UseCudaSources
 
bool m_UseCudaSplat
 
itk::Array2D< float > m_Weights
 

Detailed Description

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
class rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >

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.

ProjectionStackToFourDImageFilter implements S_theta^T R_theta^T.

dot_inline_dotgraph_23.png
Test:
rtkfourdconjugategradienttest.cxx
Author
Cyril Mory

Definition at line 105 of file rtkProjectionStackToFourDImageFilter.h.

Member Typedef Documentation

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::BackProjectionFilterType = rtk::BackProjectionImageFilter< VolumeType, VolumeType >

Definition at line 134 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::ConstantVolumeSeriesSourceType = rtk::ConstantImageSource< VolumeSeriesType >

Definition at line 137 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::ConstantVolumeSourceType = rtk::ConstantImageSource< VolumeType >

Definition at line 136 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::CPUVolumeSeriesType = typename itk::Image< typename VolumeSeriesType::PixelType, VolumeSeriesType::ImageDimension>

SFINAE type alias, depending on whether a CUDA image is used.

Definition at line 144 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::CudaConstantVolumeSeriesSourceType = ConstantVolumeSeriesSourceType

Definition at line 158 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::CudaConstantVolumeSourceType = ConstantVolumeSourceType

Definition at line 157 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::CudaSplatImageFilterType = SplatFilterType

Definition at line 156 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::ExtractFilterType = itk::ExtractImageFilter< ProjectionStackType, ProjectionStackType >

Definition at line 135 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::GeometryType = rtk::ThreeDCircularProjectionGeometry

Definition at line 140 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::Pointer = itk::SmartPointer< Self >

Definition at line 113 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::Self = ProjectionStackToFourDImageFilter

Standard class type alias.

Definition at line 111 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::SplatFilterType = rtk::SplatWithKnownWeightsImageFilter<VolumeSeriesType, VolumeType>

Definition at line 138 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::Superclass = itk::ImageToImageFilter< VolumeSeriesType, VolumeSeriesType >

Definition at line 112 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
using rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::VolumeType = ProjectionStackType

Convenient type alias

Definition at line 116 of file rtkProjectionStackToFourDImageFilter.h.

Constructor & Destructor Documentation

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::ProjectionStackToFourDImageFilter ( )
protected
template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::~ProjectionStackToFourDImageFilter ( )
overrideprotecteddefault

Member Function Documentation

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
virtual::itk::LightObject::Pointer rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::CreateAnother ( ) const
template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
void rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::GenerateData ( )
overrideprotected

Does the real work.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
void rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::GenerateInputRequestedRegion ( )
overrideprotected
template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
void rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::GenerateOutputInformation ( )
overrideprotected
template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
ProjectionStackType::ConstPointer rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::GetInputProjectionStack ( )

Set/Get the stack of projections

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
VolumeSeriesType::ConstPointer rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::GetInputVolumeSeries ( )

Set/Get the 4D image to be updated.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
virtual const char* rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::GetNameOfClass ( ) const
virtual

Run-time type information (and related methods).

Reimplemented from itk::ImageSource< TOutputImage >.

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

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
virtual bool rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::GetUseCudaSources ( )
virtual
template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
virtual bool rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::GetUseCudaSplat ( )
virtual
template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
virtual itk::Array2D<float> rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::GetWeights ( )
virtual

Macros that take care of implementing the Get and Set methods for Weights

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
void rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::InitializeConstantSource ( )
protected
template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
static Pointer rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::New ( )
static

Method for creation through the object factory.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
void rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::SetBackProjectionFilter ( const typename BackProjectionFilterType::Pointer  _arg)

Pass the backprojection filter to SingleProjectionToFourDFilter

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
virtual void rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::SetGeometry ( const GeometryType _arg)
virtual

Pass the geometry to SingleProjectionToFourDFilter

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
void rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::SetInputProjectionStack ( const ProjectionStackType *  Projections)

Set/Get the stack of projections

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
void rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::SetInputVolumeSeries ( const VolumeSeriesType *  VolumeSeries)

Set/Get the 4D image to be updated.

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

Store the phase signal in a member variable

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

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
virtual void rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::SetUseCudaSources ( bool  _arg)
virtual
template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
virtual void rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::SetUseCudaSplat ( bool  _arg)
virtual

Use CUDA splat / sources

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
virtual void rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::SetWeights ( itk::Array2D< float >  _arg)
virtual

Member Data Documentation

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
BackProjectionFilterType::Pointer rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::m_BackProjectionFilter
protected

Definition at line 195 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
ConstantVolumeSeriesSourceType::Pointer rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::m_ConstantVolumeSeriesSource
protected

Definition at line 198 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
ConstantVolumeSourceType::Pointer rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::m_ConstantVolumeSource
protected

Definition at line 197 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
ExtractFilterType::Pointer rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::m_ExtractFilter
protected

Definition at line 196 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
GeometryType::ConstPointer rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::m_Geometry
protected

Definition at line 202 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
std::vector<double> rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::m_Signal
protected

Definition at line 205 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
SplatFilterType::Pointer rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::m_SplatFilter
protected

Member pointers to the filters used internally (for convenience)

Definition at line 194 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
bool rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::m_UseCudaSources
protected

Definition at line 204 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
bool rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::m_UseCudaSplat
protected

Definition at line 203 of file rtkProjectionStackToFourDImageFilter.h.

template<typename VolumeSeriesType, typename ProjectionStackType, typename TFFTPrecision = double>
itk::Array2D<float> rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType, TFFTPrecision >::m_Weights
protected

Other member variables

Definition at line 201 of file rtkProjectionStackToFourDImageFilter.h.


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