itk::SparseKernelTransform
#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
Detailed Description
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
typedef SparseKernelTransform itk::SparseKernelTransform< TScalarType, NDimensions >::Self;
Standard class typedefs.
typedef Superclass
typedef Transform<TScalarType, NDimensions, NDimensions > itk::SparseKernelTransform< TScalarType, NDimensions >::Superclass;
typedef Pointer
typedef SmartPointer<Self> itk::SparseKernelTransform< TScalarType, NDimensions >::Pointer;
typedef ConstPointer
typedef SmartPointer<const Self> itk::SparseKernelTransform< TScalarType, NDimensions >::ConstPointer;
typedef ScalarType
typedef Superclass::ScalarType itk::SparseKernelTransform< TScalarType, NDimensions >::ScalarType;
Scalar type.
typedef ParametersType
typedef Superclass::ParametersType itk::SparseKernelTransform< TScalarType, NDimensions >::ParametersType;
Parameters type.
typedef JacobianType
typedef Superclass::JacobianType itk::SparseKernelTransform< TScalarType, NDimensions >::JacobianType;
Jacobian type.
typedef InputPointType
typedef Superclass::InputPointType itk::SparseKernelTransform< TScalarType, NDimensions >::InputPointType;
Standard coordinate point type for this class.
typedef OutputPointType
typedef Superclass::OutputPointType itk::SparseKernelTransform< TScalarType, NDimensions >::OutputPointType;
typedef InputVectorType
typedef Superclass::InputVectorType itk::SparseKernelTransform< TScalarType, NDimensions >::InputVectorType;
Standard vector type for this class.
typedef OutputVectorType
typedef Superclass::OutputVectorType itk::SparseKernelTransform< TScalarType, NDimensions >::OutputVectorType;
typedef PointSetTraitsType
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
typedef PointSet<InputPointType, NDimensions, PointSetTraitsType> itk::SparseKernelTransform< TScalarType, NDimensions >::PointSetType;
typedef PointSetPointer
typedef PointSetType::Pointer itk::SparseKernelTransform< TScalarType, NDimensions >::PointSetPointer;
typedef PointsContainer
typedef PointSetType::PointsContainer itk::SparseKernelTransform< TScalarType, NDimensions >::PointsContainer;
typedef PointsIterator
typedef PointSetType::PointsContainerIterator itk::SparseKernelTransform< TScalarType, NDimensions >::PointsIterator;
typedef PointsConstIterator
typedef PointSetType::PointsContainerConstIterator itk::SparseKernelTransform< TScalarType, NDimensions >::PointsConstIterator;
typedef VectorSetType
typedef itk::VectorContainer<unsigned long,InputVectorType> itk::SparseKernelTransform< TScalarType, NDimensions >::VectorSetType;
VectorSet typedef.
typedef VectorSetPointer
typedef VectorSetType::Pointer itk::SparseKernelTransform< TScalarType, NDimensions >::VectorSetPointer;
typedef IMatrixType
typedef Eigen::Matrix<TScalarType, NDimensions, NDimensions> itk::SparseKernelTransform< TScalarType, NDimensions >::IMatrixType;
'I' (identity) matrix typedef.
typedef TripletType
typedef Eigen::Triplet<TScalarType> itk::SparseKernelTransform< TScalarType, NDimensions >::TripletType;
triplets used to fill sparse matrices.
typedef GMatrixType
typedef Eigen::Matrix<TScalarType, NDimensions, NDimensions> itk::SparseKernelTransform< TScalarType, NDimensions >::GMatrixType;
'G' matrix typedef.
typedef LMatrixType
typedef Eigen::SparseMatrix<TScalarType> itk::SparseKernelTransform< TScalarType, NDimensions >::LMatrixType;
'L' matrix typedef.
typedef KMatrixType
typedef Eigen::SparseMatrix<TScalarType> itk::SparseKernelTransform< TScalarType, NDimensions >::KMatrixType;
'K' matrix typedef.
typedef PMatrixType
typedef Eigen::SparseMatrix<TScalarType> itk::SparseKernelTransform< TScalarType, NDimensions >::PMatrixType;
'P' matrix typedef.
typedef YMatrixType
typedef Eigen::Matrix<TScalarType, Eigen::Dynamic, Eigen::Dynamic> itk::SparseKernelTransform< TScalarType, NDimensions >::YMatrixType;
'Y' matrix typedef.
typedef WMatrixType
typedef Eigen::Matrix<TScalarType, Eigen::Dynamic, Eigen::Dynamic> itk::SparseKernelTransform< TScalarType, NDimensions >::WMatrixType;
'W' matrix typedef.
typedef DMatrixType
typedef Eigen::Matrix<TScalarType, Eigen::Dynamic, Eigen::Dynamic> itk::SparseKernelTransform< TScalarType, NDimensions >::DMatrixType;
'D' matrix typedef. Deformation component
typedef AMatrixType
typedef Eigen::Matrix<TScalarType, NDimensions, NDimensions> itk::SparseKernelTransform< TScalarType, NDimensions >::AMatrixType;
'A' matrix typedef. Rotational part of the Affine component
typedef BMatrixType
typedef Eigen::Matrix<TScalarType,NDimensions,1> itk::SparseKernelTransform< TScalarType, NDimensions >::BMatrixType;
'B' matrix typedef. Translational part of the Affine component
typedef RowMatrixType
typedef Eigen::Matrix<TScalarType,1,NDimensions> itk::SparseKernelTransform< TScalarType, NDimensions >::RowMatrixType;
Row matrix typedef.
typedef ColumnMatrixType
typedef Eigen::Matrix<TScalarType,NDimensions,1> itk::SparseKernelTransform< TScalarType, NDimensions >::ColumnMatrixType;
Column matrix typedef.
Public Functions Documentation
function itkTypeMacro
itkTypeMacro(
SparseKernelTransform ,
Transform
)
Run-time type information (and related methods).
function itkNewMacro
itkNewMacro(
Self
)
New macro for creation of through a Smart Pointer
function itkStaticConstMacro
itkStaticConstMacro(
SpaceDimension ,
unsigned int ,
NDimensions
)
Dimension of the domain space.
function itkGetObjectMacro
itkGetObjectMacro(
SourceLandmarks ,
PointSetType
)
Get the source landmarks list, which we will denote \( p \).
function SetSourceLandmarks
virtual void SetSourceLandmarks(
PointSetType *
)
Set the source landmarks list.
function itkGetObjectMacro
itkGetObjectMacro(
TargetLandmarks ,
PointSetType
)
Get the target landmarks list, which we will denote \( q \).
function SetTargetLandmarks
virtual void SetTargetLandmarks(
PointSetType *
)
Set the target landmarks list.
function itkGetObjectMacro
itkGetObjectMacro(
Displacements ,
VectorSetType
)
Get the displacements list, which we will denote \( d \), where \( d_i = q_i - p_i \).
function ComputeWMatrix
void ComputeWMatrix(
void
) const
Compute W matrix.
function TransformPoint
virtual OutputPointType TransformPoint(
const InputPointType & thisPoint
) const
Compute L matrix inverse. Compute the position of point in the new space
function SetIdentity
virtual void SetIdentity()
Compute the Jacobian Matrix of the transformation at one point Set the Transformation Parameters to be an identity transform
function SetParameters
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
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
virtual void UpdateParameters(
void
) const
Update the Parameters array from the landmarks corrdinates.
function GetParameters
virtual const ParametersType & GetParameters(
void
) const
Get the Transformation Parameters - Gets the Source Landmarks
function GetFixedParameters
virtual const ParametersType & GetFixedParameters(
void
) const
Get Transform Fixed Parameters - Gets the Target Landmarks
function ComputeJacobianWithRespectToParameters
virtual void ComputeJacobianWithRespectToParameters(
const InputPointType & in,
JacobianType & jacobian
) const
Reimplemented by: itk::CompactlySupportedRBFSparseKernelTransform::ComputeJacobianWithRespectToParameters
function SetStiffness
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
itkGetMacro(
Stiffness ,
double
)
Protected Functions Documentation
function SparseKernelTransform
SparseKernelTransform()
function ~SparseKernelTransform
virtual ~SparseKernelTransform()
function PrintSelf
void PrintSelf(
std::ostream & os,
Indent indent
) const
function ComputeG
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
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
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
void ComputeK() const
Compute K matrix.
function ComputeL
void ComputeL() const
Compute L matrix.
function ComputeP
void ComputeP() const
Compute P matrix.
function ComputeY
void ComputeY() const
Compute Y matrix.
function ComputeD
void ComputeD() const
Compute displacements \( q_i - p_i \).
function ReorganizeW
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
PointSetPointer m_SourceLandmarks;
The list of source landmarks, denoted 'p'.
variable m_TargetLandmarks
PointSetPointer m_TargetLandmarks;
The list of target landmarks, denoted 'q'.
Protected Attributes Documentation
variable m_Stiffness
double m_Stiffness;
Stiffness parameter
variable m_Displacements
VectorSetPointer m_Displacements;
The list of displacements. d[i] = q[i] - p[i];
variable m_LMatrix
LMatrixType m_LMatrix;
The L matrix.
variable m_LMatrixInverse
LMatrixType m_LMatrixInverse;
The inverse of L, which we also cache.
variable m_KMatrix
KMatrixType m_KMatrix;
The K matrix.
variable m_PMatrix
PMatrixType m_PMatrix;
The P matrix.
variable m_YMatrix
YMatrixType m_YMatrix;
The Y matrix.
variable m_WMatrix
WMatrixType m_WMatrix;
The W matrix.
variable m_DMatrix
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
AMatrixType m_AMatrix;
Rotatinoal/Shearing part of the Affine component of the Transformation
variable m_BVector
BMatrixType m_BVector;
Translational part of the Affine component of the Transformation
variable m_GMatrix
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
bool m_WMatrixComputed;
Has the W matrix been computed?
variable m_LMatrixComputed
bool m_LMatrixComputed;
Has the L matrix been computed?
variable m_LInverseComputed
bool m_LInverseComputed;
Has the L inverse matrix been computed?
variable m_I
IMatrixType m_I;
Identity matrix.
Updated on 2024-11-11 at 19:51:45 +0000