Newer
Older
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
// Copyright (C) 2010 Hauke Heibel <hauke.heibel@gmail.com>
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
#ifndef EIGEN_TRANSFORM_H
#define EIGEN_TRANSFORM_H
namespace Eigen {
namespace internal {
template<typename Transform>
struct transform_traits
{
enum
{
Dim = Transform::Dim,
HDim = Transform::HDim,
Mode = Transform::Mode,
IsProjective = (int(Mode)==int(Projective))
};
};
template< typename TransformType,
typename MatrixType,
int Case = transform_traits<TransformType>::IsProjective ? 0
: int(MatrixType::RowsAtCompileTime) == int(transform_traits<TransformType>::HDim) ? 1
struct transform_right_product_impl;
template< typename Other,
int Mode,
int Options,
int Dim,
int HDim,
int OtherRows=Other::RowsAtCompileTime,
int OtherCols=Other::ColsAtCompileTime>
struct transform_left_product_impl;
template< typename Lhs,
typename Rhs,
bool AnyProjective =
transform_traits<Lhs>::IsProjective ||
transform_traits<Rhs>::IsProjective>
struct transform_transform_product_impl;
template< typename Other,
int Mode,
int Options,
int Dim,
int HDim,
int OtherRows=Other::RowsAtCompileTime,
int OtherCols=Other::ColsAtCompileTime>
struct transform_construct_from_matrix;
template<typename TransformType> struct transform_take_affine_part;
template<typename _Scalar, int _Dim, int _Mode, int _Options>
struct traits<Transform<_Scalar,_Dim,_Mode,_Options> >
{
typedef _Scalar Scalar;
typedef Eigen::Index StorageIndex;
typedef Dense StorageKind;
enum {
Dim1 = _Dim==Dynamic ? _Dim : _Dim + 1,
RowsAtCompileTime = _Mode==Projective ? Dim1 : _Dim,
ColsAtCompileTime = Dim1,
MaxRowsAtCompileTime = RowsAtCompileTime,
MaxColsAtCompileTime = ColsAtCompileTime,
Flags = 0
};
};
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
template<int Mode> struct transform_make_affine;
} // end namespace internal
/** \geometry_module \ingroup Geometry_Module
*
* \class Transform
*
* \brief Represents an homogeneous transformation in a N dimensional space
*
* \tparam _Scalar the scalar type, i.e., the type of the coefficients
* \tparam _Dim the dimension of the space
* \tparam _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.
* \tparam _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:
*
* \code v' = T * v \endcode
*
* Therefore, an affine transformation matrix M is shaped like this:
*
* \f$ \left( \begin{array}{cc}
* linear & translation\\
* 0 ... 0 & 1
* \end{array} \right) \f$
*
* 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,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:
*
* \b Translation t (Dim)x(1):
* \f$ \left( \begin{array}{cc}
* I & t \\
* 0\,...\,0 & 1
* \end{array} \right) \f$
*
* \b Rotation R (Dim)x(Dim):
* \f$ \left( \begin{array}{cc}
* R & 0\\
* 0\,...\,0 & 1
* \end{array} \right) \f$
*<!--
* \b Linear \b Matrix L (Dim)x(Dim):
* \f$ \left( \begin{array}{cc}
* L & 0\\
* 0\,...\,0 & 1
* \end{array} \right) \f$
*
* \b Affine \b Matrix A (Dim)x(Dim+1):
* \f$ \left( \begin{array}{c}
* A\\
* 0\,...\,0\,1
* \end{array} \right) \f$
*-->
* \b Scaling \b DiagonalMatrix S (Dim)x(Dim):
* \f$ \left( \begin{array}{cc}
* S & 0\\
* 0\,...\,0 & 1
* \end{array} \right) \f$
*
* \b Column \b point v (Dim)x(1):
* \f$ \left( \begin{array}{c}
* v\\
* 1
* \end{array} \right) \f$
*
* \b Set \b of \b column \b points V1...Vn (Dim)x(n):
* \f$ \left( \begin{array}{ccc}
* v_1 & ... & v_n\\
* 1 & ... & 1
* \end{array} \right) \f$
*
* 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:
* \code
* m' = T * m.colwise().homogeneous();
* \endcode
* 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
* \ref TopicCustomizing_Plugins by defining the preprocessor symbol \c EIGEN_TRANSFORM_PLUGIN.
*
* \sa class Matrix, class Quaternion
*/
template<typename _Scalar, int _Dim, int _Mode, int _Options>
class Transform
{
public:
EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_Dim==Dynamic ? Dynamic : (_Dim+1)*(_Dim+1))
enum {
Mode = _Mode,
Options = _Options,
Dim = _Dim, ///< space dimension in which the transformation holds
HDim = _Dim+1, ///< size of a respective homogeneous vector
Rows = int(Mode)==(AffineCompact) ? Dim : HDim
};
/** the scalar type of the coefficients */
typedef _Scalar Scalar;
typedef Eigen::Index StorageIndex;
typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
/** type of the matrix used to represent the transformation */
typedef typename internal::make_proper_matrix_type<Scalar,Rows,HDim,Options>::type MatrixType;
/** constified MatrixType */
typedef const MatrixType ConstMatrixType;
/** type of the matrix used to represent the linear part of the transformation */
typedef Matrix<Scalar,Dim,Dim,Options> LinearMatrixType;
/** type of read/write reference to the linear part of the transformation */
typedef Block<MatrixType,Dim,Dim,int(Mode)==(AffineCompact) && (Options&RowMajor)==0> LinearPart;
/** type of read reference to the linear part of the transformation */
typedef const Block<ConstMatrixType,Dim,Dim,int(Mode)==(AffineCompact) && (Options&RowMajor)==0> ConstLinearPart;
/** type of read/write reference to the affine part of the transformation */
typedef typename internal::conditional<int(Mode)==int(AffineCompact),
MatrixType&,
Block<MatrixType,Dim,HDim> >::type AffinePart;
/** type of read reference to the affine part of the transformation */
typedef typename internal::conditional<int(Mode)==int(AffineCompact),
const MatrixType&,
const Block<const MatrixType,Dim,HDim> >::type ConstAffinePart;
/** type of a vector */
typedef Matrix<Scalar,Dim,1> VectorType;
/** type of a read/write reference to the translation part of the rotation */
typedef Block<MatrixType,Dim,1,!(internal::traits<MatrixType>::Flags & RowMajorBit)> TranslationPart;
/** type of a read reference to the translation part of the rotation */
typedef const Block<ConstMatrixType,Dim,1,!(internal::traits<MatrixType>::Flags & RowMajorBit)> ConstTranslationPart;
/** corresponding translation type */
typedef Translation<Scalar,Dim> TranslationType;
// this intermediate enum is needed to avoid an ICE with gcc 3.4 and 4.0
enum { TransformTimeDiagonalMode = ((Mode==int(Isometry))?Affine:int(Mode)) };
/** The return type of the product between a diagonal matrix and a transform */
typedef Transform<Scalar,Dim,TransformTimeDiagonalMode> TransformTimeDiagonalReturnType;
protected:
MatrixType m_matrix;
public:
/** Default constructor without initialization of the meaningful coefficients.
* If Mode==Affine, then the last row is set to [0 ... 0 1] */
{
check_template_params();
internal::transform_make_affine<(int(Mode)==Affine) ? Affine : AffineCompact>::run(m_matrix);
}
EIGEN_DEVICE_FUNC inline Transform(const Transform& other)
{
check_template_params();
m_matrix = other.m_matrix;
}
EIGEN_DEVICE_FUNC inline explicit Transform(const TranslationType& t)
EIGEN_DEVICE_FUNC inline explicit Transform(const UniformScaling<Scalar>& s)
{
check_template_params();
*this = s;
}
template<typename Derived>
EIGEN_DEVICE_FUNC inline explicit Transform(const RotationBase<Derived, Dim>& r)
EIGEN_DEVICE_FUNC inline Transform& operator=(const Transform& other)
{ m_matrix = other.m_matrix; return *this; }
typedef internal::transform_take_affine_part<Transform> take_affine_part;
/** Constructs and initializes a transformation from a Dim^2 or a (Dim+1)^2 matrix. */
template<typename OtherDerived>
EIGEN_DEVICE_FUNC inline explicit Transform(const EigenBase<OtherDerived>& other)
{
EIGEN_STATIC_ASSERT((internal::is_same<Scalar,typename OtherDerived::Scalar>::value),
YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY);
check_template_params();
internal::transform_construct_from_matrix<OtherDerived,Mode,Options,Dim,HDim>::run(this, other.derived());
}
/** Set \c *this from a Dim^2 or (Dim+1)^2 matrix. */
template<typename OtherDerived>
EIGEN_DEVICE_FUNC inline Transform& operator=(const EigenBase<OtherDerived>& other)
{
EIGEN_STATIC_ASSERT((internal::is_same<Scalar,typename OtherDerived::Scalar>::value),
YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY);
internal::transform_construct_from_matrix<OtherDerived,Mode,Options,Dim,HDim>::run(this, other.derived());
return *this;
}
template<int OtherOptions>
EIGEN_DEVICE_FUNC inline Transform(const Transform<Scalar,Dim,Mode,OtherOptions>& other)
{
check_template_params();
// only the options change, we can directly copy the matrices
m_matrix = other.matrix();
}
template<int OtherMode,int OtherOptions>
EIGEN_DEVICE_FUNC inline Transform(const Transform<Scalar,Dim,OtherMode,OtherOptions>& other)
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
{
check_template_params();
// prevent conversions as:
// Affine | AffineCompact | Isometry = Projective
EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(OtherMode==int(Projective), Mode==int(Projective)),
YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION)
// prevent conversions as:
// Isometry = Affine | AffineCompact
EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(OtherMode==int(Affine)||OtherMode==int(AffineCompact), Mode!=int(Isometry)),
YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION)
enum { ModeIsAffineCompact = Mode == int(AffineCompact),
OtherModeIsAffineCompact = OtherMode == int(AffineCompact)
};
if(ModeIsAffineCompact == OtherModeIsAffineCompact)
{
// We need the block expression because the code is compiled for all
// combinations of transformations and will trigger a compile time error
// if one tries to assign the matrices directly
m_matrix.template block<Dim,Dim+1>(0,0) = other.matrix().template block<Dim,Dim+1>(0,0);
makeAffine();
}
else if(OtherModeIsAffineCompact)
{
typedef typename Transform<Scalar,Dim,OtherMode,OtherOptions>::MatrixType OtherMatrixType;
internal::transform_construct_from_matrix<OtherMatrixType,Mode,Options,Dim,HDim>::run(this, other.matrix());
}
else
{
// here we know that Mode == AffineCompact and OtherMode != AffineCompact.
// if OtherMode were Projective, the static assert above would already have caught it.
// So the only possibility is that OtherMode == Affine
linear() = other.linear();
translation() = other.translation();
}
}
template<typename OtherDerived>
EIGEN_DEVICE_FUNC Transform(const ReturnByValue<OtherDerived>& other)
{
check_template_params();
other.evalTo(*this);
}
template<typename OtherDerived>
EIGEN_DEVICE_FUNC Transform& operator=(const ReturnByValue<OtherDerived>& other)
{
other.evalTo(*this);
return *this;
}
#ifdef EIGEN_QT_SUPPORT
inline Transform(const QMatrix& other);
inline Transform& operator=(const QMatrix& other);
inline QMatrix toQMatrix(void) const;
inline Transform(const QTransform& other);
inline Transform& operator=(const QTransform& other);
inline QTransform toQTransform(void) const;
#endif
EIGEN_DEVICE_FUNC Index rows() const { return int(Mode)==int(Projective) ? m_matrix.cols() : (m_matrix.cols()-1); }
EIGEN_DEVICE_FUNC Index cols() const { return m_matrix.cols(); }
/** shortcut for m_matrix(row,col);
* \sa MatrixBase::operator(Index,Index) const */
EIGEN_DEVICE_FUNC inline Scalar operator() (Index row, Index col) const { return m_matrix(row,col); }
/** shortcut for m_matrix(row,col);
* \sa MatrixBase::operator(Index,Index) */
EIGEN_DEVICE_FUNC inline Scalar& operator() (Index row, Index col) { return m_matrix(row,col); }
/** \returns a read-only expression of the transformation matrix */
EIGEN_DEVICE_FUNC inline const MatrixType& matrix() const { return m_matrix; }
/** \returns a writable expression of the transformation matrix */
EIGEN_DEVICE_FUNC inline MatrixType& matrix() { return m_matrix; }
/** \returns a read-only expression of the linear part of the transformation */
EIGEN_DEVICE_FUNC inline ConstLinearPart linear() const { return ConstLinearPart(m_matrix,0,0); }
/** \returns a writable expression of the linear part of the transformation */
EIGEN_DEVICE_FUNC inline LinearPart linear() { return LinearPart(m_matrix,0,0); }
/** \returns a read-only expression of the Dim x HDim affine part of the transformation */
EIGEN_DEVICE_FUNC inline ConstAffinePart affine() const { return take_affine_part::run(m_matrix); }
/** \returns a writable expression of the Dim x HDim affine part of the transformation */
EIGEN_DEVICE_FUNC inline AffinePart affine() { return take_affine_part::run(m_matrix); }
/** \returns a read-only expression of the translation vector of the transformation */
EIGEN_DEVICE_FUNC inline ConstTranslationPart translation() const { return ConstTranslationPart(m_matrix,0,Dim); }
/** \returns a writable expression of the translation vector of the transformation */
EIGEN_DEVICE_FUNC inline TranslationPart translation() { return TranslationPart(m_matrix,0,Dim); }
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
/** \returns an expression of the product between the transform \c *this and a matrix expression \a other.
*
* The right-hand-side \a other can be either:
* \li an homogeneous vector of size Dim+1,
* \li a set of homogeneous vectors of size Dim+1 x N,
* \li a transformation matrix of size Dim+1 x Dim+1.
*
* Moreover, if \c *this represents an affine transformation (i.e., Mode!=Projective), then \a other can also be:
* \li a point of size Dim (computes: \code this->linear() * other + this->translation()\endcode),
* \li a set of N points as a Dim x N matrix (computes: \code (this->linear() * other).colwise() + this->translation()\endcode),
*
* In all cases, the return type is a matrix or vector of same sizes as the right-hand-side \a other.
*
* If you want to interpret \a 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:
* \code
* Affine3f A;
* Vector3f v1, v2;
* v2 = A.linear() * v1;
* \endcode
*
*/
// note: this function is defined here because some compilers cannot find the respective declaration
template<typename OtherDerived>
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const typename internal::transform_right_product_impl<Transform, OtherDerived>::ResultType
operator * (const EigenBase<OtherDerived> &other) const
{ return internal::transform_right_product_impl<Transform, OtherDerived>::run(*this,other.derived()); }
/** \returns the product expression of a transformation matrix \a a times a transform \a b
*
* The left hand side \a other can be either:
* \li a linear transformation matrix of size Dim x Dim,
* \li an affine transformation matrix of size Dim x Dim+1,
* \li a general transformation matrix of size Dim+1 x Dim+1.
*/
template<typename OtherDerived> friend
EIGEN_DEVICE_FUNC inline const typename internal::transform_left_product_impl<OtherDerived,Mode,Options,_Dim,_Dim+1>::ResultType
operator * (const EigenBase<OtherDerived> &a, const Transform &b)
{ return internal::transform_left_product_impl<OtherDerived,Mode,Options,Dim,HDim>::run(a.derived(),b); }
/** \returns The product expression of a transform \a a times a diagonal matrix \a 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 DiagonalDerived>
EIGEN_DEVICE_FUNC inline const TransformTimeDiagonalReturnType
operator * (const DiagonalBase<DiagonalDerived> &b) const
{
TransformTimeDiagonalReturnType res(*this);
return res;
}
/** \returns The product expression of a diagonal matrix \a a times a transform \a 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.
*/
template<typename DiagonalDerived>
EIGEN_DEVICE_FUNC friend inline TransformTimeDiagonalReturnType
operator * (const DiagonalBase<DiagonalDerived> &a, const Transform &b)
{
TransformTimeDiagonalReturnType res;
res.linear().noalias() = a*b.linear();
res.translation().noalias() = a*b.translation();
if (Mode!=int(AffineCompact))
res.matrix().row(Dim) = b.matrix().row(Dim);
return res;
}
template<typename OtherDerived>
EIGEN_DEVICE_FUNC inline Transform& operator*=(const EigenBase<OtherDerived>& other) { return *this = *this * other; }
EIGEN_DEVICE_FUNC inline const Transform operator * (const Transform& other) const
{
return internal::transform_transform_product_impl<Transform,Transform>::run(*this,other);
}
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
private:
// this intermediate structure permits to workaround a bug in ICC 11:
// error: template instantiation resulted in unexpected function type of "Eigen::Transform<double, 3, 32, 0>
// (const Eigen::Transform<double, 3, 2, 0> &) const"
// (the meaning of a name may have changed since the template declaration -- the type of the template is:
// "Eigen::internal::transform_transform_product_impl<Eigen::Transform<double, 3, 32, 0>,
// Eigen::Transform<double, 3, Mode, Options>, <expression>>::ResultType (const Eigen::Transform<double, 3, Mode, Options> &) const")
//
template<int OtherMode,int OtherOptions> struct icc_11_workaround
{
typedef internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> > ProductType;
typedef typename ProductType::ResultType ResultType;
};
public:
/** Concatenates two different transformations */
template<int OtherMode,int OtherOptions>
inline typename icc_11_workaround<OtherMode,OtherOptions>::ResultType
operator * (const Transform<Scalar,Dim,OtherMode,OtherOptions>& other) const
{
typedef typename icc_11_workaround<OtherMode,OtherOptions>::ProductType ProductType;
return ProductType::run(*this,other);
}
#else
/** Concatenates two different transformations */
template<int OtherMode,int OtherOptions>
EIGEN_DEVICE_FUNC inline typename internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> >::ResultType
operator * (const Transform<Scalar,Dim,OtherMode,OtherOptions>& other) const
{
return internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> >::run(*this,other);
}
#endif
/** \sa MatrixBase::setIdentity() */
EIGEN_DEVICE_FUNC void setIdentity() { m_matrix.setIdentity(); }
/**
* \brief Returns an identity transformation.
* \todo In the future this function should be returning a Transform expression.
*/
{
return Transform(MatrixType::Identity());
}
template<typename OtherDerived>
inline Transform& scale(const MatrixBase<OtherDerived> &other);
template<typename OtherDerived>
inline Transform& prescale(const MatrixBase<OtherDerived> &other);
EIGEN_DEVICE_FUNC inline Transform& scale(const Scalar& s);
EIGEN_DEVICE_FUNC inline Transform& prescale(const Scalar& s);
inline Transform& translate(const MatrixBase<OtherDerived> &other);
template<typename OtherDerived>
inline Transform& pretranslate(const MatrixBase<OtherDerived> &other);
template<typename RotationType>
inline Transform& rotate(const RotationType& rotation);
template<typename RotationType>
EIGEN_DEVICE_FUNC Transform& shear(const Scalar& sx, const Scalar& sy);
EIGEN_DEVICE_FUNC Transform& preshear(const Scalar& sx, const Scalar& sy);
EIGEN_DEVICE_FUNC inline Transform& operator=(const TranslationType& t);
EIGEN_DEVICE_FUNC
inline Transform& operator*=(const TranslationType& t) { return translate(t.vector()); }
EIGEN_DEVICE_FUNC inline Transform operator*(const TranslationType& t) const;
inline Transform& operator*=(const UniformScaling<Scalar>& s) { return scale(s.factor()); }
EIGEN_DEVICE_FUNC
inline TransformTimeDiagonalReturnType operator*(const UniformScaling<Scalar>& s) const
inline Transform& operator*=(const DiagonalMatrix<Scalar,Dim>& s) { linearExt() *= s; return *this; }
EIGEN_DEVICE_FUNC inline Transform& operator=(const RotationBase<Derived,Dim>& r);
EIGEN_DEVICE_FUNC inline Transform& operator*=(const RotationBase<Derived,Dim>& r) { return rotate(r.toRotationMatrix()); }
EIGEN_DEVICE_FUNC inline Transform operator*(const RotationBase<Derived,Dim>& r) const;
EIGEN_DEVICE_FUNC const LinearMatrixType rotation() const;
void computeRotationScaling(RotationMatrixType *rotation, ScalingMatrixType *scaling) const;
template<typename ScalingMatrixType, typename RotationMatrixType>
void computeScalingRotation(ScalingMatrixType *scaling, RotationMatrixType *rotation) const;
template<typename PositionDerived, typename OrientationType, typename ScaleDerived>
Transform& fromPositionOrientationScale(const MatrixBase<PositionDerived> &position,
const OrientationType& orientation, const MatrixBase<ScaleDerived> &scale);
inline Transform inverse(TransformTraits traits = (TransformTraits)Mode) const;
/** \returns a const pointer to the column major internal matrix */
EIGEN_DEVICE_FUNC const Scalar* data() const { return m_matrix.data(); }
/** \returns a non-const pointer to the column major internal matrix */
EIGEN_DEVICE_FUNC Scalar* data() { return m_matrix.data(); }
/** \returns \c *this with scalar type casted to \a NewScalarType
*
* Note that if \a NewScalarType is equal to the current scalar type of \c *this
* then this function smartly returns a const reference to \c *this.
*/
template<typename NewScalarType>
EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<Transform,Transform<NewScalarType,Dim,Mode,Options> >::type cast() const
{ return typename internal::cast_return_type<Transform,Transform<NewScalarType,Dim,Mode,Options> >::type(*this); }
/** Copy constructor with scalar type conversion */
template<typename OtherScalarType>
EIGEN_DEVICE_FUNC inline explicit Transform(const Transform<OtherScalarType,Dim,Mode,Options>& other)
{
check_template_params();
m_matrix = other.matrix().template cast<Scalar>();
}
/** \returns \c true if \c *this is approximately equal to \a other, within the precision
* determined by \a prec.
*
* \sa MatrixBase::isApprox() */
EIGEN_DEVICE_FUNC bool isApprox(const Transform& other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const
{ return m_matrix.isApprox(other.m_matrix, prec); }
/** Sets the last row to [0 ... 0 1]
*/
{
internal::transform_make_affine<int(Mode)>::run(m_matrix);
}
/** \internal
* \returns the Dim x Dim linear part if the transformation is affine,
* and the HDim x Dim part for projective transformations.
*/
EIGEN_DEVICE_FUNC inline Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,Dim> linearExt()
{ return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,Dim>(0,0); }
/** \internal
* \returns the Dim x Dim linear part if the transformation is affine,
* and the HDim x Dim part for projective transformations.
*/
EIGEN_DEVICE_FUNC inline const Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,Dim> linearExt() const
{ return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,Dim>(0,0); }
/** \internal
* \returns the translation part if the transformation is affine,
* and the last column for projective transformations.
*/
EIGEN_DEVICE_FUNC inline Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,1> translationExt()
{ return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,1>(0,Dim); }
/** \internal
* \returns the translation part if the transformation is affine,
* and the last column for projective transformations.
*/
EIGEN_DEVICE_FUNC inline const Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,1> translationExt() const
{ return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,1>(0,Dim); }
#ifdef EIGEN_TRANSFORM_PLUGIN
#include EIGEN_TRANSFORM_PLUGIN
#endif
protected:
#ifndef EIGEN_PARSED_BY_DOXYGEN
EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE void check_template_params()
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
{
EIGEN_STATIC_ASSERT((Options & (DontAlign|RowMajor)) == Options, INVALID_MATRIX_TEMPLATE_PARAMETERS)
}
#endif
};
/** \ingroup Geometry_Module */
typedef Transform<float,2,Isometry> Isometry2f;
/** \ingroup Geometry_Module */
typedef Transform<float,3,Isometry> Isometry3f;
/** \ingroup Geometry_Module */
typedef Transform<double,2,Isometry> Isometry2d;
/** \ingroup Geometry_Module */
typedef Transform<double,3,Isometry> Isometry3d;
/** \ingroup Geometry_Module */
typedef Transform<float,2,Affine> Affine2f;
/** \ingroup Geometry_Module */
typedef Transform<float,3,Affine> Affine3f;
/** \ingroup Geometry_Module */
typedef Transform<double,2,Affine> Affine2d;
/** \ingroup Geometry_Module */
typedef Transform<double,3,Affine> Affine3d;
/** \ingroup Geometry_Module */
typedef Transform<float,2,AffineCompact> AffineCompact2f;
/** \ingroup Geometry_Module */
typedef Transform<float,3,AffineCompact> AffineCompact3f;
/** \ingroup Geometry_Module */
typedef Transform<double,2,AffineCompact> AffineCompact2d;
/** \ingroup Geometry_Module */
typedef Transform<double,3,AffineCompact> AffineCompact3d;
/** \ingroup Geometry_Module */
typedef Transform<float,2,Projective> Projective2f;
/** \ingroup Geometry_Module */
typedef Transform<float,3,Projective> Projective3f;
/** \ingroup Geometry_Module */
typedef Transform<double,2,Projective> Projective2d;
/** \ingroup Geometry_Module */
typedef Transform<double,3,Projective> Projective3d;
/**************************
*** Optional QT support ***
**************************/
#ifdef EIGEN_QT_SUPPORT
/** Initializes \c *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<Scalar,Dim,Mode,Options>::Transform(const QMatrix& other)
{
check_template_params();
*this = other;
}
/** Set \c *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<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::operator=(const QMatrix& other)
{
EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
if (Mode == int(AffineCompact))
m_matrix << other.m11(), other.m21(), other.dx(),
other.m12(), other.m22(), other.dy();
else
m_matrix << other.m11(), other.m21(), other.dx(),
other.m12(), other.m22(), other.dy(),
0, 0, 1;
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
return *this;
}
/** \returns a QMatrix from \c *this assuming the dimension is 2.
*
* \warning this conversion might loss data if \c *this is not affine
*
* This function is available only if the token EIGEN_QT_SUPPORT is defined.
*/
template<typename Scalar, int Dim, int Mode, int Options>
QMatrix Transform<Scalar,Dim,Mode,Options>::toQMatrix(void) const
{
check_template_params();
EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
return QMatrix(m_matrix.coeff(0,0), m_matrix.coeff(1,0),
m_matrix.coeff(0,1), m_matrix.coeff(1,1),
m_matrix.coeff(0,2), m_matrix.coeff(1,2));
}
/** Initializes \c *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<Scalar,Dim,Mode,Options>::Transform(const QTransform& other)
{
check_template_params();
*this = other;
}
/** Set \c *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<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::operator=(const QTransform& other)
{
check_template_params();
EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
if (Mode == int(AffineCompact))
m_matrix << other.m11(), other.m21(), other.dx(),
other.m12(), other.m22(), other.dy();
else
m_matrix << other.m11(), other.m21(), other.dx(),
other.m12(), other.m22(), other.dy(),
other.m13(), other.m23(), other.m33();
return *this;
}
/** \returns a QTransform from \c *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 Transform<Scalar,Dim,Mode,Options>::toQTransform(void) const
{
EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
if (Mode == int(AffineCompact))
return QTransform(m_matrix.coeff(0,0), m_matrix.coeff(1,0),
m_matrix.coeff(0,1), m_matrix.coeff(1,1),
m_matrix.coeff(0,2), m_matrix.coeff(1,2));
else
return QTransform(m_matrix.coeff(0,0), m_matrix.coeff(1,0), m_matrix.coeff(2,0),
m_matrix.coeff(0,1), m_matrix.coeff(1,1), m_matrix.coeff(2,1),
m_matrix.coeff(0,2), m_matrix.coeff(1,2), m_matrix.coeff(2,2));
}
#endif
/*********************
*** Procedural API ***
*********************/
/** Applies on the right the non uniform scale transformation represented
* by the vector \a other to \c *this and returns a reference to \c *this.
* \sa prescale()
*/
template<typename Scalar, int Dim, int Mode, int Options>
template<typename OtherDerived>
Transform<Scalar,Dim,Mode,Options>::scale(const MatrixBase<OtherDerived> &other)
{
EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
linearExt().noalias() = (linearExt() * other.asDiagonal());
return *this;
}
/** Applies on the right a uniform scale of a factor \a c to \c *this
* and returns a reference to \c *this.
* \sa prescale(Scalar)
*/
template<typename Scalar, int Dim, int Mode, int Options>
EIGEN_DEVICE_FUNC inline Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::scale(const Scalar& s)
{
EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
linearExt() *= s;
return *this;
}
/** Applies on the left the non uniform scale transformation represented
* by the vector \a other to \c *this and returns a reference to \c *this.
* \sa scale()
*/
template<typename Scalar, int Dim, int Mode, int Options>
template<typename OtherDerived>
Transform<Scalar,Dim,Mode,Options>::prescale(const MatrixBase<OtherDerived> &other)
{
EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
return *this;
}
/** Applies on the left a uniform scale of a factor \a c to \c *this
* and returns a reference to \c *this.
* \sa scale(Scalar)
*/
template<typename Scalar, int Dim, int Mode, int Options>
EIGEN_DEVICE_FUNC inline Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::prescale(const Scalar& s)
{
EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
m_matrix.template topRows<Dim>() *= s;
return *this;
}
/** Applies on the right the translation matrix represented by the vector \a other
* to \c *this and returns a reference to \c *this.
* \sa pretranslate()
*/
template<typename Scalar, int Dim, int Mode, int Options>
template<typename OtherDerived>
Transform<Scalar,Dim,Mode,Options>::translate(const MatrixBase<OtherDerived> &other)
{
EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
translationExt() += linearExt() * other;
return *this;
}
/** Applies on the left the translation matrix represented by the vector \a other
* to \c *this and returns a reference to \c *this.
* \sa translate()
*/
template<typename Scalar, int Dim, int Mode, int Options>
template<typename OtherDerived>
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
Transform<Scalar,Dim,Mode,Options>::pretranslate(const MatrixBase<OtherDerived> &other)
{
EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim))
if(int(Mode)==int(Projective))
affine() += other * m_matrix.row(Dim);
else
translation() += other;
return *this;
}
/** Applies on the right the rotation represented by the rotation \a rotation
* to \c *this and returns a reference to \c *this.
*
* The template parameter \a RotationType is the type of the rotation which
* must be known by internal::toRotationMatrix<>.
*
* Natively supported types includes:
* - any scalar (2D),
* - a Dim x Dim matrix expression,
* - a Quaternion (3D),
* - a AngleAxis (3D)
*
* This mechanism is easily extendable to support user types such as Euler angles,
* or a pair of Quaternion for 4D rotations.
*
* \sa rotate(Scalar), class Quaternion, class AngleAxis, prerotate(RotationType)
*/
template<typename Scalar, int Dim, int Mode, int Options>
template<typename RotationType>
Transform<Scalar,Dim,Mode,Options>::rotate(const RotationType& rotation)
{
linearExt() *= internal::toRotationMatrix<Scalar,Dim>(rotation);
return *this;
}
/** Applies on the left the rotation represented by the rotation \a rotation
* to \c *this and returns a reference to \c *this.
*
* See rotate() for further details.
*
* \sa rotate()
*/
template<typename Scalar, int Dim, int Mode, int Options>
template<typename RotationType>
Transform<Scalar,Dim,Mode,Options>::prerotate(const RotationType& rotation)
{
m_matrix.template block<Dim,HDim>(0,0) = internal::toRotationMatrix<Scalar,Dim>(rotation)
* m_matrix.template block<Dim,HDim>(0,0);
return *this;
}
/** Applies on the right the shear transformation represented
* by the vector \a other to \c *this and returns a reference to \c *this.
* \warning 2D only.
* \sa preshear()
*/
template<typename Scalar, int Dim, int Mode, int Options>
Transform<Scalar,Dim,Mode,Options>::shear(const Scalar& sx, const Scalar& sy)
{
EIGEN_STATIC_ASSERT(int(Dim)==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
VectorType tmp = linear().col(0)*sy + linear().col(1);
linear() << linear().col(0) + linear().col(1)*sx, tmp;
return *this;
}
/** Applies on the left the shear transformation represented
* by the vector \a other to \c *this and returns a reference to \c *this.
* \warning 2D only.
* \sa shear()
*/
template<typename Scalar, int Dim, int Mode, int Options>
Transform<Scalar,Dim,Mode,Options>::preshear(const Scalar& sx, const Scalar& sy)
{
EIGEN_STATIC_ASSERT(int(Dim)==2, YOU_MADE_A_PROGRAMMING_MISTAKE)
EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS)
m_matrix.template block<Dim,HDim>(0,0) = LinearMatrixType(1, sx, sy, 1) * m_matrix.template block<Dim,HDim>(0,0);
return *this;
}
/******************************************************
*** Scaling, Translation and Rotation compatibility ***
******************************************************/
template<typename Scalar, int Dim, int Mode, int Options>