ITK  6.0.0
Insight Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Attributes | List of all members
itk::Rigid3DPerspectiveTransform< TParametersValueType > Class Template Reference

#include <itkRigid3DPerspectiveTransform.h>

Detailed Description

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

Rigid3DTramsform of a vector space (e.g. space coordinates)

This transform applies a rotation and translation to the 3D space followed by a projection to 2D space along the Z axis.

Definition at line 38 of file itkRigid3DPerspectiveTransform.h.

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

Public Types

using AngleType = typename VersorType::ValueType
 
using AxisType = typename VersorType::VectorType
 
using AxisValueType = typename AxisType::ValueType
 
using ConstPointer = SmartPointer< const Self >
 
using FixedParametersValueType = typename FixedParametersType::ValueType
 
using InputCovariantVectorType = CovariantVector< TParametersValueType, VInputDimension >
 
using InputPointType = Point< TParametersValueType, Self::InputSpaceDimension >
 
using InputVectorType = Vector< TParametersValueType, Self::InputSpaceDimension >
 
using InputVnlVectorType = vnl_vector_fixed< TParametersValueType, VInputDimension >
 
using InverseJacobianPositionType = vnl_matrix_fixed< ParametersValueType, VInputDimension, VOutputDimension >
 
using JacobianPositionType = vnl_matrix_fixed< ParametersValueType, VOutputDimension, VInputDimension >
 
using JacobianType = Array2D< ParametersValueType >
 
using MatrixType = Matrix< TParametersValueType, Self::InputSpaceDimension, Self::InputSpaceDimension >
 
using OffsetType = Vector< TParametersValueType, Self::InputSpaceDimension >
 
using OffsetValueType = typename OffsetType::ValueType
 
using OutputCovariantVectorType = CovariantVector< TParametersValueType, VOutputDimension >
 
using OutputPointType = Point< TParametersValueType, Self::OutputSpaceDimension >
 
using OutputVectorType = Vector< TParametersValueType, Self::OutputSpaceDimension >
 
using OutputVnlVectorType = vnl_vector_fixed< TParametersValueType, VOutputDimension >
 
using ParametersValueType = typename ParametersType::ValueType
 
using Pointer = SmartPointer< Self >
 
using ScalarType = ParametersValueType
 
using Self = Rigid3DPerspectiveTransform
 
using Superclass = Transform< TParametersValueType, Self::InputSpaceDimension, Self::OutputSpaceDimension >
 
using VersorType = Versor< TParametersValueType >
 
using VnlQuaternionType = vnl_quaternion< TParametersValueType >
 
- Public Types inherited from itk::Transform< TParametersValueType, 3, 2 >
using ConstPointer = SmartPointer< const Self >
 
using DerivativeType = Array< ParametersValueType >
 
using DirectionChangeMatrix = Matrix< double, Self::OutputSpaceDimension, Self::InputSpaceDimension >
 
using InputCovariantVectorType = CovariantVector< TParametersValueType, VInputDimension >
 
using InputDiffusionTensor3DType = DiffusionTensor3D< TParametersValueType >
 
using InputDirectionMatrix = Matrix< double, Self::InputSpaceDimension, Self::InputSpaceDimension >
 
using InputPointType = Point< TParametersValueType, VInputDimension >
 
using InputSymmetricSecondRankTensorType = SymmetricSecondRankTensor< TParametersValueType, VInputDimension >
 
using InputVectorPixelType = VariableLengthVector< TParametersValueType >
 
using InputVectorType = Vector< TParametersValueType, VInputDimension >
 
using InputVnlVectorType = vnl_vector_fixed< TParametersValueType, VInputDimension >
 
using InverseJacobianPositionType = vnl_matrix_fixed< ParametersValueType, VInputDimension, VOutputDimension >
 
using InverseTransformBasePointer = typename InverseTransformBaseType::Pointer
 
using InverseTransformBaseType = Transform< TParametersValueType, VOutputDimension, VInputDimension >
 
using JacobianPositionType = vnl_matrix_fixed< ParametersValueType, VOutputDimension, VInputDimension >
 
using JacobianType = Array2D< ParametersValueType >
 
using MatrixType = Matrix< TParametersValueType, Self::OutputSpaceDimension, Self::InputSpaceDimension >
 
using OutputCovariantVectorType = CovariantVector< TParametersValueType, VOutputDimension >
 
using OutputDiffusionTensor3DType = DiffusionTensor3D< TParametersValueType >
 
using OutputDirectionMatrix = Matrix< double, Self::OutputSpaceDimension, Self::OutputSpaceDimension >
 
using OutputPointType = Point< TParametersValueType, VOutputDimension >
 
using OutputSymmetricSecondRankTensorType = SymmetricSecondRankTensor< TParametersValueType, VOutputDimension >
 
using OutputVectorPixelType = VariableLengthVector< TParametersValueType >
 
using OutputVectorType = Vector< TParametersValueType, VOutputDimension >
 
using OutputVnlVectorType = vnl_vector_fixed< TParametersValueType, VOutputDimension >
 
using Pointer = SmartPointer< Self >
 
using ScalarType = ParametersValueType
 
using Self = Transform
 
using Superclass = TransformBaseTemplate< TParametersValueType >
 

Public Member Functions

void ComputeJacobianWithRespectToParameters (const InputPointType &p, JacobianType &jacobian) const override
 
virtual void ComputeJacobianWithRespectToPosition (const InputPointType &, JacobianPositionType &) const
 
void ComputeJacobianWithRespectToPosition (const InputPointType &, JacobianPositionType &) const override
 
void ComputeMatrix ()
 
double GetFocalDistance () const
 
const char * GetNameOfClass () const override
 
const OffsetTypeGetOffset () const
 
const ParametersType & GetParameters () const override
 
const VersorTypeGetRotation () const
 
const MatrixTypeGetRotationMatrix () const
 
void SetFixedParameters (const FixedParametersType &) override
 
void SetFocalDistance (TParametersValueType focalDistance)
 
void SetOffset (const OffsetType &offset)
 
void SetParameters (const ParametersType &parameters) override
 
void SetRotation (const Vector< TParametersValueType, 3 > &axis, double angle)
 
void SetRotation (const VersorType &rotation)
 
virtual OutputCovariantVectorType TransformCovariantVector (const InputCovariantVectorType &) const
 
OutputCovariantVectorType TransformCovariantVector (const InputCovariantVectorType &) const override
 
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
 
OutputPointType TransformPoint (const InputPointType &point) const override
 
virtual OutputVectorPixelType TransformVector (const InputVectorPixelType &) const
 
virtual OutputVectorPixelType TransformVector (const InputVectorPixelType &vector, const InputPointType &point) const
 
virtual OutputVectorType TransformVector (const InputVectorType &) const
 
OutputVectorType TransformVector (const InputVectorType &) const override
 
virtual OutputVectorType TransformVector (const InputVectorType &vector, const InputPointType &point) const
 
virtual OutputVnlVectorType TransformVector (const InputVnlVectorType &) const
 
OutputVnlVectorType TransformVector (const InputVnlVectorType &) const override
 
virtual OutputVnlVectorType TransformVector (const InputVnlVectorType &vector, const InputPointType &point) const
 
virtual const OffsetTypeGetFixedOffset () const
 
virtual void SetFixedOffset (OffsetType _arg)
 
virtual void SetCenterOfRotation (InputPointType _arg)
 
virtual const InputPointTypeGetCenterOfRotation () const
 
- Public Member Functions inherited from itk::Transform< TParametersValueType, 3, 2 >
virtual void ComputeJacobianWithRespectToParametersCachedTemporaries (const InputPointType &p, JacobianType &jacobian, JacobianType &) const
 
void CopyInFixedParameters (const FixedParametersValueType *const begin, const FixedParametersValueType *const end) override
 
void CopyInParameters (const ParametersValueType *const begin, const ParametersValueType *const end) override
 
const FixedParametersType & GetFixedParameters () const override
 
unsigned int GetInputSpaceDimension () const override
 
bool GetInverse (Self *) const
 
virtual InverseTransformBasePointer GetInverseTransform () const
 
const char * GetNameOfClass () const override
 
virtual NumberOfParametersType GetNumberOfFixedParameters () const
 
virtual NumberOfParametersType GetNumberOfLocalParameters () const
 
NumberOfParametersType GetNumberOfParameters () const override
 
unsigned int GetOutputSpaceDimension () const override
 
const ParametersType & GetParameters () const override
 
TransformCategoryEnum GetTransformCategory () const override
 
std::string GetTransformTypeAsString () const override
 
virtual bool IsLinear () const
 
 itkCloneMacro (Self)
 
void SetParametersByValue (const ParametersType &p) override
 
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 &inputTensor, const InputPointType &point) const
 
virtual OutputVectorPixelType TransformDiffusionTensor3D (const InputVectorPixelType &) const
 
virtual OutputVectorPixelType TransformDiffusionTensor3D (const InputVectorPixelType &inputTensor, const InputPointType &point) const
 
virtual OutputPointType TransformPoint (const InputPointType &) const=0
 
virtual OutputSymmetricSecondRankTensorType TransformSymmetricSecondRankTensor (const InputSymmetricSecondRankTensorType &) const
 
virtual OutputSymmetricSecondRankTensorType TransformSymmetricSecondRankTensor (const InputSymmetricSecondRankTensorType &inputTensor, const InputPointType &point) const
 
virtual OutputVectorPixelType TransformSymmetricSecondRankTensor (const InputVectorPixelType &) const
 
virtual OutputVectorPixelType TransformSymmetricSecondRankTensor (const InputVectorPixelType &inputTensor, const InputPointType &point) const
 
virtual OutputVectorPixelType TransformVector (const InputVectorPixelType &) const
 
virtual OutputVectorPixelType TransformVector (const InputVectorPixelType &vector, 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 &vector, const InputPointType &point) const
 
virtual void UpdateTransformParameters (const DerivativeType &update, ParametersValueType factor=1.0)
 
virtual void ComputeJacobianWithRespectToParameters (const InputPointType &, JacobianType &) const=0
 
virtual void ComputeJacobianWithRespectToPosition (const InputPointType &, JacobianPositionType &) const
 
 itkLegacyMacro (virtual void ComputeJacobianWithRespectToPosition(const InputPointType &x, JacobianType &jacobian) const)
 
 itkLegacyMacro (virtual void ComputeInverseJacobianWithRespectToPosition(const InputPointType &x, JacobianType &jacobian) const)
 
virtual void ComputeInverseJacobianWithRespectToPosition (const InputPointType &pnt, InverseJacobianPositionType &jacobian) const
 
std::enable_if_t< TImage::ImageDimension==VInputDimension &&TImage::ImageDimension==VOutputDimension, void > ApplyToImageMetadata (TImage *image) const
 
std::enable_if_t< TImage::ImageDimension==VInputDimension &&TImage::ImageDimension==VOutputDimension, void > ApplyToImageMetadata (SmartPointer< TImage > image) const
 

Static Public Member Functions

static Pointer New ()
 

Static Public Attributes

static constexpr unsigned int InputSpaceDimension = 3
 
static constexpr unsigned int OutputSpaceDimension = 2
 
static constexpr unsigned int ParametersDimension = 6
 
static constexpr unsigned int SpaceDimension = 3
 
- Static Public Attributes inherited from itk::Transform< TParametersValueType, 3, 2 >
static constexpr unsigned int InputSpaceDimension
 
static constexpr unsigned int OutputSpaceDimension
 

Protected Member Functions

void PrintSelf (std::ostream &os, Indent indent) const override
 
 Rigid3DPerspectiveTransform ()
 
 ~Rigid3DPerspectiveTransform () override=default
 
- Protected Member Functions inherited from itk::Transform< TParametersValueType, 3, 2 >
LightObject::Pointer InternalClone () const override
 
OutputDiffusionTensor3DType PreservationOfPrincipalDirectionDiffusionTensor3DReorientation (const InputDiffusionTensor3DType &, const InverseJacobianPositionType &) const
 
 Transform ()=default
 
 Transform (NumberOfParametersType numberOfParameters)
 
 ~Transform () override=default
 

Private Attributes

InputPointType m_CenterOfRotation {}
 
OffsetType m_FixedOffset {}
 
TParametersValueType m_FocalDistance {}
 
OffsetType m_Offset {}
 
MatrixType m_RotationMatrix {}
 
VersorType m_Versor {}
 

Additional Inherited Members

- Static Protected Member Functions inherited from itk::Transform< TParametersValueType, 3, 2 >
static InverseTransformBasePointer InvertTransform (const TTransform &transform)
 
- Protected Attributes inherited from itk::Transform< TParametersValueType, 3, 2 >
FixedParametersType m_FixedParameters
 
ParametersType m_Parameters
 

Member Typedef Documentation

◆ AngleType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::AngleType = typename VersorType::ValueType

Definition at line 107 of file itkRigid3DPerspectiveTransform.h.

◆ AxisType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::AxisType = typename VersorType::VectorType

Definition at line 106 of file itkRigid3DPerspectiveTransform.h.

◆ AxisValueType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::AxisValueType = typename AxisType::ValueType

Definition at line 108 of file itkRigid3DPerspectiveTransform.h.

◆ ConstPointer

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

Definition at line 56 of file itkRigid3DPerspectiveTransform.h.

◆ FixedParametersValueType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::FixedParametersValueType = typename FixedParametersType::ValueType

Definition at line 69 of file itkRigid3DPerspectiveTransform.h.

◆ InputCovariantVectorType

template<typename TParametersValueType = double>
using itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::InputCovariantVectorType = CovariantVector<TParametersValueType, VInputDimension>

Standard covariant vector type for this class

Definition at line 152 of file itkTransform.h.

◆ InputPointType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::InputPointType = Point<TParametersValueType, Self::InputSpaceDimension>

Standard coordinate point type for this class.

Definition at line 94 of file itkRigid3DPerspectiveTransform.h.

◆ InputVectorType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::InputVectorType = Vector<TParametersValueType, Self::InputSpaceDimension>

Standard vector type for this class.

Definition at line 86 of file itkRigid3DPerspectiveTransform.h.

◆ InputVnlVectorType

template<typename TParametersValueType = double>
using itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::InputVnlVectorType = vnl_vector_fixed<TParametersValueType, VInputDimension>

Standard vnl_vector type for this class.

Definition at line 156 of file itkTransform.h.

◆ InverseJacobianPositionType

template<typename TParametersValueType = double>
using itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::InverseJacobianPositionType = vnl_matrix_fixed<ParametersValueType, VInputDimension, VOutputDimension>

Definition at line 132 of file itkTransform.h.

◆ JacobianPositionType

template<typename TParametersValueType = double>
using itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::JacobianPositionType = vnl_matrix_fixed<ParametersValueType, VOutputDimension, VInputDimension>

Definition at line 131 of file itkTransform.h.

◆ JacobianType

template<typename TParametersValueType = double>
using itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::JacobianType = Array2D<ParametersValueType>

Type of the Jacobian matrix.

Definition at line 130 of file itkTransform.h.

◆ MatrixType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::MatrixType = Matrix<TParametersValueType, Self::InputSpaceDimension, Self::InputSpaceDimension>

Standard matrix type for this class.

Definition at line 79 of file itkRigid3DPerspectiveTransform.h.

◆ OffsetType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::OffsetType = Vector<TParametersValueType, Self::InputSpaceDimension>

Standard vector type for this class.

Definition at line 82 of file itkRigid3DPerspectiveTransform.h.

◆ OffsetValueType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::OffsetValueType = typename OffsetType::ValueType

Definition at line 83 of file itkRigid3DPerspectiveTransform.h.

◆ OutputCovariantVectorType

template<typename TParametersValueType = double>
using itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::OutputCovariantVectorType = CovariantVector<TParametersValueType, VOutputDimension>

Definition at line 153 of file itkTransform.h.

◆ OutputPointType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::OutputPointType = Point<TParametersValueType, Self::OutputSpaceDimension>

Definition at line 95 of file itkRigid3DPerspectiveTransform.h.

◆ OutputVectorType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::OutputVectorType = Vector<TParametersValueType, Self::OutputSpaceDimension>

Definition at line 87 of file itkRigid3DPerspectiveTransform.h.

◆ OutputVnlVectorType

template<typename TParametersValueType = double>
using itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::OutputVnlVectorType = vnl_vector_fixed<TParametersValueType, VOutputDimension>

Definition at line 157 of file itkTransform.h.

◆ ParametersValueType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::ParametersValueType = typename ParametersType::ValueType

Definition at line 71 of file itkRigid3DPerspectiveTransform.h.

◆ Pointer

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

Definition at line 55 of file itkRigid3DPerspectiveTransform.h.

◆ ScalarType

template<typename TParametersValueType = double>
using itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::ScalarType = ParametersValueType

Type of the scalar representing coordinate and vector elements.

Definition at line 127 of file itkTransform.h.

◆ Self

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

Standard class type aliases.

Definition at line 52 of file itkRigid3DPerspectiveTransform.h.

◆ Superclass

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::Superclass = Transform<TParametersValueType, Self::InputSpaceDimension, Self::OutputSpaceDimension>

Definition at line 53 of file itkRigid3DPerspectiveTransform.h.

◆ VersorType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::VersorType = Versor<TParametersValueType>

Versor type.

Definition at line 105 of file itkRigid3DPerspectiveTransform.h.

◆ VnlQuaternionType

template<typename TParametersValueType = double>
using itk::Rigid3DPerspectiveTransform< TParametersValueType >::VnlQuaternionType = vnl_quaternion<TParametersValueType>

Standard vnl_quaternion type.

Definition at line 98 of file itkRigid3DPerspectiveTransform.h.

Constructor & Destructor Documentation

◆ Rigid3DPerspectiveTransform()

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

◆ ~Rigid3DPerspectiveTransform()

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

Member Function Documentation

◆ ComputeJacobianWithRespectToParameters()

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

Compute the Jacobian Matrix of the transformation at one point, allowing for thread-safety.

◆ ComputeJacobianWithRespectToPosition() [1/2]

template<typename TParametersValueType = double>
virtual void itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::ComputeJacobianWithRespectToPosition
inline

This provides the ability to get a local jacobian value in a dense/local transform, e.g. DisplacementFieldTransform. For such transforms it would be unclear what parameters would refer to. Generally, global transforms should return an identity jacobian since there is no change with respect to position.

Definition at line 528 of file itkTransform.h.

◆ ComputeJacobianWithRespectToPosition() [2/2]

template<typename TParametersValueType = double>
void itk::Rigid3DPerspectiveTransform< TParametersValueType >::ComputeJacobianWithRespectToPosition ( const InputPointType ,
JacobianPositionType  
) const
inlineoverride

Definition at line 230 of file itkRigid3DPerspectiveTransform.h.

References GetNameOfClass().

◆ ComputeMatrix()

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

Compute the matrix.

◆ GetCenterOfRotation()

template<typename TParametersValueType = double>
virtual const InputPointType& itk::Rigid3DPerspectiveTransform< TParametersValueType >::GetCenterOfRotation ( ) const
virtual

Set the center of Rotation

◆ GetFixedOffset()

template<typename TParametersValueType = double>
virtual const OffsetType& itk::Rigid3DPerspectiveTransform< TParametersValueType >::GetFixedOffset ( ) const
virtual

Set a fixed offset: this allow to center the object to be transformed

◆ GetFocalDistance()

template<typename TParametersValueType = double>
double itk::Rigid3DPerspectiveTransform< TParametersValueType >::GetFocalDistance ( ) const
inline

Return the Focal Distance

Definition at line 177 of file itkRigid3DPerspectiveTransform.h.

◆ GetNameOfClass()

template<typename TParametersValueType = double>
const char* itk::Rigid3DPerspectiveTransform< TParametersValueType >::GetNameOfClass ( ) const
override

◆ GetOffset()

template<typename TParametersValueType = double>
const OffsetType& itk::Rigid3DPerspectiveTransform< TParametersValueType >::GetOffset ( ) const
inline

Get offset of an Rigid3DPerspectiveTransform This method returns the value of the offset of the Rigid3DPerspectiveTransform.

Definition at line 114 of file itkRigid3DPerspectiveTransform.h.

◆ GetParameters()

template<typename TParametersValueType = double>
const ParametersType& itk::Rigid3DPerspectiveTransform< TParametersValueType >::GetParameters ( ) const
override

◆ GetRotation()

template<typename TParametersValueType = double>
const VersorType& itk::Rigid3DPerspectiveTransform< TParametersValueType >::GetRotation ( ) const
inline

Get rotation from an Rigid3DPerspectiveTransform. This method returns the value of the rotation of the Rigid3DPerspectiveTransform.

Definition at line 123 of file itkRigid3DPerspectiveTransform.h.

◆ GetRotationMatrix()

template<typename TParametersValueType = double>
const MatrixType& itk::Rigid3DPerspectiveTransform< TParametersValueType >::GetRotationMatrix ( ) const
inline

Return the rotation matrix

Definition at line 215 of file itkRigid3DPerspectiveTransform.h.

◆ New()

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

New macro for creation of through a Smart Pointer.

◆ PrintSelf()

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

◆ SetCenterOfRotation()

template<typename TParametersValueType = double>
virtual void itk::Rigid3DPerspectiveTransform< TParametersValueType >::SetCenterOfRotation ( InputPointType  _arg)
virtual

Set the center of Rotation

◆ SetFixedOffset()

template<typename TParametersValueType = double>
virtual void itk::Rigid3DPerspectiveTransform< TParametersValueType >::SetFixedOffset ( OffsetType  _arg)
virtual

Set a fixed offset: this allow to center the object to be transformed

◆ SetFixedParameters()

template<typename TParametersValueType = double>
void itk::Rigid3DPerspectiveTransform< TParametersValueType >::SetFixedParameters ( const FixedParametersType &  )
inlineoverridevirtual

Set the fixed parameters and update internal transformation. This transform has no fixed parameters

Implements itk::Transform< TParametersValueType, 3, 2 >.

Definition at line 142 of file itkRigid3DPerspectiveTransform.h.

◆ SetFocalDistance()

template<typename TParametersValueType = double>
void itk::Rigid3DPerspectiveTransform< TParametersValueType >::SetFocalDistance ( TParametersValueType  focalDistance)
inline

Set the Focal Distance of the projection This method sets the focal distance for the perspective projection to a value specified by the user.

Definition at line 170 of file itkRigid3DPerspectiveTransform.h.

◆ SetOffset()

template<typename TParametersValueType = double>
void itk::Rigid3DPerspectiveTransform< TParametersValueType >::SetOffset ( const OffsetType offset)
inline

This method sets the offset of an Rigid3DPerspectiveTransform to a value specified by the user.

Definition at line 148 of file itkRigid3DPerspectiveTransform.h.

◆ SetParameters()

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

Set/Get the transformation from a container of parameters. This is typically used by optimizers. There are 6 parameters. The first three represent the versor and the last three represents the offset.

Implements itk::Transform< TParametersValueType, 3, 2 >.

◆ SetRotation() [1/2]

template<typename TParametersValueType = double>
void itk::Rigid3DPerspectiveTransform< TParametersValueType >::SetRotation ( const Vector< TParametersValueType, 3 > &  axis,
double  angle 
)

Set Rotation of the Rigid transform. This method sets the rotation of an Rigid3DTransform to a value specified by the user using the axis of rotation an the angle.

◆ SetRotation() [2/2]

template<typename TParametersValueType = double>
void itk::Rigid3DPerspectiveTransform< TParametersValueType >::SetRotation ( const VersorType rotation)

This method sets the rotation of an Rigid3DPerspectiveTransform to a value specified by the user.

◆ TransformCovariantVector() [1/5]

template<typename TParametersValueType = double>
virtual OutputCovariantVectorType itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::TransformCovariantVector
inline

Method to transform a CovariantVector.

Definition at line 234 of file itkTransform.h.

◆ TransformCovariantVector() [2/5]

template<typename TParametersValueType = double>
OutputCovariantVectorType itk::Rigid3DPerspectiveTransform< TParametersValueType >::TransformCovariantVector ( const InputCovariantVectorType ) const
inlineoverridevirtual

Method to transform a CovariantVector.

Reimplemented from itk::Transform< TParametersValueType, 3, 2 >.

Definition at line 207 of file itkRigid3DPerspectiveTransform.h.

◆ TransformCovariantVector() [3/5]

template<typename TParametersValueType = double>
virtual OutputCovariantVectorType itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::TransformCovariantVector

Method to transform a CovariantVector, using a point. Global transforms can ignore the point parameter. Local transforms (e.g. deformation field transform) must override and provide required behavior. By default, point is ignored and TransformCovariantVector(vector) is called

◆ TransformCovariantVector() [4/5]

template<typename TParametersValueType = double>
virtual OutputVectorPixelType itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::TransformCovariantVector
inline

Method to transform a CovariantVector stored in a VectorImage.

Definition at line 252 of file itkTransform.h.

◆ TransformCovariantVector() [5/5]

template<typename TParametersValueType = double>
virtual OutputVectorPixelType itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::TransformCovariantVector

Method to transform a CovariantVector, using a point. Global transforms can ignore the point parameter. Local transforms (e.g. deformation field transform) must override and provide required behavior. By default, point is ignored and TransformCovariantVector(vector) is called

◆ TransformPoint()

template<typename TParametersValueType = double>
OutputPointType itk::Rigid3DPerspectiveTransform< TParametersValueType >::TransformPoint ( const InputPointType point) const
override

Transform by a Rigid3DPerspectiveTransform. This method applies the transform given by self to a given point, returning the transformed point.

◆ TransformVector() [1/8]

template<typename TParametersValueType = double>
virtual OutputVectorPixelType itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::TransformVector
inline

These vector transforms are not implemented for this transform

Definition at line 218 of file itkTransform.h.

◆ TransformVector() [2/8]

template<typename TParametersValueType = double>
virtual OutputVectorPixelType itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::TransformVector

These vector transforms are not implemented for this transform

◆ TransformVector() [3/8]

template<typename TParametersValueType = double>
virtual OutputVectorType itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::TransformVector
inline

These vector transforms are not implemented for this transform

Definition at line 186 of file itkTransform.h.

◆ TransformVector() [4/8]

template<typename TParametersValueType = double>
OutputVectorType itk::Rigid3DPerspectiveTransform< TParametersValueType >::TransformVector ( const InputVectorType ) const
inlineoverride

Definition at line 192 of file itkRigid3DPerspectiveTransform.h.

◆ TransformVector() [5/8]

template<typename TParametersValueType = double>
virtual OutputVectorType itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::TransformVector

These vector transforms are not implemented for this transform

◆ TransformVector() [6/8]

template<typename TParametersValueType = double>
virtual OutputVnlVectorType itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::TransformVector
inline

These vector transforms are not implemented for this transform

Definition at line 202 of file itkTransform.h.

◆ TransformVector() [7/8]

template<typename TParametersValueType = double>
OutputVnlVectorType itk::Rigid3DPerspectiveTransform< TParametersValueType >::TransformVector ( const InputVnlVectorType ) const
inlineoverridevirtual

Method to transform a vnl_vector.

Reimplemented from itk::Transform< TParametersValueType, 3, 2 >.

Definition at line 198 of file itkRigid3DPerspectiveTransform.h.

◆ TransformVector() [8/8]

template<typename TParametersValueType = double>
virtual OutputVnlVectorType itk::Transform< TParametersValueType, VInputDimension, VOutputDimension >::TransformVector

These vector transforms are not implemented for this transform

Member Data Documentation

◆ InputSpaceDimension

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

Dimension of the domain space.

Definition at line 44 of file itkRigid3DPerspectiveTransform.h.

◆ m_CenterOfRotation

template<typename TParametersValueType = double>
InputPointType itk::Rigid3DPerspectiveTransform< TParametersValueType >::m_CenterOfRotation {}
private

Center of rotation

Definition at line 271 of file itkRigid3DPerspectiveTransform.h.

◆ m_FixedOffset

template<typename TParametersValueType = double>
OffsetType itk::Rigid3DPerspectiveTransform< TParametersValueType >::m_FixedOffset {}
private

Fixed offset

Definition at line 268 of file itkRigid3DPerspectiveTransform.h.

◆ m_FocalDistance

template<typename TParametersValueType = double>
TParametersValueType itk::Rigid3DPerspectiveTransform< TParametersValueType >::m_FocalDistance {}
private

Set Focal distance of the projection.

Definition at line 262 of file itkRigid3DPerspectiveTransform.h.

◆ m_Offset

template<typename TParametersValueType = double>
OffsetType itk::Rigid3DPerspectiveTransform< TParametersValueType >::m_Offset {}
private

Offset of the transformation.

Definition at line 256 of file itkRigid3DPerspectiveTransform.h.

◆ m_RotationMatrix

template<typename TParametersValueType = double>
MatrixType itk::Rigid3DPerspectiveTransform< TParametersValueType >::m_RotationMatrix {}
private

Matrix representation of the rotation.

Definition at line 265 of file itkRigid3DPerspectiveTransform.h.

◆ m_Versor

template<typename TParametersValueType = double>
VersorType itk::Rigid3DPerspectiveTransform< TParametersValueType >::m_Versor {}
private

Rotation of the transformation.

Definition at line 259 of file itkRigid3DPerspectiveTransform.h.

◆ OutputSpaceDimension

template<typename TParametersValueType = double>
constexpr unsigned int itk::Rigid3DPerspectiveTransform< TParametersValueType >::OutputSpaceDimension = 2
staticconstexpr

Definition at line 45 of file itkRigid3DPerspectiveTransform.h.

◆ ParametersDimension

template<typename TParametersValueType = double>
constexpr unsigned int itk::Rigid3DPerspectiveTransform< TParametersValueType >::ParametersDimension = 6
staticconstexpr

Definition at line 49 of file itkRigid3DPerspectiveTransform.h.

◆ SpaceDimension

template<typename TParametersValueType = double>
constexpr unsigned int itk::Rigid3DPerspectiveTransform< TParametersValueType >::SpaceDimension = 3
staticconstexpr

Dimension of parameters.

Definition at line 48 of file itkRigid3DPerspectiveTransform.h.


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