![]() |
RTK
2.7.0
Reconstruction Toolkit
|
#include <rtkFourDSARTConeBeamReconstructionFilter.h>
Inheritance diagram for rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >:
Collaboration diagram for rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >:Public Types | |
| using | AddFilterType = itk::AddImageFilter< VolumeSeriesType, VolumeSeriesType > |
| using | BackProjectionFilterType = rtk::BackProjectionImageFilter< VolumeType, VolumeType > |
| using | BackProjectionType = typename Superclass::BackProjectionType |
| using | ConstantProjectionStackSourceType = rtk::ConstantImageSource< ProjectionStackType > |
| using | ConstantVolumeSeriesSourceType = rtk::ConstantImageSource< VolumeSeriesType > |
| using | ConstPointer = itk::SmartPointer< const Self > |
| using | DisplacedDetectorFilterType = rtk::DisplacedDetectorImageFilter< ProjectionStackType > |
| using | DivideFilterType = itk::DivideOrZeroOutImageFilter< ProjectionStackType, ProjectionStackType, ProjectionStackType > |
| using | ExtractFilterType = itk::ExtractImageFilter< ProjectionStackType, ProjectionStackType > |
| using | ForwardProjectionFilterType = rtk::ForwardProjectionImageFilter< ProjectionStackType, ProjectionStackType > |
| using | ForwardProjectionType = typename Superclass::ForwardProjectionType |
| using | FourDToProjectionStackFilterType = rtk::FourDToProjectionStackImageFilter< ProjectionStackType, VolumeSeriesType > |
| using | InputImageType = VolumeSeriesType |
| using | MultiplyFilterType = itk::MultiplyImageFilter< ProjectionStackType, ProjectionStackType, ProjectionStackType > |
| using | OutputImageType = VolumeSeriesType |
| using | Pointer = itk::SmartPointer< Self > |
| using | ProjectionStackToFourDFilterType = rtk::ProjectionStackToFourDImageFilter< VolumeSeriesType, ProjectionStackType > |
| using | RayBoxIntersectionFilterType = rtk::RayBoxIntersectionImageFilter< ProjectionStackType, ProjectionStackType > |
| using | Self = FourDSARTConeBeamReconstructionFilter |
| using | SubtractFilterType = itk::SubtractImageFilter< ProjectionStackType, ProjectionStackType > |
| using | Superclass = IterativeConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType > |
| using | ThresholdFilterType = itk::ThresholdImageFilter< VolumeSeriesType > |
| using | VolumeType = ProjectionStackType |
Public Types inherited from rtk::IterativeConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType > | |
| using | BackProjectionFilterType = rtk::BackProjectionImageFilter< ProjectionStackType, VolumeType > |
| using | BackProjectionPointerType = typename BackProjectionFilterType::Pointer |
| enum | BackProjectionType |
| using | ConstPointer = itk::SmartPointer< const Self > |
| using | ForwardProjectionFilterType = rtk::ForwardProjectionImageFilter< VolumeType, ProjectionStackType > |
| using | ForwardProjectionPointerType = typename ForwardProjectionFilterType::Pointer |
| enum | ForwardProjectionType |
| using | Pointer = itk::SmartPointer< Self > |
| using | Self = IterativeConeBeamReconstructionFilter |
| using | Superclass = itk::ImageToImageFilter< VolumeSeriesType, VolumeSeriesType > |
| using | TClipImageType = itk::Image< double, VolumeType::ImageDimension > |
| using | VolumeType = ProjectionStackType |
Public Member Functions | |
| virtual ::itk::LightObject::Pointer | CreateAnother () const |
| const char * | GetNameOfClass () const override |
| virtual void | SetSignal (const std::vector< double > signal) |
| void | SetWeights (const itk::Array2D< float > _arg) |
| void | SetInputVolumeSeries (const VolumeSeriesType *VolumeSeries) |
| VolumeSeriesType::ConstPointer | GetInputVolumeSeries () |
| void | SetInputProjectionStack (const ProjectionStackType *Projection) |
| ProjectionStackType::Pointer | GetInputProjectionStack () |
| virtual ThreeDCircularProjectionGeometry * | GetModifiableGeometry () |
| virtual const ThreeDCircularProjectionGeometry * | GetGeometry () const |
| virtual void | SetGeometry (ThreeDCircularProjectionGeometry *_arg) |
| virtual unsigned int | GetNumberOfIterations () |
| virtual void | SetNumberOfIterations (unsigned int _arg) |
| virtual unsigned int | GetNumberOfProjectionsPerSubset () |
| virtual void | SetNumberOfProjectionsPerSubset (unsigned int _arg) |
| virtual double | GetLambda () |
| virtual void | SetLambda (double _arg) |
| virtual bool | GetEnforcePositivity () |
| virtual void | SetEnforcePositivity (bool _arg) |
| virtual void | SetDisableDisplacedDetectorFilter (bool _arg) |
| virtual bool | GetDisableDisplacedDetectorFilter () |
Public Member Functions inherited from rtk::IterativeConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType > | |
| virtual ::itk::LightObject::Pointer | CreateAnother () const |
| const char * | GetNameOfClass () const override |
| virtual void | SetForwardProjectionFilter (ForwardProjectionType fwtype) |
| ForwardProjectionType | GetForwardProjectionFilter () |
| virtual void | SetBackProjectionFilter (BackProjectionType bptype) |
| BackProjectionType | GetBackProjectionFilter () |
| void | SetAttenuationMap (const VolumeType *attenuationMap) |
| VolumeType::ConstPointer | GetAttenuationMap () |
| void | SetInferiorClipImage (const TClipImageType *inferiorClipImage) |
| TClipImageType::ConstPointer | GetInferiorClipImage () |
| void | SetSuperiorClipImage (const TClipImageType *superiorClipImage) |
| TClipImageType::ConstPointer | GetSuperiorClipImage () |
| virtual double | GetSigmaZero () |
| virtual void | SetSigmaZero (double _arg) |
| virtual double | GetAlphaPSF () |
| virtual void | SetAlphaPSF (double _arg) |
| virtual double | GetStepSize () const |
| virtual void | SetStepSize (double _arg) |
Static Public Member Functions | |
| static Pointer | New () |
Static Public Member Functions inherited from rtk::IterativeConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType > | |
| static Pointer | New () |
Private Attributes | |
| ThreeDCircularProjectionGeometry::Pointer | m_Geometry |
| double | m_Lambda |
| unsigned int | m_NumberOfIterations |
| unsigned int | m_NumberOfProjectionsPerSubset |
Additional Inherited Members | |
Protected Types inherited from rtk::IterativeConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType > | |
| using | CPUImageType = typename itk::Image< typename ProjectionStackType::PixelType, ProjectionStackType::ImageDimension > |
| using | EnableCudaScalarAndVectorType = typename std::enable_if< !std::is_same_v< CPUImageType, ImageType > &&std::is_same_v< typename itk::PixelTraits< typename ImageType::PixelType >::ValueType, float > &&(itk::PixelTraits< typename ImageType::PixelType >::Dimension==1||itk::PixelTraits< typename ImageType::PixelType >::Dimension==2||itk::PixelTraits< typename ImageType::PixelType >::Dimension==3)>::type |
| using | DisableCudaScalarAndVectorType = typename std::enable_if< std::is_same_v< CPUImageType, ImageType >||!std::is_same_v< typename itk::PixelTraits< typename ImageType::PixelType >::ValueType, float >||(itk::PixelTraits< typename ImageType::PixelType >::Dimension !=1 &&itk::PixelTraits< typename ImageType::PixelType >::Dimension !=2 &&itk::PixelTraits< typename ImageType::PixelType >::Dimension !=3)>::type |
| using | EnableCudaScalarType = typename std::enable_if< !std::is_same_v< CPUImageType, ImageType > &&std::is_same_v< typename itk::PixelTraits< typename ImageType::PixelType >::ValueType, float > &&itk::PixelTraits< typename ImageType::PixelType >::Dimension==1 >::type |
| using | DisableCudaScalarType = typename std::enable_if< std::is_same_v< CPUImageType, ImageType >||!std::is_same_v< typename itk::PixelTraits< typename ImageType::PixelType >::ValueType, float >||itk::PixelTraits< typename ImageType::PixelType >::Dimension !=1 >::type |
| using | EnableVectorType = typename std::enable_if< itk::PixelTraits< typename ImageType::PixelType >::Dimension !=1 >::type |
| using | DisableVectorType = typename std::enable_if< itk::PixelTraits< typename ImageType::PixelType >::Dimension==1 >::type |
Implements the 4D Simultaneous Algebraic Reconstruction Technique.
FourDSARTConeBeamReconstructionFilter is a composite filter. The pipeline is essentially the same as in SARTConeBeamReconstructionFilter, with the ForwardProjectionImageFilter replaced by 4DToProjectionStackImageFilter and the BackProjectionImageFilter replaced by ProjectionStackTo4DImageFilter.
Definition at line 122 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::AddFilterType = itk::AddImageFilter<VolumeSeriesType, VolumeSeriesType> |
Definition at line 149 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::BackProjectionFilterType = rtk::BackProjectionImageFilter<VolumeType, VolumeType> |
Definition at line 150 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::BackProjectionType = typename Superclass::BackProjectionType |
Definition at line 140 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::ConstantProjectionStackSourceType = rtk::ConstantImageSource<ProjectionStackType> |
Definition at line 158 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::ConstantVolumeSeriesSourceType = rtk::ConstantImageSource<VolumeSeriesType> |
Definition at line 157 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::ConstPointer = itk::SmartPointer<const Self> |
Definition at line 132 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::DisplacedDetectorFilterType = rtk::DisplacedDetectorImageFilter<ProjectionStackType> |
Definition at line 156 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::DivideFilterType = itk::DivideOrZeroOutImageFilter<ProjectionStackType, ProjectionStackType, ProjectionStackType> |
Definition at line 155 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::ExtractFilterType = itk::ExtractImageFilter<ProjectionStackType, ProjectionStackType> |
Typedefs of each subfilter of this composite filter
Definition at line 143 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::ForwardProjectionFilterType = rtk::ForwardProjectionImageFilter<ProjectionStackType, ProjectionStackType> |
Definition at line 144 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::ForwardProjectionType = typename Superclass::ForwardProjectionType |
Definition at line 139 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::FourDToProjectionStackFilterType = rtk::FourDToProjectionStackImageFilter<ProjectionStackType, VolumeSeriesType> |
Definition at line 146 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::InputImageType = VolumeSeriesType |
Some convenient type alias.
Definition at line 135 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::MultiplyFilterType = itk::MultiplyImageFilter<ProjectionStackType, ProjectionStackType, ProjectionStackType> |
Definition at line 148 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::OutputImageType = VolumeSeriesType |
Definition at line 136 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::Pointer = itk::SmartPointer<Self> |
Definition at line 131 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::ProjectionStackToFourDFilterType = rtk::ProjectionStackToFourDImageFilter<VolumeSeriesType, ProjectionStackType> |
Definition at line 152 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::RayBoxIntersectionFilterType = rtk::RayBoxIntersectionImageFilter<ProjectionStackType, ProjectionStackType> |
Definition at line 153 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::Self = FourDSARTConeBeamReconstructionFilter |
Standard class type alias.
Definition at line 129 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::SubtractFilterType = itk::SubtractImageFilter<ProjectionStackType, ProjectionStackType> |
Definition at line 147 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::Superclass = IterativeConeBeamReconstructionFilter<VolumeSeriesType, ProjectionStackType> |
Definition at line 130 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::ThresholdFilterType = itk::ThresholdImageFilter<VolumeSeriesType> |
Definition at line 159 of file rtkFourDSARTConeBeamReconstructionFilter.h.
| using rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::VolumeType = ProjectionStackType |
Definition at line 137 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
|
overrideprotecteddefault |
| virtual::itk::LightObject::Pointer rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::CreateAnother | ( | ) | const |
|
overrideprotected |
|
overrideprotected |
|
overrideprotected |
|
virtual |
Set / Get whether the displaced detector filter should be disabled
|
virtual |
Get / Set the positivity enforcement behaviour
|
virtual |
Get / Set the object pointer to projection geometry
| ProjectionStackType::Pointer rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::GetInputProjectionStack | ( | ) |
The stack of measured projections
| VolumeSeriesType::ConstPointer rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::GetInputVolumeSeries | ( | ) |
The 4D image to be updated.
|
virtual |
Get / Set the convergence factor. Default is 0.3.
|
virtual |
Get / Set the object pointer to projection geometry
|
override |
Runtime information support.
|
virtual |
Get / Set the number of iterations. Default is 3.
|
virtual |
Get / Set the number of projections per subset. Default is 1.
|
static |
Standard New method.
|
virtual |
Set / Get whether the displaced detector filter should be disabled
|
virtual |
Get / Set the positivity enforcement behaviour
|
virtual |
Get / Set the object pointer to projection geometry
| void rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::SetInputProjectionStack | ( | const ProjectionStackType * | Projection | ) |
The stack of measured projections
| void rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::SetInputVolumeSeries | ( | const VolumeSeriesType * | VolumeSeries | ) |
The 4D image to be updated.
|
virtual |
Get / Set the convergence factor. Default is 0.3.
|
virtual |
Get / Set the number of iterations. Default is 3.
|
virtual |
Get / Set the number of projections per subset. Default is 1.
|
virtual |
Store the phase signal in a member variable
| void rtk::FourDSARTConeBeamReconstructionFilter< VolumeSeriesType, ProjectionStackType >::SetWeights | ( | const itk::Array2D< float > | _arg | ) |
Pass the interpolation weights to subfilters
|
inlineoverrideprotected |
The two inputs should not be in the same space so there is nothing to verify.
Definition at line 239 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
overrideprotected |
Checks that inputs are correctly set.
|
protected |
Definition at line 249 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 250 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 252 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 257 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 258 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 266 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 256 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 255 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 264 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Pointers to each subfilter of this composite filter
Definition at line 243 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 244 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 246 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 247 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
private |
Geometry object
Definition at line 274 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
private |
Convergence factor according to Andersen's publications which relates to the step size of the gradient descent. Default 0.3, Must be in (0,2).
Definition at line 281 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 251 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
private |
Number of iterations
Definition at line 277 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
private |
Number of projections processed before the volume is updated (1 for SART, several for OS-SART, all for SIRT)
Definition at line 271 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Miscellaneous member variables
Definition at line 262 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 263 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 253 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 254 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 265 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 248 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 259 of file rtkFourDSARTConeBeamReconstructionFilter.h.
|
protected |
Definition at line 245 of file rtkFourDSARTConeBeamReconstructionFilter.h.
1.8.14