ITK  5.3.0
Insight Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | List of all members
itk::ScaleVersor3DTransform< TParametersValueType > Class Template Reference

#include <itkScaleVersor3DTransform.h>

+ Inheritance diagram for itk::ScaleVersor3DTransform< TParametersValueType >:
+ Collaboration diagram for itk::ScaleVersor3DTransform< TParametersValueType >:

Public Types

using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using ScaleVectorType = Vector< TParametersValueType, 3 >
 
using Self = ScaleVersor3DTransform
 
using Superclass = VersorRigid3DTransform< TParametersValueType >
 
- Public Types inherited from itk::VersorRigid3DTransform< TParametersValueType >
using ConstPointer = SmartPointer< const Self >
 
using DerivativeType = Array< ParametersValueType >
 
using Pointer = SmartPointer< Self >
 
using Self = VersorRigid3DTransform
 
using Superclass = VersorTransform< TParametersValueType >
 
using VectorType = typename VersorType::VectorType
 
- Public Types inherited from itk::VersorTransform< TParametersValueType >
using AngleType = typename VersorType::ValueType
 
using AxisType = typename VersorType::VectorType
 
using AxisValueType = typename AxisType::ValueType
 
using ConstPointer = SmartPointer< const Self >
 
using ParametersValueType = typename ParametersType::ValueType
 
using Pointer = SmartPointer< Self >
 
using Self = VersorTransform
 
using Superclass = Rigid3DTransform< TParametersValueType >
 
using VersorType = Versor< TParametersValueType >
 
using VnlQuaternionType = vnl_quaternion< TParametersValueType >
 
- Public Types inherited from itk::Rigid3DTransform< TParametersValueType >
using ConstPointer = SmartPointer< const Self >
 
using InverseTransformBasePointer = typename InverseTransformBaseType::Pointer
 
using InverseTransformBaseType = typename Superclass::InverseTransformBaseType
 
using Pointer = SmartPointer< Self >
 
using Self = Rigid3DTransform
 
using Superclass = MatrixOffsetTransformBase< TParametersValueType, 3, 3 >
 
- Public Types inherited from itk::MatrixOffsetTransformBase< TParametersValueType, 3, 3 >
using CenterType = InputPointType
 
using ConstPointer = SmartPointer< const Self >
 
using InputCovariantVectorType = CovariantVector< TParametersValueType, Self::InputSpaceDimension >
 
using InputPointType = Point< TParametersValueType, Self::InputSpaceDimension >
 
using InputPointValueType = typename InputPointType::ValueType
 
using InputTensorEigenVectorType = CovariantVector< TParametersValueType, InputDiffusionTensor3DType::Dimension >
 
using InputVectorType = Vector< TParametersValueType, Self::InputSpaceDimension >
 
using InputVnlVectorType = vnl_vector_fixed< TParametersValueType, Self::InputSpaceDimension >
 
using InverseMatrixType = Matrix< TParametersValueType, Self::InputSpaceDimension, Self::OutputSpaceDimension >
 
using InverseTransformBasePointer = typename InverseTransformBaseType::Pointer
 
using InverseTransformBaseType = typename Superclass::InverseTransformBaseType
 
using MatrixType = Matrix< TParametersValueType, Self::OutputSpaceDimension, Self::InputSpaceDimension >
 
using MatrixValueType = typename MatrixType::ValueType
 
using OffsetType = OutputVectorType
 
using OffsetValueType = typename OffsetType::ValueType
 
using OutputCovariantVectorType = CovariantVector< TParametersValueType, Self::OutputSpaceDimension >
 
using OutputPointType = Point< TParametersValueType, Self::OutputSpaceDimension >
 
using OutputPointValueType = typename OutputPointType::ValueType
 
using OutputVectorType = Vector< TParametersValueType, Self::OutputSpaceDimension >
 
using OutputVectorValueType = typename OutputVectorType::ValueType
 
using OutputVnlVectorType = vnl_vector_fixed< TParametersValueType, Self::OutputSpaceDimension >
 
using Pointer = SmartPointer< Self >
 
using Self = MatrixOffsetTransformBase
 
using Superclass = Transform< TParametersValueType, NInputDimensions, NOutputDimensions >
 
using TranslationType = OutputVectorType
 
using TranslationValueType = typename TranslationType::ValueType
 
- Public Types inherited from itk::Transform< TParametersValueType, NInputDimensions, NOutputDimensions >
using ConstPointer = SmartPointer< const Self >
 
using DerivativeType = Array< ParametersValueType >
 
using DirectionChangeMatrix = Matrix< double, Self::OutputSpaceDimension, Self::InputSpaceDimension >
 
using InputCovariantVectorType = CovariantVector< TParametersValueType, NInputDimensions >
 
using InputDiffusionTensor3DType = DiffusionTensor3D< TParametersValueType >
 
using InputDirectionMatrix = Matrix< double, Self::InputSpaceDimension, Self::InputSpaceDimension >
 
using InputPointType = Point< TParametersValueType, NInputDimensions >
 
using InputSymmetricSecondRankTensorType = SymmetricSecondRankTensor< TParametersValueType, NInputDimensions >
 
using InputVectorPixelType = VariableLengthVector< TParametersValueType >
 
using InputVectorType = Vector< TParametersValueType, NInputDimensions >
 
using InputVnlVectorType = vnl_vector_fixed< TParametersValueType, NInputDimensions >
 
using InverseJacobianPositionType = vnl_matrix_fixed< ParametersValueType, NInputDimensions, NOutputDimensions >
 
using InverseTransformBasePointer = typename InverseTransformBaseType::Pointer
 
using InverseTransformBaseType = Transform< TParametersValueType, NOutputDimensions, NInputDimensions >
 
using JacobianPositionType = vnl_matrix_fixed< ParametersValueType, NOutputDimensions, NInputDimensions >
 
using JacobianType = Array2D< ParametersValueType >
 
using MatrixType = Matrix< TParametersValueType, Self::OutputSpaceDimension, Self::InputSpaceDimension >
 
using OutputCovariantVectorType = CovariantVector< TParametersValueType, NOutputDimensions >
 
using OutputDiffusionTensor3DType = DiffusionTensor3D< TParametersValueType >
 
using OutputDirectionMatrix = Matrix< double, Self::OutputSpaceDimension, Self::OutputSpaceDimension >
 
using OutputPointType = Point< TParametersValueType, NOutputDimensions >
 
using OutputSymmetricSecondRankTensorType = SymmetricSecondRankTensor< TParametersValueType, NOutputDimensions >
 
using OutputVectorPixelType = VariableLengthVector< TParametersValueType >
 
using OutputVectorType = Vector< TParametersValueType, NOutputDimensions >
 
using OutputVnlVectorType = vnl_vector_fixed< TParametersValueType, NOutputDimensions >
 
using Pointer = SmartPointer< Self >
 
using ScalarType = ParametersValueType
 
using Self = Transform
 
using Superclass = TransformBaseTemplate< TParametersValueType >
 
- Public Types inherited from itk::TransformBaseTemplate< TParametersValueType >
using ConstPointer = SmartPointer< const Self >
 
using FixedParametersType = OptimizerParameters< FixedParametersValueType >
 
using FixedParametersValueType = double
 
using NumberOfParametersType = IdentifierType
 
using ParametersType = OptimizerParameters< ParametersValueType >
 
using ParametersValueType = TParametersValueType
 
using Pointer = SmartPointer< Self >
 
using Self = TransformBaseTemplate
 
using Superclass = Object
 
using TransformCategoryEnum = TransformBaseTemplateEnums::TransformCategory
 
- Public Types inherited from itk::Object
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = Object
 
using Superclass = LightObject
 
- Public Types inherited from itk::LightObject
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = LightObject
 

Public Member Functions

virtual ::itk::LightObject::Pointer CreateAnother () const
 
virtual const char * GetNameOfClass () const
 
- Public Member Functions inherited from itk::VersorRigid3DTransform< TParametersValueType >
void ComputeJacobianWithRespectToParameters (const InputPointType &p, JacobianType &jacobian) const override
 
virtual ::itk::LightObject::Pointer CreateAnother () const
 
const ParametersTypeGetParameters () const override
 
void SetParameters (const ParametersType &parameters) override
 
void UpdateTransformParameters (const DerivativeType &update, TParametersValueType factor=1.0) override
 
- Public Member Functions inherited from itk::VersorTransform< TParametersValueType >
void ComputeJacobianWithRespectToParameters (const InputPointType &p, JacobianType &jacobian) const override
 
virtual ::itk::LightObject::Pointer CreateAnother () const
 
const ParametersTypeGetParameters () const override
 
virtual const VersorTypeGetVersor () const
 
void SetIdentity () override
 
void SetParameters (const ParametersType &parameters) override
 
void SetRotation (const AxisType &axis, AngleType angle)
 
void SetRotation (const VersorType &versor)
 
- Public Member Functions inherited from itk::Rigid3DTransform< TParametersValueType >
bool MatrixIsOrthogonal (const MatrixType &matrix, const TParametersValueType tolerance=MatrixOrthogonalityTolerance< TParametersValueType >::GetTolerance())
 
void Translate (const OffsetType &offset, bool pre=false)
 
- Public Member Functions inherited from itk::MatrixOffsetTransformBase< TParametersValueType, 3, 3 >
virtual ::itk::LightObject::Pointer CreateAnother () const
 
TransformCategoryEnum GetTransformCategory () const override
 
virtual const MatrixTypeGetMatrix () const
 
void SetOffset (const OutputVectorType &offset)
 
const OutputVectorTypeGetOffset () const
 
void SetCenter (const InputPointType &center)
 
const InputPointTypeGetCenter () const
 
void SetTranslation (const OutputVectorType &translation)
 
const OutputVectorTypeGetTranslation () const
 
void SetParameters (const ParametersType &parameters) override
 
const ParametersTypeGetParameters () const override
 
void SetFixedParameters (const FixedParametersType &) override
 
const FixedParametersTypeGetFixedParameters () const override
 
void Compose (const Self *other, bool pre=false)
 
OutputPointType TransformPoint (const InputPointType &point) const override
 
OutputVectorType TransformVector (const InputVectorType &vect) const override
 
OutputVnlVectorType TransformVector (const InputVnlVectorType &vect) const override
 
OutputVectorPixelType TransformVector (const InputVectorPixelType &vect) const override
 
OutputCovariantVectorType TransformCovariantVector (const InputCovariantVectorType &vec) const override
 
OutputVectorPixelType TransformCovariantVector (const InputVectorPixelType &vect) const override
 
OutputDiffusionTensor3DType TransformDiffusionTensor3D (const InputDiffusionTensor3DType &tensor) const override
 
OutputVectorPixelType TransformDiffusionTensor3D (const InputVectorPixelType &tensor) const override
 
OutputSymmetricSecondRankTensorType TransformSymmetricSecondRankTensor (const InputSymmetricSecondRankTensorType &inputTensor) const override
 
OutputVectorPixelType TransformSymmetricSecondRankTensor (const InputVectorPixelType &inputTensor) const override
 
void ComputeJacobianWithRespectToParameters (const InputPointType &p, JacobianType &jacobian) const override
 
void ComputeJacobianWithRespectToPosition (const InputPointType &x, JacobianPositionType &jac) const override
 
void ComputeInverseJacobianWithRespectToPosition (const InputPointType &x, InverseJacobianPositionType &jac) const override
 
bool GetInverse (Self *inverse) const
 
InverseTransformBasePointer GetInverseTransform () const override
 
bool IsLinear () const override
 
- Public Member Functions inherited from itk::Transform< TParametersValueType, NInputDimensions, NOutputDimensions >
void CopyInFixedParameters (const FixedParametersValueType *const begin, const FixedParametersValueType *const end) override
 
void CopyInParameters (const ParametersValueType *const begin, const ParametersValueType *const end) override
 
unsigned int GetInputSpaceDimension () const override
 
bool GetInverse (Self *) const
 
virtual NumberOfParametersType GetNumberOfFixedParameters () const
 
virtual NumberOfParametersType GetNumberOfLocalParameters () const
 
NumberOfParametersType GetNumberOfParameters () const override
 
unsigned int GetOutputSpaceDimension () const override
 
std::string GetTransformTypeAsString () const override
 
 itkCloneMacro (Self)
 
void SetParametersByValue (const ParametersType &p) override
 
virtual OutputCovariantVectorType TransformCovariantVector (const InputCovariantVectorType &vector, const InputPointType &point) const
 
virtual OutputVectorPixelType TransformCovariantVector (const InputVectorPixelType &vector, const InputPointType &point) const
 
virtual OutputDiffusionTensor3DType TransformDiffusionTensor3D (const InputDiffusionTensor3DType &inputTensor, const InputPointType &point) const
 
virtual OutputVectorPixelType TransformDiffusionTensor3D (const InputVectorPixelType &inputTensor, const InputPointType &point) const
 
virtual OutputSymmetricSecondRankTensorType TransformSymmetricSecondRankTensor (const InputSymmetricSecondRankTensorType &inputTensor, const InputPointType &point) const
 
virtual OutputVectorPixelType TransformSymmetricSecondRankTensor (const InputVectorPixelType &inputTensor, const InputPointType &point) const
 
virtual OutputVectorPixelType TransformVector (const InputVectorPixelType &vector, const InputPointType &point) const
 
virtual OutputVectorType TransformVector (const InputVectorType &vector, const InputPointType &point) const
 
virtual OutputVnlVectorType TransformVector (const InputVnlVectorType &vector, const InputPointType &point) const
 
virtual void ComputeJacobianWithRespectToParametersCachedTemporaries (const InputPointType &p, JacobianType &jacobian, JacobianType &) const
 
 itkLegacyMacro (virtual void ComputeJacobianWithRespectToPosition(const InputPointType &x, JacobianType &jacobian) const)
 
 itkLegacyMacro (virtual void ComputeInverseJacobianWithRespectToPosition(const InputPointType &x, JacobianType &jacobian) const)
 
template<typename TImage >
std::enable_if_t< TImage::ImageDimension==NInputDimensions &&TImage::ImageDimension==NOutputDimensions, void > ApplyToImageMetadata (TImage *image) const
 
template<typename TImage >
std::enable_if_t< TImage::ImageDimension==NInputDimensions &&TImage::ImageDimension==NOutputDimensions, void > ApplyToImageMetadata (SmartPointer< TImage > image) const
 
- Public Member Functions inherited from itk::Object
unsigned long AddObserver (const EventObject &event, Command *)
 
unsigned long AddObserver (const EventObject &event, Command *) const
 
unsigned long AddObserver (const EventObject &event, std::function< void(const EventObject &)> function) 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
 
void Register () const override
 
void RemoveAllObservers ()
 
void RemoveObserver (unsigned long tag)
 
void SetDebug (bool debugFlag) const
 
void SetReferenceCount (int) override
 
void UnRegister () const noexcept override
 
void SetMetaDataDictionary (const MetaDataDictionary &rhs)
 
void SetMetaDataDictionary (MetaDataDictionary &&rrhs)
 
virtual void SetObjectName (std::string _arg)
 
virtual const std::string & GetObjectName () const
 
- Public Member Functions inherited from itk::LightObject
Pointer Clone () const
 
virtual void Delete ()
 
virtual int GetReferenceCount () const
 
void Print (std::ostream &os, Indent indent=0) const
 

Static Public Member Functions

static Pointer New ()
 
- Static Public Member Functions inherited from itk::VersorRigid3DTransform< TParametersValueType >
static Pointer New ()
 
- Static Public Member Functions inherited from itk::VersorTransform< TParametersValueType >
static Pointer New ()
 
- Static Public Member Functions inherited from itk::Rigid3DTransform< TParametersValueType >
static Pointer New ()
 
- Static Public Member Functions inherited from itk::MatrixOffsetTransformBase< TParametersValueType, 3, 3 >
static Pointer New ()
 
- Static Public Member Functions inherited from itk::Object
static bool GetGlobalWarningDisplay ()
 
static void GlobalWarningDisplayOff ()
 
static void GlobalWarningDisplayOn ()
 
static Pointer New ()
 
static void SetGlobalWarningDisplay (bool val)
 
- Static Public Member Functions inherited from itk::LightObject
static void BreakOnError ()
 
static Pointer New ()
 

Static Public Attributes

static constexpr unsigned int InputSpaceDimension = 3
 
static constexpr unsigned int OutputSpaceDimension = 3
 
static constexpr unsigned int ParametersDimension = 9
 
- Static Public Attributes inherited from itk::VersorRigid3DTransform< TParametersValueType >
static constexpr unsigned int InputSpaceDimension = 3
 
static constexpr unsigned int OutputSpaceDimension = 3
 
static constexpr unsigned int ParametersDimension = 6
 
static constexpr unsigned int SpaceDimension = 3
 
- Static Public Attributes inherited from itk::VersorTransform< TParametersValueType >
static constexpr unsigned int InputSpaceDimension = 3
 
static constexpr unsigned int OutputSpaceDimension = 3
 
static constexpr unsigned int ParametersDimension = 3
 
static constexpr unsigned int SpaceDimension = 3
 
- Static Public Attributes inherited from itk::Rigid3DTransform< TParametersValueType >
static constexpr unsigned int InputSpaceDimension = 3
 
static constexpr unsigned int OutputSpaceDimension = 3
 
static constexpr unsigned int ParametersDimension = 12
 
static constexpr unsigned int SpaceDimension = 3
 
- Static Public Attributes inherited from itk::MatrixOffsetTransformBase< TParametersValueType, 3, 3 >
static constexpr unsigned int InputSpaceDimension
 
static constexpr unsigned int OutputSpaceDimension
 
static constexpr unsigned int ParametersDimension
 
- Static Public Attributes inherited from itk::Transform< TParametersValueType, NInputDimensions, NOutputDimensions >
static constexpr unsigned int InputSpaceDimension = NInputDimensions
 
static constexpr unsigned int OutputSpaceDimension = NOutputDimensions
 
ScaleVectorType m_Scale
 
void SetMatrix (const MatrixType &matrix) override
 
void SetMatrix (const MatrixType &matrix, const TParametersValueType tolerance) override
 
void SetParameters (const ParametersType &parameters) override
 
const ParametersTypeGetParameters () const override
 
void SetScale (const ScaleVectorType &scale)
 
virtual const ScaleVectorTypeGetScale () const
 
void SetIdentity () override
 
void ComputeJacobianWithRespectToParameters (const InputPointType &p, JacobianType &jacobian) const override
 
 ScaleVersor3DTransform ()
 
 ScaleVersor3DTransform (const MatrixType &matrix, const OutputVectorType &offset)
 
 ScaleVersor3DTransform (unsigned int parametersDimension)
 
 ~ScaleVersor3DTransform () override=default
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
void SetVarScale (const ScaleVectorType &scale)
 
void ComputeMatrix () override
 
void ComputeMatrixParameters () override
 

Additional Inherited Members

- Protected Member Functions inherited from itk::VersorRigid3DTransform< TParametersValueType >
void PrintSelf (std::ostream &os, Indent indent) const override
 
 VersorRigid3DTransform ()
 
 VersorRigid3DTransform (const MatrixType &matrix, const OutputVectorType &offset)
 
 VersorRigid3DTransform (unsigned int paramDim)
 
 ~VersorRigid3DTransform () override=default
 
- Protected Member Functions inherited from itk::VersorTransform< TParametersValueType >
 VersorTransform (const MatrixType &matrix, const OutputVectorType &offset)
 
 VersorTransform (unsigned int parametersDimension)
 
 VersorTransform ()
 
 ~VersorTransform () override=default
 
void SetVarVersor (const VersorType &newVersor)
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
void ComputeMatrix () override
 
void ComputeMatrixParameters () override
 
- Protected Member Functions inherited from itk::Rigid3DTransform< TParametersValueType >
 Rigid3DTransform ()
 
 Rigid3DTransform (const MatrixType &matrix, const OutputVectorType &offset)
 
 Rigid3DTransform (unsigned int paramDim)
 
 ~Rigid3DTransform () override=default
 
- Protected Member Functions inherited from itk::MatrixOffsetTransformBase< TParametersValueType, 3, 3 >
const InverseMatrixTypeGetInverseMatrix () const
 
 MatrixOffsetTransformBase (const MatrixType &matrix, const OutputVectorType &offset)
 
 MatrixOffsetTransformBase (unsigned int paramDims=ParametersDimension)
 
 ~MatrixOffsetTransformBase () override=default
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
const InverseMatrixTypeGetVarInverseMatrix () const
 
void SetVarInverseMatrix (const InverseMatrixType &matrix) const
 
bool InverseMatrixIsOld () const
 
void SetVarMatrix (const MatrixType &matrix)
 
virtual void ComputeTranslation ()
 
void SetVarTranslation (const OutputVectorType &translation)
 
virtual void ComputeOffset ()
 
void SetVarOffset (const OutputVectorType &offset)
 
void SetVarCenter (const InputPointType &center)
 
virtual bool GetSingular () const
 
- Protected Member Functions inherited from itk::Transform< TParametersValueType, NInputDimensions, NOutputDimensions >
LightObject::Pointer InternalClone () const override
 
 Transform ()
 
 Transform (NumberOfParametersType numberOfParameters)
 
 ~Transform () override=default
 
OutputDiffusionTensor3DType PreservationOfPrincipalDirectionDiffusionTensor3DReorientation (const InputDiffusionTensor3DType &, const InverseJacobianPositionType &) const
 
- Protected Member Functions inherited from itk::TransformBaseTemplate< TParametersValueType >
 TransformBaseTemplate ()=default
 
 ~TransformBaseTemplate () override=default
 
- Protected Member Functions inherited from itk::Object
 Object ()
 
 ~Object () override
 
bool PrintObservers (std::ostream &os, Indent indent) const
 
virtual void SetTimeStamp (const TimeStamp &timeStamp)
 
- 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 ()
 
- Protected Attributes inherited from itk::Transform< TParametersValueType, NInputDimensions, NOutputDimensions >
ParametersType m_Parameters
 
FixedParametersType m_FixedParameters
 
- Protected Attributes inherited from itk::LightObject
std::atomic< int > m_ReferenceCount
 

Detailed Description

template<typename TParametersValueType = double>
class itk::ScaleVersor3DTransform< TParametersValueType >

This transform applies a Versor rotation, translation and anisotropic scale to the space.

The transform can be described as: $ (\textbf{R}_v + \textbf{S})\textbf{x} $ where $\textbf{R}_v$ is the rotation matrix given the versor, and $S=\left( \begin{array}{ccc}s_0-1 & 0 & 0 \\ 0 & s_1-1 & 0 \\ 0 & 0 & s_2-1 \end{array} \right)\ $

Note
This transform's scale parameters are not related to the uniform scaling parameter of the Similarity3DTransform.
Author
Johnson H.J., Harris G., Williams K. University of Iowa Carver College of Medicine, Department of Psychiatry NeuroImaging Center

This implementation was taken from the Insight Journal paper: https://www.insight-journal.org/browse/publication/180

Definition at line 48 of file itkScaleVersor3DTransform.h.

Member Typedef Documentation

◆ ConstPointer

template<typename TParametersValueType = double>
using itk::ScaleVersor3DTransform< TParametersValueType >::ConstPointer = SmartPointer<const Self>

Definition at line 57 of file itkScaleVersor3DTransform.h.

◆ Pointer

template<typename TParametersValueType = double>
using itk::ScaleVersor3DTransform< TParametersValueType >::Pointer = SmartPointer<Self>

Definition at line 56 of file itkScaleVersor3DTransform.h.

◆ ScaleVectorType

template<typename TParametersValueType = double>
using itk::ScaleVersor3DTransform< TParametersValueType >::ScaleVectorType = Vector<TParametersValueType, 3>

Scale Vector Type.

Definition at line 96 of file itkScaleVersor3DTransform.h.

◆ Self

template<typename TParametersValueType = double>
using itk::ScaleVersor3DTransform< TParametersValueType >::Self = ScaleVersor3DTransform

Standard class type aliases.

Definition at line 54 of file itkScaleVersor3DTransform.h.

◆ Superclass

template<typename TParametersValueType = double>
using itk::ScaleVersor3DTransform< TParametersValueType >::Superclass = VersorRigid3DTransform<TParametersValueType>

Definition at line 55 of file itkScaleVersor3DTransform.h.

Constructor & Destructor Documentation

◆ ScaleVersor3DTransform() [1/3]

template<typename TParametersValueType = double>
itk::ScaleVersor3DTransform< TParametersValueType >::ScaleVersor3DTransform ( )
protected

Vector containing the scale.

◆ ScaleVersor3DTransform() [2/3]

template<typename TParametersValueType = double>
itk::ScaleVersor3DTransform< TParametersValueType >::ScaleVersor3DTransform ( const MatrixType matrix,
const OutputVectorType offset 
)
protected

Vector containing the scale.

◆ ScaleVersor3DTransform() [3/3]

template<typename TParametersValueType = double>
itk::ScaleVersor3DTransform< TParametersValueType >::ScaleVersor3DTransform ( unsigned int  parametersDimension)
protected

Vector containing the scale.

◆ ~ScaleVersor3DTransform()

template<typename TParametersValueType = double>
itk::ScaleVersor3DTransform< TParametersValueType >::~ScaleVersor3DTransform ( )
overrideprotecteddefault

Vector containing the scale.

Member Function Documentation

◆ ComputeJacobianWithRespectToParameters()

template<typename TParametersValueType = double>
void itk::ScaleVersor3DTransform< TParametersValueType >::ComputeJacobianWithRespectToParameters ( const InputPointType p,
JacobianType jacobian 
) const
overridevirtual

This method computes the Jacobian matrix of the transformation. given point or vector, returning the transformed point or vector. The rank of the Jacobian will also indicate if the transform is invertible at this point.

Implements itk::Transform< TParametersValueType, NInputDimensions, NOutputDimensions >.

◆ ComputeMatrix()

template<typename TParametersValueType = double>
void itk::ScaleVersor3DTransform< TParametersValueType >::ComputeMatrix ( )
overrideprotectedvirtual

Compute the components of the rotation matrix in the superclass.

Reimplemented from itk::MatrixOffsetTransformBase< TParametersValueType, 3, 3 >.

◆ ComputeMatrixParameters()

template<typename TParametersValueType = double>
void itk::ScaleVersor3DTransform< TParametersValueType >::ComputeMatrixParameters ( )
overrideprotectedvirtual

Vector containing the scale.

Reimplemented from itk::MatrixOffsetTransformBase< TParametersValueType, 3, 3 >.

◆ CreateAnother()

template<typename TParametersValueType = double>
virtual::itk::LightObject::Pointer itk::ScaleVersor3DTransform< TParametersValueType >::CreateAnother ( ) const
virtual

Create an object from an instance, potentially deferring to a factory. This method allows you to create an instance of an object that is exactly the same type as the referring object. This is useful in cases where an object has been cast back to a base class.

Reimplemented from itk::Rigid3DTransform< TParametersValueType >.

◆ GetNameOfClass()

template<typename TParametersValueType = double>
virtual const char* itk::ScaleVersor3DTransform< TParametersValueType >::GetNameOfClass ( ) const
virtual

Run-time type information (and related methods).

Reimplemented from itk::VersorRigid3DTransform< TParametersValueType >.

◆ GetParameters()

template<typename TParametersValueType = double>
const ParametersType& itk::ScaleVersor3DTransform< TParametersValueType >::GetParameters ( ) const
overridevirtual

◆ GetScale()

template<typename TParametersValueType = double>
virtual const ScaleVectorType& itk::ScaleVersor3DTransform< TParametersValueType >::GetScale ( ) const
virtual

Vector containing the scale.

◆ New()

template<typename TParametersValueType = double>
static Pointer itk::ScaleVersor3DTransform< TParametersValueType >::New ( )
static

New macro for creation of through a Smart Pointer.

◆ PrintSelf()

template<typename TParametersValueType = double>
void itk::ScaleVersor3DTransform< TParametersValueType >::PrintSelf ( std::ostream &  os,
Indent  indent 
) const
overrideprotectedvirtual

Vector containing the scale.

Reimplemented from itk::Rigid3DTransform< TParametersValueType >.

◆ SetIdentity()

template<typename TParametersValueType = double>
void itk::ScaleVersor3DTransform< TParametersValueType >::SetIdentity ( )
overridevirtual

Set the internal parameters of the transform in order to represent an Identity transform.

Reimplemented from itk::MatrixOffsetTransformBase< TParametersValueType, 3, 3 >.

◆ SetMatrix() [1/2]

template<typename TParametersValueType = double>
void itk::ScaleVersor3DTransform< TParametersValueType >::SetMatrix ( const MatrixType matrix)
overridevirtual

Directly set the matrix of the transform.

Orthogonality testing is bypassed in this case.

\sa MatrixOffsetTransformBase::SetMatrix() 

Reimplemented from itk::Rigid3DTransform< TParametersValueType >.

◆ SetMatrix() [2/2]

template<typename TParametersValueType = double>
void itk::ScaleVersor3DTransform< TParametersValueType >::SetMatrix ( const MatrixType matrix,
const TParametersValueType  tolerance 
)
overridevirtual

Vector containing the scale.

Reimplemented from itk::Rigid3DTransform< TParametersValueType >.

◆ SetParameters()

template<typename TParametersValueType = double>
void itk::ScaleVersor3DTransform< TParametersValueType >::SetParameters ( const ParametersType parameters)
overridevirtual

Set the transformation from a container of parameters This is typically used by optimizers. There are 9 parameters: 0-2 versor 3-5 translation 6-8 Scale

Reimplemented from itk::Rigid3DTransform< TParametersValueType >.

◆ SetScale()

template<typename TParametersValueType = double>
void itk::ScaleVersor3DTransform< TParametersValueType >::SetScale ( const ScaleVectorType scale)

Set/Get the scale vector. These scale factors are associated to the axis of coordinates.

◆ SetVarScale()

template<typename TParametersValueType = double>
void itk::ScaleVersor3DTransform< TParametersValueType >::SetVarScale ( const ScaleVectorType scale)
inlineprotected

Vector containing the scale.

Definition at line 151 of file itkScaleVersor3DTransform.h.

Member Data Documentation

◆ InputSpaceDimension

template<typename TParametersValueType = double>
constexpr unsigned int itk::ScaleVersor3DTransform< TParametersValueType >::InputSpaceDimension = 3
staticconstexpr

Dimension of parameters.

Definition at line 66 of file itkScaleVersor3DTransform.h.

◆ m_Scale

template<typename TParametersValueType = double>
ScaleVectorType itk::ScaleVersor3DTransform< TParametersValueType >::m_Scale
private

Vector containing the scale.

Definition at line 165 of file itkScaleVersor3DTransform.h.

◆ OutputSpaceDimension

template<typename TParametersValueType = double>
constexpr unsigned int itk::ScaleVersor3DTransform< TParametersValueType >::OutputSpaceDimension = 3
staticconstexpr

Definition at line 67 of file itkScaleVersor3DTransform.h.

◆ ParametersDimension

template<typename TParametersValueType = double>
constexpr unsigned int itk::ScaleVersor3DTransform< TParametersValueType >::ParametersDimension = 9
staticconstexpr

Definition at line 68 of file itkScaleVersor3DTransform.h.


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