ITK  4.8.0
Insight Segmentation and Registration Toolkit
Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder > Class Template Reference

#include <itkBSplineDeformableTransform.h>

+ Inheritance diagram for itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >:
+ Collaboration diagram for itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >:

Detailed Description

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
class itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >

Deformable transform using a BSpline representation.

Note
BSplineTransform is a newer version of this class, and it is preferred.

This class encapsulates a deformable transform of points from one N-dimensional space to another N-dimensional space. The deformation field is modelled using B-splines. A deformation is defined on a sparse regular grid of control points $ \vec{\lambda}_j $ and is varied by defining a deformation $ \vec{g}(\vec{\lambda}_j) $ of each control point. The deformation $ D(\vec{x}) $ at any point $ \vec{x} $ is obtained by using a B-spline interpolation kernel.

The deformation field grid is defined by a user specified GridRegion, GridSpacing and GridOrigin. Each grid/control point has associated with it N deformation coefficients $ \vec{\delta}_j $, representing the N directional components of the deformation. Deformation outside the grid plus support region for the BSpline interpolation is assumed to be zero.

Additionally, the user can specified an addition bulk transform $ B $ such that the transformed point is given by:

\[ \vec{y} = B(\vec{x}) + D(\vec{x}) \]

The parameters for this transform is an N x N-D grid of spline coefficients. The user specifies the parameters as one flat array: each N-D grid is represented by an array in the same way an N-D image is represented in the buffer; the N arrays are then concatentated together on form a single array.

For efficiency, this transform does not make a copy of the parameters. It only keeps a pointer to the input parameters and assumes that the memory is managed by the caller.

The following illustrates the typical usage of this class:

* typedef BSplineDeformableTransform<double,2,3> TransformType;
* TransformType::Pointer transform = TransformType::New();
*
* transform->SetGridRegion( region );
* transform->SetGridSpacing( spacing );
* transform->SetGridOrigin( origin );
*
* // NB: the region must be set first before setting the parameters
*
* TransformType::ParametersType parameters(
*                                       transform->GetNumberOfParameters() );
*
* // Fill the parameters with values
*
* transform->SetParameters( parameters )
*
* outputPoint = transform->TransformPoint( inputPoint );
*
* 

An alternative way to set the B-spline coefficients is via array of images. The grid region, spacing and origin information is taken directly from the first image. It is assumed that the subsequent images are the same buffered region. The following illustrates the API:

*
* TransformType::ImageConstPointer images[2];
*
* // Fill the images up with values
*
* transform->SetCoefficientImages( images );
* outputPoint = transform->TransformPoint( inputPoint );
*
* 

Warning: use either the SetParameters() or SetCoefficientImages() API. Mixing the two modes may results in unexpected results.

The class is templated coordinate representation type (float or double), the space dimension and the spline order.

See Also
BSplineTransform
Wiki Examples:
Examples:
WikiExamples/Registration/ImageRegistrationMethodBSpline.cxx.

Definition at line 116 of file itkBSplineDeformableTransform.h.

Public Types

typedef
BulkTransformType::ConstPointer 
BulkTransformPointer
 
typedef Transform
< TParametersValueType,
itkGetStaticConstMacro(SpaceDimension),
itkGetStaticConstMacro(SpaceDimension)> 
BulkTransformType
 
typedef
Superclass::CoefficientImageArray 
CoefficientImageArray
 
typedef SmartPointer< const SelfConstPointer
 
typedef
Superclass::ContinuousIndexType 
ContinuousIndexType
 
typedef Superclass::DirectionType DirectionType
 
typedef
Superclass::FixedParametersType 
FixedParametersType
 
typedef Superclass::ImagePointer ImagePointer
 
typedef Superclass::ImageType ImageType
 
typedef Superclass::IndexType IndexType
 
typedef
Superclass::InputCovariantVectorType 
InputCovariantVectorType
 
typedef Superclass::InputVectorType InputVectorType
 
typedef
Superclass::InputVnlVectorType 
InputVnlVectorType
 
typedef Superclass::JacobianType JacobianType
 
typedef Superclass::MeshSizeType MeshSizeType
 
typedef
Superclass::NumberOfParametersType 
NumberOfParametersType
 
typedef Superclass::OriginType OriginType
 
typedef
Superclass::OutputCovariantVectorType 
OutputCovariantVectorType
 
typedef
Superclass::OutputVectorType 
OutputVectorType
 
typedef
Superclass::OutputVnlVectorType 
OutputVnlVectorType
 
typedef
Superclass::ParameterIndexArrayType 
ParameterIndexArrayType
 
typedef Superclass::ParametersType ParametersType
 
typedef
Superclass::ParametersValueType 
ParametersValueType
 
typedef Superclass::SpacingType PhysicalDimensionsType
 
typedef Superclass::PixelType PixelType
 
typedef SmartPointer< SelfPointer
 
typedef Superclass::RegionType RegionType
 
typedef TParametersValueType ScalarType
 
typedef BSplineDeformableTransform Self
 
typedef Superclass::SizeType SizeType
 
typedef Superclass::SpacingType SpacingType
 
typedef BSplineBaseTransform
< TParametersValueType,
NDimensions, VSplineOrder > 
Superclass
 
typedef
Superclass::WeightsFunctionType 
WeightsFunctionType
 
typedef Superclass::WeightsType WeightsType
 
typedef Point
< TParametersValueType,
itkGetStaticConstMacro(SpaceDimension)> 
InputPointType
 
typedef Point
< TParametersValueType,
itkGetStaticConstMacro(SpaceDimension)> 
OutputPointType
 
- Public Types inherited from itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >
typedef FixedArray
< ImagePointer, NDimensions > 
CoefficientImageArray
 
typedef SmartPointer< const SelfConstPointer
 
typedef
WeightsFunctionType::ContinuousIndexType 
ContinuousIndexType
 
typedef Superclass::DerivativeType DerivativeType
 
typedef ImageType::DirectionType DirectionType
 
typedef
Superclass::FixedParametersType 
FixedParametersType
 
typedef ImageType::Pointer ImagePointer
 
typedef Image
< ParametersValueType,
itkGetStaticConstMacro(SpaceDimension)> 
ImageType
 
typedef RegionType::IndexType IndexType
 
typedef vnl_vector_fixed
< TParametersValueType,
SpaceDimension
InputVnlVectorType
 
typedef Superclass::JacobianType JacobianType
 
typedef SizeType MeshSizeType
 
typedef
Superclass::NumberOfParametersType 
NumberOfParametersType
 
typedef ImageType::PointType OriginType
 
typedef vnl_vector_fixed
< TParametersValueType,
SpaceDimension
OutputVnlVectorType
 
typedef Array< unsigned long > ParameterIndexArrayType
 
typedef Superclass::ParametersType ParametersType
 
typedef ParametersType::ValueType ParametersValueType
 
typedef ImageType::SpacingType PhysicalDimensionsType
 
typedef ImageType::PixelType PixelType
 
typedef SmartPointer< SelfPointer
 
typedef ImageRegion
< itkGetStaticConstMacro(SpaceDimension)> 
RegionType
 
typedef Superclass::ScalarType ScalarType
 
typedef BSplineBaseTransform Self
 
typedef RegionType::SizeType SizeType
 
typedef ImageType::SpacingType SpacingType
 
typedef Transform
< TParametersValueType,
NDimensions, NDimensions > 
Superclass
 
typedef
Superclass::TransformCategoryType 
TransformCategoryType
 
typedef
BSplineInterpolationWeightFunction
< ScalarType,
itkGetStaticConstMacro(SpaceDimension),
itkGetStaticConstMacro(SplineOrder)> 
WeightsFunctionType
 
typedef
WeightsFunctionType::WeightsType 
WeightsType
 
typedef Vector
< TParametersValueType,
itkGetStaticConstMacro(SpaceDimension)> 
InputVectorType
 
typedef Vector
< TParametersValueType,
itkGetStaticConstMacro(SpaceDimension)> 
OutputVectorType
 
typedef CovariantVector
< TParametersValueType,
itkGetStaticConstMacro(SpaceDimension)> 
InputCovariantVectorType
 
typedef CovariantVector
< TParametersValueType,
itkGetStaticConstMacro(SpaceDimension)> 
OutputCovariantVectorType
 
typedef Point
< TParametersValueType,
itkGetStaticConstMacro(SpaceDimension)> 
InputPointType
 
typedef Point
< TParametersValueType,
itkGetStaticConstMacro(SpaceDimension)> 
OutputPointType
 
- Public Types inherited from itk::Transform< TParametersValueType, NDimensions, NDimensions >
typedef SmartPointer< const SelfConstPointer
 
typedef Array
< ParametersValueType
DerivativeType
 
typedef Matrix< double,
itkGetStaticConstMacro(OutputSpaceDimension),
itkGetStaticConstMacro(InputSpaceDimension)> 
DirectionChangeMatrix
 
typedef
Superclass::FixedParametersType 
FixedParametersType
 
typedef
Superclass::FixedParametersValueType 
FixedParametersValueType
 
typedef CovariantVector
< TParametersValueType,
NInputDimensions > 
InputCovariantVectorType
 
typedef DiffusionTensor3D
< TParametersValueType > 
InputDiffusionTensor3DType
 
typedef Matrix< double,
itkGetStaticConstMacro(InputSpaceDimension),
itkGetStaticConstMacro(InputSpaceDimension)> 
InputDirectionMatrix
 
typedef Point
< TParametersValueType,
NInputDimensions > 
InputPointType
 
typedef
SymmetricSecondRankTensor
< TParametersValueType,
NInputDimensions > 
InputSymmetricSecondRankTensorType
 
typedef VariableLengthVector
< TParametersValueType > 
InputVectorPixelType
 
typedef Vector
< TParametersValueType,
NInputDimensions > 
InputVectorType
 
typedef vnl_vector_fixed
< TParametersValueType,
NInputDimensions > 
InputVnlVectorType
 
typedef
InverseTransformBaseType::Pointer 
InverseTransformBasePointer
 
typedef Transform
< TParametersValueType,
NOutputDimensions,
NInputDimensions > 
InverseTransformBaseType
 
typedef Array2D
< ParametersValueType
JacobianType
 
typedef Matrix
< TParametersValueType,
itkGetStaticConstMacro(OutputSpaceDimension),
itkGetStaticConstMacro(InputSpaceDimension)> 
MatrixType
 
typedef
Superclass::NumberOfParametersType 
NumberOfParametersType
 
typedef CovariantVector
< TParametersValueType,
NOutputDimensions > 
OutputCovariantVectorType
 
typedef DiffusionTensor3D
< TParametersValueType > 
OutputDiffusionTensor3DType
 
typedef Matrix< double,
itkGetStaticConstMacro(OutputSpaceDimension),
itkGetStaticConstMacro(OutputSpaceDimension)> 
OutputDirectionMatrix
 
typedef Point
< TParametersValueType,
NOutputDimensions > 
OutputPointType
 
typedef
SymmetricSecondRankTensor
< TParametersValueType,
NOutputDimensions > 
OutputSymmetricSecondRankTensorType
 
typedef VariableLengthVector
< TParametersValueType > 
OutputVectorPixelType
 
typedef Vector
< TParametersValueType,
NOutputDimensions > 
OutputVectorType
 
typedef vnl_vector_fixed
< TParametersValueType,
NOutputDimensions > 
OutputVnlVectorType
 
typedef Superclass::ParametersType ParametersType
 
typedef
Superclass::ParametersValueType 
ParametersValueType
 
typedef SmartPointer< SelfPointer
 
typedef ParametersValueType ScalarType
 
typedef Transform Self
 
typedef TransformBaseTemplate
< TParametersValueType > 
Superclass
 
typedef
Superclass::TransformCategoryType 
TransformCategoryType
 
- Public Types inherited from itk::TransformBaseTemplate< TParametersValueType >
typedef SmartPointer< const SelfConstPointer
 
typedef OptimizerParameters
< FixedParametersValueType
FixedParametersType
 
typedef ParametersValueType FixedParametersValueType
 
typedef IdentifierType NumberOfParametersType
 
typedef OptimizerParameters
< ParametersValueType
ParametersType
 
typedef TParametersValueType ParametersValueType
 
typedef SmartPointer< SelfPointer
 
typedef TransformBaseTemplate Self
 
typedef Object Superclass
 
enum  TransformCategoryType {
  UnknownTransformCategory =0,
  Linear =1,
  BSpline =2,
  Spline =3,
  DisplacementField =4,
  VelocityField =5
}
 
- Public Types inherited from itk::Object
typedef SmartPointer< const SelfConstPointer
 
typedef SmartPointer< SelfPointer
 
typedef Object Self
 
typedef LightObject Superclass
 
- Public Types inherited from itk::LightObject
typedef SmartPointer< const SelfConstPointer
 
typedef SmartPointer< SelfPointer
 
typedef LightObject Self
 

Public Member Functions

virtual void ComputeJacobianWithRespectToParameters (const InputPointType &, JacobianType &) const override
 
virtual DirectionType GetGridDirection () const
 
virtual OriginType GetGridOrigin () const
 
virtual RegionType GetGridRegion () const
 
virtual SpacingType GetGridSpacing () const
 
virtual const char * GetNameOfClass () const
 
virtual NumberOfParametersType GetNumberOfParameters () const override
 
NumberOfParametersType GetNumberOfParametersPerDimension () const override
 
virtual const RegionTypeGetValidRegion () const
 
 itkCloneMacro (Self)
 
virtual void SetCoefficientImages (const CoefficientImageArray &images) override
 
virtual void SetGridDirection (const DirectionType &)
 
virtual void SetGridOrigin (const OriginType &)
 
virtual void SetGridRegion (const RegionType &)
 
virtual void SetGridSpacing (const SpacingType &)
 
virtual void SetFixedParameters (const FixedParametersType &parameters) override
 
virtual void TransformPoint (const InputPointType &inputPoint, OutputPointType &outputPoint, WeightsType &weights, ParameterIndexArrayType &indices, bool &inside) const override
 
virtual void SetBulkTransform (const BulkTransformType *_arg)
 
virtual const BulkTransformTypeGetBulkTransform () const
 
- Public Member Functions inherited from itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >
void ComputeJacobianFromBSplineWeightsWithRespectToPosition (const InputPointType &, WeightsType &, ParameterIndexArrayType &) const
 
virtual void ComputeJacobianWithRespectToPosition (const InputPointType &, JacobianType &) const override
 
const CoefficientImageArray GetCoefficientImages () const
 
virtual const FixedParametersTypeGetFixedParameters () const override
 
unsigned int GetNumberOfAffectedWeights () const
 
virtual NumberOfParametersType GetNumberOfLocalParameters () const override
 
unsigned long GetNumberOfWeights () const
 
virtual const ParametersTypeGetParameters () const override
 
virtual TransformCategoryType GetTransformCategory () const override
 
 itkCloneMacro (Self)
 
void SetIdentity ()
 
void SetParameters (const ParametersType &parameters) override
 
void SetParametersByValue (const ParametersType &parameters) override
 
OutputPointType TransformPoint (const InputPointType &point) const override
 
virtual OutputVnlVectorType TransformVector (const InputVnlVectorType &) const override
 
virtual void UpdateTransformParameters (const DerivativeType &update, TParametersValueType factor=1.0) override
 
virtual OutputVectorType TransformVector (const InputVectorType &) const override
 
virtual OutputCovariantVectorType TransformCovariantVector (const InputCovariantVectorType &) const override
 
- Public Member Functions inherited from itk::Transform< TParametersValueType, NDimensions, NDimensions >
virtual void ComputeInverseJacobianWithRespectToPosition (const InputPointType &x, JacobianType &jacobian) const
 
virtual void ComputeJacobianWithRespectToParameters (const InputPointType &, JacobianType &) const =0
 
virtual void ComputeJacobianWithRespectToParametersCachedTemporaries (const InputPointType &p, JacobianType &jacobian, JacobianType &) const
 
virtual void ComputeJacobianWithRespectToPosition (const InputPointType &, JacobianType &) const
 
virtual void CopyInFixedParameters (const FixedParametersValueType *const begin, const FixedParametersValueType *const end) override
 
virtual void CopyInParameters (const ParametersValueType *const begin, const ParametersValueType *const end) override
 
unsigned int GetInputSpaceDimension (void) const override
 
bool GetInverse (Self *) const
 
virtual InverseTransformBasePointer GetInverseTransform () const
 
virtual NumberOfParametersType GetNumberOfFixedParameters () const
 
unsigned int GetOutputSpaceDimension (void) const override
 
virtual std::string GetTransformTypeAsString () const override
 
virtual bool IsLinear () const
 
 itkCloneMacro (Self)
 
virtual void SetFixedParameters (const FixedParametersType &) override=0
 
virtual void SetParameters (const ParametersType &) override=0
 
virtual void SetParametersByValue (const ParametersType &p) override
 
virtual OutputCovariantVectorType TransformCovariantVector (const InputCovariantVectorType &) const
 
virtual OutputCovariantVectorType TransformCovariantVector (const InputCovariantVectorType &vector, const InputPointType &point) const
 
virtual OutputVectorPixelType TransformCovariantVector (const InputVectorPixelType &) const
 
virtual OutputVectorPixelType TransformCovariantVector (const InputVectorPixelType &vector, const InputPointType &point) const
 
virtual OutputDiffusionTensor3DType TransformDiffusionTensor3D (const InputDiffusionTensor3DType &) const
 
virtual OutputDiffusionTensor3DType TransformDiffusionTensor3D (const InputDiffusionTensor3DType &tensor, const InputPointType &point) const
 
virtual OutputVectorPixelType TransformDiffusionTensor3D (const InputVectorPixelType &) const
 
virtual OutputVectorPixelType TransformDiffusionTensor3D (const InputVectorPixelType &tensor, const InputPointType &point) const
 
virtual OutputPointType TransformPoint (const InputPointType &) const =0
 
virtual
OutputSymmetricSecondRankTensorType 
TransformSymmetricSecondRankTensor (const InputSymmetricSecondRankTensorType &tensor, const InputPointType &point) const
 
virtual
OutputSymmetricSecondRankTensorType 
TransformSymmetricSecondRankTensor (const InputSymmetricSecondRankTensorType &) const
 
virtual OutputVectorPixelType TransformSymmetricSecondRankTensor (const InputVectorPixelType &) const
 
virtual OutputVectorPixelType TransformSymmetricSecondRankTensor (const InputVectorPixelType &tensor, const InputPointType &point) const
 
virtual OutputVectorType TransformVector (const InputVectorType &) const
 
virtual OutputVectorType TransformVector (const InputVectorType &vector, const InputPointType &point) const
 
virtual OutputVnlVectorType TransformVector (const InputVnlVectorType &) const
 
virtual OutputVnlVectorType TransformVector (const InputVnlVectorType &vector, const InputPointType &point) const
 
virtual OutputVectorPixelType TransformVector (const InputVectorPixelType &) const
 
virtual OutputVectorPixelType TransformVector (const InputVectorPixelType &vector, const InputPointType &point) const
 
virtual void UpdateTransformParameters (const DerivativeType &update, ParametersValueType factor=1.0)
 
- Public Member Functions inherited from itk::Object
unsigned long AddObserver (const EventObject &event, Command *)
 
unsigned long AddObserver (const EventObject &event, Command *) const
 
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
 
virtual void Register () const override
 
void RemoveAllObservers ()
 
void RemoveObserver (unsigned long tag)
 
void SetDebug (bool debugFlag) const
 
void SetMetaDataDictionary (const MetaDataDictionary &rhs)
 
virtual void SetReferenceCount (int) override
 
virtual void UnRegister () const noexceptoverride
 
virtual void SetObjectName (std::string _arg)
 
virtual const std::string & GetObjectName () const
 
- Public Member Functions inherited from itk::LightObject
virtual void Delete ()
 
virtual int GetReferenceCount () const
 
 itkCloneMacro (Self)
 
void Print (std::ostream &os, Indent indent=0) const
 

Static Public Attributes

static const unsigned int SpaceDimension = NDimensions
 
static const unsigned int SplineOrder = VSplineOrder
 
- Static Public Attributes inherited from itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >
static const unsigned int SpaceDimension = NDimensions
 
static const unsigned int SplineOrder = VSplineOrder
 
- Static Public Attributes inherited from itk::Transform< TParametersValueType, NDimensions, NDimensions >
static const unsigned int InputSpaceDimension
 
static const unsigned int OutputSpaceDimension
 

Protected Member Functions

 BSplineDeformableTransform ()
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
virtual ~BSplineDeformableTransform ()
 
- Protected Member Functions inherited from itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >
 BSplineBaseTransform ()
 
void SetFixedParametersFromTransformDomainInformation () const
 
void WrapAsImages ()
 
virtual ~BSplineBaseTransform ()
 
virtual void SetWeightsFunction (WeightsFunctionType *_arg)
 
virtual WeightsFunctionTypeGetModifiableWeightsFunction ()
 
virtual const WeightsFunctionTypeGetWeightsFunction () const
 
- Protected Member Functions inherited from itk::Transform< TParametersValueType, NDimensions, NDimensions >
virtual LightObject::Pointer InternalClone () const override
 
OutputDiffusionTensor3DType PreservationOfPrincipalDirectionDiffusionTensor3DReorientation (const InputDiffusionTensor3DType, const JacobianType) const
 
 Transform ()
 
 Transform (NumberOfParametersType NumberOfParameters)
 
virtual ~Transform ()
 
- Protected Member Functions inherited from itk::TransformBaseTemplate< TParametersValueType >
 TransformBaseTemplate ()
 
virtual ~TransformBaseTemplate ()
 
- Protected Member Functions inherited from itk::Object
 Object ()
 
bool PrintObservers (std::ostream &os, Indent indent) const
 
virtual void SetTimeStamp (const TimeStamp &time)
 
virtual ~Object ()
 
- Protected Member Functions inherited from itk::LightObject
 LightObject ()
 
virtual void PrintHeader (std::ostream &os, Indent indent) const
 
virtual void PrintTrailer (std::ostream &os, Indent indent) const
 
virtual ~LightObject ()
 

Private Member Functions

 BSplineDeformableTransform (const Self &)
 
virtual bool InsideValidRegion (ContinuousIndexType &) const override
 
void operator= (const Self &)
 
virtual void SetCoefficientImageInformationFromFixedParameters () override
 
virtual void SetFixedParametersGridDirectionFromTransformDomainInformation () const override
 
virtual void SetFixedParametersGridOriginFromTransformDomainInformation () const override
 
virtual void SetFixedParametersGridSizeFromTransformDomainInformation () const override
 
virtual void SetFixedParametersGridSpacingFromTransformDomainInformation () const override
 
void UpdateValidGridRegion ()
 

Private Attributes

BulkTransformPointer m_BulkTransform
 
const DirectionTypem_GridDirection
 
const OriginTypem_GridOrigin
 
const RegionTypem_GridRegion
 
const SpacingTypem_GridSpacing
 
unsigned long m_Offset
 
bool m_SplineOrderOdd
 
RegionType m_ValidRegion
 
IndexType m_ValidRegionFirst
 
IndexType m_ValidRegionLast
 
static Pointer New ()
 
virtual ::itk::LightObject::Pointer CreateAnother (void) const override
 

Additional Inherited Members

- Static Public Member Functions inherited from itk::Object
static bool GetGlobalWarningDisplay ()
 
static void GlobalWarningDisplayOff ()
 
static void GlobalWarningDisplayOn ()
 
static Pointer New ()
 
static void SetGlobalWarningDisplay (bool flag)
 
- Static Public Member Functions inherited from itk::LightObject
static void BreakOnError ()
 
static Pointer New ()
 
- Protected Attributes inherited from itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >
CoefficientImageArray m_CoefficientImages
 
ParametersType m_InternalParametersBuffer
 
WeightsFunctionType::Pointer m_WeightsFunction
 
- Protected Attributes inherited from itk::Transform< TParametersValueType, NDimensions, NDimensions >
DirectionChangeMatrix m_DirectionChange
 
FixedParametersType m_FixedParameters
 
ParametersType m_Parameters
 
- Protected Attributes inherited from itk::LightObject
AtomicInt< int > m_ReferenceCount
 

Member Typedef Documentation

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef BulkTransformType::ConstPointer itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::BulkTransformPointer

Definition at line 307 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Transform<TParametersValueType, itkGetStaticConstMacro(SpaceDimension), itkGetStaticConstMacro(SpaceDimension)> itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::BulkTransformType

Definition at line 302 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::CoefficientImageArray itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::CoefficientImageArray

Definition at line 209 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef SmartPointer<const Self> itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::ConstPointer

Definition at line 124 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::ContinuousIndexType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::ContinuousIndexType

Definition at line 249 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::DirectionType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::DirectionType

Definition at line 242 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::FixedParametersType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::FixedParametersType

Standard parameters container.

Definition at line 159 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::ImagePointer itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::ImagePointer

Definition at line 208 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::ImageType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::ImageType

Definition at line 207 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::IndexType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::IndexType

Definition at line 239 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::InputCovariantVectorType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::InputCovariantVectorType

Standard covariant vector type for this class.

Definition at line 173 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Point<TParametersValueType, itkGetStaticConstMacro( SpaceDimension )> itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::InputPointType

Standard coordinate point type for this class.

Definition at line 181 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::InputVectorType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::InputVectorType

Standard vector type for this class.

Definition at line 169 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::InputVnlVectorType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::InputVnlVectorType

Standard vnl_vector type for this class.

Definition at line 177 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::JacobianType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::JacobianType

Standard Jacobian container.

Definition at line 163 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::MeshSizeType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::MeshSizeType

Definition at line 278 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::NumberOfParametersType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::NumberOfParametersType

The number of parameters defininig this transform.

Definition at line 166 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::OriginType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::OriginType

Definition at line 243 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::OutputCovariantVectorType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::OutputCovariantVectorType

Definition at line 174 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Point<TParametersValueType, itkGetStaticConstMacro( SpaceDimension )> itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::OutputPointType

Standard coordinate point type for this class.

Definition at line 182 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::OutputVectorType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::OutputVectorType

Definition at line 170 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::OutputVnlVectorType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::OutputVnlVectorType

Definition at line 178 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::ParameterIndexArrayType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::ParameterIndexArrayType

Parameter index array type.

Definition at line 252 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::ParametersType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::ParametersType

Definition at line 160 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::ParametersValueType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::ParametersValueType

Parameters as SpaceDimension number of images.

Definition at line 206 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::SpacingType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::PhysicalDimensionsType

Definition at line 275 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::PixelType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::PixelType

Definition at line 276 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef SmartPointer<Self> itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::Pointer

Definition at line 123 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::RegionType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::RegionType

Typedefs for specifying the extent of the grid.

Definition at line 237 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef TParametersValueType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::ScalarType

Standard scalar type for this class.

Definition at line 156 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef BSplineDeformableTransform itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::Self

Standard class typedefs.

Definition at line 121 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::SizeType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SizeType

Definition at line 240 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::SpacingType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SpacingType

Definition at line 241 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef BSplineBaseTransform<TParametersValueType,NDimensions,VSplineOrder> itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::Superclass

Definition at line 122 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::WeightsFunctionType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::WeightsFunctionType

Interpolation weights function type.

Definition at line 246 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
typedef Superclass::WeightsType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::WeightsType

Definition at line 248 of file itkBSplineDeformableTransform.h.

Constructor & Destructor Documentation

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::BSplineDeformableTransform ( )
protected
template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::~BSplineDeformableTransform ( )
protectedvirtual
template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::BSplineDeformableTransform ( const Self )
private

Member Function Documentation

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::ComputeJacobianWithRespectToParameters ( const InputPointType ,
JacobianType  
) const
overridevirtual
template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual ::itk::LightObject::Pointer itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::CreateAnother ( void  ) const
inlineoverridevirtual
template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual const BulkTransformType* itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::GetBulkTransform ( ) const
virtual

This method specifies the bulk transform to be applied. The default is the identity transform.

Referenced by itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::CreateAnother().

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual DirectionType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::GetGridDirection ( ) const
virtual

Function to retrieve the transform domain direction.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual OriginType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::GetGridOrigin ( ) const
virtual

Function to retrieve the transform domain origin.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual RegionType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::GetGridRegion ( ) const
virtual

Function to retrieve the transform domain mesh size.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual SpacingType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::GetGridSpacing ( ) const
virtual

This method retrieve the grid spacing or resolution.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual const char* itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::GetNameOfClass ( ) const
virtual

Run-time type information (and related methods).

Reimplemented from itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual NumberOfParametersType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::GetNumberOfParameters ( ) const
overridevirtual

Return the number of parameters that completely define the Transfom

Implements itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
NumberOfParametersType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::GetNumberOfParametersPerDimension ( ) const
overridevirtual

Return the number of parameters per dimension

Implements itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual const RegionType& itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::GetValidRegion ( ) const
virtual

Return the region of the grid wholly within the support region

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual bool itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::InsideValidRegion ( ContinuousIndexType ) const
overrideprivatevirtual

Check if a continuous index is inside the valid region.

Implements itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::itkCloneMacro ( Self  )

implement type-specific clone method

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
static Pointer itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::New ( )
static

New macro for creation of through the object factory.

Referenced by itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::CreateAnother().

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::operator= ( const Self )
private
template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::PrintSelf ( std::ostream &  os,
Indent  indent 
) const
overrideprotectedvirtual
template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetBulkTransform ( const BulkTransformType _arg)
virtual

This method specifies the bulk transform to be applied. The default is the identity transform.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetCoefficientImageInformationFromFixedParameters ( )
overrideprivatevirtual

Construct control point grid size from transform domain information

Implements itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetCoefficientImages ( const CoefficientImageArray images)
overridevirtual

Set the array of coefficient images.

This is an alternative API for setting the BSpline coefficients as an array of SpaceDimension images. The fixed parameters are taken from the first image. It is assumed that the buffered region of all the subsequent images are the same as the first image. Note that no error checking is done.

Warning: use either the SetParameters() or SetCoefficientImages() API. Mixing the two modes may results in unexpected results.

Implements itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetFixedParameters ( const FixedParametersType parameters)
overridevirtual

This method sets the fixed parameters of the transform. For a BSpline deformation transform, the parameters are the following: Grid Size, Grid Origin, and Grid Spacing

The fixed parameters are the three times the size of the templated dimensions. This function has the effect of make the following calls: transform->SetGridSpacing( spacing ); transform->SetGridOrigin( origin ); transform->SetGridDirection( direction ); transform->SetGridRegion( bsplineRegion );

This function was added to allow the transform to work with the itkTransformReader/Writer I/O filters.

Implements itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetFixedParametersGridDirectionFromTransformDomainInformation ( ) const
overrideprivatevirtual

Construct control point grid direction from transform domain information

Implements itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetFixedParametersGridOriginFromTransformDomainInformation ( ) const
overrideprivatevirtual

Construct control point grid origin from transform domain information

Implements itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetFixedParametersGridSizeFromTransformDomainInformation ( ) const
overrideprivatevirtual

Construct control point grid size from transform domain information

Implements itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetFixedParametersGridSpacingFromTransformDomainInformation ( ) const
overrideprivatevirtual

Construct control point grid spacing from transform domain information

Implements itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetGridDirection ( const DirectionType )
virtual

Function to specify the transform domain direction.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetGridOrigin ( const OriginType )
virtual

Function to specify the transform domain origin.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetGridRegion ( const RegionType )
virtual

Function to specify the transform domain mesh size.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SetGridSpacing ( const SpacingType )
virtual

This method specifies the grid spacing or resolution.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
virtual void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::TransformPoint ( const InputPointType inputPoint,
OutputPointType outputPoint,
WeightsType weights,
ParameterIndexArrayType indices,
bool &  inside 
) const
overridevirtual

Transform points by a BSpline deformable transformation. On return, weights contains the interpolation weights used to compute the deformation and indices of the x (zeroth) dimension coefficient parameters in the support region used to compute the deformation. Parameter indices for the i-th dimension can be obtained by adding ( i * this->GetNumberOfParametersPerDimension() ) to the indices array.

Implements itk::BSplineBaseTransform< TParametersValueType, NDimensions, VSplineOrder >.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
void itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::UpdateValidGridRegion ( )
private

Member Data Documentation

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
BulkTransformPointer itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::m_BulkTransform
private

The bulk transform.

Definition at line 362 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
const DirectionType& itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::m_GridDirection
private

Definition at line 359 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
const OriginType& itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::m_GridOrigin
private

Definition at line 357 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
const RegionType& itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::m_GridRegion
private

The variables defining the coefficient grid domain for the InternalParametersBuffer are taken from the m_CoefficientImages[0] image, and must be kept in sync with them. by using references to that instance, this is more naturally enforced and does not introduce a speed penalty of dereferencing through the pointers (although it does enforce some internal class synchronization).

Definition at line 356 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
const SpacingType& itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::m_GridSpacing
private

Definition at line 358 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
unsigned long itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::m_Offset
private

Variables defining the interpolation support region.

Definition at line 367 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
bool itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::m_SplineOrderOdd
private

Definition at line 368 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
RegionType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::m_ValidRegion
private

Definition at line 364 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
IndexType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::m_ValidRegionFirst
private

Definition at line 370 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
IndexType itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::m_ValidRegionLast
private

Definition at line 369 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
const unsigned int itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SpaceDimension = NDimensions
static

Dimension of the domain space.

Definition at line 150 of file itkBSplineDeformableTransform.h.

template<typename TParametersValueType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
const unsigned int itk::BSplineDeformableTransform< TParametersValueType, NDimensions, VSplineOrder >::SplineOrder = VSplineOrder
static

The BSpline order.

Definition at line 153 of file itkBSplineDeformableTransform.h.


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