ITK  6.0.0
Insight Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Attributes | List of all members
itk::BSplineControlPointImageFunction< TInputImage, TCoordinate > Class Template Reference

#include <itkBSplineControlPointImageFunction.h>

Detailed Description

template<typename TInputImage, typename TCoordinate = double>
class itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >

Evaluate a B-spline object given a grid of control points.

The output of the class itkBSplineScatteredDataPointSetToImageFilter
is a control point grid defining a B-spline object. This class is used to hold various routines meant to operate on that control point grid. In addition to specifying the control point grid as the input, the user must also supply the spline order and the parametric domain (i.e. size, domain, origin, spacing).

Operations include

  1. Evaluation of the B-spline object at any point in the domain.
  2. Evaluation of the gradient of the B-spline object at any point in the domain.
  3. Evaluation of the Hessian of the B-spline object at any point in the domain.
Author
Nicholas J. Tustison

Definition at line 58 of file itkBSplineControlPointImageFunction.h.

+ Inheritance diagram for itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >:
+ Collaboration diagram for itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >:

Public Types

using ArrayType = FixedArray< unsigned int, ImageDimension >
 
using ConstPointer = SmartPointer< const Self >
 
using ControlPointLatticeType = TInputImage
 
using CoordinateType = TCoordinate
 
using GradientType = VariableSizeMatrix< CoordinateType >
 
using HessianComponentType = VariableSizeMatrix< CoordinateType >
 
using IndexType = typename InputImageType::IndexType
 
using InputImageRegionType = typename InputImageType::RegionType
 
using InputImageType = TInputImage
 
using KernelOrder0Type = BSplineKernelFunction< 0 >
 
using KernelOrder1Type = BSplineKernelFunction< 1 >
 
using KernelOrder2Type = BSplineKernelFunction< 2 >
 
using KernelOrder3Type = BSplineKernelFunction< 3 >
 
using KernelType = CoxDeBoorBSplineKernelFunction< 3 >
 
using OriginType = typename InputImageType::PointType
 
using OutputType = PixelType
 
using PixelType = typename InputImageType::PixelType
 
using Pointer = SmartPointer< Self >
 
using RealImagePointer = typename RealImageType::Pointer
 
using RealImageType = Image< CoordinateType, ImageDimension >
 
using RealType = float
 
using RegionType = typename InputImageType::RegionType
 
using Self = BSplineControlPointImageFunction
 
using SizeType = typename InputImageType::SizeType
 
using SpacingType = typename InputImageType::SpacingType
 
using Superclass = ImageFunction< TInputImage, typename TInputImage::PixelType, TCoordinate >
 
- Public Types inherited from itk::ImageFunction< TInputImage, TInputImage::PixelType, TCoordinate >
using ConstPointer = SmartPointer< const Self >
 
using ContinuousIndexType = ContinuousIndex< TCoordinate, Self::ImageDimension >
 
using CoordinateType = TCoordinate
 
using IndexType = typename InputImageType::IndexType
 
using IndexValueType = typename InputImageType::IndexValueType
 
using InputImageConstPointer = typename InputImageType::ConstPointer
 
using InputImageType = TInputImage
 
using InputPixelType = typename InputImageType::PixelType
 
using OutputType = TInputImage::PixelType
 
using Pointer = SmartPointer< Self >
 
using PointType = Point< TCoordinate, Self::ImageDimension >
 
using Self = ImageFunction
 
using Superclass = FunctionBase< Point< TCoordinate, Self::ImageDimension >, TInputImage::PixelType >
 
- Public Types inherited from itk::FunctionBase< Point< TCoordinate, TInputImage::ImageDimension >, TInputImage::PixelType >
using ConstPointer = SmartPointer< const Self >
 
using InputType = Point< TCoordinate, TInputImage::ImageDimension >
 
using OutputType = TInputImage::PixelType
 
using Pointer = SmartPointer< Self >
 
using Self = FunctionBase
 
using Superclass = Object
 
- Public Types inherited from itk::Object
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = Object
 
using Superclass = LightObject
 
- Public Types inherited from itk::LightObject
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = LightObject
 

Public Member Functions

OutputType Evaluate (const PointType &) const override
 
OutputType EvaluateAtContinuousIndex (const ContinuousIndexType &) const override
 
OutputType EvaluateAtIndex (const IndexType &) const override
 
OutputType EvaluateAtParametricPoint (const PointType &) const
 
GradientType EvaluateGradient (const PointType &) const
 
GradientType EvaluateGradientAtContinuousIndex (const ContinuousIndexType &) const
 
GradientType EvaluateGradientAtIndex (const IndexType &) const
 
GradientType EvaluateGradientAtParametricPoint (const PointType &) const
 
HessianComponentType EvaluateHessian (const PointType &, const unsigned int) const
 
HessianComponentType EvaluateHessianAtContinuousIndex (const ContinuousIndexType &, const unsigned int) const
 
HessianComponentType EvaluateHessianAtIndex (const IndexType &, const unsigned int) const
 
HessianComponentType EvaluateHessianAtParametricPoint (const PointType &, const unsigned int) const
 
virtual const ArrayTypeGetCloseDimension () const
 
const char * GetNameOfClass () const override
 
virtual const ArrayTypeGetSplineOrder () const
 
virtual void SetCloseDimension (ArrayType _arg)
 
void SetInputImage (const InputImageType *) override
 
void SetSplineOrder (const ArrayType &)
 
void SetSplineOrder (const unsigned int)
 
virtual void SetSpacing (SpacingType _arg)
 
virtual SpacingType GetSpacing () const
 
virtual void SetOrigin (OriginType _arg)
 
virtual OriginType GetOrigin () const
 
virtual void SetSize (SizeType _arg)
 
virtual SizeType GetSize () const
 
virtual void SetBSplineEpsilon (RealType _arg)
 
virtual RealType GetBSplineEpsilon () const
 
- Public Member Functions inherited from itk::ImageFunction< TInputImage, TInputImage::PixelType, TCoordinate >
void ConvertContinuousIndexToNearestIndex (const ContinuousIndexType &cindex, IndexType &index) const
 
void ConvertPointToContinuousIndex (const PointType &point, ContinuousIndexType &cindex) const
 
virtual const ContinuousIndexTypeGetEndContinuousIndex () const
 
virtual const IndexTypeGetEndIndex () const
 
const InputImageTypeGetInputImage () const
 
const char * GetNameOfClass () const override
 
virtual const ContinuousIndexTypeGetStartContinuousIndex () const
 
virtual const IndexTypeGetStartIndex () const
 
virtual bool IsInsideBuffer (const IndexType &index) const
 
virtual bool IsInsideBuffer (const ContinuousIndexType &index) const
 
virtual bool IsInsideBuffer (const PointType &point) const
 
void ConvertPointToNearestIndex (const PointType &point, IndexType &index) const
 
- Public Member Functions inherited from itk::Object
unsigned long AddObserver (const EventObject &event, Command *cmd) const
 
unsigned long AddObserver (const EventObject &event, std::function< void(const EventObject &)> function) const
 
LightObject::Pointer CreateAnother () const override
 
virtual void DebugOff () const
 
virtual void DebugOn () const
 
CommandGetCommand (unsigned long tag)
 
bool GetDebug () const
 
MetaDataDictionaryGetMetaDataDictionary ()
 
const MetaDataDictionaryGetMetaDataDictionary () const
 
virtual ModifiedTimeType GetMTime () const
 
virtual const TimeStampGetTimeStamp () const
 
bool HasObserver (const EventObject &event) const
 
void InvokeEvent (const EventObject &)
 
void InvokeEvent (const EventObject &) const
 
virtual void Modified () const
 
void Register () const override
 
void RemoveAllObservers ()
 
void RemoveObserver (unsigned long tag) const
 
void SetDebug (bool debugFlag) const
 
void SetReferenceCount (int) override
 
void UnRegister () const noexcept override
 
void SetMetaDataDictionary (const MetaDataDictionary &rhs)
 
void SetMetaDataDictionary (MetaDataDictionary &&rrhs)
 
virtual void SetObjectName (std::string _arg)
 
virtual const std::string & GetObjectName () const
 
- Public Member Functions inherited from itk::LightObject
Pointer Clone () const
 
virtual void Delete ()
 
virtual int GetReferenceCount () const
 
void Print (std::ostream &os, Indent indent=0) const
 

Static Public Member Functions

static Pointer New ()
 
- Static Public Member Functions inherited from itk::Object
static bool GetGlobalWarningDisplay ()
 
static void GlobalWarningDisplayOff ()
 
static void GlobalWarningDisplayOn ()
 
static Pointer New ()
 
static void SetGlobalWarningDisplay (bool val)
 
- Static Public Member Functions inherited from itk::LightObject
static void BreakOnError ()
 
static Pointer New ()
 

Static Public Attributes

static constexpr unsigned int ImageDimension = TInputImage::ImageDimension
 
- Static Public Attributes inherited from itk::ImageFunction< TInputImage, TInputImage::PixelType, TCoordinate >
static constexpr unsigned int ImageDimension
 

Protected Member Functions

 BSplineControlPointImageFunction ()
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
 ~BSplineControlPointImageFunction () override=default
 
- Protected Member Functions inherited from itk::ImageFunction< TInputImage, TInputImage::PixelType, TCoordinate >
 ImageFunction ()
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
 ~ImageFunction () override=default
 
- Protected Member Functions inherited from itk::FunctionBase< Point< TCoordinate, TInputImage::ImageDimension >, TInputImage::PixelType >
 FunctionBase ()=default
 
 ~FunctionBase () override=default
 
- Protected Member Functions inherited from itk::Object
 Object ()
 
bool PrintObservers (std::ostream &os, Indent indent) const
 
virtual void SetTimeStamp (const TimeStamp &timeStamp)
 
 ~Object () override
 
- Protected Member Functions inherited from itk::LightObject
virtual LightObject::Pointer InternalClone () const
 
 LightObject ()
 
virtual void PrintHeader (std::ostream &os, Indent indent) const
 
virtual void PrintTrailer (std::ostream &os, Indent indent) const
 
virtual ~LightObject ()
 

Private Attributes

CoordinateType m_BSplineEpsilon {}
 
ArrayType m_CloseDimension {}
 
KernelType::Pointer m_Kernel [ImageDimension] {}
 
KernelOrder0Type::Pointer m_KernelOrder0 {}
 
KernelOrder1Type::Pointer m_KernelOrder1 {}
 
KernelOrder2Type::Pointer m_KernelOrder2 {}
 
KernelOrder3Type::Pointer m_KernelOrder3 {}
 
RealImagePointer m_NeighborhoodWeightImage {}
 
ArrayType m_NumberOfControlPoints {}
 
OriginType m_Origin {}
 
SizeType m_Size {}
 
SpacingType m_Spacing {}
 
ArrayType m_SplineOrder {}
 

Additional Inherited Members

- Protected Attributes inherited from itk::ImageFunction< TInputImage, TInputImage::PixelType, TCoordinate >
ContinuousIndexType m_EndContinuousIndex
 
IndexType m_EndIndex
 
InputImageConstPointer m_Image
 
ContinuousIndexType m_StartContinuousIndex
 
IndexType m_StartIndex
 
- Protected Attributes inherited from itk::LightObject
std::atomic< int > m_ReferenceCount {}
 

Member Typedef Documentation

◆ ArrayType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::ArrayType = FixedArray<unsigned int, ImageDimension>

Other type alias

Definition at line 102 of file itkBSplineControlPointImageFunction.h.

◆ ConstPointer

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::ConstPointer = SmartPointer<const Self>

Definition at line 67 of file itkBSplineControlPointImageFunction.h.

◆ ControlPointLatticeType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::ControlPointLatticeType = TInputImage

Image type alias support

Definition at line 79 of file itkBSplineControlPointImageFunction.h.

◆ CoordinateType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::CoordinateType = TCoordinate

Definition at line 81 of file itkBSplineControlPointImageFunction.h.

◆ GradientType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::GradientType = VariableSizeMatrix<CoordinateType>

Definition at line 98 of file itkBSplineControlPointImageFunction.h.

◆ HessianComponentType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::HessianComponentType = VariableSizeMatrix<CoordinateType>

Definition at line 99 of file itkBSplineControlPointImageFunction.h.

◆ IndexType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::IndexType = typename InputImageType::IndexType

Definition at line 88 of file itkBSplineControlPointImageFunction.h.

◆ InputImageRegionType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::InputImageRegionType = typename InputImageType::RegionType

Definition at line 90 of file itkBSplineControlPointImageFunction.h.

◆ InputImageType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::InputImageType = TInputImage

Definition at line 80 of file itkBSplineControlPointImageFunction.h.

◆ KernelOrder0Type

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::KernelOrder0Type = BSplineKernelFunction<0>

Definition at line 110 of file itkBSplineControlPointImageFunction.h.

◆ KernelOrder1Type

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::KernelOrder1Type = BSplineKernelFunction<1>

Definition at line 111 of file itkBSplineControlPointImageFunction.h.

◆ KernelOrder2Type

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::KernelOrder2Type = BSplineKernelFunction<2>

Definition at line 112 of file itkBSplineControlPointImageFunction.h.

◆ KernelOrder3Type

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::KernelOrder3Type = BSplineKernelFunction<3>

Definition at line 113 of file itkBSplineControlPointImageFunction.h.

◆ KernelType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::KernelType = CoxDeBoorBSplineKernelFunction<3>

Interpolation kernel type (default spline order = 3)

Definition at line 109 of file itkBSplineControlPointImageFunction.h.

◆ OriginType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::OriginType = typename InputImageType::PointType

Definition at line 93 of file itkBSplineControlPointImageFunction.h.

◆ OutputType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::OutputType = PixelType

Output type alias support

Definition at line 97 of file itkBSplineControlPointImageFunction.h.

◆ PixelType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::PixelType = typename InputImageType::PixelType

Definition at line 86 of file itkBSplineControlPointImageFunction.h.

◆ Pointer

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::Pointer = SmartPointer<Self>

Definition at line 66 of file itkBSplineControlPointImageFunction.h.

◆ RealImagePointer

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::RealImagePointer = typename RealImageType::Pointer

Definition at line 104 of file itkBSplineControlPointImageFunction.h.

◆ RealImageType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::RealImageType = Image<CoordinateType, ImageDimension>

Definition at line 103 of file itkBSplineControlPointImageFunction.h.

◆ RealType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::RealType = float

Definition at line 106 of file itkBSplineControlPointImageFunction.h.

◆ RegionType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::RegionType = typename InputImageType::RegionType

Definition at line 87 of file itkBSplineControlPointImageFunction.h.

◆ Self

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::Self = BSplineControlPointImageFunction

Definition at line 64 of file itkBSplineControlPointImageFunction.h.

◆ SizeType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::SizeType = typename InputImageType::SizeType

Definition at line 94 of file itkBSplineControlPointImageFunction.h.

◆ SpacingType

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::SpacingType = typename InputImageType::SpacingType

Definition at line 92 of file itkBSplineControlPointImageFunction.h.

◆ Superclass

template<typename TInputImage , typename TCoordinate = double>
using itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::Superclass = ImageFunction<TInputImage, typename TInputImage::PixelType, TCoordinate>

Definition at line 65 of file itkBSplineControlPointImageFunction.h.

Constructor & Destructor Documentation

◆ BSplineControlPointImageFunction()

template<typename TInputImage , typename TCoordinate = double>
itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::BSplineControlPointImageFunction ( )
protected

◆ ~BSplineControlPointImageFunction()

template<typename TInputImage , typename TCoordinate = double>
itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::~BSplineControlPointImageFunction ( )
overrideprotecteddefault

Member Function Documentation

◆ Evaluate()

template<typename TInputImage , typename TCoordinate = double>
OutputType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::Evaluate ( const PointType ) const
overridevirtual

Evaluate the resulting B-spline object at a specified internal parametric point. Note that the internal parameterization over each dimension of the B-spline object is [0, 1).

Implements itk::ImageFunction< TInputImage, TInputImage::PixelType, TCoordinate >.

◆ EvaluateAtContinuousIndex()

template<typename TInputImage , typename TCoordinate = double>
OutputType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::EvaluateAtContinuousIndex ( const ContinuousIndexType ) const
overridevirtual

Evaluate the resulting B-spline object at a specified continuous index in the parametric domain.

Implements itk::ImageFunction< TInputImage, TInputImage::PixelType, TCoordinate >.

◆ EvaluateAtIndex()

template<typename TInputImage , typename TCoordinate = double>
OutputType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::EvaluateAtIndex ( const IndexType ) const
overridevirtual

Evaluate the resulting B-spline object at a specified index in the parametric domain.

Implements itk::ImageFunction< TInputImage, TInputImage::PixelType, TCoordinate >.

◆ EvaluateAtParametricPoint()

template<typename TInputImage , typename TCoordinate = double>
OutputType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::EvaluateAtParametricPoint ( const PointType ) const

Evaluate the resulting B-spline object at a specified point in the parametric domain.

◆ EvaluateGradient()

template<typename TInputImage , typename TCoordinate = double>
GradientType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::EvaluateGradient ( const PointType ) const

Evaluate the gradient of the resulting B-spline object at a specified internal parametric point. Note that the internal parameterization over each dimension of the B-spline object is [0, 1).

◆ EvaluateGradientAtContinuousIndex()

template<typename TInputImage , typename TCoordinate = double>
GradientType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::EvaluateGradientAtContinuousIndex ( const ContinuousIndexType ) const

Evaluate the gradient of the resulting B-spline object at a specified continuous index in the parametric domain.

◆ EvaluateGradientAtIndex()

template<typename TInputImage , typename TCoordinate = double>
GradientType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::EvaluateGradientAtIndex ( const IndexType ) const

Evaluate the gradient of the resulting B-spline object at a specified index in the parametric domain.

◆ EvaluateGradientAtParametricPoint()

template<typename TInputImage , typename TCoordinate = double>
GradientType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::EvaluateGradientAtParametricPoint ( const PointType ) const

Evaluate the gradient of the resulting B-spline object at a specified point in the parametric domain.

◆ EvaluateHessian()

template<typename TInputImage , typename TCoordinate = double>
HessianComponentType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::EvaluateHessian ( const PointType ,
const unsigned int   
) const

Evaluate the hessian of the resulting B-spline object at a specified internal parametric point. Note that the internal parameterization over each dimension of the B-spline object is [0, 1).

◆ EvaluateHessianAtContinuousIndex()

template<typename TInputImage , typename TCoordinate = double>
HessianComponentType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::EvaluateHessianAtContinuousIndex ( const ContinuousIndexType ,
const unsigned int   
) const

Evaluate the Hessian of the resulting B-spline object at a specified con- tinuous index within the parametric domain. Since the Hessian for a vector function is a 3-tensor, one must specify the component.

◆ EvaluateHessianAtIndex()

template<typename TInputImage , typename TCoordinate = double>
HessianComponentType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::EvaluateHessianAtIndex ( const IndexType ,
const unsigned int   
) const

Evaluate the Hessian of the resulting B-spline object at a specified index within the parametric domain. Since the Hessian for a vector function is a 3-tensor, one must specify the component.

◆ EvaluateHessianAtParametricPoint()

template<typename TInputImage , typename TCoordinate = double>
HessianComponentType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::EvaluateHessianAtParametricPoint ( const PointType ,
const unsigned int   
) const

Evaluate the Hessian of the resulting B-spline object at a specified point within the parametric domain. Since the Hessian for a vector function is a 3-tensor, one must specify the component.

◆ GetBSplineEpsilon()

template<typename TInputImage , typename TCoordinate = double>
virtual RealType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::GetBSplineEpsilon ( ) const
virtual

Set/Get the epsilon used for B-splines. The B-spline parametric domain in 1-D is defined on the half-closed interval [a,b). Extension to n-D is defined similarly. This presents some difficulty for defining the the image domain to be co-extensive with the parametric domain. We use the B-spline epsilon to push the edge of the image boundary inside the B-spline parametric domain.

◆ GetCloseDimension()

template<typename TInputImage , typename TCoordinate = double>
virtual const ArrayType& itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::GetCloseDimension ( ) const
virtual

Get the boolean array indicating which dimensions are closed.

◆ GetNameOfClass()

template<typename TInputImage , typename TCoordinate = double>
const char* itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::GetNameOfClass ( ) const
overridevirtual
See also
LightObject::GetNameOfClass()

Reimplemented from itk::Object.

◆ GetOrigin()

template<typename TInputImage , typename TCoordinate = double>
virtual OriginType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::GetOrigin ( ) const
virtual

Set/Get the parametric origin of the B-spline object domain.

◆ GetSize()

template<typename TInputImage , typename TCoordinate = double>
virtual SizeType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::GetSize ( ) const
virtual

Set/Get the parametric size of the B-spline object domain.

◆ GetSpacing()

template<typename TInputImage , typename TCoordinate = double>
virtual SpacingType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::GetSpacing ( ) const
virtual

Set/Get the parametric spacing of the B-spline object domain.

◆ GetSplineOrder()

template<typename TInputImage , typename TCoordinate = double>
virtual const ArrayType& itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::GetSplineOrder ( ) const
virtual

Get the spline order array of the B-spline object. Default = 3.

◆ New()

template<typename TInputImage , typename TCoordinate = double>
static Pointer itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::New ( )
static

Method for creation through the object factory.

◆ PrintSelf()

template<typename TInputImage , typename TCoordinate = double>
void itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::PrintSelf ( std::ostream &  os,
Indent  indent 
) const
overrideprotectedvirtual

Methods invoked by Print() to print information about the object including superclasses. Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.

Reimplemented from itk::Object.

◆ SetBSplineEpsilon()

template<typename TInputImage , typename TCoordinate = double>
virtual void itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::SetBSplineEpsilon ( RealType  _arg)
virtual

Set/Get the epsilon used for B-splines. The B-spline parametric domain in 1-D is defined on the half-closed interval [a,b). Extension to n-D is defined similarly. This presents some difficulty for defining the the image domain to be co-extensive with the parametric domain. We use the B-spline epsilon to push the edge of the image boundary inside the B-spline parametric domain.

◆ SetCloseDimension()

template<typename TInputImage , typename TCoordinate = double>
virtual void itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::SetCloseDimension ( ArrayType  _arg)
virtual

Set the boolean array indicating the periodicity of the B-spline object. This array of 0/1 values defines whether a particular dimension of the parametric space is to be considered periodic or not. For example, if you are using interpolating along a 1D closed curve, the array type will have size 1, and you should set the first element of this array to the value "1". In the case that you were interpolating in a planar surface with cylindrical topology, the array type will have two components, and you should set to "1" the component that goes around the cylinder, and set to "0" the component that goes from the top of the cylinder to the bottom. This will indicate the periodicity of that parameter to the filter. Internally, in order to make periodic the domain of the parameter, the filter will reuse some of the points at the beginning of the domain as if they were also located at the end of the domain. The number of points to be reused will depend on the spline order. As a user, you don't need to replicate the points, the filter will do this for you.

◆ SetInputImage()

template<typename TInputImage , typename TCoordinate = double>
void itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::SetInputImage ( const InputImageType )
overridevirtual

Set the input image. Note that the size, spacing, origin, and spline order must be called prior to setting the input image.

Reimplemented from itk::ImageFunction< TInputImage, TInputImage::PixelType, TCoordinate >.

◆ SetOrigin()

template<typename TInputImage , typename TCoordinate = double>
virtual void itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::SetOrigin ( OriginType  _arg)
virtual

Set/Get the parametric origin of the B-spline object domain.

◆ SetSize()

template<typename TInputImage , typename TCoordinate = double>
virtual void itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::SetSize ( SizeType  _arg)
virtual

Set/Get the parametric size of the B-spline object domain.

◆ SetSpacing()

template<typename TInputImage , typename TCoordinate = double>
virtual void itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::SetSpacing ( SpacingType  _arg)
virtual

Set/Get the parametric spacing of the B-spline object domain.

◆ SetSplineOrder() [1/2]

template<typename TInputImage , typename TCoordinate = double>
void itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::SetSplineOrder ( const ArrayType )

Set the spline order array where each element of the array corresponds to a single parametric dimension of the B-spline object. Default = 3.

◆ SetSplineOrder() [2/2]

template<typename TInputImage , typename TCoordinate = double>
void itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::SetSplineOrder ( const unsigned int  )

Set the spline order of the B-spline object for all parametric dimensions. Default = 3.

Member Data Documentation

◆ ImageDimension

template<typename TInputImage , typename TCoordinate = double>
constexpr unsigned int itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::ImageDimension = TInputImage::ImageDimension
staticconstexpr

Extract dimension from input image.

Definition at line 76 of file itkBSplineControlPointImageFunction.h.

◆ m_BSplineEpsilon

template<typename TInputImage , typename TCoordinate = double>
CoordinateType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_BSplineEpsilon {}
private

Definition at line 311 of file itkBSplineControlPointImageFunction.h.

◆ m_CloseDimension

template<typename TInputImage , typename TCoordinate = double>
ArrayType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_CloseDimension {}
private

Definition at line 300 of file itkBSplineControlPointImageFunction.h.

◆ m_Kernel

template<typename TInputImage , typename TCoordinate = double>
KernelType::Pointer itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_Kernel[ImageDimension] {}
private

Definition at line 305 of file itkBSplineControlPointImageFunction.h.

◆ m_KernelOrder0

template<typename TInputImage , typename TCoordinate = double>
KernelOrder0Type::Pointer itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_KernelOrder0 {}
private

Definition at line 306 of file itkBSplineControlPointImageFunction.h.

◆ m_KernelOrder1

template<typename TInputImage , typename TCoordinate = double>
KernelOrder1Type::Pointer itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_KernelOrder1 {}
private

Definition at line 307 of file itkBSplineControlPointImageFunction.h.

◆ m_KernelOrder2

template<typename TInputImage , typename TCoordinate = double>
KernelOrder2Type::Pointer itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_KernelOrder2 {}
private

Definition at line 308 of file itkBSplineControlPointImageFunction.h.

◆ m_KernelOrder3

template<typename TInputImage , typename TCoordinate = double>
KernelOrder3Type::Pointer itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_KernelOrder3 {}
private

Definition at line 309 of file itkBSplineControlPointImageFunction.h.

◆ m_NeighborhoodWeightImage

template<typename TInputImage , typename TCoordinate = double>
RealImagePointer itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_NeighborhoodWeightImage {}
private

Definition at line 303 of file itkBSplineControlPointImageFunction.h.

◆ m_NumberOfControlPoints

template<typename TInputImage , typename TCoordinate = double>
ArrayType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_NumberOfControlPoints {}
private

Definition at line 299 of file itkBSplineControlPointImageFunction.h.

◆ m_Origin

template<typename TInputImage , typename TCoordinate = double>
OriginType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_Origin {}
private

Definition at line 297 of file itkBSplineControlPointImageFunction.h.

◆ m_Size

template<typename TInputImage , typename TCoordinate = double>
SizeType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_Size {}
private

Parameters for the B-spline object domain

Definition at line 295 of file itkBSplineControlPointImageFunction.h.

◆ m_Spacing

template<typename TInputImage , typename TCoordinate = double>
SpacingType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_Spacing {}
private

Definition at line 296 of file itkBSplineControlPointImageFunction.h.

◆ m_SplineOrder

template<typename TInputImage , typename TCoordinate = double>
ArrayType itk::BSplineControlPointImageFunction< TInputImage, TCoordinate >::m_SplineOrder {}
private

Definition at line 301 of file itkBSplineControlPointImageFunction.h.


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