template<typename _Scalar, int _Dim, int _Mode, int _Options>
class Eigen::Transform< _Scalar, _Dim, _Mode, _Options >
Represents an homogeneous transformation in a N dimensional space.
This is defined in the Geometry module.
#include <Eigen/Geometry>
- 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. |
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
|
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 const Block< ConstMatrixType, Dim, Dim, int(Mode)==(AffineCompact) &&(int(Options)&RowMajor)==0 > |
ConstLinearPart |
|
typedef const MatrixType |
ConstMatrixType |
|
typedef const Block< ConstMatrixType, Dim, 1,!(internal::traits< MatrixType >::Flags &RowMajorBit)> |
ConstTranslationPart |
|
typedef Eigen::Index |
Index |
|
typedef Matrix< Scalar, Dim, Dim, Options > |
LinearMatrixType |
|
typedef Block< MatrixType, Dim, Dim, int(Mode)==(AffineCompact) &&(int(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 |
|
|
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 |
|
Scalar * |
data () |
|
const Scalar * |
data () 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 () |
|
MatrixType & |
matrix () |
|
const MatrixType & |
matrix () const |
|
Scalar & |
operator() (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 > |
Transform & |
operator= (const EigenBase< OtherDerived > &other) |
|
Transform & |
operator= (const QMatrix &other) |
|
Transform & |
operator= (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) |
|
Transform & |
prescale (const Scalar &s) |
|
Transform & |
preshear (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) |
|
Transform & |
scale (const Scalar &s) |
|
void |
setIdentity () |
|
Transform & |
shear (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 |
|
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 | | 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:
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;