itk::KernelTransform2
#include <itkKernelTransform2.h>
Inherits from itk::AdvancedTransform< TScalarType, NDimensions, NDimensions >, Transform< TScalarType, 3, 3 >
Public Types
Name | |
---|---|
typedef KernelTransform2 | Self |
typedef AdvancedTransform< TScalarType, NDimensions, NDimensions > | Superclass |
typedef SmartPointer< Self > | Pointer |
typedef SmartPointer< const Self > | ConstPointer |
typedef Superclass::ScalarType | ScalarType |
typedef Superclass::ParametersType | ParametersType |
typedef Superclass::NumberOfParametersType | NumberOfParametersType |
typedef Superclass::JacobianType | JacobianType |
typedef Superclass::InputPointType | InputPointType |
typedef Superclass::OutputPointType | OutputPointType |
typedef Superclass::InputVectorType | InputVectorType |
typedef Superclass::OutputVectorType | OutputVectorType |
typedef Superclass::InputCovariantVectorType | InputCovariantVectorType |
typedef Superclass::OutputCovariantVectorType | OutputCovariantVectorType |
typedef Superclass::InputVnlVectorType | InputVnlVectorType |
typedef Superclass::OutputVnlVectorType | OutputVnlVectorType |
typedef Superclass ::NonZeroJacobianIndicesType | NonZeroJacobianIndicesType |
typedef Superclass::SpatialJacobianType | SpatialJacobianType |
typedef Superclass ::JacobianOfSpatialJacobianType | JacobianOfSpatialJacobianType |
typedef Superclass::SpatialHessianType | SpatialHessianType |
typedef Superclass ::JacobianOfSpatialHessianType | JacobianOfSpatialHessianType |
typedef Superclass::InternalMatrixType | InternalMatrixType |
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 VectorContainer< unsigned long, InputVectorType > | VectorSetType |
typedef VectorSetType::Pointer | VectorSetPointer |
typedef vnl_matrix_fixed< TScalarType, NDimensions, NDimensions > | IMatrixType |
typedef vnl_matrix_fixed< TScalarType, NDimensions, NDimensions > | GMatrixType |
typedef vnl_matrix< TScalarType > | LMatrixType |
typedef vnl_matrix< TScalarType > | KMatrixType |
typedef vnl_matrix< TScalarType > | PMatrixType |
typedef vnl_matrix< TScalarType > | YMatrixType |
typedef vnl_matrix< TScalarType > | WMatrixType |
typedef vnl_matrix< TScalarType > | DMatrixType |
typedef vnl_matrix_fixed< TScalarType, NDimensions, NDimensions > | AMatrixType |
typedef vnl_vector_fixed< TScalarType, NDimensions > | BMatrixType |
typedef vnl_matrix_fixed< TScalarType, 1, NDimensions > | RowMatrixType |
typedef vnl_matrix_fixed< TScalarType, NDimensions, 1 > | ColumnMatrixType |
Protected Types
Name | |
---|---|
typedef vnl_svd< ScalarType > | SVDDecompositionType |
typedef vnl_qr< ScalarType > | QRDecompositionType |
Public Functions
Name | |
---|---|
itkTypeMacro(KernelTransform2 , AdvancedTransform ) | |
itkNewMacro(Self ) | |
itkStaticConstMacro(SpaceDimension , unsigned int , NDimensions ) | |
virtual NumberOfParametersType | GetNumberOfParameters(void ) const |
itkGetObjectMacro(SourceLandmarks , PointSetType ) | |
virtual void | SetSourceLandmarks(PointSetType * ) |
itkGetObjectMacro(TargetLandmarks , PointSetType ) | |
virtual void | SetTargetLandmarks(PointSetType * ) |
itkGetObjectMacro(Displacements , VectorSetType ) | |
void | ComputeWMatrix(void ) |
void | ComputeLInverse(void ) |
virtual OutputPointType | TransformPoint(const InputPointType & thisPoint) const |
virtual OutputVectorType | TransformVector(const InputVectorType & ) const |
virtual OutputVnlVectorType | TransformVector(const InputVnlVectorType & ) const |
virtual OutputCovariantVectorType | TransformCovariantVector(const InputCovariantVectorType & ) const |
virtual void | GetJacobian(const InputPointType & , JacobianType & , NonZeroJacobianIndicesType & ) const |
virtual void | SetIdentity(void ) |
virtual void | SetParameters(const ParametersType & ) |
virtual void | SetFixedParameters(const ParametersType & ) |
virtual void | UpdateParameters(void ) |
virtual const ParametersType & | GetParameters(void ) const |
virtual const ParametersType & | GetFixedParameters(void ) const |
virtual void | SetStiffness(double stiffness) |
itkGetMacro(Stiffness , double ) | |
virtual void | SetAlpha(TScalarType itkNotUsedAlpha) |
virtual TScalarType | GetAlpha(void ) const |
itkSetMacro(PoissonRatio , TScalarType ) | |
virtual const TScalarType | GetPoissonRatio(void ) const |
itkSetMacro(MatrixInversionMethod , std::string ) | |
itkGetConstReferenceMacro(MatrixInversionMethod , std::string ) | |
virtual void | GetSpatialJacobian(const InputPointType & ipp, SpatialJacobianType & sj) const |
virtual void | GetSpatialHessian(const InputPointType & ipp, SpatialHessianType & sh) const |
virtual void | GetJacobianOfSpatialJacobian(const InputPointType & ipp, JacobianOfSpatialJacobianType & jsj, NonZeroJacobianIndicesType & nonZeroJacobianIndices) const |
virtual void | GetJacobianOfSpatialJacobian(const InputPointType & ipp, SpatialJacobianType & sj, JacobianOfSpatialJacobianType & jsj, NonZeroJacobianIndicesType & nonZeroJacobianIndices) const |
virtual void | GetJacobianOfSpatialHessian(const InputPointType & ipp, JacobianOfSpatialHessianType & jsh, NonZeroJacobianIndicesType & nonZeroJacobianIndices) const |
virtual void | GetJacobianOfSpatialHessian(const InputPointType & ipp, SpatialHessianType & sh, JacobianOfSpatialHessianType & jsh, NonZeroJacobianIndicesType & nonZeroJacobianIndices) const |
Protected Functions
Name | |
---|---|
KernelTransform2() | |
virtual | ~KernelTransform2() |
void | PrintSelf(std::ostream & os, Indent indent) const |
virtual void | ComputeG(const InputVectorType & landmarkVector, GMatrixType & GMatrix) const |
virtual void | ComputeReflexiveG(PointsIterator , GMatrixType & GMatrix) const |
virtual void | ComputeDeformationContribution(const InputPointType & inputPoint, OutputPointType & result) const |
void | ComputeK(void ) |
void | ComputeL(void ) |
void | ComputeP(void ) |
void | ComputeY(void ) |
void | ComputeD(void ) |
void | ReorganizeW(void ) |
Public Attributes
Name | |
---|---|
PointSetPointer | m_SourceLandmarks |
PointSetPointer | m_TargetLandmarks |
Protected Attributes
Additional inherited members
Public Types inherited from itk::AdvancedTransform< TScalarType, NDimensions, NDimensions >
Name | |
---|---|
typedef Superclass::ParametersValueType | ParametersValueType |
typedef Superclass::DerivativeType | DerivativeType |
typedef Superclass::InverseTransformBaseType | InverseTransformBaseType |
typedef Superclass::InverseTransformBasePointer | InverseTransformBasePointer |
typedef Transform< TScalarType, NInputDimensions, NOutputDimensions > | TransformType |
typedef TransformType::Pointer | TransformTypePointer |
typedef TransformType::ConstPointer | TransformTypeConstPointer |
typedef OutputCovariantVectorType | MovingImageGradientType |
typedef MovingImageGradientType::ValueType | MovingImageGradientValueType |
Public Functions inherited from itk::AdvancedTransform< TScalarType, NDimensions, NDimensions >
Name | |
---|---|
virtual NumberOfParametersType | GetNumberOfNonZeroJacobianIndices(void ) const |
itkGetConstMacro(HasNonZeroSpatialHessian , bool ) | |
itkGetConstMacro(HasNonZeroJacobianOfSpatialHessian , bool ) | |
virtual void | EvaluateJacobianWithImageGradientProduct(const InputPointType & ipp, const MovingImageGradientType & movingImageGradient, DerivativeType & imageJacobian, NonZeroJacobianIndicesType & nonZeroJacobianIndices) const |
virtual void | ComputeJacobianWithRespectToParameters(const InputPointType & itkNotUsedp, JacobianType & itkNotUsedj) const |
virtual void | ComputeJacobianWithRespectToPosition(const InputPointType & itkNotUsedp, JacobianType & itkNotUsedj) const |
Protected Functions inherited from itk::AdvancedTransform< TScalarType, NDimensions, NDimensions >
Name | |
---|---|
AdvancedTransform() | |
AdvancedTransform(NumberOfParametersType numberOfParameters) | |
virtual | ~AdvancedTransform() |
Protected Attributes inherited from itk::AdvancedTransform< TScalarType, NDimensions, NDimensions >
Name | |
---|---|
bool | m_HasNonZeroSpatialHessian |
bool | m_HasNonZeroJacobianOfSpatialHessian |
Detailed Description
template <class TScalarType ,
unsigned int NDimensions>
class itk::KernelTransform2;
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.
KernelTransform2: 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."
Modified to include it in elastix:
- style
- make it inherit from AdvancedTransform
- make it threadsafe, like was done in the itk as well.
- Support for matrix inversion by QR decomposition, instead of SVD. QR is much faster. Used in SetParameters() and SetFixedParameters().
- Much faster Jacobian computation for some of the derived kernel transforms.
Public Types Documentation
typedef Self
typedef KernelTransform2 itk::KernelTransform2< TScalarType, NDimensions >::Self;
Standard class typedefs.
typedef Superclass
typedef AdvancedTransform< TScalarType, NDimensions, NDimensions > itk::KernelTransform2< TScalarType, NDimensions >::Superclass;
typedef Pointer
typedef SmartPointer< Self > itk::KernelTransform2< TScalarType, NDimensions >::Pointer;
typedef ConstPointer
typedef SmartPointer< const Self > itk::KernelTransform2< TScalarType, NDimensions >::ConstPointer;
typedef ScalarType
typedef Superclass::ScalarType itk::KernelTransform2< TScalarType, NDimensions >::ScalarType;
Typedefs.
typedef ParametersType
typedef Superclass::ParametersType itk::KernelTransform2< TScalarType, NDimensions >::ParametersType;
typedef NumberOfParametersType
typedef Superclass::NumberOfParametersType itk::KernelTransform2< TScalarType, NDimensions >::NumberOfParametersType;
typedef JacobianType
typedef Superclass::JacobianType itk::KernelTransform2< TScalarType, NDimensions >::JacobianType;
typedef InputPointType
typedef Superclass::InputPointType itk::KernelTransform2< TScalarType, NDimensions >::InputPointType;
typedef OutputPointType
typedef Superclass::OutputPointType itk::KernelTransform2< TScalarType, NDimensions >::OutputPointType;
typedef InputVectorType
typedef Superclass::InputVectorType itk::KernelTransform2< TScalarType, NDimensions >::InputVectorType;
typedef OutputVectorType
typedef Superclass::OutputVectorType itk::KernelTransform2< TScalarType, NDimensions >::OutputVectorType;
typedef InputCovariantVectorType
typedef Superclass::InputCovariantVectorType itk::KernelTransform2< TScalarType, NDimensions >::InputCovariantVectorType;
typedef OutputCovariantVectorType
typedef Superclass::OutputCovariantVectorType itk::KernelTransform2< TScalarType, NDimensions >::OutputCovariantVectorType;
typedef InputVnlVectorType
typedef Superclass::InputVnlVectorType itk::KernelTransform2< TScalarType, NDimensions >::InputVnlVectorType;
typedef OutputVnlVectorType
typedef Superclass::OutputVnlVectorType itk::KernelTransform2< TScalarType, NDimensions >::OutputVnlVectorType;
typedef NonZeroJacobianIndicesType
typedef Superclass ::NonZeroJacobianIndicesType itk::KernelTransform2< TScalarType, NDimensions >::NonZeroJacobianIndicesType;
AdvancedTransform typedefs.
typedef SpatialJacobianType
typedef Superclass::SpatialJacobianType itk::KernelTransform2< TScalarType, NDimensions >::SpatialJacobianType;
typedef JacobianOfSpatialJacobianType
typedef Superclass ::JacobianOfSpatialJacobianType itk::KernelTransform2< TScalarType, NDimensions >::JacobianOfSpatialJacobianType;
typedef SpatialHessianType
typedef Superclass::SpatialHessianType itk::KernelTransform2< TScalarType, NDimensions >::SpatialHessianType;
typedef JacobianOfSpatialHessianType
typedef Superclass ::JacobianOfSpatialHessianType itk::KernelTransform2< TScalarType, NDimensions >::JacobianOfSpatialHessianType;
typedef InternalMatrixType
typedef Superclass::InternalMatrixType itk::KernelTransform2< TScalarType, NDimensions >::InternalMatrixType;
typedef PointSetTraitsType
typedef DefaultStaticMeshTraits< TScalarType, NDimensions, NDimensions, TScalarType, TScalarType > itk::KernelTransform2< 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::KernelTransform2< TScalarType, NDimensions >::PointSetType;
typedef PointSetPointer
typedef PointSetType::Pointer itk::KernelTransform2< TScalarType, NDimensions >::PointSetPointer;
typedef PointsContainer
typedef PointSetType::PointsContainer itk::KernelTransform2< TScalarType, NDimensions >::PointsContainer;
typedef PointsIterator
typedef PointSetType::PointsContainerIterator itk::KernelTransform2< TScalarType, NDimensions >::PointsIterator;
typedef PointsConstIterator
typedef PointSetType::PointsContainerConstIterator itk::KernelTransform2< TScalarType, NDimensions >::PointsConstIterator;
typedef VectorSetType
typedef VectorContainer< unsigned long, InputVectorType > itk::KernelTransform2< TScalarType, NDimensions >::VectorSetType;
VectorSet typedef.
typedef VectorSetPointer
typedef VectorSetType::Pointer itk::KernelTransform2< TScalarType, NDimensions >::VectorSetPointer;
typedef IMatrixType
typedef vnl_matrix_fixed< TScalarType, NDimensions, NDimensions > itk::KernelTransform2< TScalarType, NDimensions >::IMatrixType;
'I' (identity) matrix typedef.
typedef GMatrixType
typedef vnl_matrix_fixed< TScalarType, NDimensions, NDimensions > itk::KernelTransform2< TScalarType, NDimensions >::GMatrixType;
'G' matrix typedef.
typedef LMatrixType
typedef vnl_matrix< TScalarType > itk::KernelTransform2< TScalarType, NDimensions >::LMatrixType;
'L' matrix typedef.
typedef KMatrixType
typedef vnl_matrix< TScalarType > itk::KernelTransform2< TScalarType, NDimensions >::KMatrixType;
'K' matrix typedef.
typedef PMatrixType
typedef vnl_matrix< TScalarType > itk::KernelTransform2< TScalarType, NDimensions >::PMatrixType;
'P' matrix typedef.
typedef YMatrixType
typedef vnl_matrix< TScalarType > itk::KernelTransform2< TScalarType, NDimensions >::YMatrixType;
'Y' matrix typedef.
typedef WMatrixType
typedef vnl_matrix< TScalarType > itk::KernelTransform2< TScalarType, NDimensions >::WMatrixType;
'W' matrix typedef.
typedef DMatrixType
typedef vnl_matrix< TScalarType > itk::KernelTransform2< TScalarType, NDimensions >::DMatrixType;
'D' matrix typedef. Deformation component
typedef AMatrixType
typedef vnl_matrix_fixed< TScalarType, NDimensions, NDimensions > itk::KernelTransform2< TScalarType, NDimensions >::AMatrixType;
'A' matrix typedef. Rotational part of the Affine component
typedef BMatrixType
typedef vnl_vector_fixed< TScalarType, NDimensions > itk::KernelTransform2< TScalarType, NDimensions >::BMatrixType;
'B' matrix typedef. Translational part of the Affine component
typedef RowMatrixType
typedef vnl_matrix_fixed< TScalarType, 1, NDimensions > itk::KernelTransform2< TScalarType, NDimensions >::RowMatrixType;
Row matrix typedef.
typedef ColumnMatrixType
typedef vnl_matrix_fixed< TScalarType, NDimensions, 1 > itk::KernelTransform2< TScalarType, NDimensions >::ColumnMatrixType;
Column matrix typedef.
Protected Types Documentation
typedef SVDDecompositionType
typedef vnl_svd< ScalarType > itk::KernelTransform2< TScalarType, NDimensions >::SVDDecompositionType;
Decompositions, needed for the L matrix. These decompositions are cached for performance reasons during registration. During registration, in every iteration SetParameters() is called, which in turn calls ComputeWMatrix(). The L matrix is not changed however, and therefore it is not needed to redo the decomposition.
typedef QRDecompositionType
typedef vnl_qr< ScalarType > itk::KernelTransform2< TScalarType, NDimensions >::QRDecompositionType;
Public Functions Documentation
function itkTypeMacro
itkTypeMacro(
KernelTransform2 ,
AdvancedTransform
)
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 GetNumberOfParameters
inline virtual NumberOfParametersType GetNumberOfParameters(
void
) const
Return the number of parameters that completely define the Transform.
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
)
Compute W matrix.
function ComputeLInverse
void ComputeLInverse(
void
)
Compute L matrix inverse.
function TransformPoint
virtual OutputPointType TransformPoint(
const InputPointType & thisPoint
) const
Compute the position of point in the new space
function TransformVector
inline virtual OutputVectorType TransformVector(
const InputVectorType &
) const
These vector transforms are not implemented for this transform.
function TransformVector
inline virtual OutputVnlVectorType TransformVector(
const InputVnlVectorType &
) const
function TransformCovariantVector
inline virtual OutputCovariantVectorType TransformCovariantVector(
const InputCovariantVectorType &
) const
function GetJacobian
virtual void GetJacobian(
const InputPointType & ,
JacobianType & ,
NonZeroJacobianIndicesType &
) const
Compute the Jacobian of the transformation.
function SetIdentity
virtual void SetIdentity(
void
)
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
)
Update the Parameters array from the landmarks coordinates.
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 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
)
function SetAlpha
inline virtual void SetAlpha(
TScalarType itkNotUsedAlpha
)
This method makes only sense for the ElasticBody splines. Declare here, so that you can always call it if you don't know the type of kernel beforehand. It will be overridden in the ElasticBodySplineKernelTransform and in the ElasticBodyReciprocalSplineKernelTransform.
function GetAlpha
inline virtual TScalarType GetAlpha(
void
) const
function itkSetMacro
itkSetMacro(
PoissonRatio ,
TScalarType
)
This method makes only sense for the ElasticBody splines. Declare here, so that you can always call it if you don't know the type of kernel beforehand. It will be overridden in the ElasticBodySplineKernelTransform and in the ElasticBodyReciprocalSplineKernelTransform.
function GetPoissonRatio
inline virtual const TScalarType GetPoissonRatio(
void
) const
function itkSetMacro
itkSetMacro(
MatrixInversionMethod ,
std::string
)
Matrix inversion by SVD or QR decomposition.
function itkGetConstReferenceMacro
itkGetConstReferenceMacro(
MatrixInversionMethod ,
std::string
)
function GetSpatialJacobian
inline virtual void GetSpatialJacobian(
const InputPointType & ipp,
SpatialJacobianType & sj
) const
Reimplements: itk::AdvancedTransform::GetSpatialJacobian
Must be provided.
function GetSpatialHessian
inline virtual void GetSpatialHessian(
const InputPointType & ipp,
SpatialHessianType & sh
) const
Reimplements: itk::AdvancedTransform::GetSpatialHessian
Compute the spatial Hessian of the transformation.
The spatial Hessian is the vector of matrices of partial second order derivatives of the transformation components with respect to the spatial position \(x\), evaluated at a point \(p\).
with i the i-th component of the transformation.
function GetJacobianOfSpatialJacobian
inline virtual void GetJacobianOfSpatialJacobian(
const InputPointType & ipp,
JacobianOfSpatialJacobianType & jsj,
NonZeroJacobianIndicesType & nonZeroJacobianIndices
) const
function GetJacobianOfSpatialJacobian
inline virtual void GetJacobianOfSpatialJacobian(
const InputPointType & ipp,
SpatialJacobianType & sj,
JacobianOfSpatialJacobianType & jsj,
NonZeroJacobianIndicesType & nonZeroJacobianIndices
) const
function GetJacobianOfSpatialHessian
inline virtual void GetJacobianOfSpatialHessian(
const InputPointType & ipp,
JacobianOfSpatialHessianType & jsh,
NonZeroJacobianIndicesType & nonZeroJacobianIndices
) const
function GetJacobianOfSpatialHessian
inline virtual void GetJacobianOfSpatialHessian(
const InputPointType & ipp,
SpatialHessianType & sh,
JacobianOfSpatialHessianType & jsh,
NonZeroJacobianIndicesType & nonZeroJacobianIndices
) const
Protected Functions Documentation
function KernelTransform2
KernelTransform2()
function ~KernelTransform2
virtual ~KernelTransform2()
function PrintSelf
void PrintSelf(
std::ostream & os,
Indent indent
) const
function ComputeG
virtual void ComputeG(
const InputVectorType & landmarkVector,
GMatrixType & GMatrix
) const
Reimplemented by: itk::ThinPlateSplineKernelTransform2::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 void ComputeReflexiveG(
PointsIterator ,
GMatrixType & GMatrix
) 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::ThinPlateSplineKernelTransform2::ComputeDeformationContribution
Compute the contribution of the landmarks weighted by the kernel function to the global deformation of the space.
function ComputeK
void ComputeK(
void
)
Compute K matrix.
function ComputeL
void ComputeL(
void
)
Compute L matrix.
function ComputeP
void ComputeP(
void
)
Compute P matrix.
function ComputeY
void ComputeY(
void
)
Compute Y matrix.
function ComputeD
void ComputeD(
void
)
Compute displacements \( q_i - p_i \).
function ReorganizeW
void ReorganizeW(
void
)
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;
Rotational/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_WMatrixComputed
bool m_WMatrixComputed;
The G matrix. It used to be mutable because m_GMatrix was made an ivar only to avoid copying the matrix at return time but this is not necessary. SK: we don't need this matrix anymore as a member. 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_LMatrixDecompositionComputed
bool m_LMatrixDecompositionComputed;
Has the L matrix decomposition been computed?
variable m_LMatrixDecompositionSVD
SVDDecompositionType * m_LMatrixDecompositionSVD;
variable m_LMatrixDecompositionQR
QRDecompositionType * m_LMatrixDecompositionQR;
variable m_I
IMatrixType m_I;
Identity matrix.
variable m_NonZeroJacobianIndices
NonZeroJacobianIndicesType m_NonZeroJacobianIndices;
Precomputed nonzero Jacobian indices (simply all params)
variable m_NonZeroJacobianIndicesTemp
NonZeroJacobianIndicesType m_NonZeroJacobianIndicesTemp;
for old GetJacobian() method:
variable m_FastComputationPossible
bool m_FastComputationPossible;
The Jacobian can be computed much faster for some of the derived kerbel transforms, most notably the TPS.
Updated on 2022-07-23 at 16:40:05 -0600