itk::AdvancedTransform
Transform maps points, vectors and covariant vectors from an input space to an output space. More...
#include <itkAdvancedTransform.h>
Inherits from Transform< TScalarType, 3, 3 >
Public Types
Name | |
---|---|
typedef AdvancedTransform | Self |
typedef Transform< TScalarType, NInputDimensions, NOutputDimensions > | Superclass |
typedef SmartPointer< Self > | Pointer |
typedef SmartPointer< const Self > | ConstPointer |
typedef Superclass::ScalarType | ScalarType |
typedef Superclass::ParametersType | ParametersType |
typedef Superclass::ParametersValueType | ParametersValueType |
typedef Superclass::NumberOfParametersType | NumberOfParametersType |
typedef Superclass::DerivativeType | DerivativeType |
typedef Superclass::JacobianType | JacobianType |
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::InputPointType | InputPointType |
typedef Superclass::OutputPointType | OutputPointType |
typedef Superclass::InverseTransformBaseType | InverseTransformBaseType |
typedef Superclass::InverseTransformBasePointer | InverseTransformBasePointer |
typedef Transform< TScalarType, NInputDimensions, NOutputDimensions > | TransformType |
typedef TransformType::Pointer | TransformTypePointer |
typedef TransformType::ConstPointer | TransformTypeConstPointer |
typedef std::vector< unsigned long > | NonZeroJacobianIndicesType |
typedef Matrix< ScalarType, OutputSpaceDimension, InputSpaceDimension > | SpatialJacobianType |
typedef std::vector< SpatialJacobianType > | JacobianOfSpatialJacobianType |
typedef FixedArray< Matrix< ScalarType, InputSpaceDimension, InputSpaceDimension >, OutputSpaceDimension > | SpatialHessianType |
typedef std::vector< SpatialHessianType > | JacobianOfSpatialHessianType |
typedef SpatialJacobianType::InternalMatrixType | InternalMatrixType |
typedef OutputCovariantVectorType | MovingImageGradientType |
typedef MovingImageGradientType::ValueType | MovingImageGradientValueType |
Public Functions
Name | |
---|---|
itkTypeMacro(AdvancedTransform , Transform ) | |
itkStaticConstMacro(InputSpaceDimension , unsigned int , NInputDimensions ) | |
itkStaticConstMacro(OutputSpaceDimension , unsigned int , NOutputDimensions ) | |
virtual NumberOfParametersType | GetNumberOfNonZeroJacobianIndices(void ) const |
itkGetConstMacro(HasNonZeroSpatialHessian , bool ) | |
itkGetConstMacro(HasNonZeroJacobianOfSpatialHessian , bool ) | |
virtual void | GetJacobian(const InputPointType & ipp, JacobianType & j, NonZeroJacobianIndicesType & nonZeroJacobianIndices) const =0 |
virtual void | EvaluateJacobianWithImageGradientProduct(const InputPointType & ipp, const MovingImageGradientType & movingImageGradient, DerivativeType & imageJacobian, NonZeroJacobianIndicesType & nonZeroJacobianIndices) const |
virtual void | GetSpatialJacobian(const InputPointType & ipp, SpatialJacobianType & sj) const =0 |
virtual void | ComputeJacobianWithRespectToParameters(const InputPointType & itkNotUsedp, JacobianType & itkNotUsedj) const |
virtual void | ComputeJacobianWithRespectToPosition(const InputPointType & itkNotUsedp, JacobianType & itkNotUsedj) const |
virtual void | GetSpatialHessian(const InputPointType & ipp, SpatialHessianType & sh) const =0 |
virtual void | GetJacobianOfSpatialJacobian(const InputPointType & ipp, JacobianOfSpatialJacobianType & jsj, NonZeroJacobianIndicesType & nonZeroJacobianIndices) const =0 |
virtual void | GetJacobianOfSpatialJacobian(const InputPointType & ipp, SpatialJacobianType & sj, JacobianOfSpatialJacobianType & jsj, NonZeroJacobianIndicesType & nonZeroJacobianIndices) const =0 |
virtual void | GetJacobianOfSpatialHessian(const InputPointType & ipp, JacobianOfSpatialHessianType & jsh, NonZeroJacobianIndicesType & nonZeroJacobianIndices) const =0 |
virtual void | GetJacobianOfSpatialHessian(const InputPointType & ipp, SpatialHessianType & sh, JacobianOfSpatialHessianType & jsh, NonZeroJacobianIndicesType & nonZeroJacobianIndices) const =0 |
Protected Functions
Name | |
---|---|
AdvancedTransform() | |
AdvancedTransform(NumberOfParametersType numberOfParameters) | |
virtual | ~AdvancedTransform() |
Protected Attributes
Name | |
---|---|
bool | m_HasNonZeroSpatialHessian |
bool | m_HasNonZeroJacobianOfSpatialHessian |
Detailed Description
template <class TScalarType ,
unsigned int NInputDimensions =3,
unsigned int NOutputDimensions =3>
class itk::AdvancedTransform;
Transform maps points, vectors and covariant vectors from an input space to an output space.
Par: Registration Framework Support
Typically a Transform class has several methods for setting its parameters. For use in the registration framework, the parameters must also be represented by an array of doubles to allow communication with generic optimizers. The Array of transformation parameters is set using the SetParameters() method.
This abstract class define the generic interface for a geometrical transformation from one space to another. The class provides methods for mapping points, vectors and covariant vectors from the input space to the output space.
Given that transformation are not necessarily invertible, this basic class does not provide the methods for back transformation. Back transform methods are implemented in derived classes where appropriate.
Another requirement of the registration framework is the computation of the Jacobian of the transform T. In general, an ImageToImageMetric requires the knowledge of this Jacobian in order to compute the metric derivatives. The Jacobian is a matrix whose element are the partial derivatives of the transformation with respect to the array of parameters mu that defines the transform, evaluated at a point p: dT/dmu(p).
If penalty terms are included in the registration, the transforms also need to implement other derivatives of T. Often, penalty terms are functions of the spatial derivatives of T. Therefore, e.g. the SpatialJacobian dT/dx and the SpatialHessian d^2T/dx_idx_j require implementation. The GetValueAndDerivative() requires the d/dmu of those terms. Therefore, we additionally define GetJacobianOfSpatialJacobian() and GetJacobianOfSpatialHessian().
Public Types Documentation
typedef Self
typedef AdvancedTransform itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::Self;
Standard class typedefs.
typedef Superclass
typedef Transform< TScalarType, NInputDimensions, NOutputDimensions > itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::Superclass;
typedef Pointer
typedef SmartPointer< Self > itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::Pointer;
typedef ConstPointer
typedef SmartPointer< const Self > itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::ConstPointer;
typedef ScalarType
typedef Superclass::ScalarType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::ScalarType;
Typedefs from the Superclass.
typedef ParametersType
typedef Superclass::ParametersType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::ParametersType;
typedef ParametersValueType
typedef Superclass::ParametersValueType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::ParametersValueType;
typedef NumberOfParametersType
typedef Superclass::NumberOfParametersType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::NumberOfParametersType;
typedef DerivativeType
typedef Superclass::DerivativeType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::DerivativeType;
typedef JacobianType
typedef Superclass::JacobianType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::JacobianType;
typedef InputVectorType
typedef Superclass::InputVectorType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::InputVectorType;
typedef OutputVectorType
typedef Superclass::OutputVectorType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::OutputVectorType;
typedef InputCovariantVectorType
typedef Superclass::InputCovariantVectorType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::InputCovariantVectorType;
typedef OutputCovariantVectorType
typedef Superclass::OutputCovariantVectorType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::OutputCovariantVectorType;
typedef InputVnlVectorType
typedef Superclass::InputVnlVectorType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::InputVnlVectorType;
typedef OutputVnlVectorType
typedef Superclass::OutputVnlVectorType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::OutputVnlVectorType;
typedef InputPointType
typedef Superclass::InputPointType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::InputPointType;
typedef OutputPointType
typedef Superclass::OutputPointType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::OutputPointType;
typedef InverseTransformBaseType
typedef Superclass::InverseTransformBaseType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::InverseTransformBaseType;
typedef InverseTransformBasePointer
typedef Superclass::InverseTransformBasePointer itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::InverseTransformBasePointer;
typedef TransformType
typedef Transform< TScalarType, NInputDimensions, NOutputDimensions > itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::TransformType;
Transform typedefs for the from Superclass.
typedef TransformTypePointer
typedef TransformType::Pointer itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::TransformTypePointer;
typedef TransformTypeConstPointer
typedef TransformType::ConstPointer itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::TransformTypeConstPointer;
typedef NonZeroJacobianIndicesType
typedef std::vector< unsigned long > itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::NonZeroJacobianIndicesType;
Types for the (Spatial)Jacobian/Hessian. Using an itk::FixedArray instead of an std::vector gives a performance gain for the SpatialHessianType.
typedef SpatialJacobianType
typedef Matrix< ScalarType, OutputSpaceDimension, InputSpaceDimension > itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::SpatialJacobianType;
typedef JacobianOfSpatialJacobianType
typedef std::vector< SpatialJacobianType > itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::JacobianOfSpatialJacobianType;
typedef SpatialHessianType
typedef FixedArray< Matrix< ScalarType, InputSpaceDimension, InputSpaceDimension >, OutputSpaceDimension > itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::SpatialHessianType;
typedef JacobianOfSpatialHessianType
typedef std::vector< SpatialHessianType > itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::JacobianOfSpatialHessianType;
typedef InternalMatrixType
typedef SpatialJacobianType::InternalMatrixType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::InternalMatrixType;
typedef MovingImageGradientType
typedef OutputCovariantVectorType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::MovingImageGradientType;
Typedef for the moving image gradient type. This type is defined by the B-spline interpolator as typedef CovariantVector< RealType, ImageDimension > As we cannot access this type we simply re-construct it to be identical.
typedef MovingImageGradientValueType
typedef MovingImageGradientType::ValueType itk::AdvancedTransform< TScalarType, NInputDimensions, NOutputDimensions >::MovingImageGradientValueType;
Public Functions Documentation
function itkTypeMacro
itkTypeMacro(
AdvancedTransform ,
Transform
)
New method for creating an object using a factory. Run-time type information (and related methods).
function itkStaticConstMacro
itkStaticConstMacro(
InputSpaceDimension ,
unsigned int ,
NInputDimensions
)
Dimension of the domain space.
function itkStaticConstMacro
itkStaticConstMacro(
OutputSpaceDimension ,
unsigned int ,
NOutputDimensions
)
function GetNumberOfNonZeroJacobianIndices
virtual NumberOfParametersType GetNumberOfNonZeroJacobianIndices(
void
) const
Get the number of nonzero Jacobian indices. By default all.
function itkGetConstMacro
itkGetConstMacro(
HasNonZeroSpatialHessian ,
bool
)
Whether the advanced transform has nonzero matrices.
function itkGetConstMacro
itkGetConstMacro(
HasNonZeroJacobianOfSpatialHessian ,
bool
)
function GetJacobian
virtual void GetJacobian(
const InputPointType & ipp,
JacobianType & j,
NonZeroJacobianIndicesType & nonZeroJacobianIndices
) const =0
Reimplemented by: itk::KernelTransform2::GetJacobian, itk::KernelTransform2::GetJacobian
This returns a sparse version of the Jacobian of the transformation.
The Jacobian is expressed as a vector of partial derivatives of the transformation components with respect to the parameters \(\mu\) that define the transformation \(T\), evaluated at a point \(p\).
with \(m\) the number of parameters, i.e. the size of \(\mu\), and \(d\) the dimension of the image.
function EvaluateJacobianWithImageGradientProduct
virtual void EvaluateJacobianWithImageGradientProduct(
const InputPointType & ipp,
const MovingImageGradientType & movingImageGradient,
DerivativeType & imageJacobian,
NonZeroJacobianIndicesType & nonZeroJacobianIndices
) const
Compute the inner product of the Jacobian with the moving image gradient. The Jacobian is (partially) constructed inside this function, but not returned.
function GetSpatialJacobian
virtual void GetSpatialJacobian(
const InputPointType & ipp,
SpatialJacobianType & sj
) const =0
Reimplemented by: itk::KernelTransform2::GetSpatialJacobian, itk::KernelTransform2::GetSpatialJacobian
Compute the spatial Jacobian of the transformation.
The spatial Jacobian is expressed as a vector of partial derivatives of the transformation components with respect to the spatial position \(x\), evaluated at a point \(p\).
with \(m\) the number of parameters, i.e. the size of \(\mu\), and \(d\) the dimension of the image.
function ComputeJacobianWithRespectToParameters
inline virtual void ComputeJacobianWithRespectToParameters(
const InputPointType & itkNotUsedp,
JacobianType & itkNotUsedj
) const
Override some pure virtual ITK4 functions.
function ComputeJacobianWithRespectToPosition
inline virtual void ComputeJacobianWithRespectToPosition(
const InputPointType & itkNotUsedp,
JacobianType & itkNotUsedj
) const
function GetSpatialHessian
virtual void GetSpatialHessian(
const InputPointType & ipp,
SpatialHessianType & sh
) const =0
Reimplemented by: itk::KernelTransform2::GetSpatialHessian, itk::KernelTransform2::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
virtual void GetJacobianOfSpatialJacobian(
const InputPointType & ipp,
JacobianOfSpatialJacobianType & jsj,
NonZeroJacobianIndicesType & nonZeroJacobianIndices
) const =0
Reimplemented by: itk::KernelTransform2::GetJacobianOfSpatialJacobian, itk::KernelTransform2::GetJacobianOfSpatialJacobian
Compute the Jacobian of the spatial Jacobian of the transformation.
The Jacobian of the spatial Jacobian is the derivative of the spatial Jacobian to the transformation parameters \(\mu\), evaluated at a point \(p\).
function GetJacobianOfSpatialJacobian
virtual void GetJacobianOfSpatialJacobian(
const InputPointType & ipp,
SpatialJacobianType & sj,
JacobianOfSpatialJacobianType & jsj,
NonZeroJacobianIndicesType & nonZeroJacobianIndices
) const =0
Reimplemented by: itk::KernelTransform2::GetJacobianOfSpatialJacobian, itk::KernelTransform2::GetJacobianOfSpatialJacobian
Compute both the spatial Jacobian and the Jacobian of the spatial Jacobian of the transformation.
function GetJacobianOfSpatialHessian
virtual void GetJacobianOfSpatialHessian(
const InputPointType & ipp,
JacobianOfSpatialHessianType & jsh,
NonZeroJacobianIndicesType & nonZeroJacobianIndices
) const =0
Reimplemented by: itk::KernelTransform2::GetJacobianOfSpatialHessian, itk::KernelTransform2::GetJacobianOfSpatialHessian
Compute the Jacobian of the spatial Hessian of the transformation.
The Jacobian of the spatial Hessian is the derivative of the spatial Hessian to the transformation parameters \(\mu\), evaluated at a point \(p\).
function GetJacobianOfSpatialHessian
virtual void GetJacobianOfSpatialHessian(
const InputPointType & ipp,
SpatialHessianType & sh,
JacobianOfSpatialHessianType & jsh,
NonZeroJacobianIndicesType & nonZeroJacobianIndices
) const =0
Reimplemented by: itk::KernelTransform2::GetJacobianOfSpatialHessian, itk::KernelTransform2::GetJacobianOfSpatialHessian
Compute both the spatial Hessian and the Jacobian of the spatial Hessian of the transformation.
Protected Functions Documentation
function AdvancedTransform
AdvancedTransform()
function AdvancedTransform
AdvancedTransform(
NumberOfParametersType numberOfParameters
)
function ~AdvancedTransform
inline virtual ~AdvancedTransform()
Protected Attributes Documentation
variable m_HasNonZeroSpatialHessian
bool m_HasNonZeroSpatialHessian;
variable m_HasNonZeroJacobianOfSpatialHessian
bool m_HasNonZeroJacobianOfSpatialHessian;
Updated on 2024-03-17 at 12:58:44 -0600