template<typename _Scalar, int _Dim, int _Mode, int _Options>
Eigen::Transform class

Represents an homogeneous transformation in a N dimensional space.

Template parameters
_Scalar the scalar type, i.e., the type of the coefficients
_Dim the dimension of the space
_Mode

the 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().
_Options has 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.

This is defined in the Geometry module. #include <Eigen/Geometry>

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.

Public types

enum (anonymous) { TransformTimeDiagonalMode = ((Mode==int(Isometry))?Affine:int(Mode)) }
using AffinePart = internal::conditional<int(Mode)==int(AffineCompact), MatrixType&, Block<MatrixType, Dim, HDim>>::type
using ConstAffinePart = internal::conditional<int(Mode)==int(AffineCompact), const MatrixType&, const Block<const MatrixType, Dim, HDim>>::type
using ConstLinearPart = const Block<ConstMatrixType, Dim, Dim, int(Mode)==(AffineCompact) && (Options&RowMajor)==0>
using ConstMatrixType = const MatrixType
using ConstTranslationPart = const Block<ConstMatrixType, Dim, 1,!(internal::traits<MatrixType>::Flags&RowMajorBit)>
using Index = Eigen::Index deprecated
using LinearMatrixType = Matrix<Scalar, Dim, Dim, Options>
using LinearPart = Block<MatrixType, Dim, Dim, int(Mode)==(AffineCompact) && (Options&RowMajor)==0>
using MatrixType = internal::make_proper_matrix_type<Scalar, Rows, HDim, Options>::type
using RotationReturnType = internal::conditional<int(Mode)==Isometry, ConstLinearPart, const LinearMatrixType>::type
using Scalar = _Scalar
using StorageIndex = Eigen::Index
using take_affine_part = internal::transform_take_affine_part<Transform>
using TransformTimeDiagonalReturnType = Transform<Scalar, Dim, TransformTimeDiagonalMode>
using TranslationPart = Block<MatrixType, Dim, 1,!(internal::traits<MatrixType>::Flags&RowMajorBit)>
using TranslationType = Translation<Scalar, Dim>
using VectorType = Matrix<Scalar, Dim, 1>

Public static functions

static auto Identity() -> const Transform
Returns an identity transformation.

Constructors, destructors, conversion operators

EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar, _Dim = =Dynamic ? Dynamic :(_Dim+1)*(_Dim+1)) enum
Transform()
Transform(const Transform& other)
Transform(const TranslationType& t) explicit
Transform(const UniformScaling<Scalar>& s) explicit
template<typename Derived>
Transform(const RotationBase<Derived, Dim>& r) explicit
template<typename OtherDerived>
Transform(const EigenBase<OtherDerived>& other) explicit
template<int OtherOptions>
Transform(const Transform<Scalar, Dim, Mode, OtherOptions>& other)
template<int OtherMode, int OtherOptions>
Transform(const Transform<Scalar, Dim, OtherMode, OtherOptions>& other)
template<typename OtherDerived>
Transform(const ReturnByValue<OtherDerived>& other)
Transform(const QMatrix& other)
Transform(const QTransform& other)
template<typename OtherScalarType>
Transform(const Transform<OtherScalarType, Dim, Mode, Options>& other) explicit

Public functions

auto affine() const -> ConstAffinePart
auto affine() -> AffinePart
template<typename NewScalarType>
auto cast() const -> internal::cast_return_type<Transform, Transform<NewScalarType, Dim, Mode, Options>>::type
auto cols() const -> Index
template<typename RotationMatrixType, typename ScalingMatrixType>
void computeRotationScaling(RotationMatrixType* rotation, ScalingMatrixType* scaling) const
template<typename ScalingMatrixType, typename RotationMatrixType>
void computeScalingRotation(ScalingMatrixType* scaling, RotationMatrixType* rotation) const
auto data() const -> const Scalar*
auto data() -> Scalar*
template<typename PositionDerived, typename OrientationType, typename ScaleDerived>
auto fromPositionOrientationScale(const MatrixBase<PositionDerived>& position, const OrientationType& orientation, const MatrixBase<ScaleDerived>& scale) -> Transform&
template<typename PositionDerived, typename OrientationType, typename ScaleDerived>
auto fromPositionOrientationScale(const MatrixBase<PositionDerived>& position, const OrientationType& orientation, const MatrixBase<ScaleDerived>& scale) -> Transform<Scalar, Dim, Mode, Options>&
auto inverse(TransformTraits traits = (TransformTraits) Mode) const -> Transform
auto isApprox(const Transform& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
auto linear() const -> ConstLinearPart
auto linear() -> LinearPart
auto linearExt() -> Block<MatrixType, int(Mode)==int(Projective)?HDim:Dim, Dim>
auto linearExt() const -> const Block<MatrixType, int(Mode)==int(Projective)?HDim:Dim, Dim>
void makeAffine()
auto matrix() const -> const MatrixType&
auto matrix() -> MatrixType&
template<typename OtherDerived>
auto operator*(const EigenBase<OtherDerived>& other) const -> const internal::transform_right_product_impl<Transform, OtherDerived>::ResultType
template<typename DiagonalDerived>
auto operator*(const DiagonalBase<DiagonalDerived>& b) const -> const TransformTimeDiagonalReturnType
auto operator*(const Transform& other) const -> const Transform
template<int OtherMode, int OtherOptions>
auto operator*(const Transform<Scalar, Dim, OtherMode, OtherOptions>& other) const -> internal::transform_transform_product_impl<Transform, Transform<Scalar, Dim, OtherMode, OtherOptions>>::ResultType
auto operator*(const TranslationType& t) const -> Transform
auto operator*(const UniformScaling<Scalar>& s) const -> TransformTimeDiagonalReturnType
template<typename Derived>
auto operator*(const RotationBase<Derived, Dim>& r) const -> Transform
template<typename Derived>
auto operator*(const RotationBase<Derived, Dim>& r) const -> Transform<Scalar, Dim, Mode, Options>
template<typename OtherDerived>
auto operator*=(const EigenBase<OtherDerived>& other) -> Transform&
auto operator*=(const TranslationType& t) -> Transform&
auto operator*=(const UniformScaling<Scalar>& s) -> Transform&
auto operator*=(const DiagonalMatrix<Scalar, Dim>& s) -> Transform&
template<typename Derived>
auto operator*=(const RotationBase<Derived, Dim>& r) -> Transform&
auto operator()(Index row, Index col) const -> Scalar
auto operator()(Index row, Index col) -> Scalar&
auto operator=(const Transform& other) -> Transform&
template<typename OtherDerived>
auto operator=(const EigenBase<OtherDerived>& other) -> Transform&
template<typename OtherDerived>
auto operator=(const ReturnByValue<OtherDerived>& other) -> Transform&
auto operator=(const QMatrix& other) -> Transform&
auto operator=(const QTransform& other) -> Transform&
auto operator=(const TranslationType& t) -> Transform&
auto operator=(const UniformScaling<Scalar>& t) -> Transform&
template<typename Derived>
auto operator=(const RotationBase<Derived, Dim>& r) -> Transform&
template<typename Derived>
auto operator=(const RotationBase<Derived, Dim>& r) -> Transform<Scalar, Dim, Mode, Options>&
template<typename RotationType>
auto prerotate(const RotationType& rotation) -> Transform&
template<typename RotationType>
auto prerotate(const RotationType& rotation) -> Transform<Scalar, Dim, Mode, Options>&
template<typename OtherDerived>
auto prescale(const MatrixBase<OtherDerived>& other) -> Transform&
auto prescale(const Scalar& s) -> Transform&
template<typename OtherDerived>
auto prescale(const MatrixBase<OtherDerived>& other) -> Transform<Scalar, Dim, Mode, Options>&
auto preshear(const Scalar& sx, const Scalar& sy) -> Transform&
template<typename OtherDerived>
auto pretranslate(const MatrixBase<OtherDerived>& other) -> Transform&
template<typename OtherDerived>
auto pretranslate(const MatrixBase<OtherDerived>& other) -> Transform<Scalar, Dim, Mode, Options>&
template<typename RotationType>
auto rotate(const RotationType& rotation) -> Transform&
template<typename RotationType>
auto rotate(const RotationType& rotation) -> Transform<Scalar, Dim, Mode, Options>&
auto rotation() const -> RotationReturnType
auto rows() const -> Index
template<typename OtherDerived>
auto scale(const MatrixBase<OtherDerived>& other) -> Transform&
auto scale(const Scalar& s) -> Transform&
template<typename OtherDerived>
auto scale(const MatrixBase<OtherDerived>& other) -> Transform<Scalar, Dim, Mode, Options>&
void setIdentity()
auto shear(const Scalar& sx, const Scalar& sy) -> Transform&
auto toQMatrix(void) const -> QMatrix
auto toQTransform(void) const -> QTransform
template<typename OtherDerived>
auto translate(const MatrixBase<OtherDerived>& other) -> Transform&
template<typename OtherDerived>
auto translate(const MatrixBase<OtherDerived>& other) -> Transform<Scalar, Dim, Mode, Options>&
auto translation() const -> ConstTranslationPart
auto translation() -> TranslationPart
auto translationExt() -> Block<MatrixType, int(Mode)==int(Projective)?HDim:Dim, 1>
auto translationExt() const -> const Block<MatrixType, int(Mode)==int(Projective)?HDim:Dim, 1>

Protected variables

MatrixType m_matrix

Friends

template<typename OtherDerived>
auto operator*(const EigenBase<OtherDerived>& a, const Transform& b) -> const internal::transform_left_product_impl<OtherDerived, Mode, Options, _Dim, _Dim+1>::ResultType
template<typename DiagonalDerived>
auto operator*(const DiagonalBase<DiagonalDerived>& a, const Transform& b) -> TransformTimeDiagonalReturnType

Enum documentation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
enum Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::(anonymous)

Enumerators
TransformTimeDiagonalMode

Typedef documentation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef internal::conditional<int(Mode)==int(AffineCompact), MatrixType&, Block<MatrixType, Dim, HDim>>::type Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::AffinePart

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef internal::conditional<int(Mode)==int(AffineCompact), const MatrixType&, const Block<const MatrixType, Dim, HDim>>::type Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::ConstAffinePart

type of read reference to the affine part of the transformation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef const Block<ConstMatrixType, Dim, Dim, int(Mode)==(AffineCompact) && (Options&RowMajor)==0> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::ConstLinearPart

type of read reference to the linear part of the transformation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef const MatrixType Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::ConstMatrixType

constified MatrixType

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef const Block<ConstMatrixType, Dim, 1,!(internal::traits<MatrixType>::Flags&RowMajorBit)> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::ConstTranslationPart

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef Eigen::Index Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Index

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef Matrix<Scalar, Dim, Dim, Options> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::LinearMatrixType

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef Block<MatrixType, Dim, Dim, int(Mode)==(AffineCompact) && (Options&RowMajor)==0> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::LinearPart

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef internal::make_proper_matrix_type<Scalar, Rows, HDim, Options>::type Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::MatrixType

type of the matrix used to represent the transformation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef internal::conditional<int(Mode)==Isometry, ConstLinearPart, const LinearMatrixType>::type Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::RotationReturnType

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef _Scalar Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Scalar

the scalar type of the coefficients

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef Eigen::Index Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::StorageIndex

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef internal::transform_take_affine_part<Transform> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::take_affine_part

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef Transform<Scalar, Dim, TransformTimeDiagonalMode> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::TransformTimeDiagonalReturnType

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef Block<MatrixType, Dim, 1,!(internal::traits<MatrixType>::Flags&RowMajorBit)> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::TranslationPart

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef Translation<Scalar, Dim> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::TranslationType

corresponding translation type

template<typename _Scalar, int _Dim, int _Mode, int _Options>
typedef Matrix<Scalar, Dim, 1> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::VectorType

type of a vector

Function documentation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
static const Transform Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Identity()

Returns an identity transformation.

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar, _Dim = =Dynamic ? Dynamic :(_Dim+1)*(_Dim+1)) enum

< space dimension in which the transformation holds

< size of a respective homogeneous vector

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform()

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform(const Transform& other)

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform(const TranslationType& t) explicit

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform(const UniformScaling<Scalar>& s) explicit

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename Derived>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform(const RotationBase<Derived, Dim>& r) explicit

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform(const EigenBase<OtherDerived>& other) explicit

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

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<int OtherOptions>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform(const Transform<Scalar, Dim, Mode, OtherOptions>& other)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<int OtherMode, int OtherOptions>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform(const Transform<Scalar, Dim, OtherMode, OtherOptions>& other)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform(const ReturnByValue<OtherDerived>& other)

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform(const QMatrix& other)

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

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform(const QTransform& other)

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

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

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherScalarType>
Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::Transform(const Transform<OtherScalarType, Dim, Mode, Options>& other) explicit

Copy constructor with scalar type conversion

template<typename _Scalar, int _Dim, int _Mode, int _Options>
ConstAffinePart Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::affine() const

Returns a read-only expression of the Dim x HDim affine part of the transformation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
AffinePart Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::affine()

Returns a writable expression of the Dim x HDim affine part of the transformation

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename NewScalarType>
internal::cast_return_type<Transform, Transform<NewScalarType, Dim, Mode, Options>>::type Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::cast() const

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.

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Index Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::cols() const

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename RotationMatrixType, typename ScalingMatrixType>
void Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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>

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename ScalingMatrixType, typename RotationMatrixType>
void Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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>

template<typename _Scalar, int _Dim, int _Mode, int _Options>
const Scalar* Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::data() const

Returns a const pointer to the column major internal matrix

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Scalar* Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::data()

Returns a non-const pointer to the column major internal matrix

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename PositionDerived, typename OrientationType, typename ScaleDerived>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::fromPositionOrientationScale(const MatrixBase<PositionDerived>& position, const OrientationType& orientation, const MatrixBase<ScaleDerived>& scale)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename PositionDerived, typename OrientationType, typename ScaleDerived>
Transform<Scalar, Dim, Mode, Options>& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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.

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::inverse(TransformTraits traits = (TransformTraits) Mode) const

Returns the inverse transformation according to some given knowledge on *this.

template<typename _Scalar, int _Dim, int _Mode, int _Options>
bool Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::isApprox(const Transform& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const

Returns true if *this is approximately equal to other, within the precision determined by prec.

template<typename _Scalar, int _Dim, int _Mode, int _Options>
ConstLinearPart Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::linear() const

Returns a read-only expression of the linear part of the transformation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
LinearPart Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::linear()

Returns a writable expression of the linear part of the transformation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Block<MatrixType, int(Mode)==int(Projective)?HDim:Dim, Dim> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::linearExt()

template<typename _Scalar, int _Dim, int _Mode, int _Options>
const Block<MatrixType, int(Mode)==int(Projective)?HDim:Dim, Dim> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::linearExt() const

template<typename _Scalar, int _Dim, int _Mode, int _Options>
void Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::makeAffine()

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
const MatrixType& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::matrix() const

Returns a read-only expression of the transformation matrix

template<typename _Scalar, int _Dim, int _Mode, int _Options>
MatrixType& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::matrix()

Returns a writable expression of the transformation matrix

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
const internal::transform_right_product_impl<Transform, OtherDerived>::ResultType Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*(const EigenBase<OtherDerived>& other) const

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;

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename DiagonalDerived>
const TransformTimeDiagonalReturnType Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*(const DiagonalBase<DiagonalDerived>& b) const

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.

template<typename _Scalar, int _Dim, int _Mode, int _Options>
const Transform Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*(const Transform& other) const

Concatenates two transformations

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<int OtherMode, int OtherOptions>
internal::transform_transform_product_impl<Transform, Transform<Scalar, Dim, OtherMode, OtherOptions>>::ResultType Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*(const Transform<Scalar, Dim, OtherMode, OtherOptions>& other) const

Concatenates two different transformations

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*(const TranslationType& t) const

template<typename _Scalar, int _Dim, int _Mode, int _Options>
TransformTimeDiagonalReturnType Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*(const UniformScaling<Scalar>& s) const

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename Derived>
Transform Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*(const RotationBase<Derived, Dim>& r) const

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename Derived>
Transform<Scalar, Dim, Mode, Options> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*(const RotationBase<Derived, Dim>& r) const

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*=(const EigenBase<OtherDerived>& other)

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*=(const TranslationType& t)

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*=(const UniformScaling<Scalar>& s)

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*=(const DiagonalMatrix<Scalar, Dim>& s)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename Derived>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator*=(const RotationBase<Derived, Dim>& r)

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Scalar Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator()(Index row, Index col) const

shortcut for m_matrix(row,col);

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Scalar& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator()(Index row, Index col)

shortcut for m_matrix(row,col);

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator=(const Transform& other)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator=(const EigenBase<OtherDerived>& other)

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

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator=(const ReturnByValue<OtherDerived>& other)

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator=(const QMatrix& other)

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

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator=(const QTransform& other)

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

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator=(const TranslationType& t)

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator=(const UniformScaling<Scalar>& t)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename Derived>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator=(const RotationBase<Derived, Dim>& r)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename Derived>
Transform<Scalar, Dim, Mode, Options>& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::operator=(const RotationBase<Derived, Dim>& r)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename RotationType>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::prerotate(const RotationType& rotation)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename RotationType>
Transform<Scalar, Dim, Mode, Options>& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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.

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::prescale(const MatrixBase<OtherDerived>& other)

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::prescale(const Scalar& s)

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

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Transform<Scalar, Dim, Mode, Options>& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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.

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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.

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::pretranslate(const MatrixBase<OtherDerived>& other)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Transform<Scalar, Dim, Mode, Options>& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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.

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename RotationType>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::rotate(const RotationType& rotation)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename RotationType>
Transform<Scalar, Dim, Mode, Options>& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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.

template<typename _Scalar, int _Dim, int _Mode, int _Options>
RotationReturnType Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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>

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Index Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::rows() const

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::scale(const MatrixBase<OtherDerived>& other)

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::scale(const Scalar& s)

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

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Transform<Scalar, Dim, Mode, Options>& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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.

template<typename _Scalar, int _Dim, int _Mode, int _Options>
void Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::setIdentity()

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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.

template<typename _Scalar, int _Dim, int _Mode, int _Options>
QMatrix Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::toQMatrix(void) const

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

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

template<typename _Scalar, int _Dim, int _Mode, int _Options>
QTransform Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::toQTransform(void) const

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

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

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Transform& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::translate(const MatrixBase<OtherDerived>& other)

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
Transform<Scalar, Dim, Mode, Options>& Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::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.

template<typename _Scalar, int _Dim, int _Mode, int _Options>
ConstTranslationPart Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::translation() const

Returns a read-only expression of the translation vector of the transformation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
TranslationPart Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::translation()

Returns a writable expression of the translation vector of the transformation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
Block<MatrixType, int(Mode)==int(Projective)?HDim:Dim, 1> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::translationExt()

template<typename _Scalar, int _Dim, int _Mode, int _Options>
const Block<MatrixType, int(Mode)==int(Projective)?HDim:Dim, 1> Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::translationExt() const

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename OtherDerived>
const internal::transform_left_product_impl<OtherDerived, Mode, Options, _Dim, _Dim+1>::ResultType operator*(const EigenBase<OtherDerived>& a, const Transform& b)

Returns the product expression of a transformation matrix a times a transform b

The left hand side other can be either:

  • a linear transformation matrix of size Dim x Dim,
  • an affine transformation matrix of size Dim x Dim+1,
  • a general transformation matrix of size Dim+1 x Dim+1.

template<typename _Scalar, int _Dim, int _Mode, int _Options> template<typename DiagonalDerived>
TransformTimeDiagonalReturnType operator*(const DiagonalBase<DiagonalDerived>& a, const Transform& b)

Returns The product expression of a diagonal matrix a times a transform b

The lhs 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.

Variable documentation

template<typename _Scalar, int _Dim, int _Mode, int _Options>
MatrixType Eigen::Transform<_Scalar, _Dim, _Mode, _Options>::m_matrix protected