RTK  2.5.0
Reconstruction Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Attributes | List of all members
rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision > Class Template Reference

#include <rtkFFTRampImageFilter.h>

+ Inheritance diagram for rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >:
+ Collaboration diagram for rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >:

Public Types

using ConstPointer = itk::SmartPointer< const Self >
 
using FFTInputImagePointer = typename FFTInputImageType::Pointer
 
using FFTInputImageType = typename Superclass::FFTInputImageType
 
using FFTOutputImagePointer = typename FFTOutputImageType::Pointer
 
using FFTOutputImageType = typename Superclass::FFTOutputImageType
 
using FFTPrecisionType = TFFTPrecision
 
using IndexType = typename InputImageType::IndexType
 
using InputImageType = TInputImage
 
using OutputImageType = TOutputImage
 
using Pointer = itk::SmartPointer< Self >
 
using Self = FFTRampImageFilter
 
using SizeType = typename InputImageType::SizeType
 
using Superclass = rtk::FFTProjectionsConvolutionImageFilter< TInputImage, TOutputImage, TFFTPrecision >
 
- Public Types inherited from rtk::FFTProjectionsConvolutionImageFilter< TInputImage, TOutputImage, TFFTPrecision >
using ConstPointer = itk::SmartPointer< const Self >
 
using FFTInputImagePointer = typename FFTInputImageType::Pointer
 
using FFTInputImageType = typename itk::Image< TFFTPrecision, TInputImage::ImageDimension >
 
using FFTOutputImagePointer = typename FFTOutputImageType::Pointer
 
using FFTOutputImageType = typename itk::Image< std::complex< TFFTPrecision >, TInputImage::ImageDimension >
 
using IndexType = typename InputImageType::IndexType
 
using InputImageType = TInputImage
 
using OutputImageType = TOutputImage
 
using Pointer = itk::SmartPointer< Self >
 
using RegionType = typename InputImageType::RegionType
 
using Self = FFTProjectionsConvolutionImageFilter
 
using SizeType = typename InputImageType::SizeType
 
using Superclass = itk::ImageToImageFilter< TInputImage, TOutputImage >
 
using ZeroPadFactorsType = itk::Vector< int, 2 >
 

Public Member Functions

virtual ::itk::LightObject::Pointer CreateAnother () const
 
virtual const char * GetNameOfClass () const
 
virtual double GetHannCutFrequency () const
 
 itkSetMacro (HannCutFrequency, double)
 
virtual double GetCosineCutFrequency () const
 
 itkSetMacro (CosineCutFrequency, double)
 
virtual double GetHammingFrequency () const
 
 itkSetMacro (HammingFrequency, double)
 
virtual double GetHannCutFrequencyY () const
 
 itkSetMacro (HannCutFrequencyY, double)
 
virtual double GetRamLakCutFrequency () const
 
 itkSetMacro (RamLakCutFrequency, double)
 
virtual double GetSheppLoganCutFrequency () const
 
 itkSetMacro (SheppLoganCutFrequency, double)
 
- Public Member Functions inherited from rtk::FFTProjectionsConvolutionImageFilter< TInputImage, TOutputImage, TFFTPrecision >
virtual int GetGreatestPrimeFactor () const
 
virtual void SetGreatestPrimeFactor (int _arg)
 
virtual double GetTruncationCorrection () const
 
virtual void SetTruncationCorrection (double _arg)
 
virtual ZeroPadFactorsType GetZeroPadFactors () const
 
virtual void SetZeroPadFactors (ZeroPadFactorsType _arg)
 

Static Public Member Functions

static Pointer New ()
 

Protected Member Functions

 FFTRampImageFilter ()
 
void GenerateInputRequestedRegion () override
 
void UpdateFFTProjectionsConvolutionKernel (const SizeType s) override
 
 ~FFTRampImageFilter () override=default
 
- Protected Member Functions inherited from rtk::FFTProjectionsConvolutionImageFilter< TInputImage, TOutputImage, TFFTPrecision >
void AfterThreadedGenerateData () override
 
void BeforeThreadedGenerateData () override
 
 FFTProjectionsConvolutionImageFilter ()
 
int GreatestPrimeFactor (int n) const
 
bool IsPrime (int n) const
 
void PrintSelf (std::ostream &os, itk::Indent indent) const override
 
void ThreadedGenerateData (const RegionType &outputRegionForThread, ThreadIdType threadId) override
 
virtual void UpdateTruncationMirrorWeights ()
 
 ~FFTProjectionsConvolutionImageFilter () override=default
 
virtual FFTInputImagePointer PadInputImageRegion (const RegionType &inputRegion)
 
RegionType GetPaddedImageRegion (const RegionType &inputRegion)
 

Private Attributes

double m_CosineCutFrequency { 0. }
 
double m_HammingFrequency { 0. }
 
double m_HannCutFrequency { 0. }
 
double m_HannCutFrequencyY { 0. }
 
SizeType m_PreviousKernelUpdateSize
 
double m_RamLakCutFrequency { 0. }
 
double m_SheppLoganCutFrequency { 0. }
 

Additional Inherited Members

- Static Public Attributes inherited from rtk::FFTProjectionsConvolutionImageFilter< TInputImage, TOutputImage, TFFTPrecision >
static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension
 
- Protected Attributes inherited from rtk::FFTProjectionsConvolutionImageFilter< TInputImage, TOutputImage, TFFTPrecision >
int m_KernelDimension { 1 }
 
FFTOutputImagePointer m_KernelFFT
 
std::vector< TFFTPrecision > m_TruncationMirrorWeights
 

Detailed Description

template<class TInputImage, class TOutputImage = TInputImage, class TFFTPrecision = double>
class rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >

Implements the ramp image filter of the filtered backprojection algorithm.

The filter code is based on FFTProjectionsConvolutionImageFilter by Gaetan Lehmann (see https://hdl.handle.net/10380/3154)

Test:
rtkrampfiltertest.cxx
Author
Simon Rit

Definition at line 65 of file rtkFFTRampImageFilter.h.

Member Typedef Documentation

◆ ConstPointer

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::ConstPointer = itk::SmartPointer<const Self>

Definition at line 75 of file rtkFFTRampImageFilter.h.

◆ FFTInputImagePointer

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::FFTInputImagePointer = typename FFTInputImageType::Pointer

Definition at line 85 of file rtkFFTRampImageFilter.h.

◆ FFTInputImageType

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::FFTInputImageType = typename Superclass::FFTInputImageType

Definition at line 84 of file rtkFFTRampImageFilter.h.

◆ FFTOutputImagePointer

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::FFTOutputImagePointer = typename FFTOutputImageType::Pointer

Definition at line 87 of file rtkFFTRampImageFilter.h.

◆ FFTOutputImageType

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::FFTOutputImageType = typename Superclass::FFTOutputImageType

Definition at line 86 of file rtkFFTRampImageFilter.h.

◆ FFTPrecisionType

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::FFTPrecisionType = TFFTPrecision

Definition at line 80 of file rtkFFTRampImageFilter.h.

◆ IndexType

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::IndexType = typename InputImageType::IndexType

Definition at line 81 of file rtkFFTRampImageFilter.h.

◆ InputImageType

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::InputImageType = TInputImage

Some convenient type alias.

Definition at line 78 of file rtkFFTRampImageFilter.h.

◆ OutputImageType

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::OutputImageType = TOutputImage

Definition at line 79 of file rtkFFTRampImageFilter.h.

◆ Pointer

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::Pointer = itk::SmartPointer<Self>

Definition at line 74 of file rtkFFTRampImageFilter.h.

◆ Self

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::Self = FFTRampImageFilter

Standard class type alias.

Definition at line 72 of file rtkFFTRampImageFilter.h.

◆ SizeType

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::SizeType = typename InputImageType::SizeType

Definition at line 82 of file rtkFFTRampImageFilter.h.

◆ Superclass

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
using rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::Superclass = rtk::FFTProjectionsConvolutionImageFilter<TInputImage, TOutputImage, TFFTPrecision>

Definition at line 73 of file rtkFFTRampImageFilter.h.

Constructor & Destructor Documentation

◆ FFTRampImageFilter()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::FFTRampImageFilter ( )
protected

◆ ~FFTRampImageFilter()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::~FFTRampImageFilter ( )
overrideprotecteddefault

Member Function Documentation

◆ CreateAnother()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
virtual::itk::LightObject::Pointer rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::CreateAnother ( ) const
virtual

Reimplemented from itk::Object.

◆ GenerateInputRequestedRegion()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
void rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::GenerateInputRequestedRegion ( )
overrideprotectedvirtual

◆ GetCosineCutFrequency()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
virtual double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::GetCosineCutFrequency ( ) const
virtual

Set/Get the Cosine Cut window frequency. 0 (default) disables it

◆ GetHammingFrequency()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
virtual double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::GetHammingFrequency ( ) const
virtual

Set/Get the Hamming window frequency. 0 (default) disables it

◆ GetHannCutFrequency()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
virtual double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::GetHannCutFrequency ( ) const
virtual

Set/Get the Hann window frequency. 0 (default) disables it

◆ GetHannCutFrequencyY()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
virtual double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::GetHannCutFrequencyY ( ) const
virtual

Set/Get the Hann window frequency in Y direction. 0 (default) disables it

◆ GetNameOfClass()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
virtual const char* rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::GetNameOfClass ( ) const
virtual

◆ GetRamLakCutFrequency()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
virtual double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::GetRamLakCutFrequency ( ) const
virtual

Set/Get the Ram-Lak window frequency (0...1). 0 (default) disable it. Equation and further explanation about Ram-Lak filter could be found in:

  1. Fundamentals of 2D and 3D reconstruction (by Dr. Gunter Lauritsch). https://campar.in.tum.de/twiki/pub/Chair/TeachingWs04IOIV/08CTReconstruction.pdf
  2. Reconstruction. http://oftankonyv.reak.bme.hu/tiki-index.php?page=Reconstruction

◆ GetSheppLoganCutFrequency()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
virtual double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::GetSheppLoganCutFrequency ( ) const
virtual

Set/Get the Shepp-Logan window frequency (0...1). 0 (default) disable it. Equation and further explanation about Shepp-Logan filter could be found in:

  1. Fundamentals of 2D and 3D reconstruction (by Dr. Gunter Lauritsch). https://campar.in.tum.de/twiki/pub/Chair/TeachingWs04IOIV/08CTReconstruction.pdf
  2. Reconstruction. http://oftankonyv.reak.bme.hu/tiki-index.php?page=Reconstruction

◆ itkSetMacro() [1/6]

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::itkSetMacro ( HannCutFrequency  ,
double   
)

Set/Get the Hann window frequency. 0 (default) disables it

◆ itkSetMacro() [2/6]

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::itkSetMacro ( CosineCutFrequency  ,
double   
)

Set/Get the Cosine Cut window frequency. 0 (default) disables it

◆ itkSetMacro() [3/6]

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::itkSetMacro ( HammingFrequency  ,
double   
)

Set/Get the Hamming window frequency. 0 (default) disables it

◆ itkSetMacro() [4/6]

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::itkSetMacro ( HannCutFrequencyY  ,
double   
)

Set/Get the Hann window frequency in Y direction. 0 (default) disables it

◆ itkSetMacro() [5/6]

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::itkSetMacro ( RamLakCutFrequency  ,
double   
)

Set/Get the Ram-Lak window frequency (0...1). 0 (default) disable it. Equation and further explanation about Ram-Lak filter could be found in:

  1. Fundamentals of 2D and 3D reconstruction (by Dr. Gunter Lauritsch). https://campar.in.tum.de/twiki/pub/Chair/TeachingWs04IOIV/08CTReconstruction.pdf
  2. Reconstruction. http://oftankonyv.reak.bme.hu/tiki-index.php?page=Reconstruction

◆ itkSetMacro() [6/6]

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::itkSetMacro ( SheppLoganCutFrequency  ,
double   
)

Set/Get the Shepp-Logan window frequency (0...1). 0 (default) disable it. Equation and further explanation about Shepp-Logan filter could be found in:

  1. Fundamentals of 2D and 3D reconstruction (by Dr. Gunter Lauritsch). https://campar.in.tum.de/twiki/pub/Chair/TeachingWs04IOIV/08CTReconstruction.pdf
  2. Reconstruction. http://oftankonyv.reak.bme.hu/tiki-index.php?page=Reconstruction

◆ New()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
static Pointer rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::New ( )
static

Standard New method.

◆ UpdateFFTProjectionsConvolutionKernel()

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
void rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::UpdateFFTProjectionsConvolutionKernel ( const SizeType  s)
overrideprotectedvirtual

Creates and return a pointer to one line of the ramp kernel in Fourier space. Used in generate data functions.

Implements rtk::FFTProjectionsConvolutionImageFilter< TInputImage, TOutputImage, TFFTPrecision >.

Member Data Documentation

◆ m_CosineCutFrequency

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::m_CosineCutFrequency { 0. }
private

Definition at line 160 of file rtkFFTRampImageFilter.h.

◆ m_HammingFrequency

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::m_HammingFrequency { 0. }
private

Definition at line 161 of file rtkFFTRampImageFilter.h.

◆ m_HannCutFrequency

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::m_HannCutFrequency { 0. }
private

Cut frequency of Hann, Cosine and Hamming windows. The first one which is non-zero is used.

Definition at line 159 of file rtkFFTRampImageFilter.h.

◆ m_HannCutFrequencyY

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::m_HannCutFrequencyY { 0. }
private

Definition at line 162 of file rtkFFTRampImageFilter.h.

◆ m_PreviousKernelUpdateSize

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
SizeType rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::m_PreviousKernelUpdateSize
private

Definition at line 169 of file rtkFFTRampImageFilter.h.

◆ m_RamLakCutFrequency

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::m_RamLakCutFrequency { 0. }
private

Cut frequency of Ram-Lak and Shepp-Logan

Definition at line 166 of file rtkFFTRampImageFilter.h.

◆ m_SheppLoganCutFrequency

template<class TInputImage , class TOutputImage = TInputImage, class TFFTPrecision = double>
double rtk::FFTRampImageFilter< TInputImage, TOutputImage, TFFTPrecision >::m_SheppLoganCutFrequency { 0. }
private

Definition at line 167 of file rtkFFTRampImageFilter.h.


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