Eigen::Transform< Scalar_, Dim_, Mode_, Options_ > Class Template Reference

Represents an homogeneous transformation in a N dimensional space. More...

Public Types

enum  { TransformTimeDiagonalMode }
 
typedef std::conditional_t< int(Mode)==int(AffineCompact), MatrixType &, Block< MatrixType, Dim, HDim > > AffinePart
 
typedef std::conditional_t< int(Mode)==int(AffineCompact), const MatrixType &, const Block< const MatrixType, Dim, HDim > > 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 std::conditional_t< int(Mode)==Isometry, ConstLinearPart, const LinearMatrixTypeRotationReturnType
 
typedef Scalar_ Scalar
 
typedef Eigen::Index StorageIndex
 
typedef internal::transform_take_affine_part< Transformtake_affine_part
 
typedef Transform< Scalar, Dim, TransformTimeDiagonalModeTransformTimeDiagonalReturnType
 
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
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
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 >
TransformfromPositionOrientationScale (const MatrixBase< PositionDerived > &position, const OrientationType &orientation, const MatrixBase< ScaleDerived > &scale)
 
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
 
Block< MatrixType, int(Mode)==int(Projective)?HDim:Dim, Dim > linearExt ()
 
const Block< MatrixType, int(Mode)==int(Projective)?HDim:Dim, Dim > linearExt () 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
 
template<typename Derived >
Transform operator* (const RotationBase< Derived, Dim > &r) const
 
template<typename Derived >
Transform< Scalar, Dim, Mode, Options > operator* (const RotationBase< Derived, Dim > &r) 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
 
Transform operator* (const TranslationType &t) const
 
TransformTimeDiagonalReturnType operator* (const UniformScaling< Scalar > &s) const
 
Transformoperator*= (const DiagonalMatrix< Scalar, Dim > &s)
 
template<typename OtherDerived >
Transformoperator*= (const EigenBase< OtherDerived > &other)
 
template<typename Derived >
Transformoperator*= (const RotationBase< Derived, Dim > &r)
 
Transformoperator*= (const TranslationType &t)
 
Transformoperator*= (const UniformScaling< Scalar > &s)
 
template<typename OtherDerived >
Transformoperator= (const EigenBase< OtherDerived > &other)
 
Transformoperator= (const QTransform &other)
 
template<typename OtherDerived >
Transformoperator= (const ReturnByValue< OtherDerived > &other)
 
template<typename Derived >
Transformoperator= (const RotationBase< Derived, Dim > &r)
 
template<typename Derived >
Transform< Scalar, Dim, Mode, Options > & operator= (const RotationBase< Derived, Dim > &r)
 
Transformoperator= (const TranslationType &t)
 
Transformoperator= (const UniformScaling< Scalar > &t)
 
template<typename RotationType >
Transformprerotate (const RotationType &rotation)
 
template<typename RotationType >
Transform< Scalar, Dim, Mode, Options > & prerotate (const RotationType &rotation)
 
template<typename OtherDerived >
Transformprescale (const MatrixBase< OtherDerived > &other)
 
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 >
Transformpretranslate (const MatrixBase< OtherDerived > &other)
 
template<typename OtherDerived >
Transform< Scalar, Dim, Mode, Options > & pretranslate (const MatrixBase< OtherDerived > &other)
 
template<typename RotationType >
Transformrotate (const RotationType &rotation)
 
template<typename RotationType >
Transform< Scalar, Dim, Mode, Options > & rotate (const RotationType &rotation)
 
RotationReturnType rotation () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
template<typename OtherDerived >
Transformscale (const MatrixBase< OtherDerived > &other)
 
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)
 
QTransform toQTransform (void) const
 
 Transform ()
 
template<typename OtherDerived >
 Transform (const EigenBase< OtherDerived > &other)
 
 Transform (const QTransform &other)
 
template<typename OtherDerived >
 Transform (const ReturnByValue< OtherDerived > &other)
 
template<typename Derived >
 Transform (const RotationBase< Derived, Dim > &r)
 
template<typename OtherScalarType >
 Transform (const Transform< OtherScalarType, Dim, Mode, Options > &other)
 
template<int OtherOptions>
 Transform (const Transform< Scalar, Dim, Mode, OtherOptions > &other)
 
template<int OtherMode, int OtherOptions>
 Transform (const Transform< Scalar, Dim, OtherMode, OtherOptions > &other)
 
 Transform (const TranslationType &t)
 
 Transform (const UniformScaling< Scalar > &s)
 
template<typename OtherDerived >
Transformtranslate (const MatrixBase< OtherDerived > &other)
 
template<typename OtherDerived >
Transform< Scalar, Dim, Mode, Options > & translate (const MatrixBase< OtherDerived > &other)
 
TranslationPart translation ()
 
ConstTranslationPart translation () const
 
Block< MatrixType, int(Mode)==int(Projective)?HDim:Dim, 1 > translationExt ()
 
const Block< MatrixType, int(Mode)==int(Projective)?HDim:Dim, 1 > translationExt () const
 

Static Public Member Functions

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

Protected Attributes

MatrixType m_matrix
 

Detailed Description

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.

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
Array< int, Dynamic, 1 > 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();
Matrix3f m

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

Definition at line 206 of file Transform.h.

Member Typedef Documentation

◆ AffinePart

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
typedef std::conditional_t<int(Mode)==int(AffineCompact), MatrixType&, Block<MatrixType,Dim,HDim> > Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::AffinePart

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

Definition at line 234 of file Transform.h.

◆ ConstAffinePart

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
typedef std::conditional_t<int(Mode)==int(AffineCompact), const MatrixType&, const Block<const MatrixType,Dim,HDim> > Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::ConstAffinePart

type of read reference to the affine part of the transformation

Definition at line 238 of file Transform.h.

◆ ConstLinearPart

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
typedef const Block<ConstMatrixType,Dim,Dim,int(Mode)==(AffineCompact) && (int(Options)&RowMajor)==0> Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::ConstLinearPart

type of read reference to the linear part of the transformation

Definition at line 230 of file Transform.h.

◆ ConstMatrixType

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
typedef const MatrixType Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::ConstMatrixType

constified MatrixType

Definition at line 224 of file Transform.h.

◆ ConstTranslationPart

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

Definition at line 244 of file Transform.h.

◆ Index

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
typedef Eigen::Index Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::Index
Deprecated:
since Eigen 3.3

Definition at line 220 of file Transform.h.

◆ LinearMatrixType

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

Definition at line 226 of file Transform.h.

◆ LinearPart

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
typedef Block<MatrixType,Dim,Dim,int(Mode)==(AffineCompact) && (int(Options)&RowMajor)==0> Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::LinearPart

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

Definition at line 228 of file Transform.h.

◆ MatrixType

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

Definition at line 222 of file Transform.h.

◆ RotationReturnType

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
typedef std::conditional_t<int(Mode)==Isometry,ConstLinearPart,const LinearMatrixType> Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::RotationReturnType

Definition at line 603 of file Transform.h.

◆ Scalar

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
typedef Scalar_ Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::Scalar

the scalar type of the coefficients

Definition at line 218 of file Transform.h.

◆ StorageIndex

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
typedef Eigen::Index Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::StorageIndex

Definition at line 219 of file Transform.h.

◆ take_affine_part

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

Definition at line 284 of file Transform.h.

◆ TransformTimeDiagonalReturnType

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

Definition at line 251 of file Transform.h.

◆ TranslationPart

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

Definition at line 242 of file Transform.h.

◆ TranslationType

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
typedef Translation<Scalar,Dim> Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::TranslationType

corresponding translation type

Definition at line 246 of file Transform.h.

◆ VectorType

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

Definition at line 240 of file Transform.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
anonymous enum
Enumerator
TransformTimeDiagonalMode 

Definition at line 249 of file Transform.h.

249 { TransformTimeDiagonalMode = ((Mode==int(Isometry))?Affine:int(Mode)) };
@ TransformTimeDiagonalMode
Definition: Transform.h:249
@ Affine
Definition: Constants.h:464
@ Isometry
Definition: Constants.h:461

Constructor & Destructor Documentation

◆ Transform() [1/10]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::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]

Definition at line 261 of file Transform.h.

262  {
263  check_template_params();
264  internal::transform_make_affine<(int(Mode)==Affine || int(Mode)==Isometry) ? Affine : AffineCompact>::run(m_matrix);
265  }
MatrixType m_matrix
Definition: Transform.h:255
@ AffineCompact
Definition: Constants.h:466

◆ Transform() [2/10]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::Transform ( const TranslationType t)
inlineexplicit

Definition at line 267 of file Transform.h.

268  {
269  check_template_params();
270  *this = t;
271  }

◆ Transform() [3/10]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::Transform ( const UniformScaling< Scalar > &  s)
inlineexplicit

Definition at line 272 of file Transform.h.

273  {
274  check_template_params();
275  *this = s;
276  }

◆ Transform() [4/10]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename Derived >
Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::Transform ( const RotationBase< Derived, Dim > &  r)
inlineexplicit

Definition at line 278 of file Transform.h.

279  {
280  check_template_params();
281  *this = r;
282  }

◆ Transform() [5/10]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename OtherDerived >
Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::Transform ( const EigenBase< OtherDerived > &  other)
inlineexplicit

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

Definition at line 288 of file Transform.h.

289  {
290  EIGEN_STATIC_ASSERT((internal::is_same<Scalar,typename OtherDerived::Scalar>::value),
291  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY);
292 
293  check_template_params();
294  internal::transform_construct_from_matrix<OtherDerived,Mode,Options,Dim,HDim>::run(this, other.derived());
295  }
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26

◆ Transform() [6/10]

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)
inline

Definition at line 309 of file Transform.h.

310  {
311  check_template_params();
312  // only the options change, we can directly copy the matrices
313  m_matrix = other.matrix();
314  }
MatrixBase< Derived > & matrix()
Definition: MatrixBase.h:319

◆ Transform() [7/10]

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)
inline

Definition at line 317 of file Transform.h.

318  {
319  check_template_params();
320  // prevent conversions as:
321  // Affine | AffineCompact | Isometry = Projective
323  YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION)
324 
325  // prevent conversions as:
326  // Isometry = Affine | AffineCompact
327  EIGEN_STATIC_ASSERT(internal::check_implication(OtherMode==int(Affine)||OtherMode==int(AffineCompact), Mode!=int(Isometry)),
328  YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION)
329 
330  enum { ModeIsAffineCompact = Mode == int(AffineCompact),
331  OtherModeIsAffineCompact = OtherMode == int(AffineCompact)
332  };
333 
334  if(EIGEN_CONST_CONDITIONAL(ModeIsAffineCompact == OtherModeIsAffineCompact))
335  {
336  // We need the block expression because the code is compiled for all
337  // combinations of transformations and will trigger a compile time error
338  // if one tries to assign the matrices directly
339  m_matrix.template block<Dim,Dim+1>(0,0) = other.matrix().template block<Dim,Dim+1>(0,0);
340  makeAffine();
341  }
342  else if(EIGEN_CONST_CONDITIONAL(OtherModeIsAffineCompact))
343  {
344  typedef typename Transform<Scalar,Dim,OtherMode,OtherOptions>::MatrixType OtherMatrixType;
345  internal::transform_construct_from_matrix<OtherMatrixType,Mode,Options,Dim,HDim>::run(this, other.matrix());
346  }
347  else
348  {
349  // here we know that Mode == AffineCompact and OtherMode != AffineCompact.
350  // if OtherMode were Projective, the static assert above would already have caught it.
351  // So the only possibility is that OtherMode == Affine
352  linear() = other.linear();
353  translation() = other.translation();
354  }
355  }
#define EIGEN_CONST_CONDITIONAL(cond)
Definition: Macros.h:1048
ConstLinearPart linear() const
Definition: Transform.h:398
internal::make_proper_matrix_type< Scalar, Rows, HDim, Options >::type MatrixType
Definition: Transform.h:222
ConstTranslationPart translation() const
Definition: Transform.h:408
void makeAffine()
Definition: Transform.h:652
@ Projective
Definition: Constants.h:468
constexpr bool check_implication(bool a, bool b)
Definition: Meta.h:579

◆ Transform() [8/10]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename OtherDerived >
Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::Transform ( const ReturnByValue< OtherDerived > &  other)
inline

Definition at line 358 of file Transform.h.

359  {
360  check_template_params();
361  other.evalTo(*this);
362  }

◆ Transform() [9/10]

template<typename Scalar , int Dim, int Mode, int Options>
Eigen::Transform< Scalar, Dim, Mode, Options >::Transform ( const QTransform< Scalar_, Dim_, Mode_, Options_ > &  other)
inline

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

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

Definition at line 792 of file Transform.h.

793 {
794  check_template_params();
795  *this = other;
796 }

◆ Transform() [10/10]

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)
inlineexplicit

Copy constructor with scalar type conversion

Definition at line 637 of file Transform.h.

638  {
639  check_template_params();
640  m_matrix = other.matrix().template cast<Scalar>();
641  }

Member Function Documentation

◆ affine() [1/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
AffinePart Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::affine ( )
inline
Returns
a writable expression of the Dim x HDim affine part of the transformation

Definition at line 405 of file Transform.h.

405 { return take_affine_part::run(m_matrix); }

◆ affine() [2/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
ConstAffinePart Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::affine ( ) const
inline
Returns
a read-only expression of the Dim x HDim affine part of the transformation

Definition at line 403 of file Transform.h.

403 { return take_affine_part::run(m_matrix); }

◆ cast()

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

Definition at line 632 of file Transform.h.

633  { return typename internal::cast_return_type<Transform,Transform<NewScalarType,Dim,Mode,Options> >::type(*this); }

◆ cols()

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
EIGEN_CONSTEXPR Index Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::cols ( void  ) const
inline

Definition at line 383 of file Transform.h.

383 { return m_matrix.cols(); }
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT

◆ computeRotationScaling()

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>
See also
computeScalingRotation(), rotation(), class SVD

Definition at line 1105 of file Transform.h.

1106 {
1107  // Note that JacobiSVD is faster than BDCSVD for small matrices.
1108  JacobiSVD<LinearMatrixType, ComputeFullU | ComputeFullV> svd(linear());
1109 
1110  Scalar x = (svd.matrixU() * svd.matrixV().adjoint()).determinant() < Scalar(0) ? Scalar(-1) : Scalar(1); // so x has absolute value 1
1111  VectorType sv(svd.singularValues());
1112  sv.coeffRef(Dim-1) *= x;
1113  if(scaling) *scaling = svd.matrixV() * sv.asDiagonal() * svd.matrixV().adjoint();
1114  if(rotation)
1115  {
1116  LinearMatrixType m(svd.matrixU());
1117  m.col(Dim-1) *= x;
1118  *rotation = m * svd.matrixV().adjoint();
1119  }
1120 }
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf, ComputeThinU|ComputeThinV > svd(m)
Scalar_ Scalar
Definition: Transform.h:216
Matrix< Scalar, Dim, 1 > VectorType
Definition: Transform.h:240
RotationReturnType rotation() const
Definition: Transform.h:1086
Matrix< Scalar, Dim, Dim, Options > LinearMatrixType
Definition: Transform.h:226

◆ computeScalingRotation()

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>
See also
computeRotationScaling(), rotation(), class SVD

Definition at line 1135 of file Transform.h.

1136 {
1137  // Note that JacobiSVD is faster than BDCSVD for small matrices.
1138  JacobiSVD<LinearMatrixType, ComputeFullU | ComputeFullV> svd(linear());
1139 
1140  Scalar x = (svd.matrixU() * svd.matrixV().adjoint()).determinant() < Scalar(0) ? Scalar(-1) : Scalar(1); // so x has absolute value 1
1141  VectorType sv(svd.singularValues());
1142  sv.coeffRef(Dim-1) *= x;
1143  if(scaling) *scaling = svd.matrixU() * sv.asDiagonal() * svd.matrixU().adjoint();
1144  if(rotation)
1145  {
1146  LinearMatrixType m(svd.matrixU());
1147  m.col(Dim-1) *= x;
1148  *rotation = m * svd.matrixV().adjoint();
1149  }
1150 }

◆ data() [1/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
Scalar* Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::data ( )
inline
Returns
a non-const pointer to the column major internal matrix

Definition at line 624 of file Transform.h.

624 { return m_matrix.data(); }
const Scalar * data() const

◆ data() [2/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
const Scalar* Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::data ( ) const
inline
Returns
a const pointer to the column major internal matrix

Definition at line 622 of file Transform.h.

622 { return m_matrix.data(); }

◆ EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE()

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) 
)
inline

< space dimension in which the transformation holds

< size of a respective homogeneous vector

Definition at line 209 of file Transform.h.

209  : (Dim_+1)*(Dim_+1))
210  enum {
211  Mode = Mode_,
212  Options = Options_,
213  Dim = Dim_,
214  HDim = Dim_+1,
215  Rows = int(Mode)==(AffineCompact) ? Dim : HDim
216  };

◆ fromPositionOrientationScale() [1/2]

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 
)

◆ fromPositionOrientationScale() [2/2]

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.

Definition at line 1158 of file Transform.h.

1160 {
1161  linear() = internal::toRotationMatrix<Scalar,Dim>(orientation);
1162  linear() *= scale.asDiagonal();
1163  translation() = position;
1164  makeAffine();
1165  return *this;
1166 }
Transform & scale(const MatrixBase< OtherDerived > &other)

◆ Identity()

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
static const Transform Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::Identity ( )
inlinestatic

Returns an identity transformation.

Definition at line 537 of file Transform.h.

538  {
540  }
static const IdentityReturnType Identity()

◆ inverse()

template<typename Scalar , int Dim, int Mode, int Options>
Transform< Scalar, Dim, Mode, Options > Eigen::Transform< Scalar, Dim, Mode, Options >::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()

Definition at line 1230 of file Transform.h.

1231 {
1232  Transform res;
1233  if (hint == Projective)
1234  {
1235  internal::projective_transform_inverse<Transform>::run(*this, res);
1236  }
1237  else
1238  {
1239  if (hint == Isometry)
1240  {
1241  res.matrix().template topLeftCorner<Dim,Dim>() = linear().transpose();
1242  }
1243  else if(hint&Affine)
1244  {
1245  res.matrix().template topLeftCorner<Dim,Dim>() = linear().inverse();
1246  }
1247  else
1248  {
1249  eigen_assert(false && "Invalid transform traits in Transform::Inverse");
1250  }
1251  // translation and remaining parts
1252  res.matrix().template topRightCorner<Dim,1>()
1253  = - res.matrix().template topLeftCorner<Dim,Dim>() * translation();
1254  res.makeAffine(); // we do need this, because in the beginning res is uninitialized
1255  }
1256  return res;
1257 }
#define eigen_assert(x)
Definition: Macros.h:902
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res

◆ isApprox()

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
bool Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::isApprox ( const Transform< Scalar_, Dim_, Mode_, Options_ > &  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()

Definition at line 647 of file Transform.h.

648  { return m_matrix.isApprox(other.m_matrix, prec); }
bool isApprox(const DenseBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Fuzzy.h:105

◆ linear() [1/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
LinearPart Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::linear ( )
inline
Returns
a writable expression of the linear part of the transformation

Definition at line 400 of file Transform.h.

400 { return LinearPart(m_matrix,0,0); }
Block< MatrixType, Dim, Dim, int(Mode)==(AffineCompact) &&(int(Options)&RowMajor)==0 > LinearPart
Definition: Transform.h:228

◆ linear() [2/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
ConstLinearPart Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::linear ( ) const
inline
Returns
a read-only expression of the linear part of the transformation

Definition at line 398 of file Transform.h.

398 { return ConstLinearPart(m_matrix,0,0); }
const Block< ConstMatrixType, Dim, Dim, int(Mode)==(AffineCompact) &&(int(Options)&RowMajor)==0 > ConstLinearPart
Definition: Transform.h:230

◆ linearExt() [1/2]

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 ( )
inline

Definition at line 661 of file Transform.h.

661  :Dim,Dim> linearExt()
662  { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,Dim>(0,0); }
Block< MatrixType, int(Mode)==int(Projective)?HDim:Dim, Dim > linearExt()
Definition: Transform.h:661

◆ linearExt() [2/2]

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
inline

Definition at line 667 of file Transform.h.

667  :Dim,Dim> linearExt() const
668  { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,Dim>(0,0); }

◆ makeAffine()

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
void Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::makeAffine ( )
inline

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

Definition at line 652 of file Transform.h.

653  {
654  internal::transform_make_affine<int(Mode)>::run(m_matrix);
655  }

◆ matrix() [1/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
MatrixType& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::matrix ( )
inline
Returns
a writable expression of the transformation matrix

Definition at line 395 of file Transform.h.

395 { return m_matrix; }

◆ matrix() [2/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
const MatrixType& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::matrix ( ) const
inline
Returns
a read-only expression of the transformation matrix

Definition at line 393 of file Transform.h.

393 { return m_matrix; }

◆ operator()() [1/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
Scalar& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::operator() ( Index  row,
Index  col 
)
inline

shortcut for m_matrix(row,col);

See also
MatrixBase::operator(Index,Index)

Definition at line 390 of file Transform.h.

390 { return m_matrix(row,col); }
RowXpr row(Index i)
This is the const version of row(). *‍/.
ColXpr col(Index i)
This is the const version of col().

◆ operator()() [2/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
Scalar Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::operator() ( Index  row,
Index  col 
) const
inline

shortcut for m_matrix(row,col);

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

Definition at line 387 of file Transform.h.

387 { return m_matrix(row,col); }

◆ operator*() [1/8]

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

Definition at line 462 of file Transform.h.

463  {
465  res.linearExt() *= b;
466  return res;
467  }
Array< int, 3, 1 > b
Transform< Scalar, Dim, TransformTimeDiagonalMode > TransformTimeDiagonalReturnType
Definition: Transform.h:251

◆ operator*() [2/8]

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:

  • 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:

v2 = A.linear() * v1;
MatrixXcf A
Map< RowVectorXf > v2(M2.data(), M2.size())
M1<< 1, 2, 3, 4, 5, 6, 7, 8, 9;Map< RowVectorXf > v1(M1.data(), M1.size())
Transform< float, 3, Affine > Affine3f
Definition: Transform.h:710
Matrix< float, 3, 1 > Vector3f
3×1 vector of type float.
Definition: Matrix.h:501

Definition at line 439 of file Transform.h.

440  { return internal::transform_right_product_impl<Transform, OtherDerived>::run(*this,other.derived()); }

◆ operator*() [3/8]

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
inline

◆ operator*() [4/8]

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
inline

Definition at line 1041 of file Transform.h.

1042 {
1043  Transform res = *this;
1044  res.rotate(r.derived());
1045  return res;
1046 }

◆ operator*() [5/8]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
const Transform Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::operator* ( const Transform< Scalar_, Dim_, Mode_, Options_ > &  other) const
inline

Concatenates two transformations

Definition at line 491 of file Transform.h.

492  {
493  return internal::transform_transform_product_impl<Transform,Transform>::run(*this,other);
494  }

◆ operator*() [6/8]

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
inline

Concatenates two different transformations

Definition at line 524 of file Transform.h.

525  {
526  return internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> >::run(*this,other);
527  }

◆ operator*() [7/8]

template<typename Scalar , int Dim, int Mode, int Options>
Transform< Scalar, Dim, Mode, Options > Eigen::Transform< Scalar, Dim, Mode, Options >::operator* ( const TranslationType t) const
inline

Definition at line 1013 of file Transform.h.

1014 {
1015  Transform res = *this;
1016  res.translate(t.vector());
1017  return res;
1018 }

◆ operator*() [8/8]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
TransformTimeDiagonalReturnType Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::operator* ( const UniformScaling< Scalar > &  s) const
inline

Definition at line 586 of file Transform.h.

587  {
589  res.scale(s.factor());
590  return res;
591  }

◆ operator*=() [1/5]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::operator*= ( const DiagonalMatrix< Scalar, Dim > &  s)
inline

Definition at line 594 of file Transform.h.

594 { linearExt() *= s; return *this; }

◆ operator*=() [2/5]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename OtherDerived >
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::operator*= ( const EigenBase< OtherDerived > &  other)
inline

Definition at line 488 of file Transform.h.

488 { return *this = *this * other; }

◆ operator*=() [3/5]

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)
inline

Definition at line 599 of file Transform.h.

599 { return rotate(r.toRotationMatrix()); }
Transform & rotate(const RotationType &rotation)

◆ operator*=() [4/5]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::operator*= ( const TranslationType t)
inline

Definition at line 575 of file Transform.h.

575 { return translate(t.vector()); }
Transform & translate(const MatrixBase< OtherDerived > &other)

◆ operator*=() [5/5]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::operator*= ( const UniformScaling< Scalar > &  s)
inline

Definition at line 583 of file Transform.h.

583 { return scale(s.factor()); }

◆ operator=() [1/7]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename OtherDerived >
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::operator= ( const EigenBase< OtherDerived > &  other)
inline

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

Definition at line 299 of file Transform.h.

300  {
301  EIGEN_STATIC_ASSERT((internal::is_same<Scalar,typename OtherDerived::Scalar>::value),
302  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY);
303 
304  internal::transform_construct_from_matrix<OtherDerived,Mode,Options,Dim,HDim>::run(this, other.derived());
305  return *this;
306  }

◆ operator=() [2/7]

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

Definition at line 803 of file Transform.h.

804 {
805  check_template_params();
806  EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
807  if (EIGEN_CONST_CONDITIONAL(Mode == int(AffineCompact)))
808  m_matrix << other.m11(), other.m21(), other.dx(),
809  other.m12(), other.m22(), other.dy();
810  else
811  m_matrix << other.m11(), other.m21(), other.dx(),
812  other.m12(), other.m22(), other.dy(),
813  other.m13(), other.m23(), other.m33();
814  return *this;
815 }

◆ operator=() [3/7]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename OtherDerived >
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::operator= ( const ReturnByValue< OtherDerived > &  other)
inline

Definition at line 365 of file Transform.h.

366  {
367  other.evalTo(*this);
368  return *this;
369  }

◆ operator=() [4/7]

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)
inline

◆ operator=() [5/7]

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)
inline

Definition at line 1031 of file Transform.h.

1032 {
1033  linear() = internal::toRotationMatrix<Scalar,Dim>(r);
1034  translation().setZero();
1035  makeAffine();
1036  return *this;
1037 }

◆ operator=() [6/7]

template<typename Scalar , int Dim, int Mode, int Options>
Transform< Scalar, Dim, Mode, Options > & Eigen::Transform< Scalar, Dim, Mode, Options >::operator= ( const TranslationType t)
inline

Scaling, Translation and Rotation compatibility ***

Definition at line 1004 of file Transform.h.

1005 {
1006  linear().setIdentity();
1007  translation() = t.vector();
1008  makeAffine();
1009  return *this;
1010 }

◆ operator=() [7/7]

template<typename Scalar , int Dim, int Mode, int Options>
Transform< Scalar, Dim, Mode, Options > & Eigen::Transform< Scalar, Dim, Mode, Options >::operator= ( const UniformScaling< Scalar > &  t)
inline

Definition at line 1021 of file Transform.h.

1022 {
1023  m_matrix.setZero();
1024  linear().diagonal().fill(s.factor());
1025  makeAffine();
1026  return *this;
1027 }
Derived & setZero(Index size)

◆ prerotate() [1/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename RotationType >
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::prerotate ( const RotationType &  rotation)
inline

◆ prerotate() [2/2]

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.

See also
rotate()

Definition at line 961 of file Transform.h.

962 {
963  m_matrix.template block<Dim,HDim>(0,0) = internal::toRotationMatrix<Scalar,Dim>(rotation)
964  * m_matrix.template block<Dim,HDim>(0,0);
965  return *this;
966 }

◆ prescale() [1/3]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename OtherDerived >
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::prescale ( const MatrixBase< OtherDerived > &  other)
inline

◆ prescale() [2/3]

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.

See also
scale()

Definition at line 874 of file Transform.h.

875 {
876  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
877  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
878  affine().noalias() = (other.asDiagonal() * affine());
879  return *this;
880 }
#define EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(TYPE, SIZE)
Definition: StaticAssert.h:51
ConstAffinePart affine() const
Definition: Transform.h:403

◆ prescale() [3/3]

template<typename Scalar , int Dim, int Mode, int Options>
Transform< Scalar, Dim, Mode, Options > & Eigen::Transform< Scalar, Dim, Mode, Options >::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)

Definition at line 887 of file Transform.h.

888 {
889  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
890  m_matrix.template topRows<Dim>() *= s;
891  return *this;
892 }

◆ preshear()

template<typename Scalar , int Dim, int Mode, int Options>
Transform< Scalar, Dim, Mode, Options > & 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.

Warning
2D only.
See also
shear()

Definition at line 991 of file Transform.h.

992 {
993  EIGEN_STATIC_ASSERT(int(Dim)==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
994  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
995  m_matrix.template block<Dim,HDim>(0,0) = LinearMatrixType(1, sx, sy, 1) * m_matrix.template block<Dim,HDim>(0,0);
996  return *this;
997 }

◆ pretranslate() [1/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename OtherDerived >
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::pretranslate ( const MatrixBase< OtherDerived > &  other)
inline

◆ pretranslate() [2/2]

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.

See also
translate()

Definition at line 915 of file Transform.h.

916 {
917  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
918  if(EIGEN_CONST_CONDITIONAL(int(Mode)==int(Projective)))
919  affine() += other * m_matrix.row(Dim);
920  else
921  translation() += other;
922  return *this;
923 }

◆ rotate() [1/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename RotationType >
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::rotate ( const RotationType &  rotation)
inline

◆ rotate() [2/2]

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.

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

Definition at line 945 of file Transform.h.

946 {
947  linearExt() *= internal::toRotationMatrix<Scalar,Dim>(rotation);
948  return *this;
949 }

◆ rotation()

template<typename Scalar , int Dim, int Mode, int Options>
Transform< Scalar, Dim, Mode, Options >::RotationReturnType Eigen::Transform< Scalar, Dim, Mode, Options >::rotation
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

Definition at line 1086 of file Transform.h.

1087 {
1088  return internal::transform_rotation_impl<Mode>::run(*this);
1089 }

◆ rows()

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
EIGEN_CONSTEXPR Index Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::rows ( void  ) const
inline

Definition at line 382 of file Transform.h.

382 { return int(Mode)==int(Projective) ? m_matrix.cols() : (m_matrix.cols()-1); }

◆ scale() [1/3]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename OtherDerived >
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::scale ( const MatrixBase< OtherDerived > &  other)
inline

◆ scale() [2/3]

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)

Procedural API *** 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()

Definition at line 847 of file Transform.h.

848 {
849  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
850  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
851  linearExt().noalias() = (linearExt() * other.asDiagonal());
852  return *this;
853 }

◆ scale() [3/3]

template<typename Scalar , int Dim, int Mode, int Options>
Transform< Scalar, Dim, Mode, Options > & Eigen::Transform< Scalar, Dim, Mode, Options >::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)

Definition at line 860 of file Transform.h.

861 {
862  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
863  linearExt() *= s;
864  return *this;
865 }

◆ setIdentity()

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
void Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::setIdentity ( )
inline
See also
MatrixBase::setIdentity()

Definition at line 531 of file Transform.h.

531 { m_matrix.setIdentity(); }
Derived & setIdentity()

◆ shear()

template<typename Scalar , int Dim, int Mode, int Options>
Transform< Scalar, Dim, Mode, Options > & 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.

Warning
2D only.
See also
preshear()

Definition at line 975 of file Transform.h.

976 {
977  EIGEN_STATIC_ASSERT(int(Dim)==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
978  EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
979  VectorType tmp = linear().col(0)*sy + linear().col(1);
980  linear() << linear().col(0) + linear().col(1)*sx, tmp;
981  return *this;
982 }

◆ toQTransform()

template<typename Scalar , int Dim, int Mode, int Options>
QTransform Eigen::Transform< Scalar, Dim, Mode, Options >::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.

Definition at line 822 of file Transform.h.

823 {
824  EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
825  if (EIGEN_CONST_CONDITIONAL(Mode == int(AffineCompact)))
826  return QTransform(m_matrix.coeff(0,0), m_matrix.coeff(1,0),
827  m_matrix.coeff(0,1), m_matrix.coeff(1,1),
828  m_matrix.coeff(0,2), m_matrix.coeff(1,2));
829  else
830  return QTransform(m_matrix.coeff(0,0), m_matrix.coeff(1,0), m_matrix.coeff(2,0),
831  m_matrix.coeff(0,1), m_matrix.coeff(1,1), m_matrix.coeff(2,1),
832  m_matrix.coeff(0,2), m_matrix.coeff(1,2), m_matrix.coeff(2,2));
833 }
constexpr const Scalar & coeff(Index rowId, Index colId) const

◆ translate() [1/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
template<typename OtherDerived >
Transform& Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::translate ( const MatrixBase< OtherDerived > &  other)
inline

◆ translate() [2/2]

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.

See also
pretranslate()

Definition at line 901 of file Transform.h.

902 {
903  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
904  translationExt() += linearExt() * other;
905  return *this;
906 }
Block< MatrixType, int(Mode)==int(Projective)?HDim:Dim, 1 > translationExt()
Definition: Transform.h:674

◆ translation() [1/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
TranslationPart Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::translation ( )
inline
Returns
a writable expression of the translation vector of the transformation

Definition at line 410 of file Transform.h.

410 { return TranslationPart(m_matrix,0,Dim); }
Block< MatrixType, Dim, 1,!(internal::traits< MatrixType >::Flags &RowMajorBit)> TranslationPart
Definition: Transform.h:242

◆ translation() [2/2]

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
ConstTranslationPart Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::translation ( ) const
inline
Returns
a read-only expression of the translation vector of the transformation

Definition at line 408 of file Transform.h.

408 { return ConstTranslationPart(m_matrix,0,Dim); }
const Block< ConstMatrixType, Dim, 1,!(internal::traits< MatrixType >::Flags &RowMajorBit)> ConstTranslationPart
Definition: Transform.h:244

◆ translationExt() [1/2]

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 ( )
inline

Definition at line 674 of file Transform.h.

674  :Dim,1> translationExt()
675  { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,1>(0,Dim); }

◆ translationExt() [2/2]

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
inline

Definition at line 680 of file Transform.h.

680  :Dim,1> translationExt() const
681  { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,1>(0,Dim); }

Member Data Documentation

◆ m_matrix

template<typename Scalar_ , int Dim_, int Mode_, int Options_>
MatrixType Eigen::Transform< Scalar_, Dim_, Mode_, Options_ >::m_matrix
protected

Definition at line 255 of file Transform.h.


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