Please, help us to better know about our user community by answering the following short survey: https://forms.gle/wpyrxWi18ox9Z5ae9
Eigen  3.3.90 (git rev 30960d485ec7e45b095d3ad206b2dbcc8bc835ba)
Eigen::Transform Class Reference

Detailed Description

Represents an homogeneous transformation in a N dimensional space.

This is defined in the Geometry module.

#include <Eigen/Geometry>
Template Parameters
_Scalarthe scalar type, i.e., the type of the coefficients
_Dimthe dimension of the space
_Modethe type of the transformation. Can be:
  • Affine: the transformation is stored as a (Dim+1)^2 matrix, where the last row is assumed to be [0 ... 0 1].
  • AffineCompact: the transformation is stored as a (Dim)x(Dim+1) matrix.
  • Projective: the transformation is stored as a (Dim+1)^2 matrix without any assumption.
  • Isometry: same as Affine with the additional assumption that the linear part represents a rotation. This assumption is exploited to speed up some functions such as inverse() and rotation().
_Optionshas the same meaning as in class Matrix. It allows to specify DontAlign and/or RowMajor. These Options are passed directly to the underlying matrix type.

The homography is internally represented and stored by a matrix which is available through the matrix() method. To understand the behavior of this class you have to think a Transform object as its internal matrix representation. The chosen convention is right multiply:

v' = T * v

Therefore, an affine transformation matrix M is shaped like this:

$ \left( \begin{array}{cc} linear & translation\\ 0 ... 0 & 1 \end{array} \right) $

Note that for a projective transformation the last row can be anything, and then the interpretation of different parts might be slightly different.

However, unlike a plain matrix, the Transform class provides many features simplifying both its assembly and usage. In particular, it can be composed with any other transformations (Transform,Translation,RotationBase,DiagonalMatrix) and can be directly used to transform implicit homogeneous vectors. All these operations are handled via the operator*. For the composition of transformations, its principle consists to first convert the right/left hand sides of the product to a compatible (Dim+1)^2 matrix and then perform a pure matrix product. Of course, internally, operator* tries to perform the minimal number of operations according to the nature of each terms. Likewise, when applying the transform to points, the latters are automatically promoted to homogeneous vectors before doing the matrix product. The conventions to homogeneous representations are performed as follow:

Translation t (Dim)x(1): $ \left( \begin{array}{cc} I & t \\ 0\,...\,0 & 1 \end{array} \right) $

Rotation R (Dim)x(Dim): $ \left( \begin{array}{cc} R & 0\\ 0\,...\,0 & 1 \end{array} \right) $

Scaling DiagonalMatrix S (Dim)x(Dim): $ \left( \begin{array}{cc} S & 0\\ 0\,...\,0 & 1 \end{array} \right) $

Column point v (Dim)x(1): $ \left( \begin{array}{c} v\\ 1 \end{array} \right) $

Set of column points V1...Vn (Dim)x(n): $ \left( \begin{array}{ccc} v_1 & ... & v_n\\ 1 & ... & 1 \end{array} \right) $

The concatenation of a Transform object with any kind of other transformation always returns a Transform object.

A little exception to the "as pure matrix product" rule is the case of the transformation of non homogeneous vectors by an affine transformation. In that case the last matrix row can be ignored, and the product returns non homogeneous vectors.

Since, for instance, a Dim x Dim matrix is interpreted as a linear transformation, it is not possible to directly transform Dim vectors stored in a Dim x Dim matrix. The solution is either to use a Dim x Dynamic matrix or explicitly request a vector transformation by making the vector homogeneous:

m' = T * m.colwise().homogeneous();

Note that there is zero overhead.

Conversion methods from/to Qt's QMatrix and QTransform are available if the preprocessor token EIGEN_QT_SUPPORT is defined.

This class can be extended with the help of the plugin mechanism described on the page Extending MatrixBase (and other classes) by defining the preprocessor symbol EIGEN_TRANSFORM_PLUGIN.

See also
class Matrix, class Quaternion

Public Types

typedef internal::conditional< int(Mode)==int(AffineCompact), MatrixType &, Block< MatrixType, Dim, HDim > >::type AffinePart
 
typedef internal::conditional< int(Mode)==int(AffineCompact), const MatrixType &, const Block< const MatrixType, Dim, HDim > >::type ConstAffinePart
 
typedef Eigen::Index Index
 
typedef Matrix< Scalar, Dim, Dim, Options > LinearMatrixType
 
typedef Block< MatrixType, Dim, Dim, int(Mode)==(AffineCompact) &&(Options &RowMajor)==0 > LinearPart
 
typedef internal::make_proper_matrix_type< Scalar, Rows, HDim, Options >::type MatrixType
 
typedef _Scalar Scalar
 
typedef Transform< Scalar, Dim, TransformTimeDiagonalMode > TransformTimeDiagonalReturnType
 
typedef Block< MatrixType, Dim, 1,!(internal::traits< MatrixType >::Flags &RowMajorBit)> TranslationPart
 
typedef Translation< Scalar, Dim > TranslationType
 
typedef Matrix< Scalar, Dim, 1 > VectorType
 

Public Member Functions

AffinePart affine ()
 
ConstAffinePart affine () const
 
template<typename NewScalarType >
internal::cast_return_type< Transform, Transform< NewScalarType, Dim, Mode, Options > >::type cast () const
 
template<typename RotationMatrixType , typename ScalingMatrixType >
void computeRotationScaling (RotationMatrixType *rotation, ScalingMatrixType *scaling) const
 
template<typename ScalingMatrixType , typename RotationMatrixType >
void computeScalingRotation (ScalingMatrixType *scaling, RotationMatrixType *rotation) const
 
Scalardata ()
 
const Scalardata () const
 
 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE (_Scalar, _Dim==Dynamic ? Dynamic :(_Dim+1) *(_Dim+1)) enum
 
template<typename PositionDerived , typename OrientationType , typename ScaleDerived >
Transform< Scalar, Dim, Mode, Options > & fromPositionOrientationScale (const MatrixBase< PositionDerived > &position, const OrientationType &orientation, const MatrixBase< ScaleDerived > &scale)
 
Transform inverse (TransformTraits traits=(TransformTraits) Mode) const
 
bool isApprox (const Transform &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const
 
LinearPart linear ()
 
ConstLinearPart linear () const
 
void makeAffine ()
 
MatrixTypematrix ()
 
const MatrixTypematrix () const
 
Scalaroperator() (Index row, Index col)
 
Scalar operator() (Index row, Index col) const
 
template<typename DiagonalDerived >
const TransformTimeDiagonalReturnType operator* (const DiagonalBase< DiagonalDerived > &b) const
 
template<typename OtherDerived >
const internal::transform_right_product_impl< Transform, OtherDerived >::ResultType operator* (const EigenBase< OtherDerived > &other) const
 
const Transform operator* (const Transform &other) const
 
template<int OtherMode, int OtherOptions>
internal::transform_transform_product_impl< Transform, Transform< Scalar, Dim, OtherMode, OtherOptions > >::ResultType operator* (const Transform< Scalar, Dim, OtherMode, OtherOptions > &other) const
 
template<typename OtherDerived >
Transformoperator= (const EigenBase< OtherDerived > &other)
 
Transformoperator= (const QMatrix &other)
 
Transformoperator= (const QTransform &other)
 
template<typename RotationType >
Transform< Scalar, Dim, Mode, Options > & prerotate (const RotationType &rotation)
 
template<typename OtherDerived >
Transform< Scalar, Dim, Mode, Options > & prescale (const MatrixBase< OtherDerived > &other)
 
Transformprescale (const Scalar &s)
 
Transformpreshear (const Scalar &sx, const Scalar &sy)
 
template<typename OtherDerived >
Transform< Scalar, Dim, Mode, Options > & pretranslate (const MatrixBase< OtherDerived > &other)
 
template<typename RotationType >
Transform< Scalar, Dim, Mode, Options > & rotate (const RotationType &rotation)
 
RotationReturnType rotation () const
 
template<typename OtherDerived >
Transform< Scalar, Dim, Mode, Options > & scale (const MatrixBase< OtherDerived > &other)
 
Transformscale (const Scalar &s)
 
void setIdentity ()
 
Transformshear (const Scalar &sx, const Scalar &sy)
 
QMatrix toQMatrix (void) const
 
QTransform toQTransform (void) const
 
 Transform ()
 
template<typename OtherDerived >
 Transform (const EigenBase< OtherDerived > &other)
 
 Transform (const QMatrix &other)
 
 Transform (const QTransform &other)
 
template<typename OtherScalarType >
 Transform (const Transform< OtherScalarType, Dim, Mode, Options > &other)
 
template<typename OtherDerived >
Transform< Scalar, Dim, Mode, Options > & translate (const MatrixBase< OtherDerived > &other)
 
TranslationPart translation ()
 
ConstTranslationPart translation () const
 

Static Public Member Functions

static const Transform Identity ()
 Returns an identity transformation. More...
 

Public Attributes

const typedef Block< ConstMatrixType, Dim, Dim, int(Mode)==(AffineCompact) &&(Options &RowMajor)==0 > ConstLinearPart
 
const typedef MatrixType ConstMatrixType
 
const typedef Block< ConstMatrixType, Dim, 1,!(internal::traits< MatrixType >::Flags &RowMajorBit)> ConstTranslationPart
 

Member Typedef Documentation

◆ AffinePart

typedef internal::conditional<int(Mode)==int(AffineCompact), MatrixType&, Block<MatrixType,Dim,HDim> >::type Eigen::Transform::AffinePart

type of read/write reference to the affine part of the transformation

◆ ConstAffinePart

typedef internal::conditional<int(Mode)==int(AffineCompact), const MatrixType&, const Block<const MatrixType,Dim,HDim> >::type Eigen::Transform::ConstAffinePart

type of read reference to the affine part of the transformation

◆ Index

◆ LinearMatrixType

type of the matrix used to represent the linear part of the transformation

◆ LinearPart

typedef Block<MatrixType,Dim,Dim,int(Mode)==(AffineCompact) && (Options&RowMajor)==0> Eigen::Transform::LinearPart

type of read/write reference to the linear part of the transformation

◆ MatrixType

typedef internal::make_proper_matrix_type<Scalar,Rows,HDim,Options>::type Eigen::Transform::MatrixType

type of the matrix used to represent the transformation

◆ Scalar

typedef _Scalar Eigen::Transform::Scalar

the scalar type of the coefficients

◆ TransformTimeDiagonalReturnType

typedef Transform<Scalar,Dim,TransformTimeDiagonalMode> Eigen::Transform::TransformTimeDiagonalReturnType

The return type of the product between a diagonal matrix and a transform

◆ TranslationPart

typedef Block<MatrixType,Dim,1,!(internal::traits<MatrixType>::Flags & RowMajorBit)> Eigen::Transform::TranslationPart

type of a read/write reference to the translation part of the rotation

◆ TranslationType

corresponding translation type

◆ VectorType

type of a vector

Constructor & Destructor Documentation

◆ Transform() [1/5]

Eigen::Transform::Transform ( )
inline

Default constructor without initialization of the meaningful coefficients. If Mode==Affine or Mode==Isometry, then the last row is set to [0 ... 0 1]

◆ Transform() [2/5]

template<typename OtherDerived >
Eigen::Transform::Transform ( const EigenBase< OtherDerived > &  other)
inlineexplicit

Constructs and initializes a transformation from a Dim^2 or a (Dim+1)^2 matrix.

◆ Transform() [3/5]

Eigen::Transform::Transform ( const QMatrix &  other)
inline

Initializes *this from a QMatrix assuming the dimension is 2.

This function is available only if the token EIGEN_QT_SUPPORT is defined.

◆ Transform() [4/5]

Eigen::Transform::Transform ( const QTransform &  other)
inline

Initializes *this from a QTransform assuming the dimension is 2.

This function is available only if the token EIGEN_QT_SUPPORT is defined.

◆ Transform() [5/5]

template<typename OtherScalarType >
Eigen::Transform::Transform ( const Transform< OtherScalarType, Dim, Mode, Options > &  other)
inlineexplicit

Copy constructor with scalar type conversion

Member Function Documentation

◆ affine() [1/2]

AffinePart Eigen::Transform::affine ( )
inline
Returns
a writable expression of the Dim x HDim affine part of the transformation

◆ affine() [2/2]

ConstAffinePart Eigen::Transform::affine ( ) const
inline
Returns
a read-only expression of the Dim x HDim affine part of the transformation

◆ cast()

template<typename NewScalarType >
internal::cast_return_type<Transform,Transform<NewScalarType,Dim,Mode,Options> >::type Eigen::Transform::cast ( ) const
inline
Returns
*this with scalar type casted to NewScalarType

Note that if NewScalarType is equal to the current scalar type of *this then this function smartly returns a const reference to *this.

◆ computeRotationScaling()

template<typename RotationMatrixType , typename ScalingMatrixType >
void Eigen::Transform::computeRotationScaling ( RotationMatrixType *  rotation,
ScalingMatrixType *  scaling 
) const

decomposes the linear part of the transformation as a product rotation x scaling, the scaling being not necessarily positive.

If either pointer is zero, the corresponding computation is skipped.

This is defined in the SVD module.

#include <Eigen/SVD>
See also
computeScalingRotation(), rotation(), class SVD

◆ computeScalingRotation()

template<typename ScalingMatrixType , typename RotationMatrixType >
void Eigen::Transform::computeScalingRotation ( ScalingMatrixType *  scaling,
RotationMatrixType *  rotation 
) const

decomposes the linear part of the transformation as a product scaling x rotation, the scaling being not necessarily positive.

If either pointer is zero, the corresponding computation is skipped.

This is defined in the SVD module.

#include <Eigen/SVD>
See also
computeRotationScaling(), rotation(), class SVD

◆ data() [1/2]

Scalar* Eigen::Transform::data ( )
inline
Returns
a non-const pointer to the column major internal matrix

◆ data() [2/2]

const Scalar* Eigen::Transform::data ( ) const
inline
Returns
a const pointer to the column major internal matrix

◆ EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE()

Eigen::Transform::EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE ( _Scalar  ,
_Dim  = =Dynamic ? Dynamic : (_Dim+1)*(_Dim+1) 
)
inline

< space dimension in which the transformation holds

< size of a respective homogeneous vector

◆ fromPositionOrientationScale()

template<typename PositionDerived , typename OrientationType , typename ScaleDerived >
Transform<Scalar,Dim,Mode,Options>& Eigen::Transform::fromPositionOrientationScale ( const MatrixBase< PositionDerived > &  position,
const OrientationType &  orientation,
const MatrixBase< ScaleDerived > &  scale 
)

Convenient method to set *this from a position, orientation and scale of a 3D object.

◆ Identity()

static const Transform Eigen::Transform::Identity ( )
inlinestatic

Returns an identity transformation.

◆ inverse()

Transform< Scalar, Dim, Mode, Options > Eigen::Transform::inverse ( TransformTraits  hint = (TransformTraits)Mode) const
inline
Returns
the inverse transformation according to some given knowledge on *this.
Parameters
hintallows to optimize the inversion process when the transformation is known to be not a general transformation (optional). The possible values are:
  • Projective if the transformation is not necessarily affine, i.e., if the last row is not guaranteed to be [0 ... 0 1]
  • Affine if the last row can be assumed to be [0 ... 0 1]
  • Isometry if the transformation is only a concatenations of translations and rotations. The default is the template class parameter Mode.
Warning
unless traits is always set to NoShear or NoScaling, this function requires the generic inverse method of MatrixBase defined in the LU module. If you forget to include this module, then you will get hard to debug linking errors.
See also
MatrixBase::inverse()

◆ isApprox()

bool Eigen::Transform::isApprox ( const Transform other,
const typename NumTraits< Scalar >::Real &  prec = NumTraits<Scalar>::dummy_precision() 
) const
inline
Returns
true if *this is approximately equal to other, within the precision determined by prec.
See also
MatrixBase::isApprox()

◆ linear() [1/2]

LinearPart Eigen::Transform::linear ( )
inline
Returns
a writable expression of the linear part of the transformation

◆ linear() [2/2]

ConstLinearPart Eigen::Transform::linear ( ) const
inline
Returns
a read-only expression of the linear part of the transformation

◆ makeAffine()

void Eigen::Transform::makeAffine ( )
inline

Sets the last row to [0 ... 0 1]

◆ matrix() [1/2]

MatrixType& Eigen::Transform::matrix ( )
inline
Returns
a writable expression of the transformation matrix

◆ matrix() [2/2]

const MatrixType& Eigen::Transform::matrix ( ) const
inline
Returns
a read-only expression of the transformation matrix

◆ operator()() [1/2]

Scalar& Eigen::Transform::operator() ( Index  row,
Index  col 
)
inline

shortcut for m_matrix(row,col);

See also
MatrixBase::operator(Index,Index)

◆ operator()() [2/2]

Scalar Eigen::Transform::operator() ( Index  row,
Index  col 
) const
inline

shortcut for m_matrix(row,col);

See also
MatrixBase::operator(Index,Index) const

◆ operator*() [1/4]

template<typename DiagonalDerived >
const TransformTimeDiagonalReturnType Eigen::Transform::operator* ( const DiagonalBase< DiagonalDerived > &  b) const
inline
Returns
The product expression of a transform a times a diagonal matrix b

The rhs diagonal matrix is interpreted as an affine scaling transformation. The product results in a Transform of the same type (mode) as the lhs only if the lhs mode is no isometry. In that case, the returned transform is an affinity.

◆ operator*() [2/4]

template<typename OtherDerived >
const internal::transform_right_product_impl<Transform, OtherDerived>::ResultType Eigen::Transform::operator* ( const EigenBase< OtherDerived > &  other) const
inline
Returns
an expression of the product between the transform *this and a matrix expression other.

The right-hand-side other can be either:

  • an homogeneous vector of size Dim+1,
  • a set of homogeneous vectors of size Dim+1 x N,
  • a transformation matrix of size Dim+1 x Dim+1.

Moreover, if *this represents an affine transformation (i.e., Mode!=Projective), then other can also be:

  • a point of size Dim (computes:
    this->linear() * other + this->translation()
    ),
  • a set of N points as a Dim x N matrix (computes:
    (this->linear() * other).colwise() + this->translation()
    ),

In all cases, the return type is a matrix or vector of same sizes as the right-hand-side other.

If you want to interpret other as a linear or affine transformation, then first convert it to a Transform<> type, or do your own cooking.

Finally, if you want to apply Affine transformations to vectors, then explicitly apply the linear part only:

Affine3f A;
Vector3f v1, v2;
v2 = A.linear() * v1;

◆ operator*() [3/4]

const Transform Eigen::Transform::operator* ( const Transform other) const
inline

Concatenates two transformations

◆ operator*() [4/4]

template<int OtherMode, int OtherOptions>
internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> >::ResultType Eigen::Transform::operator* ( const Transform< Scalar, Dim, OtherMode, OtherOptions > &  other) const
inline

Concatenates two different transformations

◆ operator=() [1/3]

template<typename OtherDerived >
Transform& Eigen::Transform::operator= ( const EigenBase< OtherDerived > &  other)
inline

Set *this from a Dim^2 or (Dim+1)^2 matrix.

◆ operator=() [2/3]

Transform< Scalar, Dim, Mode, Options > & Eigen::Transform::operator= ( const QMatrix &  other)
inline

Set *this from a QMatrix assuming the dimension is 2.

This function is available only if the token EIGEN_QT_SUPPORT is defined.

◆ operator=() [3/3]

Transform< Scalar, Dim, Mode, Options > & Eigen::Transform::operator= ( const QTransform &  other)
inline

Set *this from a QTransform assuming the dimension is 2.

This function is available only if the token EIGEN_QT_SUPPORT is defined.

◆ prerotate()

template<typename RotationType >
Transform<Scalar,Dim,Mode,Options>& Eigen::Transform::prerotate ( const RotationType &  rotation)

Applies on the left the rotation represented by the rotation rotation to *this and returns a reference to *this.

See rotate() for further details.

See also
rotate()

◆ prescale() [1/2]

template<typename OtherDerived >
Transform<Scalar,Dim,Mode,Options>& Eigen::Transform::prescale ( const MatrixBase< OtherDerived > &  other)

Applies on the left the non uniform scale transformation represented by the vector other to *this and returns a reference to *this.

See also
scale()

◆ prescale() [2/2]

Transform< Scalar, Dim, Mode, Options > & Eigen::Transform::prescale ( const Scalar s)
inline

Applies on the left a uniform scale of a factor c to *this and returns a reference to *this.

See also
scale(Scalar)

◆ preshear()

Transform< Scalar, Dim, Mode, Options > & Eigen::Transform::preshear ( const Scalar sx,
const Scalar sy 
)

Applies on the left the shear transformation represented by the vector other to *this and returns a reference to *this.

Warning
2D only.
See also
shear()

◆ pretranslate()

template<typename OtherDerived >
Transform<Scalar,Dim,Mode,Options>& Eigen::Transform::pretranslate ( const MatrixBase< OtherDerived > &  other)

Applies on the left the translation matrix represented by the vector other to *this and returns a reference to *this.

See also
translate()

◆ rotate()

template<typename RotationType >
Transform<Scalar,Dim,Mode,Options>& Eigen::Transform::rotate ( const RotationType &  rotation)

Applies on the right the rotation represented by the rotation rotation to *this and returns a reference to *this.

The template parameter RotationType is the type of the rotation which must be known by internal::toRotationMatrix<>.

Natively supported types includes:

This mechanism is easily extendable to support user types such as Euler angles, or a pair of Quaternion for 4D rotations.

See also
rotate(Scalar), class Quaternion, class AngleAxis, prerotate(RotationType)

◆ rotation()

Transform< Scalar, Dim, Mode, Options >::RotationReturnType Eigen::Transform::rotation ( ) const
Returns
the rotation part of the transformation

If Mode==Isometry, then this method is an alias for linear(), otherwise it calls computeRotationScaling() to extract the rotation through a SVD decomposition.

This is defined in the SVD module.

#include <Eigen/SVD>
See also
computeRotationScaling(), computeScalingRotation(), class SVD

◆ scale() [1/2]

template<typename OtherDerived >
Transform<Scalar,Dim,Mode,Options>& Eigen::Transform::scale ( const MatrixBase< OtherDerived > &  other)

Applies on the right the non uniform scale transformation represented by the vector other to *this and returns a reference to *this.

See also
prescale()

◆ scale() [2/2]

Transform< Scalar, Dim, Mode, Options > & Eigen::Transform::scale ( const Scalar s)
inline

Applies on the right a uniform scale of a factor c to *this and returns a reference to *this.

See also
prescale(Scalar)

◆ setIdentity()

void Eigen::Transform::setIdentity ( )
inline

◆ shear()

Transform< Scalar, Dim, Mode, Options > & Eigen::Transform::shear ( const Scalar sx,
const Scalar sy 
)

Applies on the right the shear transformation represented by the vector other to *this and returns a reference to *this.

Warning
2D only.
See also
preshear()

◆ toQMatrix()

QMatrix Eigen::Transform::toQMatrix ( void  ) const
inline
Returns
a QMatrix from *this assuming the dimension is 2.
Warning
this conversion might loss data if *this is not affine

This function is available only if the token EIGEN_QT_SUPPORT is defined.

◆ toQTransform()

QTransform Eigen::Transform::toQTransform ( void  ) const
inline
Returns
a QTransform from *this assuming the dimension is 2.

This function is available only if the token EIGEN_QT_SUPPORT is defined.

◆ translate()

template<typename OtherDerived >
Transform<Scalar,Dim,Mode,Options>& Eigen::Transform::translate ( const MatrixBase< OtherDerived > &  other)

Applies on the right the translation matrix represented by the vector other to *this and returns a reference to *this.

See also
pretranslate()

◆ translation() [1/2]

TranslationPart Eigen::Transform::translation ( )
inline
Returns
a writable expression of the translation vector of the transformation

◆ translation() [2/2]

ConstTranslationPart Eigen::Transform::translation ( ) const
inline
Returns
a read-only expression of the translation vector of the transformation

Member Data Documentation

◆ ConstLinearPart

const typedef Block<ConstMatrixType,Dim,Dim,int(Mode)==(AffineCompact) && (Options&RowMajor)==0> Eigen::Transform::ConstLinearPart

type of read reference to the linear part of the transformation

◆ ConstMatrixType

const typedef MatrixType Eigen::Transform::ConstMatrixType

constified MatrixType

◆ ConstTranslationPart

const typedef Block<ConstMatrixType,Dim,1,!(internal::traits<MatrixType>::Flags & RowMajorBit)> Eigen::Transform::ConstTranslationPart

type of a read reference to the translation part of the rotation


The documentation for this class was generated from the following files:
Eigen::Transform::linear
ConstLinearPart linear() const
Definition: Transform.h:394
Eigen::Transform::translation
ConstTranslationPart translation() const
Definition: Transform.h:404