Represents an homogeneous transformation in a N dimensional space. More...
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 const Block < ConstMatrixType, Dim, Dim > | ConstLinearPart |
typedef const MatrixType | ConstMatrixType |
typedef const Block < ConstMatrixType, Dim, 1 > | ConstTranslationPart |
typedef Matrix< Scalar, Dim, Dim, Options > | LinearMatrixType |
typedef Block< MatrixType, Dim, Dim > | 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 > | 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 |
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 & | fromPositionOrientationScale (const MatrixBase< PositionDerived > &position, const OrientationType &orientation, const MatrixBase< ScaleDerived > &scale) |
Transform | inverse (TransformTraits traits=(TransformTraits) Mode) const |
bool | isApprox (const Transform &other, 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<int OtherMode, int OtherOptions> | |
const internal::transform_transform_product_impl < Transform, Transform< Scalar, Dim, OtherMode, OtherOptions > >::ResultType | operator* (const Transform< Scalar, Dim, OtherMode, OtherOptions > &other) const |
const Transform | operator* (const Transform &other) 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 |
Transform & | operator= (const QTransform &other) |
Transform & | operator= (const QMatrix &other) |
template<typename OtherDerived > | |
Transform & | operator= (const EigenBase< OtherDerived > &other) |
template<typename RotationType > | |
Transform & | prerotate (const RotationType &rotation) |
Transform & | prescale (Scalar s) |
template<typename OtherDerived > | |
Transform & | prescale (const MatrixBase< OtherDerived > &other) |
Transform & | preshear (Scalar sx, Scalar sy) |
template<typename OtherDerived > | |
Transform & | pretranslate (const MatrixBase< OtherDerived > &other) |
template<typename RotationType > | |
Transform & | rotate (const RotationType &rotation) |
const LinearMatrixType | rotation () const |
Transform & | scale (Scalar s) |
template<typename OtherDerived > | |
Transform & | scale (const MatrixBase< OtherDerived > &other) |
void | setIdentity () |
Transform & | shear (Scalar sx, Scalar sy) |
QMatrix | toQMatrix (void) const |
QTransform | toQTransform (void) const |
template<typename OtherScalarType > | |
Transform (const Transform< OtherScalarType, Dim, Mode, Options > &other) | |
Transform (const QTransform &other) | |
Transform (const QMatrix &other) | |
template<typename OtherDerived > | |
Transform (const EigenBase< OtherDerived > &other) | |
Transform () | |
template<typename OtherDerived > | |
Transform & | translate (const MatrixBase< OtherDerived > &other) |
TranslationPart | translation () |
ConstTranslationPart | translation () const |
Static Public Member Functions | |
static const Transform | Identity () |
Returns an identity transformation. | |
Friends | |
template<typename DiagonalDerived > | |
TransformTimeDiagonalReturnType | operator* (const DiagonalBase< DiagonalDerived > &a, const Transform &b) |
template<typename OtherDerived > | |
const internal::transform_left_product_impl < OtherDerived, Mode, Options, _Dim, _Dim+1 >::ResultType | operator* (const EigenBase< OtherDerived > &a, const Transform &b) |
Represents an homogeneous transformation in a N dimensional space.
This is defined in the Geometry module.
#include <Eigen/Geometry>
_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:
| |
_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:
Note that for a projective transformation the last row can be anything, and then the interpretation of different parts might be sightly 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,Matrix) 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 non homogeneous vectors, the latters are automatically promoted to homogeneous one before doing the matrix product. The convertions to homogeneous representations are performed as follow:
Translation t (Dim)x(1):
Rotation R (Dim)x(Dim):
Linear Matrix L (Dim)x(Dim):
Affine Matrix A (Dim)x(Dim+1):
Column vector v (Dim)x(1):
Set of column vectors V1...Vn (Dim)x(n):
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 Customizing/Extending Eigen by defining the preprocessor symbol EIGEN_TRANSFORM_PLUGIN
.
typedef internal::conditional<int(Mode)==int(AffineCompact), MatrixType&, Block<MatrixType,Dim,HDim> >::type AffinePart |
type of read/write reference to the affine part of the transformation
typedef internal::conditional<int(Mode)==int(AffineCompact), const MatrixType&, const Block<const MatrixType,Dim,HDim> >::type ConstAffinePart |
type of read reference to the affine part of the transformation
typedef const Block<ConstMatrixType,Dim,Dim> ConstLinearPart |
type of read reference to the linear part of the transformation
typedef const MatrixType ConstMatrixType |
constified MatrixType
typedef const Block<ConstMatrixType,Dim,1> ConstTranslationPart |
type of a read reference to the translation part of the rotation
typedef Matrix<Scalar,Dim,Dim,Options> LinearMatrixType |
type of the matrix used to represent the linear part of the transformation
typedef Block<MatrixType,Dim,Dim> LinearPart |
type of read/write reference to the linear part of the transformation
typedef internal::make_proper_matrix_type<Scalar,Rows,HDim,Options>::type MatrixType |
type of the matrix used to represent the transformation
typedef _Scalar Scalar |
the scalar type of the coefficients
typedef Transform<Scalar,Dim,TransformTimeDiagonalMode> TransformTimeDiagonalReturnType |
The return type of the product between a diagonal matrix and a transform
typedef Block<MatrixType,Dim,1> TranslationPart |
type of a read/write reference to the translation part of the rotation
typedef Translation<Scalar,Dim> TranslationType |
corresponding translation type
typedef Matrix<Scalar,Dim,1> VectorType |
type of a vector
Transform | ( | ) | [inline] |
Default constructor without initialization of the meaningful coefficients. If Mode==Affine, then the last row is set to [0 ... 0 1]
Constructs and initializes a transformation from a Dim^2 or a (Dim+1)^2 matrix.
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 | ( | const QTransform< _Scalar, _Dim, _Mode, _Options > & | 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.
Copy constructor with scalar type conversion
AffinePart affine | ( | ) | [inline] |
ConstAffinePart affine | ( | ) | const [inline] |
internal::cast_return_type<Transform,Transform<NewScalarType,Dim,Mode,Options> >::type cast | ( | ) | const [inline] |
*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
.
void computeRotationScaling | ( | RotationMatrixType * | rotation, | |
ScalingMatrixType * | scaling | |||
) | const [inline] |
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>
void computeScalingRotation | ( | ScalingMatrixType * | scaling, | |
RotationMatrixType * | rotation | |||
) | const [inline] |
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>
Scalar* data | ( | ) | [inline] |
const Scalar* data | ( | ) | const [inline] |
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
Transform< Scalar, Dim, Mode, Options > & fromPositionOrientationScale | ( | const MatrixBase< PositionDerived > & | position, | |
const OrientationType & | orientation, | |||
const MatrixBase< ScaleDerived > & | scale | |||
) | [inline] |
Convenient method to set *this
from a position, orientation and scale of a 3D object.
static const Transform Identity | ( | ) | [inline, static] |
Returns an identity transformation.
Transform< Scalar, Dim, Mode, Options > inverse | ( | TransformTraits | hint = (TransformTraits)Mode |
) | const [inline] |
*this
.hint | allows to optimize the inversion process when the transformation is known to be not a general transformation (optional). The possible values are:
|
bool isApprox | ( | const Transform< _Scalar, _Dim, _Mode, _Options > & | other, | |
typename NumTraits< Scalar >::Real | prec = NumTraits<Scalar>::dummy_precision() | |||
) | const [inline] |
true
if *this
is approximately equal to other, within the precision determined by prec.LinearPart linear | ( | ) | [inline] |
ConstLinearPart linear | ( | ) | const [inline] |
void makeAffine | ( | ) | [inline] |
Sets the last row to [0 ... 0 1]
MatrixType& matrix | ( | ) | [inline] |
const MatrixType& matrix | ( | ) | const [inline] |
Scalar& operator() | ( | Index | row, | |
Index | col | |||
) | [inline] |
shortcut for m_matrix(row,col);
Scalar operator() | ( | Index | row, | |
Index | col | |||
) | const [inline] |
shortcut for m_matrix(row,col);
const internal::transform_transform_product_impl< Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> >::ResultType operator* | ( | const Transform< Scalar, Dim, OtherMode, OtherOptions > & | other | ) | const [inline] |
Concatenates two different transformations
const Transform operator* | ( | const Transform< _Scalar, _Dim, _Mode, _Options > & | other | ) | const [inline] |
Concatenates two transformations
const TransformTimeDiagonalReturnType operator* | ( | const DiagonalBase< DiagonalDerived > & | b | ) | const [inline] |
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.
const internal::transform_right_product_impl<Transform, OtherDerived>::ResultType operator* | ( | const EigenBase< OtherDerived > & | other | ) | const [inline] |
*this
and a matrix expression other The right hand side other might be either:
Transform< Scalar, Dim, Mode, Options > & operator= | ( | const QTransform< _Scalar, _Dim, _Mode, _Options > & | 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.
Set *this
from a QMatrix assuming the dimension is 2.
This function is available only if the token EIGEN_QT_SUPPORT is defined.
Set *this
from a Dim^2 or (Dim+1)^2 matrix.
Applies on the left a uniform scale of a factor c to *this
and returns a reference to *this
.
Transform< Scalar, Dim, Mode, Options > & prescale | ( | const MatrixBase< OtherDerived > & | other | ) | [inline] |
Applies on the left the non uniform scale transformation represented by the vector other to *this
and returns a reference to *this
.
Applies on the left the shear transformation represented by the vector other to *this
and returns a reference to *this
.
Transform< Scalar, Dim, Mode, Options > & pretranslate | ( | const MatrixBase< OtherDerived > & | other | ) | [inline] |
Applies on the left the translation matrix represented by the vector other to *this
and returns a reference to *this
.
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.
const Transform< Scalar, Dim, Mode, Options >::LinearMatrixType rotation | ( | ) | const [inline] |
This is defined in the SVD module.
#include <Eigen/SVD>
Applies on the right a uniform scale of a factor c to *this
and returns a reference to *this
.
Transform< Scalar, Dim, Mode, Options > & scale | ( | const MatrixBase< OtherDerived > & | other | ) | [inline] |
Applies on the right the non uniform scale transformation represented by the vector other to *this
and returns a reference to *this
.
void setIdentity | ( | ) | [inline] |
Applies on the right the shear transformation represented by the vector other to *this
and returns a reference to *this
.
QMatrix toQMatrix | ( | void | ) | const [inline] |
*this
assuming the dimension is 2.*this
is not affineThis function is available only if the token EIGEN_QT_SUPPORT is defined.
QTransform toQTransform | ( | void | ) | const [inline] |
*this
assuming the dimension is 2.This function is available only if the token EIGEN_QT_SUPPORT is defined.
Transform< Scalar, Dim, Mode, Options > & translate | ( | const MatrixBase< OtherDerived > & | other | ) | [inline] |
Applies on the right the translation matrix represented by the vector other to *this
and returns a reference to *this
.
TranslationPart translation | ( | ) | [inline] |
ConstTranslationPart translation | ( | ) | const [inline] |
TransformTimeDiagonalReturnType operator* | ( | const DiagonalBase< DiagonalDerived > & | a, | |
const Transform< _Scalar, _Dim, _Mode, _Options > & | b | |||
) | [friend] |
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.
const internal::transform_left_product_impl<OtherDerived,Mode,Options,_Dim,_Dim+1>::ResultType operator* | ( | const EigenBase< OtherDerived > & | a, | |
const Transform< _Scalar, _Dim, _Mode, _Options > & | b | |||
) | [friend] |
The left hand side other might be either: