Skip to content

itk::SparseKernelTransform

More...

#include <itkSparseKernelTransform.h>

Inherits from Transform< TScalarType, NDimensions, NDimensions >

Public Types

Name
typedef SparseKernelTransform Self
typedef Transform< TScalarType, NDimensions, NDimensions > Superclass
typedef SmartPointer< Self > Pointer
typedef SmartPointer< const Self > ConstPointer
typedef Superclass::ScalarType ScalarType
typedef Superclass::ParametersType ParametersType
typedef Superclass::JacobianType JacobianType
typedef Superclass::InputPointType InputPointType
typedef Superclass::OutputPointType OutputPointType
typedef Superclass::InputVectorType InputVectorType
typedef Superclass::OutputVectorType OutputVectorType
typedef DefaultStaticMeshTraits< TScalarType, NDimensions, NDimensions, TScalarType, TScalarType > PointSetTraitsType
typedef PointSet< InputPointType, NDimensions, PointSetTraitsType > PointSetType
typedef PointSetType::Pointer PointSetPointer
typedef PointSetType::PointsContainer PointsContainer
typedef PointSetType::PointsContainerIterator PointsIterator
typedef PointSetType::PointsContainerConstIterator PointsConstIterator
typedef itk::VectorContainer< unsigned long, InputVectorType > VectorSetType
typedef VectorSetType::Pointer VectorSetPointer
typedef Eigen::Matrix< TScalarType, NDimensions, NDimensions > IMatrixType
typedef Eigen::Triplet< TScalarType > TripletType
typedef Eigen::Matrix< TScalarType, NDimensions, NDimensions > GMatrixType
typedef Eigen::SparseMatrix< TScalarType > LMatrixType
typedef Eigen::SparseMatrix< TScalarType > KMatrixType
typedef Eigen::SparseMatrix< TScalarType > PMatrixType
typedef Eigen::Matrix< TScalarType, Eigen::Dynamic, Eigen::Dynamic > YMatrixType
typedef Eigen::Matrix< TScalarType, Eigen::Dynamic, Eigen::Dynamic > WMatrixType
typedef Eigen::Matrix< TScalarType, Eigen::Dynamic, Eigen::Dynamic > DMatrixType
typedef Eigen::Matrix< TScalarType, NDimensions, NDimensions > AMatrixType
typedef Eigen::Matrix< TScalarType, NDimensions, 1 > BMatrixType
typedef Eigen::Matrix< TScalarType, 1, NDimensions > RowMatrixType
typedef Eigen::Matrix< TScalarType, NDimensions, 1 > ColumnMatrixType

Public Functions

Name
itkTypeMacro(SparseKernelTransform , Transform )
itkNewMacro(Self )
itkStaticConstMacro(SpaceDimension , unsigned int , NDimensions )
itkGetObjectMacro(SourceLandmarks , PointSetType )
virtual void SetSourceLandmarks(PointSetType * )
itkGetObjectMacro(TargetLandmarks , PointSetType )
virtual void SetTargetLandmarks(PointSetType * )
itkGetObjectMacro(Displacements , VectorSetType )
void ComputeWMatrix(void ) const
virtual OutputPointType TransformPoint(const InputPointType & thisPoint) const
virtual void SetIdentity()
virtual void SetParameters(const ParametersType & )
virtual void SetFixedParameters(const ParametersType & )
virtual void UpdateParameters(void ) const
virtual const ParametersType & GetParameters(void ) const
virtual const ParametersType & GetFixedParameters(void ) const
virtual void ComputeJacobianWithRespectToParameters(const InputPointType & in, JacobianType & jacobian) const
virtual void SetStiffness(double stiffness)
itkGetMacro(Stiffness , double )

Protected Functions

Name
SparseKernelTransform()
virtual ~SparseKernelTransform()
void PrintSelf(std::ostream & os, Indent indent) const
virtual const GMatrixType & ComputeG(const InputVectorType & landmarkVector) const
virtual const GMatrixType & ComputeReflexiveG(PointsIterator ) const
virtual void ComputeDeformationContribution(const InputPointType & inputPoint, OutputPointType & result) const
void ComputeK() const
void ComputeL() const
void ComputeP() const
void ComputeY() const
void ComputeD() const
void ReorganizeW(void ) const

Public Attributes

Name
PointSetPointer m_SourceLandmarks
PointSetPointer m_TargetLandmarks

Protected Attributes

Name
double m_Stiffness
VectorSetPointer m_Displacements
LMatrixType m_LMatrix
LMatrixType m_LMatrixInverse
KMatrixType m_KMatrix
PMatrixType m_PMatrix
YMatrixType m_YMatrix
WMatrixType m_WMatrix
DMatrixType m_DMatrix
AMatrixType m_AMatrix
BMatrixType m_BVector
GMatrixType m_GMatrix
bool m_WMatrixComputed
bool m_LMatrixComputed
bool m_LInverseComputed
IMatrixType m_I

Detailed Description

cpp template <class TScalarType , unsigned int NDimensions> class itk::SparseKernelTransform;

Intended to be a base class for elastic body spline and thin plate spline. This is implemented in as straightforward a manner as possible from the IEEE TMI paper by Davis, Khotanzad, Flamig, and Harms, Vol. 16, No. 3 June 1997. Notation closely follows their paper, so if you have it in front of you, this code will make a lot more sense.

SparseKernelTransform: Provides support for defining source and target landmarks Defines a number of data types used in the computations Defines the mathematical framework used to compute all splines, so that subclasses need only provide a kernel specific to that spline

This formulation allows the stiffness of the spline to be adjusted, allowing the spline to vary from interpolating the landmarks to approximating the landmarks. This part of the formulation is based on the short paper by R. Sprengel, K. Rohr, H. Stiehl. "Thin-Plate Spline Approximation for Image Registration". In 18th International Conference of the IEEE Engineering in Medicine and Biology Society. 1996.

This class was modified to support its use in the ITK registration framework by Rupert Brooks, McGill Centre for Intelligent Machines, Montreal, Canada March 2007. See the Insight Journal Paper by Brooks, R., Arbel, T. "Improvements to the itk::KernelTransform and its subclasses."

Public Types Documentation

typedef Self

cpp typedef SparseKernelTransform itk::SparseKernelTransform< TScalarType, NDimensions >::Self;

Standard class typedefs.

typedef Superclass

cpp typedef Transform<TScalarType, NDimensions, NDimensions > itk::SparseKernelTransform< TScalarType, NDimensions >::Superclass;

typedef Pointer

cpp typedef SmartPointer<Self> itk::SparseKernelTransform< TScalarType, NDimensions >::Pointer;

typedef ConstPointer

cpp typedef SmartPointer<const Self> itk::SparseKernelTransform< TScalarType, NDimensions >::ConstPointer;

typedef ScalarType

cpp typedef Superclass::ScalarType itk::SparseKernelTransform< TScalarType, NDimensions >::ScalarType;

Scalar type.

typedef ParametersType

cpp typedef Superclass::ParametersType itk::SparseKernelTransform< TScalarType, NDimensions >::ParametersType;

Parameters type.

typedef JacobianType

cpp typedef Superclass::JacobianType itk::SparseKernelTransform< TScalarType, NDimensions >::JacobianType;

Jacobian type.

typedef InputPointType

cpp typedef Superclass::InputPointType itk::SparseKernelTransform< TScalarType, NDimensions >::InputPointType;

Standard coordinate point type for this class.

typedef OutputPointType

cpp typedef Superclass::OutputPointType itk::SparseKernelTransform< TScalarType, NDimensions >::OutputPointType;

typedef InputVectorType

cpp typedef Superclass::InputVectorType itk::SparseKernelTransform< TScalarType, NDimensions >::InputVectorType;

Standard vector type for this class.

typedef OutputVectorType

cpp typedef Superclass::OutputVectorType itk::SparseKernelTransform< TScalarType, NDimensions >::OutputVectorType;

typedef PointSetTraitsType

cpp typedef DefaultStaticMeshTraits<TScalarType, NDimensions, NDimensions, TScalarType, TScalarType> itk::SparseKernelTransform< TScalarType, NDimensions >::PointSetTraitsType;

PointList typedef. This type is used for maintaining lists of points, specifically, the source and target landmark lists.

typedef PointSetType

cpp typedef PointSet<InputPointType, NDimensions, PointSetTraitsType> itk::SparseKernelTransform< TScalarType, NDimensions >::PointSetType;

typedef PointSetPointer

cpp typedef PointSetType::Pointer itk::SparseKernelTransform< TScalarType, NDimensions >::PointSetPointer;

typedef PointsContainer

cpp typedef PointSetType::PointsContainer itk::SparseKernelTransform< TScalarType, NDimensions >::PointsContainer;

typedef PointsIterator

cpp typedef PointSetType::PointsContainerIterator itk::SparseKernelTransform< TScalarType, NDimensions >::PointsIterator;

typedef PointsConstIterator

cpp typedef PointSetType::PointsContainerConstIterator itk::SparseKernelTransform< TScalarType, NDimensions >::PointsConstIterator;

typedef VectorSetType

cpp typedef itk::VectorContainer<unsigned long,InputVectorType> itk::SparseKernelTransform< TScalarType, NDimensions >::VectorSetType;

VectorSet typedef.

typedef VectorSetPointer

cpp typedef VectorSetType::Pointer itk::SparseKernelTransform< TScalarType, NDimensions >::VectorSetPointer;

typedef IMatrixType

cpp typedef Eigen::Matrix<TScalarType, NDimensions, NDimensions> itk::SparseKernelTransform< TScalarType, NDimensions >::IMatrixType;

'I' (identity) matrix typedef.

typedef TripletType

cpp typedef Eigen::Triplet<TScalarType> itk::SparseKernelTransform< TScalarType, NDimensions >::TripletType;

triplets used to fill sparse matrices.

typedef GMatrixType

cpp typedef Eigen::Matrix<TScalarType, NDimensions, NDimensions> itk::SparseKernelTransform< TScalarType, NDimensions >::GMatrixType;

'G' matrix typedef.

typedef LMatrixType

cpp typedef Eigen::SparseMatrix<TScalarType> itk::SparseKernelTransform< TScalarType, NDimensions >::LMatrixType;

'L' matrix typedef.

typedef KMatrixType

cpp typedef Eigen::SparseMatrix<TScalarType> itk::SparseKernelTransform< TScalarType, NDimensions >::KMatrixType;

'K' matrix typedef.

typedef PMatrixType

cpp typedef Eigen::SparseMatrix<TScalarType> itk::SparseKernelTransform< TScalarType, NDimensions >::PMatrixType;

'P' matrix typedef.

typedef YMatrixType

cpp typedef Eigen::Matrix<TScalarType, Eigen::Dynamic, Eigen::Dynamic> itk::SparseKernelTransform< TScalarType, NDimensions >::YMatrixType;

'Y' matrix typedef.

typedef WMatrixType

cpp typedef Eigen::Matrix<TScalarType, Eigen::Dynamic, Eigen::Dynamic> itk::SparseKernelTransform< TScalarType, NDimensions >::WMatrixType;

'W' matrix typedef.

typedef DMatrixType

cpp typedef Eigen::Matrix<TScalarType, Eigen::Dynamic, Eigen::Dynamic> itk::SparseKernelTransform< TScalarType, NDimensions >::DMatrixType;

'D' matrix typedef. Deformation component

typedef AMatrixType

cpp typedef Eigen::Matrix<TScalarType, NDimensions, NDimensions> itk::SparseKernelTransform< TScalarType, NDimensions >::AMatrixType;

'A' matrix typedef. Rotational part of the Affine component

typedef BMatrixType

cpp typedef Eigen::Matrix<TScalarType,NDimensions,1> itk::SparseKernelTransform< TScalarType, NDimensions >::BMatrixType;

'B' matrix typedef. Translational part of the Affine component

typedef RowMatrixType

cpp typedef Eigen::Matrix<TScalarType,1,NDimensions> itk::SparseKernelTransform< TScalarType, NDimensions >::RowMatrixType;

Row matrix typedef.

typedef ColumnMatrixType

cpp typedef Eigen::Matrix<TScalarType,NDimensions,1> itk::SparseKernelTransform< TScalarType, NDimensions >::ColumnMatrixType;

Column matrix typedef.

Public Functions Documentation

function itkTypeMacro

cpp itkTypeMacro( SparseKernelTransform , Transform )

Run-time type information (and related methods).

function itkNewMacro

cpp itkNewMacro( Self )

New macro for creation of through a Smart Pointer

function itkStaticConstMacro

cpp itkStaticConstMacro( SpaceDimension , unsigned int , NDimensions )

Dimension of the domain space.

function itkGetObjectMacro

cpp itkGetObjectMacro( SourceLandmarks , PointSetType )

Get the source landmarks list, which we will denote \( p \).

function SetSourceLandmarks

cpp virtual void SetSourceLandmarks( PointSetType * )

Set the source landmarks list.

function itkGetObjectMacro

cpp itkGetObjectMacro( TargetLandmarks , PointSetType )

Get the target landmarks list, which we will denote \( q \).

function SetTargetLandmarks

cpp virtual void SetTargetLandmarks( PointSetType * )

Set the target landmarks list.

function itkGetObjectMacro

cpp itkGetObjectMacro( Displacements , VectorSetType )

Get the displacements list, which we will denote \( d \), where \( d_i = q_i - p_i \).

function ComputeWMatrix

cpp void ComputeWMatrix( void ) const

Compute W matrix.

function TransformPoint

cpp virtual OutputPointType TransformPoint( const InputPointType & thisPoint ) const

Compute L matrix inverse. Compute the position of point in the new space

function SetIdentity

cpp virtual void SetIdentity()

Compute the Jacobian Matrix of the transformation at one point Set the Transformation Parameters to be an identity transform

function SetParameters

cpp virtual void SetParameters( const ParametersType & )

Set the Transformation Parameters and update the internal transformation. The parameters represent the source landmarks. Each landmark point is represented by NDimensions doubles. All the landmarks are concatenated to form one flat Array.

function SetFixedParameters

cpp virtual void SetFixedParameters( const ParametersType & )

Set Transform Fixed Parameters: To support the transform file writer this function was added to set the target landmarks similar to the SetParameters function setting the source landmarks

function UpdateParameters

cpp virtual void UpdateParameters( void ) const

Update the Parameters array from the landmarks corrdinates.

function GetParameters

cpp virtual const ParametersType & GetParameters( void ) const

Get the Transformation Parameters - Gets the Source Landmarks

function GetFixedParameters

cpp virtual const ParametersType & GetFixedParameters( void ) const

Get Transform Fixed Parameters - Gets the Target Landmarks

function ComputeJacobianWithRespectToParameters

cpp virtual void ComputeJacobianWithRespectToParameters( const InputPointType & in, JacobianType & jacobian ) const

Reimplemented by: itk::CompactlySupportedRBFSparseKernelTransform::ComputeJacobianWithRespectToParameters

function SetStiffness

cpp inline virtual void SetStiffness( double stiffness )

Stiffness of the spline. A stiffness of zero results in the standard interpolating spline. A non-zero stiffness allows the spline to approximate rather than interpolate the landmarks. Stiffness values are usually rather small, typically in the range of 0.001 to 0.1. The approximating spline formulation is based on the short paper by R. Sprengel, K. Rohr, H. Stiehl. "Thin-Plate Spline Approximation for Image Registration". In 18th International Conference of the IEEE Engineering in Medicine and Biology Society. 1996.

function itkGetMacro

cpp itkGetMacro( Stiffness , double )

Protected Functions Documentation

function SparseKernelTransform

cpp SparseKernelTransform()

function ~SparseKernelTransform

cpp virtual ~SparseKernelTransform()

function PrintSelf

cpp void PrintSelf( std::ostream & os, Indent indent ) const

function ComputeG

cpp virtual const GMatrixType & ComputeG( const InputVectorType & landmarkVector ) const

Reimplemented by: itk::CompactlySupportedRBFSparseKernelTransform::ComputeG

Compute G(x) This is essentially the kernel of the transform. By overriding this method, we can obtain (among others): Elastic body spline Thin plate spline Volume spline

function ComputeReflexiveG

cpp virtual const GMatrixType & ComputeReflexiveG( PointsIterator ) const

Compute a G(x) for a point to itself (i.e. for the block diagonal elements of the matrix K. Parameter indicates for which landmark the reflexive G is to be computed. The default implementation for the reflexive contribution is a diagonal matrix where the diagonal elements are the stiffness of the spline.

function ComputeDeformationContribution

cpp virtual void ComputeDeformationContribution( const InputPointType & inputPoint, OutputPointType & result ) const

Reimplemented by: itk::CompactlySupportedRBFSparseKernelTransform::ComputeDeformationContribution

Compute the contribution of the landmarks weighted by the kernel funcion to the global deformation of the space

function ComputeK

cpp void ComputeK() const

Compute K matrix.

function ComputeL

cpp void ComputeL() const

Compute L matrix.

function ComputeP

cpp void ComputeP() const

Compute P matrix.

function ComputeY

cpp void ComputeY() const

Compute Y matrix.

function ComputeD

cpp void ComputeD() const

Compute displacements \( q_i - p_i \).

function ReorganizeW

cpp void ReorganizeW( void ) const

Warning: This method release the memory of the W Matrix

Reorganize the components of W into D (deformable), A (rotation part of affine) and B (translational part of affine ) components.

Public Attributes Documentation

variable m_SourceLandmarks

cpp PointSetPointer m_SourceLandmarks;

The list of source landmarks, denoted 'p'.

variable m_TargetLandmarks

cpp PointSetPointer m_TargetLandmarks;

The list of target landmarks, denoted 'q'.

Protected Attributes Documentation

variable m_Stiffness

cpp double m_Stiffness;

Stiffness parameter

variable m_Displacements

cpp VectorSetPointer m_Displacements;

The list of displacements. d[i] = q[i] - p[i];

variable m_LMatrix

cpp LMatrixType m_LMatrix;

The L matrix.

variable m_LMatrixInverse

cpp LMatrixType m_LMatrixInverse;

The inverse of L, which we also cache.

variable m_KMatrix

cpp KMatrixType m_KMatrix;

The K matrix.

variable m_PMatrix

cpp PMatrixType m_PMatrix;

The P matrix.

variable m_YMatrix

cpp YMatrixType m_YMatrix;

The Y matrix.

variable m_WMatrix

cpp WMatrixType m_WMatrix;

The W matrix.

variable m_DMatrix

cpp DMatrixType m_DMatrix;

The Deformation matrix. This is an auxiliary matrix that will hold the Deformation (non-affine) part of the transform. Those are the coefficients that will multiply the Kernel function

variable m_AMatrix

cpp AMatrixType m_AMatrix;

Rotatinoal/Shearing part of the Affine component of the Transformation

variable m_BVector

cpp BMatrixType m_BVector;

Translational part of the Affine component of the Transformation

variable m_GMatrix

cpp GMatrixType m_GMatrix;

The G matrix. It is made mutable because m_GMatrix was made an ivar only to avoid copying the matrix at return time

variable m_WMatrixComputed

cpp bool m_WMatrixComputed;

Has the W matrix been computed?

variable m_LMatrixComputed

cpp bool m_LMatrixComputed;

Has the L matrix been computed?

variable m_LInverseComputed

cpp bool m_LInverseComputed;

Has the L inverse matrix been computed?

variable m_I

cpp IMatrixType m_I;

Identity matrix.


Updated on 2026-03-31 at 16:02:09 +0000