Eigne Matrix Class Library
Dependents: Eigen_test Odometry_test AttitudeEstimation_usingTicker MPU9250_Quaternion_Binary_Serial ... more
Transform.h
00001 // This file is part of Eigen, a lightweight C++ template library 00002 // for linear algebra. 00003 // 00004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr> 00005 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com> 00006 // Copyright (C) 2010 Hauke Heibel <hauke.heibel@gmail.com> 00007 // 00008 // This Source Code Form is subject to the terms of the Mozilla 00009 // Public License v. 2.0. If a copy of the MPL was not distributed 00010 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 00011 00012 #ifndef EIGEN_TRANSFORM_H 00013 #define EIGEN_TRANSFORM_H 00014 00015 namespace Eigen { 00016 00017 namespace internal { 00018 00019 template<typename Transform> 00020 struct transform_traits 00021 { 00022 enum 00023 { 00024 Dim = Transform::Dim, 00025 HDim = Transform::HDim, 00026 Mode = Transform::Mode, 00027 IsProjective = (int(Mode)==int(Projective)) 00028 }; 00029 }; 00030 00031 template< typename TransformType, 00032 typename MatrixType, 00033 int Case = transform_traits<TransformType>::IsProjective ? 0 00034 : int(MatrixType::RowsAtCompileTime) == int(transform_traits<TransformType>::HDim) ? 1 00035 : 2> 00036 struct transform_right_product_impl; 00037 00038 template< typename Other, 00039 int Mode, 00040 int Options, 00041 int Dim, 00042 int HDim, 00043 int OtherRows=Other::RowsAtCompileTime, 00044 int OtherCols=Other::ColsAtCompileTime> 00045 struct transform_left_product_impl; 00046 00047 template< typename Lhs, 00048 typename Rhs, 00049 bool AnyProjective = 00050 transform_traits<Lhs>::IsProjective || 00051 transform_traits<Rhs>::IsProjective> 00052 struct transform_transform_product_impl; 00053 00054 template< typename Other, 00055 int Mode, 00056 int Options, 00057 int Dim, 00058 int HDim, 00059 int OtherRows=Other::RowsAtCompileTime, 00060 int OtherCols=Other::ColsAtCompileTime> 00061 struct transform_construct_from_matrix; 00062 00063 template<typename TransformType> struct transform_take_affine_part; 00064 00065 template<int Mode> struct transform_make_affine; 00066 00067 } // end namespace internal 00068 00069 /** \geometry_module \ingroup Geometry_Module 00070 * 00071 * \class Transform 00072 * 00073 * \brief Represents an homogeneous transformation in a N dimensional space 00074 * 00075 * \tparam _Scalar the scalar type, i.e., the type of the coefficients 00076 * \tparam _Dim the dimension of the space 00077 * \tparam _Mode the type of the transformation. Can be: 00078 * - #Affine: the transformation is stored as a (Dim+1)^2 matrix, 00079 * where the last row is assumed to be [0 ... 0 1]. 00080 * - #AffineCompact: the transformation is stored as a (Dim)x(Dim+1) matrix. 00081 * - #Projective: the transformation is stored as a (Dim+1)^2 matrix 00082 * without any assumption. 00083 * \tparam _Options has the same meaning as in class Matrix. It allows to specify DontAlign and/or RowMajor. 00084 * These Options are passed directly to the underlying matrix type. 00085 * 00086 * The homography is internally represented and stored by a matrix which 00087 * is available through the matrix() method. To understand the behavior of 00088 * this class you have to think a Transform object as its internal 00089 * matrix representation. The chosen convention is right multiply: 00090 * 00091 * \code v' = T * v \endcode 00092 * 00093 * Therefore, an affine transformation matrix M is shaped like this: 00094 * 00095 * \f$ \left( \begin{array}{cc} 00096 * linear & translation\\ 00097 * 0 ... 0 & 1 00098 * \end{array} \right) \f$ 00099 * 00100 * Note that for a projective transformation the last row can be anything, 00101 * and then the interpretation of different parts might be sightly different. 00102 * 00103 * However, unlike a plain matrix, the Transform class provides many features 00104 * simplifying both its assembly and usage. In particular, it can be composed 00105 * with any other transformations (Transform,Translation,RotationBase,DiagonalMatrix) 00106 * and can be directly used to transform implicit homogeneous vectors. All these 00107 * operations are handled via the operator*. For the composition of transformations, 00108 * its principle consists to first convert the right/left hand sides of the product 00109 * to a compatible (Dim+1)^2 matrix and then perform a pure matrix product. 00110 * Of course, internally, operator* tries to perform the minimal number of operations 00111 * according to the nature of each terms. Likewise, when applying the transform 00112 * to points, the latters are automatically promoted to homogeneous vectors 00113 * before doing the matrix product. The conventions to homogeneous representations 00114 * are performed as follow: 00115 * 00116 * \b Translation t (Dim)x(1): 00117 * \f$ \left( \begin{array}{cc} 00118 * I & t \\ 00119 * 0\,...\,0 & 1 00120 * \end{array} \right) \f$ 00121 * 00122 * \b Rotation R (Dim)x(Dim): 00123 * \f$ \left( \begin{array}{cc} 00124 * R & 0\\ 00125 * 0\,...\,0 & 1 00126 * \end{array} \right) \f$ 00127 *<!-- 00128 * \b Linear \b Matrix L (Dim)x(Dim): 00129 * \f$ \left( \begin{array}{cc} 00130 * L & 0\\ 00131 * 0\,...\,0 & 1 00132 * \end{array} \right) \f$ 00133 * 00134 * \b Affine \b Matrix A (Dim)x(Dim+1): 00135 * \f$ \left( \begin{array}{c} 00136 * A\\ 00137 * 0\,...\,0\,1 00138 * \end{array} \right) \f$ 00139 *--> 00140 * \b Scaling \b DiagonalMatrix S (Dim)x(Dim): 00141 * \f$ \left( \begin{array}{cc} 00142 * S & 0\\ 00143 * 0\,...\,0 & 1 00144 * \end{array} \right) \f$ 00145 * 00146 * \b Column \b point v (Dim)x(1): 00147 * \f$ \left( \begin{array}{c} 00148 * v\\ 00149 * 1 00150 * \end{array} \right) \f$ 00151 * 00152 * \b Set \b of \b column \b points V1...Vn (Dim)x(n): 00153 * \f$ \left( \begin{array}{ccc} 00154 * v_1 & ... & v_n\\ 00155 * 1 & ... & 1 00156 * \end{array} \right) \f$ 00157 * 00158 * The concatenation of a Transform object with any kind of other transformation 00159 * always returns a Transform object. 00160 * 00161 * A little exception to the "as pure matrix product" rule is the case of the 00162 * transformation of non homogeneous vectors by an affine transformation. In 00163 * that case the last matrix row can be ignored, and the product returns non 00164 * homogeneous vectors. 00165 * 00166 * Since, for instance, a Dim x Dim matrix is interpreted as a linear transformation, 00167 * it is not possible to directly transform Dim vectors stored in a Dim x Dim matrix. 00168 * The solution is either to use a Dim x Dynamic matrix or explicitly request a 00169 * vector transformation by making the vector homogeneous: 00170 * \code 00171 * m' = T * m.colwise().homogeneous(); 00172 * \endcode 00173 * Note that there is zero overhead. 00174 * 00175 * Conversion methods from/to Qt's QMatrix and QTransform are available if the 00176 * preprocessor token EIGEN_QT_SUPPORT is defined. 00177 * 00178 * This class can be extended with the help of the plugin mechanism described on the page 00179 * \ref TopicCustomizingEigen by defining the preprocessor symbol \c EIGEN_TRANSFORM_PLUGIN. 00180 * 00181 * \sa class Matrix, class Quaternion 00182 */ 00183 template<typename _Scalar, int _Dim, int _Mode, int _Options> 00184 class Transform 00185 { 00186 public: 00187 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_Dim==Dynamic ? Dynamic : (_Dim+1)*(_Dim+1)) 00188 enum { 00189 Mode = _Mode, 00190 Options = _Options, 00191 Dim = _Dim, ///< space dimension in which the transformation holds 00192 HDim = _Dim+1, ///< size of a respective homogeneous vector 00193 Rows = int(Mode)==(AffineCompact) ? Dim : HDim 00194 }; 00195 /** the scalar type of the coefficients */ 00196 typedef _Scalar Scalar; 00197 typedef DenseIndex Index; 00198 /** type of the matrix used to represent the transformation */ 00199 typedef typename internal::make_proper_matrix_type<Scalar,Rows,HDim,Options>::type MatrixType; 00200 /** constified MatrixType */ 00201 typedef const MatrixType ConstMatrixType; 00202 /** type of the matrix used to represent the linear part of the transformation */ 00203 typedef Matrix<Scalar,Dim,Dim,Options> LinearMatrixType; 00204 /** type of read/write reference to the linear part of the transformation */ 00205 typedef Block<MatrixType,Dim,Dim,int(Mode)==(AffineCompact) && (Options&RowMajor)==0> LinearPart; 00206 /** type of read reference to the linear part of the transformation */ 00207 typedef const Block<ConstMatrixType,Dim,Dim,int(Mode)==(AffineCompact) && (Options&RowMajor)==0> ConstLinearPart; 00208 /** type of read/write reference to the affine part of the transformation */ 00209 typedef typename internal::conditional<int(Mode)==int(AffineCompact), 00210 MatrixType&, 00211 Block<MatrixType,Dim,HDim> >::type AffinePart; 00212 /** type of read reference to the affine part of the transformation */ 00213 typedef typename internal::conditional<int(Mode)==int(AffineCompact), 00214 const MatrixType&, 00215 const Block<const MatrixType,Dim,HDim> >::type ConstAffinePart; 00216 /** type of a vector */ 00217 typedef Matrix<Scalar,Dim,1> VectorType; 00218 /** type of a read/write reference to the translation part of the rotation */ 00219 typedef Block<MatrixType,Dim,1,int(Mode)==(AffineCompact)> TranslationPart; 00220 /** type of a read reference to the translation part of the rotation */ 00221 typedef const Block<ConstMatrixType,Dim,1,int(Mode)==(AffineCompact)> ConstTranslationPart; 00222 /** corresponding translation type */ 00223 typedef Translation<Scalar,Dim> TranslationType; 00224 00225 // this intermediate enum is needed to avoid an ICE with gcc 3.4 and 4.0 00226 enum { TransformTimeDiagonalMode = ((Mode==int(Isometry))?Affine:int(Mode)) }; 00227 /** The return type of the product between a diagonal matrix and a transform */ 00228 typedef Transform<Scalar,Dim,TransformTimeDiagonalMode> TransformTimeDiagonalReturnType; 00229 00230 protected: 00231 00232 MatrixType m_matrix; 00233 00234 public: 00235 00236 /** Default constructor without initialization of the meaningful coefficients. 00237 * If Mode==Affine, then the last row is set to [0 ... 0 1] */ 00238 inline Transform () 00239 { 00240 check_template_params(); 00241 internal::transform_make_affine<(int(Mode)==Affine) ? Affine : AffineCompact>::run(m_matrix); 00242 } 00243 00244 inline Transform (const Transform & other) 00245 { 00246 check_template_params(); 00247 m_matrix = other.m_matrix; 00248 } 00249 00250 inline explicit Transform(const TranslationType& t) 00251 { 00252 check_template_params(); 00253 *this = t; 00254 } 00255 inline explicit Transform(const UniformScaling<Scalar>& s) 00256 { 00257 check_template_params(); 00258 *this = s; 00259 } 00260 template<typename Derived> 00261 inline explicit Transform(const RotationBase<Derived, Dim>& r) 00262 { 00263 check_template_params(); 00264 *this = r; 00265 } 00266 00267 inline Transform& operator=(const Transform& other) 00268 { m_matrix = other.m_matrix; return *this; } 00269 00270 typedef internal::transform_take_affine_part<Transform> take_affine_part; 00271 00272 /** Constructs and initializes a transformation from a Dim^2 or a (Dim+1)^2 matrix. */ 00273 template<typename OtherDerived> 00274 inline explicit Transform (const EigenBase<OtherDerived>& other) 00275 { 00276 EIGEN_STATIC_ASSERT((internal::is_same<Scalar,typename OtherDerived::Scalar>::value), 00277 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY); 00278 00279 check_template_params(); 00280 internal::transform_construct_from_matrix<OtherDerived,Mode,Options,Dim,HDim>::run(this, other.derived ()); 00281 } 00282 00283 /** Set \c *this from a Dim^2 or (Dim+1)^2 matrix. */ 00284 template<typename OtherDerived> 00285 inline Transform & operator=(const EigenBase<OtherDerived>& other) 00286 { 00287 EIGEN_STATIC_ASSERT((internal::is_same<Scalar,typename OtherDerived::Scalar>::value), 00288 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY); 00289 00290 internal::transform_construct_from_matrix<OtherDerived,Mode,Options,Dim,HDim>::run(this, other.derived ()); 00291 return *this; 00292 } 00293 00294 template<int OtherOptions> 00295 inline Transform (const Transform<Scalar,Dim,Mode,OtherOptions> & other) 00296 { 00297 check_template_params(); 00298 // only the options change, we can directly copy the matrices 00299 m_matrix = other.matrix (); 00300 } 00301 00302 template<int OtherMode,int OtherOptions> 00303 inline Transform(const Transform<Scalar,Dim,OtherMode,OtherOptions>& other) 00304 { 00305 check_template_params(); 00306 // prevent conversions as: 00307 // Affine | AffineCompact | Isometry = Projective 00308 EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(OtherMode==int(Projective), Mode==int(Projective)), 00309 YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION) 00310 00311 // prevent conversions as: 00312 // Isometry = Affine | AffineCompact 00313 EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(OtherMode==int(Affine)||OtherMode==int(AffineCompact), Mode!=int(Isometry)), 00314 YOU_PERFORMED_AN_INVALID_TRANSFORMATION_CONVERSION) 00315 00316 enum { ModeIsAffineCompact = Mode == int(AffineCompact), 00317 OtherModeIsAffineCompact = OtherMode == int(AffineCompact) 00318 }; 00319 00320 if(ModeIsAffineCompact == OtherModeIsAffineCompact) 00321 { 00322 // We need the block expression because the code is compiled for all 00323 // combinations of transformations and will trigger a compile time error 00324 // if one tries to assign the matrices directly 00325 m_matrix.template block<Dim,Dim+1>(0,0) = other.matrix().template block<Dim,Dim+1>(0,0); 00326 makeAffine(); 00327 } 00328 else if(OtherModeIsAffineCompact) 00329 { 00330 typedef typename Transform<Scalar,Dim,OtherMode,OtherOptions>::MatrixType OtherMatrixType; 00331 internal::transform_construct_from_matrix<OtherMatrixType,Mode,Options,Dim,HDim>::run(this, other.matrix()); 00332 } 00333 else 00334 { 00335 // here we know that Mode == AffineCompact and OtherMode != AffineCompact. 00336 // if OtherMode were Projective, the static assert above would already have caught it. 00337 // So the only possibility is that OtherMode == Affine 00338 linear() = other.linear(); 00339 translation() = other.translation(); 00340 } 00341 } 00342 00343 template<typename OtherDerived> 00344 Transform(const ReturnByValue<OtherDerived>& other) 00345 { 00346 check_template_params(); 00347 other.evalTo(*this); 00348 } 00349 00350 template<typename OtherDerived> 00351 Transform& operator=(const ReturnByValue<OtherDerived>& other) 00352 { 00353 other.evalTo(*this); 00354 return *this; 00355 } 00356 00357 #ifdef EIGEN_QT_SUPPORT 00358 inline Transform(const QMatrix& other); 00359 inline Transform& operator=(const QMatrix& other); 00360 inline QMatrix toQMatrix(void) const; 00361 inline Transform(const QTransform& other); 00362 inline Transform& operator=(const QTransform& other); 00363 inline QTransform toQTransform(void) const; 00364 #endif 00365 00366 /** shortcut for m_matrix(row,col); 00367 * \sa MatrixBase::operator(Index,Index) const */ 00368 inline Scalar operator() (Index row, Index col) const { return m_matrix(row,col); } 00369 /** shortcut for m_matrix(row,col); 00370 * \sa MatrixBase::operator(Index,Index) */ 00371 inline Scalar& operator() (Index row, Index col) { return m_matrix(row,col); } 00372 00373 /** \returns a read-only expression of the transformation matrix */ 00374 inline const MatrixType& matrix () const { return m_matrix; } 00375 /** \returns a writable expression of the transformation matrix */ 00376 inline MatrixType& matrix () { return m_matrix; } 00377 00378 /** \returns a read-only expression of the linear part of the transformation */ 00379 inline ConstLinearPart linear () const { return ConstLinearPart(m_matrix,0,0); } 00380 /** \returns a writable expression of the linear part of the transformation */ 00381 inline LinearPart linear () { return LinearPart(m_matrix,0,0); } 00382 00383 /** \returns a read-only expression of the Dim x HDim affine part of the transformation */ 00384 inline ConstAffinePart affine () const { return take_affine_part::run(m_matrix); } 00385 /** \returns a writable expression of the Dim x HDim affine part of the transformation */ 00386 inline AffinePart affine () { return take_affine_part::run(m_matrix); } 00387 00388 /** \returns a read-only expression of the translation vector of the transformation */ 00389 inline ConstTranslationPart translation () const { return ConstTranslationPart(m_matrix,0,Dim); } 00390 /** \returns a writable expression of the translation vector of the transformation */ 00391 inline TranslationPart translation () { return TranslationPart(m_matrix,0,Dim); } 00392 00393 /** \returns an expression of the product between the transform \c *this and a matrix expression \a other. 00394 * 00395 * The right-hand-side \a other can be either: 00396 * \li an homogeneous vector of size Dim+1, 00397 * \li a set of homogeneous vectors of size Dim+1 x N, 00398 * \li a transformation matrix of size Dim+1 x Dim+1. 00399 * 00400 * Moreover, if \c *this represents an affine transformation (i.e., Mode!=Projective), then \a other can also be: 00401 * \li a point of size Dim (computes: \code this->linear() * other + this->translation()\endcode), 00402 * \li a set of N points as a Dim x N matrix (computes: \code (this->linear() * other).colwise() + this->translation()\endcode), 00403 * 00404 * In all cases, the return type is a matrix or vector of same sizes as the right-hand-side \a other. 00405 * 00406 * If you want to interpret \a other as a linear or affine transformation, then first convert it to a Transform<> type, 00407 * or do your own cooking. 00408 * 00409 * Finally, if you want to apply Affine transformations to vectors, then explicitly apply the linear part only: 00410 * \code 00411 * Affine3f A; 00412 * Vector3f v1, v2; 00413 * v2 = A.linear() * v1; 00414 * \endcode 00415 * 00416 */ 00417 // note: this function is defined here because some compilers cannot find the respective declaration 00418 template<typename OtherDerived> 00419 EIGEN_STRONG_INLINE const typename OtherDerived::PlainObject 00420 operator * (const EigenBase<OtherDerived> &other) const 00421 { return internal::transform_right_product_impl<Transform, OtherDerived>::run(*this,other.derived ()); } 00422 00423 /** \returns the product expression of a transformation matrix \a a times a transform \a b 00424 * 00425 * The left hand side \a other can be either: 00426 * \li a linear transformation matrix of size Dim x Dim, 00427 * \li an affine transformation matrix of size Dim x Dim+1, 00428 * \li a general transformation matrix of size Dim+1 x Dim+1. 00429 */ 00430 template<typename OtherDerived> friend 00431 inline const typename internal::transform_left_product_impl<OtherDerived,Mode,Options,_Dim,_Dim+1>::ResultType 00432 operator * (const EigenBase<OtherDerived> &a, const Transform &b) 00433 { return internal::transform_left_product_impl<OtherDerived,Mode,Options,Dim,HDim>::run(a.derived (),b); } 00434 00435 /** \returns The product expression of a transform \a a times a diagonal matrix \a b 00436 * 00437 * The rhs diagonal matrix is interpreted as an affine scaling transformation. The 00438 * product results in a Transform of the same type (mode) as the lhs only if the lhs 00439 * mode is no isometry. In that case, the returned transform is an affinity. 00440 */ 00441 template<typename DiagonalDerived> 00442 inline const TransformTimeDiagonalReturnType 00443 operator * (const DiagonalBase<DiagonalDerived> &b) const 00444 { 00445 TransformTimeDiagonalReturnType res(*this); 00446 res.linear () *= b; 00447 return res; 00448 } 00449 00450 /** \returns The product expression of a diagonal matrix \a a times a transform \a b 00451 * 00452 * The lhs diagonal matrix is interpreted as an affine scaling transformation. The 00453 * product results in a Transform of the same type (mode) as the lhs only if the lhs 00454 * mode is no isometry. In that case, the returned transform is an affinity. 00455 */ 00456 template<typename DiagonalDerived> 00457 friend inline TransformTimeDiagonalReturnType 00458 operator * (const DiagonalBase<DiagonalDerived> &a, const Transform &b) 00459 { 00460 TransformTimeDiagonalReturnType res; 00461 res.linear ().noalias() = a*b.linear (); 00462 res.translation ().noalias() = a*b.translation (); 00463 if (Mode!=int(AffineCompact)) 00464 res.matrix ().row(Dim) = b.matrix ().row(Dim); 00465 return res; 00466 } 00467 00468 template<typename OtherDerived> 00469 inline Transform & operator*=(const EigenBase<OtherDerived>& other) { return *this = *this * other; } 00470 00471 /** Concatenates two transformations */ 00472 inline const Transform operator * (const Transform & other) const 00473 { 00474 return internal::transform_transform_product_impl<Transform,Transform>::run(*this,other); 00475 } 00476 00477 #ifdef __INTEL_COMPILER 00478 private: 00479 // this intermediate structure permits to workaround a bug in ICC 11: 00480 // error: template instantiation resulted in unexpected function type of "Eigen::Transform<double, 3, 32, 0> 00481 // (const Eigen::Transform<double, 3, 2, 0> &) const" 00482 // (the meaning of a name may have changed since the template declaration -- the type of the template is: 00483 // "Eigen::internal::transform_transform_product_impl<Eigen::Transform<double, 3, 32, 0>, 00484 // Eigen::Transform<double, 3, Mode, Options>, <expression>>::ResultType (const Eigen::Transform<double, 3, Mode, Options> &) const") 00485 // 00486 template<int OtherMode,int OtherOptions> struct icc_11_workaround 00487 { 00488 typedef internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> > ProductType; 00489 typedef typename ProductType::ResultType ResultType; 00490 }; 00491 00492 public: 00493 /** Concatenates two different transformations */ 00494 template<int OtherMode,int OtherOptions> 00495 inline typename icc_11_workaround<OtherMode,OtherOptions>::ResultType 00496 operator * (const Transform<Scalar,Dim,OtherMode,OtherOptions> & other) const 00497 { 00498 typedef typename icc_11_workaround<OtherMode,OtherOptions>::ProductType ProductType; 00499 return ProductType::run(*this,other); 00500 } 00501 #else 00502 /** Concatenates two different transformations */ 00503 template<int OtherMode,int OtherOptions> 00504 inline typename internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> >::ResultType 00505 operator * (const Transform<Scalar,Dim,OtherMode,OtherOptions> & other) const 00506 { 00507 return internal::transform_transform_product_impl<Transform,Transform<Scalar,Dim,OtherMode,OtherOptions> >::run(*this,other); 00508 } 00509 #endif 00510 00511 /** \sa MatrixBase::setIdentity() */ 00512 void setIdentity () { m_matrix.setIdentity(); } 00513 00514 /** 00515 * \brief Returns an identity transformation. 00516 * \todo In the future this function should be returning a Transform expression. 00517 */ 00518 static const Transform Identity() 00519 { 00520 return Transform (MatrixType::Identity()); 00521 } 00522 00523 template<typename OtherDerived> 00524 inline Transform & scale(const MatrixBase<OtherDerived> &other); 00525 00526 template<typename OtherDerived> 00527 inline Transform & prescale(const MatrixBase<OtherDerived> &other); 00528 00529 inline Transform & scale(const Scalar& s); 00530 inline Transform & prescale(const Scalar& s); 00531 00532 template<typename OtherDerived> 00533 inline Transform & translate(const MatrixBase<OtherDerived> &other); 00534 00535 template<typename OtherDerived> 00536 inline Transform & pretranslate(const MatrixBase<OtherDerived> &other); 00537 00538 template<typename RotationType> 00539 inline Transform & rotate(const RotationType& rotation); 00540 00541 template<typename RotationType> 00542 inline Transform & prerotate(const RotationType& rotation); 00543 00544 Transform & shear(const Scalar& sx, const Scalar& sy); 00545 Transform & preshear(const Scalar& sx, const Scalar& sy); 00546 00547 inline Transform & operator=(const TranslationType& t); 00548 inline Transform & operator*=(const TranslationType& t) { return translate(t.vector()); } 00549 inline Transform operator* (const TranslationType& t) const; 00550 00551 inline Transform& operator=(const UniformScaling<Scalar>& t); 00552 inline Transform& operator*=(const UniformScaling<Scalar>& s) { return scale(s.factor()); } 00553 inline Transform<Scalar,Dim,(int(Mode)==int(Isometry)?int(Affine):int(Mode))> operator* (const UniformScaling<Scalar>& s) const 00554 { 00555 Transform<Scalar,Dim,(int(Mode)==int(Isometry)?int(Affine):int(Mode)),Options> res = *this; 00556 res.scale(s.factor()); 00557 return res; 00558 } 00559 00560 inline Transform& operator*=(const DiagonalMatrix<Scalar,Dim>& s) { linear() *= s; return *this; } 00561 00562 template<typename Derived> 00563 inline Transform& operator=(const RotationBase<Derived,Dim>& r); 00564 template<typename Derived> 00565 inline Transform& operator*=(const RotationBase<Derived,Dim>& r) { return rotate(r.toRotationMatrix()); } 00566 template<typename Derived> 00567 inline Transform operator* (const RotationBase<Derived,Dim>& r) const; 00568 00569 const LinearMatrixType rotation() const; 00570 template<typename RotationMatrixType, typename ScalingMatrixType> 00571 void computeRotationScaling(RotationMatrixType *rotation, ScalingMatrixType *scaling) const; 00572 template<typename ScalingMatrixType, typename RotationMatrixType> 00573 void computeScalingRotation(ScalingMatrixType *scaling, RotationMatrixType *rotation) const; 00574 00575 template<typename PositionDerived, typename OrientationType, typename ScaleDerived> 00576 Transform& fromPositionOrientationScale(const MatrixBase<PositionDerived> &position, 00577 const OrientationType& orientation, const MatrixBase<ScaleDerived> &scale); 00578 00579 inline Transform inverse(TransformTraits traits = (TransformTraits)Mode) const; 00580 00581 /** \returns a const pointer to the column major internal matrix */ 00582 const Scalar* data () const { return m_matrix.data (); } 00583 /** \returns a non-const pointer to the column major internal matrix */ 00584 Scalar* data () { return m_matrix.data (); } 00585 00586 /** \returns \c *this with scalar type casted to \a NewScalarType 00587 * 00588 * Note that if \a NewScalarType is equal to the current scalar type of \c *this 00589 * then this function smartly returns a const reference to \c *this. 00590 */ 00591 template<typename NewScalarType> 00592 inline typename internal::cast_return_type<Transform,Transform<NewScalarType,Dim,Mode,Options> >::type cast() const 00593 { return typename internal::cast_return_type<Transform,Transform<NewScalarType,Dim,Mode,Options> >::type(*this); } 00594 00595 /** Copy constructor with scalar type conversion */ 00596 template<typename OtherScalarType> 00597 inline explicit Transform (const Transform<OtherScalarType,Dim,Mode,Options> & other) 00598 { 00599 check_template_params(); 00600 m_matrix = other.matrix ().template cast<Scalar>(); 00601 } 00602 00603 /** \returns \c true if \c *this is approximately equal to \a other, within the precision 00604 * determined by \a prec. 00605 * 00606 * \sa MatrixBase::isApprox() */ 00607 bool isApprox(const Transform & other, const typename NumTraits<Scalar>::Real& prec = NumTraits<Scalar>::dummy_precision()) const 00608 { return m_matrix.isApprox(other.m_matrix, prec); } 00609 00610 /** Sets the last row to [0 ... 0 1] 00611 */ 00612 void makeAffine() 00613 { 00614 internal::transform_make_affine<int(Mode)>::run(m_matrix); 00615 } 00616 00617 /** \internal 00618 * \returns the Dim x Dim linear part if the transformation is affine, 00619 * and the HDim x Dim part for projective transformations. 00620 */ 00621 inline Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,Dim> linearExt() 00622 { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,Dim>(0,0); } 00623 /** \internal 00624 * \returns the Dim x Dim linear part if the transformation is affine, 00625 * and the HDim x Dim part for projective transformations. 00626 */ 00627 inline const Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,Dim> linearExt() const 00628 { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,Dim>(0,0); } 00629 00630 /** \internal 00631 * \returns the translation part if the transformation is affine, 00632 * and the last column for projective transformations. 00633 */ 00634 inline Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,1> translationExt() 00635 { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,1>(0,Dim); } 00636 /** \internal 00637 * \returns the translation part if the transformation is affine, 00638 * and the last column for projective transformations. 00639 */ 00640 inline const Block<MatrixType,int(Mode)==int(Projective)?HDim:Dim,1> translationExt() const 00641 { return m_matrix.template block<int(Mode)==int(Projective)?HDim:Dim,1>(0,Dim); } 00642 00643 00644 #ifdef EIGEN_TRANSFORM_PLUGIN 00645 #include EIGEN_TRANSFORM_PLUGIN 00646 #endif 00647 00648 protected: 00649 #ifndef EIGEN_PARSED_BY_DOXYGEN 00650 static EIGEN_STRONG_INLINE void check_template_params() 00651 { 00652 EIGEN_STATIC_ASSERT((Options & (DontAlign|RowMajor)) == Options, INVALID_MATRIX_TEMPLATE_PARAMETERS) 00653 } 00654 #endif 00655 00656 }; 00657 00658 /** \ingroup Geometry_Module */ 00659 typedef Transform<float,2,Isometry> Isometry2f; 00660 /** \ingroup Geometry_Module */ 00661 typedef Transform<float,3,Isometry> Isometry3f; 00662 /** \ingroup Geometry_Module */ 00663 typedef Transform<double,2,Isometry> Isometry2d; 00664 /** \ingroup Geometry_Module */ 00665 typedef Transform<double,3,Isometry> Isometry3d; 00666 00667 /** \ingroup Geometry_Module */ 00668 typedef Transform<float,2,Affine> Affine2f; 00669 /** \ingroup Geometry_Module */ 00670 typedef Transform<float,3,Affine> Affine3f; 00671 /** \ingroup Geometry_Module */ 00672 typedef Transform<double,2,Affine> Affine2d; 00673 /** \ingroup Geometry_Module */ 00674 typedef Transform<double,3,Affine> Affine3d; 00675 00676 /** \ingroup Geometry_Module */ 00677 typedef Transform<float,2,AffineCompact> AffineCompact2f; 00678 /** \ingroup Geometry_Module */ 00679 typedef Transform<float,3,AffineCompact> AffineCompact3f; 00680 /** \ingroup Geometry_Module */ 00681 typedef Transform<double,2,AffineCompact> AffineCompact2d; 00682 /** \ingroup Geometry_Module */ 00683 typedef Transform<double,3,AffineCompact> AffineCompact3d; 00684 00685 /** \ingroup Geometry_Module */ 00686 typedef Transform<float,2,Projective> Projective2f; 00687 /** \ingroup Geometry_Module */ 00688 typedef Transform<float,3,Projective> Projective3f; 00689 /** \ingroup Geometry_Module */ 00690 typedef Transform<double,2,Projective> Projective2d; 00691 /** \ingroup Geometry_Module */ 00692 typedef Transform<double,3,Projective> Projective3d; 00693 00694 /************************** 00695 *** Optional QT support *** 00696 **************************/ 00697 00698 #ifdef EIGEN_QT_SUPPORT 00699 /** Initializes \c *this from a QMatrix assuming the dimension is 2. 00700 * 00701 * This function is available only if the token EIGEN_QT_SUPPORT is defined. 00702 */ 00703 template<typename Scalar, int Dim, int Mode,int Options> 00704 Transform<Scalar,Dim,Mode,Options>::Transform (const QMatrix& other) 00705 { 00706 check_template_params(); 00707 *this = other; 00708 } 00709 00710 /** Set \c *this from a QMatrix assuming the dimension is 2. 00711 * 00712 * This function is available only if the token EIGEN_QT_SUPPORT is defined. 00713 */ 00714 template<typename Scalar, int Dim, int Mode,int Options> 00715 Transform<Scalar,Dim,Mode,Options> & Transform<Scalar,Dim,Mode,Options>::operator= (const QMatrix& other) 00716 { 00717 EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE) 00718 m_matrix << other.m11(), other.m21(), other.dx(), 00719 other.m12(), other.m22(), other.dy(), 00720 0, 0, 1; 00721 return *this; 00722 } 00723 00724 /** \returns a QMatrix from \c *this assuming the dimension is 2. 00725 * 00726 * \warning this conversion might loss data if \c *this is not affine 00727 * 00728 * This function is available only if the token EIGEN_QT_SUPPORT is defined. 00729 */ 00730 template<typename Scalar, int Dim, int Mode, int Options> 00731 QMatrix Transform<Scalar,Dim,Mode,Options>::toQMatrix (void) const 00732 { 00733 check_template_params(); 00734 EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE) 00735 return QMatrix(m_matrix.coeff(0,0), m_matrix.coeff(1,0), 00736 m_matrix.coeff(0,1), m_matrix.coeff(1,1), 00737 m_matrix.coeff(0,2), m_matrix.coeff(1,2)); 00738 } 00739 00740 /** Initializes \c *this from a QTransform assuming the dimension is 2. 00741 * 00742 * This function is available only if the token EIGEN_QT_SUPPORT is defined. 00743 */ 00744 template<typename Scalar, int Dim, int Mode,int Options> 00745 Transform<Scalar,Dim,Mode,Options>::Transform (const QTransform& other) 00746 { 00747 check_template_params(); 00748 *this = other; 00749 } 00750 00751 /** Set \c *this from a QTransform assuming the dimension is 2. 00752 * 00753 * This function is available only if the token EIGEN_QT_SUPPORT is defined. 00754 */ 00755 template<typename Scalar, int Dim, int Mode, int Options> 00756 Transform<Scalar,Dim,Mode,Options> & Transform<Scalar,Dim,Mode,Options>::operator= (const QTransform& other) 00757 { 00758 check_template_params(); 00759 EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE) 00760 if (Mode == int(AffineCompact)) 00761 m_matrix << other.m11(), other.m21(), other.dx(), 00762 other.m12(), other.m22(), other.dy(); 00763 else 00764 m_matrix << other.m11(), other.m21(), other.dx(), 00765 other.m12(), other.m22(), other.dy(), 00766 other.m13(), other.m23(), other.m33(); 00767 return *this; 00768 } 00769 00770 /** \returns a QTransform from \c *this assuming the dimension is 2. 00771 * 00772 * This function is available only if the token EIGEN_QT_SUPPORT is defined. 00773 */ 00774 template<typename Scalar, int Dim, int Mode, int Options> 00775 QTransform Transform<Scalar,Dim,Mode,Options>::toQTransform (void) const 00776 { 00777 EIGEN_STATIC_ASSERT(Dim==2, YOU_MADE_A_PROGRAMMING_MISTAKE) 00778 if (Mode == int(AffineCompact)) 00779 return QTransform(m_matrix.coeff(0,0), m_matrix.coeff(1,0), 00780 m_matrix.coeff(0,1), m_matrix.coeff(1,1), 00781 m_matrix.coeff(0,2), m_matrix.coeff(1,2)); 00782 else 00783 return QTransform(m_matrix.coeff(0,0), m_matrix.coeff(1,0), m_matrix.coeff(2,0), 00784 m_matrix.coeff(0,1), m_matrix.coeff(1,1), m_matrix.coeff(2,1), 00785 m_matrix.coeff(0,2), m_matrix.coeff(1,2), m_matrix.coeff(2,2)); 00786 } 00787 #endif 00788 00789 /********************* 00790 *** Procedural API *** 00791 *********************/ 00792 00793 /** Applies on the right the non uniform scale transformation represented 00794 * by the vector \a other to \c *this and returns a reference to \c *this. 00795 * \sa prescale() 00796 */ 00797 template<typename Scalar, int Dim, int Mode, int Options> 00798 template<typename OtherDerived> 00799 Transform<Scalar,Dim,Mode,Options> & 00800 Transform<Scalar,Dim,Mode,Options>::scale (const MatrixBase<OtherDerived> &other) 00801 { 00802 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim)) 00803 EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS) 00804 linearExt().noalias() = (linearExt() * other.asDiagonal ()); 00805 return *this; 00806 } 00807 00808 /** Applies on the right a uniform scale of a factor \a c to \c *this 00809 * and returns a reference to \c *this. 00810 * \sa prescale(Scalar) 00811 */ 00812 template<typename Scalar, int Dim, int Mode, int Options> 00813 inline Transform<Scalar,Dim,Mode,Options> & Transform<Scalar,Dim,Mode,Options>::scale (const Scalar& s) 00814 { 00815 EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS) 00816 linearExt() *= s; 00817 return *this; 00818 } 00819 00820 /** Applies on the left the non uniform scale transformation represented 00821 * by the vector \a other to \c *this and returns a reference to \c *this. 00822 * \sa scale() 00823 */ 00824 template<typename Scalar, int Dim, int Mode, int Options> 00825 template<typename OtherDerived> 00826 Transform<Scalar,Dim,Mode,Options> & 00827 Transform<Scalar,Dim,Mode,Options>::prescale (const MatrixBase<OtherDerived> &other) 00828 { 00829 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim)) 00830 EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS) 00831 m_matrix.template block<Dim,HDim>(0,0).noalias() = (other.asDiagonal () * m_matrix.template block<Dim,HDim>(0,0)); 00832 return *this; 00833 } 00834 00835 /** Applies on the left a uniform scale of a factor \a c to \c *this 00836 * and returns a reference to \c *this. 00837 * \sa scale(Scalar) 00838 */ 00839 template<typename Scalar, int Dim, int Mode, int Options> 00840 inline Transform<Scalar,Dim,Mode,Options> & Transform<Scalar,Dim,Mode,Options>::prescale (const Scalar& s) 00841 { 00842 EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS) 00843 m_matrix.template topRows<Dim>() *= s; 00844 return *this; 00845 } 00846 00847 /** Applies on the right the translation matrix represented by the vector \a other 00848 * to \c *this and returns a reference to \c *this. 00849 * \sa pretranslate() 00850 */ 00851 template<typename Scalar, int Dim, int Mode, int Options> 00852 template<typename OtherDerived> 00853 Transform<Scalar,Dim,Mode,Options> & 00854 Transform<Scalar,Dim,Mode,Options>::translate (const MatrixBase<OtherDerived> &other) 00855 { 00856 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim)) 00857 translationExt() += linearExt() * other; 00858 return *this; 00859 } 00860 00861 /** Applies on the left the translation matrix represented by the vector \a other 00862 * to \c *this and returns a reference to \c *this. 00863 * \sa translate() 00864 */ 00865 template<typename Scalar, int Dim, int Mode, int Options> 00866 template<typename OtherDerived> 00867 Transform<Scalar,Dim,Mode,Options> & 00868 Transform<Scalar,Dim,Mode,Options>::pretranslate (const MatrixBase<OtherDerived> &other) 00869 { 00870 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,int(Dim)) 00871 if(int(Mode)==int(Projective)) 00872 affine() += other * m_matrix.row(Dim); 00873 else 00874 translation() += other; 00875 return *this; 00876 } 00877 00878 /** Applies on the right the rotation represented by the rotation \a rotation 00879 * to \c *this and returns a reference to \c *this. 00880 * 00881 * The template parameter \a RotationType is the type of the rotation which 00882 * must be known by internal::toRotationMatrix<>. 00883 * 00884 * Natively supported types includes: 00885 * - any scalar (2D), 00886 * - a Dim x Dim matrix expression, 00887 * - a Quaternion (3D), 00888 * - a AngleAxis (3D) 00889 * 00890 * This mechanism is easily extendable to support user types such as Euler angles, 00891 * or a pair of Quaternion for 4D rotations. 00892 * 00893 * \sa rotate(Scalar), class Quaternion, class AngleAxis, prerotate(RotationType) 00894 */ 00895 template<typename Scalar, int Dim, int Mode, int Options> 00896 template<typename RotationType> 00897 Transform<Scalar,Dim,Mode,Options> & 00898 Transform<Scalar,Dim,Mode,Options>::rotate (const RotationType& rotation) 00899 { 00900 linearExt() *= internal::toRotationMatrix<Scalar,Dim>(rotation); 00901 return *this; 00902 } 00903 00904 /** Applies on the left the rotation represented by the rotation \a rotation 00905 * to \c *this and returns a reference to \c *this. 00906 * 00907 * See rotate() for further details. 00908 * 00909 * \sa rotate() 00910 */ 00911 template<typename Scalar, int Dim, int Mode, int Options> 00912 template<typename RotationType> 00913 Transform<Scalar,Dim,Mode,Options> & 00914 Transform<Scalar,Dim,Mode,Options>::prerotate (const RotationType& rotation) 00915 { 00916 m_matrix.template block<Dim,HDim>(0,0) = internal::toRotationMatrix<Scalar,Dim>(rotation) 00917 * m_matrix.template block<Dim,HDim>(0,0); 00918 return *this; 00919 } 00920 00921 /** Applies on the right the shear transformation represented 00922 * by the vector \a other to \c *this and returns a reference to \c *this. 00923 * \warning 2D only. 00924 * \sa preshear() 00925 */ 00926 template<typename Scalar, int Dim, int Mode, int Options> 00927 Transform<Scalar,Dim,Mode,Options> & 00928 Transform<Scalar,Dim,Mode,Options>::shear (const Scalar& sx, const Scalar& sy) 00929 { 00930 EIGEN_STATIC_ASSERT(int(Dim)==2, YOU_MADE_A_PROGRAMMING_MISTAKE) 00931 EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS) 00932 VectorType tmp = linear().col(0)*sy + linear().col(1); 00933 linear() << linear().col(0) + linear().col(1)*sx, tmp; 00934 return *this; 00935 } 00936 00937 /** Applies on the left the shear transformation represented 00938 * by the vector \a other to \c *this and returns a reference to \c *this. 00939 * \warning 2D only. 00940 * \sa shear() 00941 */ 00942 template<typename Scalar, int Dim, int Mode, int Options> 00943 Transform<Scalar,Dim,Mode,Options> & 00944 Transform<Scalar,Dim,Mode,Options>::preshear (const Scalar& sx, const Scalar& sy) 00945 { 00946 EIGEN_STATIC_ASSERT(int(Dim)==2, YOU_MADE_A_PROGRAMMING_MISTAKE) 00947 EIGEN_STATIC_ASSERT(Mode!=int(Isometry), THIS_METHOD_IS_ONLY_FOR_SPECIFIC_TRANSFORMATIONS) 00948 m_matrix.template block<Dim,HDim>(0,0) = LinearMatrixType(1, sx, sy, 1) * m_matrix.template block<Dim,HDim>(0,0); 00949 return *this; 00950 } 00951 00952 /****************************************************** 00953 *** Scaling, Translation and Rotation compatibility *** 00954 ******************************************************/ 00955 00956 template<typename Scalar, int Dim, int Mode, int Options> 00957 inline Transform<Scalar,Dim,Mode,Options> & Transform<Scalar,Dim,Mode,Options>::operator= (const TranslationType& t) 00958 { 00959 linear().setIdentity (); 00960 translation() = t.vector(); 00961 makeAffine(); 00962 return *this; 00963 } 00964 00965 template<typename Scalar, int Dim, int Mode, int Options> 00966 inline Transform<Scalar,Dim,Mode,Options> Transform<Scalar,Dim,Mode,Options>::operator* (const TranslationType& t) const 00967 { 00968 Transform res = *this; 00969 res.translate(t.vector()); 00970 return res; 00971 } 00972 00973 template<typename Scalar, int Dim, int Mode, int Options> 00974 inline Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::operator=(const UniformScaling<Scalar>& s) 00975 { 00976 m_matrix.setZero(); 00977 linear().diagonal().fill(s.factor()); 00978 makeAffine(); 00979 return *this; 00980 } 00981 00982 template<typename Scalar, int Dim, int Mode, int Options> 00983 template<typename Derived> 00984 inline Transform<Scalar,Dim,Mode,Options>& Transform<Scalar,Dim,Mode,Options>::operator=(const RotationBase<Derived,Dim>& r) 00985 { 00986 linear() = internal::toRotationMatrix<Scalar,Dim>(r); 00987 translation().setZero(); 00988 makeAffine(); 00989 return *this; 00990 } 00991 00992 template<typename Scalar, int Dim, int Mode, int Options> 00993 template<typename Derived> 00994 inline Transform<Scalar,Dim,Mode,Options> Transform<Scalar,Dim,Mode,Options>::operator* (const RotationBase<Derived,Dim>& r) const 00995 { 00996 Transform res = *this; 00997 res.rotate(r.derived()); 00998 return res; 00999 } 01000 01001 /************************ 01002 *** Special functions *** 01003 ************************/ 01004 01005 /** \returns the rotation part of the transformation 01006 * 01007 * 01008 * \svd_module 01009 * 01010 * \sa computeRotationScaling(), computeScalingRotation(), class SVD 01011 */ 01012 template<typename Scalar, int Dim, int Mode, int Options> 01013 const typename Transform<Scalar,Dim,Mode,Options>::LinearMatrixType 01014 Transform<Scalar,Dim,Mode,Options>::rotation () const 01015 { 01016 LinearMatrixType result; 01017 computeRotationScaling(&result, (LinearMatrixType*)0); 01018 return result; 01019 } 01020 01021 01022 /** decomposes the linear part of the transformation as a product rotation x scaling, the scaling being 01023 * not necessarily positive. 01024 * 01025 * If either pointer is zero, the corresponding computation is skipped. 01026 * 01027 * 01028 * 01029 * \svd_module 01030 * 01031 * \sa computeScalingRotation(), rotation(), class SVD 01032 */ 01033 template<typename Scalar, int Dim, int Mode, int Options> 01034 template<typename RotationMatrixType, typename ScalingMatrixType> 01035 void Transform<Scalar,Dim,Mode,Options>::computeRotationScaling (RotationMatrixType *rotation, ScalingMatrixType *scaling) const 01036 { 01037 JacobiSVD<LinearMatrixType> svd(linear(), ComputeFullU | ComputeFullV); 01038 01039 Scalar x = (svd.matrixU () * svd.matrixV ().adjoint()).determinant(); // so x has absolute value 1 01040 VectorType sv(svd.singularValues ()); 01041 sv.coeffRef(0) *= x; 01042 if(scaling) scaling->lazyAssign(svd.matrixV () * sv.asDiagonal() * svd.matrixV ().adjoint()); 01043 if(rotation) 01044 { 01045 LinearMatrixType m(svd.matrixU ()); 01046 m.col(0) /= x; 01047 rotation->lazyAssign(m * svd.matrixV ().adjoint()); 01048 } 01049 } 01050 01051 /** decomposes the linear part of the transformation as a product rotation x scaling, the scaling being 01052 * not necessarily positive. 01053 * 01054 * If either pointer is zero, the corresponding computation is skipped. 01055 * 01056 * 01057 * 01058 * \svd_module 01059 * 01060 * \sa computeRotationScaling(), rotation(), class SVD 01061 */ 01062 template<typename Scalar, int Dim, int Mode, int Options> 01063 template<typename ScalingMatrixType, typename RotationMatrixType> 01064 void Transform<Scalar,Dim,Mode,Options>::computeScalingRotation (ScalingMatrixType *scaling, RotationMatrixType *rotation) const 01065 { 01066 JacobiSVD<LinearMatrixType> svd(linear(), ComputeFullU | ComputeFullV); 01067 01068 Scalar x = (svd.matrixU () * svd.matrixV ().adjoint()).determinant(); // so x has absolute value 1 01069 VectorType sv(svd.singularValues ()); 01070 sv.coeffRef(0) *= x; 01071 if(scaling) scaling->lazyAssign(svd.matrixU () * sv.asDiagonal() * svd.matrixU ().adjoint()); 01072 if(rotation) 01073 { 01074 LinearMatrixType m(svd.matrixU ()); 01075 m.col(0) /= x; 01076 rotation->lazyAssign(m * svd.matrixV ().adjoint()); 01077 } 01078 } 01079 01080 /** Convenient method to set \c *this from a position, orientation and scale 01081 * of a 3D object. 01082 */ 01083 template<typename Scalar, int Dim, int Mode, int Options> 01084 template<typename PositionDerived, typename OrientationType, typename ScaleDerived> 01085 Transform<Scalar,Dim,Mode,Options> & 01086 Transform<Scalar,Dim,Mode,Options>::fromPositionOrientationScale (const MatrixBase<PositionDerived> &position, 01087 const OrientationType& orientation, const MatrixBase<ScaleDerived> &scale) 01088 { 01089 linear() = internal::toRotationMatrix<Scalar,Dim>(orientation); 01090 linear() *= scale.asDiagonal (); 01091 translation() = position; 01092 makeAffine(); 01093 return *this; 01094 } 01095 01096 namespace internal { 01097 01098 template<int Mode> 01099 struct transform_make_affine 01100 { 01101 template<typename MatrixType> 01102 static void run(MatrixType &mat) 01103 { 01104 static const int Dim = MatrixType::ColsAtCompileTime-1; 01105 mat.template block<1,Dim>(Dim,0).setZero(); 01106 mat.coeffRef(Dim,Dim) = typename MatrixType::Scalar(1); 01107 } 01108 }; 01109 01110 template<> 01111 struct transform_make_affine<AffineCompact> 01112 { 01113 template<typename MatrixType> static void run(MatrixType &) { } 01114 }; 01115 01116 // selector needed to avoid taking the inverse of a 3x4 matrix 01117 template<typename TransformType, int Mode=TransformType::Mode> 01118 struct projective_transform_inverse 01119 { 01120 static inline void run(const TransformType&, TransformType&) 01121 {} 01122 }; 01123 01124 template<typename TransformType> 01125 struct projective_transform_inverse<TransformType, Projective> 01126 { 01127 static inline void run(const TransformType& m, TransformType& res) 01128 { 01129 res.matrix() = m.matrix().inverse(); 01130 } 01131 }; 01132 01133 } // end namespace internal 01134 01135 01136 /** 01137 * 01138 * \returns the inverse transformation according to some given knowledge 01139 * on \c *this. 01140 * 01141 * \param hint allows to optimize the inversion process when the transformation 01142 * is known to be not a general transformation (optional). The possible values are: 01143 * - #Projective if the transformation is not necessarily affine, i.e., if the 01144 * last row is not guaranteed to be [0 ... 0 1] 01145 * - #Affine if the last row can be assumed to be [0 ... 0 1] 01146 * - #Isometry if the transformation is only a concatenations of translations 01147 * and rotations. 01148 * The default is the template class parameter \c Mode. 01149 * 01150 * \warning unless \a traits is always set to NoShear or NoScaling, this function 01151 * requires the generic inverse method of MatrixBase defined in the LU module. If 01152 * you forget to include this module, then you will get hard to debug linking errors. 01153 * 01154 * \sa MatrixBase::inverse() 01155 */ 01156 template<typename Scalar, int Dim, int Mode, int Options> 01157 Transform<Scalar,Dim,Mode,Options> 01158 Transform<Scalar,Dim,Mode,Options>::inverse (TransformTraits hint) const 01159 { 01160 Transform res; 01161 if (hint == Projective) 01162 { 01163 internal::projective_transform_inverse<Transform>::run(*this, res); 01164 } 01165 else 01166 { 01167 if (hint == Isometry) 01168 { 01169 res.matrix ().template topLeftCorner<Dim,Dim>() = linear().transpose(); 01170 } 01171 else if(hint&Affine) 01172 { 01173 res.matrix ().template topLeftCorner<Dim,Dim>() = linear().inverse(); 01174 } 01175 else 01176 { 01177 eigen_assert(false && "Invalid transform traits in Transform::Inverse"); 01178 } 01179 // translation and remaining parts 01180 res.matrix ().template topRightCorner<Dim,1>() 01181 = - res.matrix ().template topLeftCorner<Dim,Dim>() * translation(); 01182 res.makeAffine(); // we do need this, because in the beginning res is uninitialized 01183 } 01184 return res; 01185 } 01186 01187 namespace internal { 01188 01189 /***************************************************** 01190 *** Specializations of take affine part *** 01191 *****************************************************/ 01192 01193 template<typename TransformType> struct transform_take_affine_part { 01194 typedef typename TransformType::MatrixType MatrixType; 01195 typedef typename TransformType::AffinePart AffinePart; 01196 typedef typename TransformType::ConstAffinePart ConstAffinePart; 01197 static inline AffinePart run(MatrixType& m) 01198 { return m.template block<TransformType::Dim,TransformType::HDim>(0,0); } 01199 static inline ConstAffinePart run(const MatrixType& m) 01200 { return m.template block<TransformType::Dim,TransformType::HDim>(0,0); } 01201 }; 01202 01203 template<typename Scalar, int Dim, int Options> 01204 struct transform_take_affine_part<Transform<Scalar,Dim,AffineCompact, Options> > { 01205 typedef typename Transform<Scalar,Dim,AffineCompact,Options>::MatrixType MatrixType; 01206 static inline MatrixType& run(MatrixType& m) { return m; } 01207 static inline const MatrixType& run(const MatrixType& m) { return m; } 01208 }; 01209 01210 /***************************************************** 01211 *** Specializations of construct from matrix *** 01212 *****************************************************/ 01213 01214 template<typename Other, int Mode, int Options, int Dim, int HDim> 01215 struct transform_construct_from_matrix<Other, Mode,Options,Dim,HDim, Dim,Dim> 01216 { 01217 static inline void run(Transform<typename Other::Scalar,Dim,Mode,Options> *transform, const Other& other) 01218 { 01219 transform->linear() = other; 01220 transform->translation().setZero(); 01221 transform->makeAffine(); 01222 } 01223 }; 01224 01225 template<typename Other, int Mode, int Options, int Dim, int HDim> 01226 struct transform_construct_from_matrix<Other, Mode,Options,Dim,HDim, Dim,HDim> 01227 { 01228 static inline void run(Transform<typename Other::Scalar,Dim,Mode,Options> *transform, const Other& other) 01229 { 01230 transform->affine() = other; 01231 transform->makeAffine(); 01232 } 01233 }; 01234 01235 template<typename Other, int Mode, int Options, int Dim, int HDim> 01236 struct transform_construct_from_matrix<Other, Mode,Options,Dim,HDim, HDim,HDim> 01237 { 01238 static inline void run(Transform<typename Other::Scalar,Dim,Mode,Options> *transform, const Other& other) 01239 { transform->matrix() = other; } 01240 }; 01241 01242 template<typename Other, int Options, int Dim, int HDim> 01243 struct transform_construct_from_matrix<Other, AffineCompact,Options,Dim,HDim, HDim,HDim> 01244 { 01245 static inline void run(Transform<typename Other::Scalar,Dim,AffineCompact,Options> *transform, const Other& other) 01246 { transform->matrix() = other.template block<Dim,HDim>(0,0); } 01247 }; 01248 01249 /********************************************************** 01250 *** Specializations of operator* with rhs EigenBase *** 01251 **********************************************************/ 01252 01253 template<int LhsMode,int RhsMode> 01254 struct transform_product_result 01255 { 01256 enum 01257 { 01258 Mode = 01259 (LhsMode == (int)Projective || RhsMode == (int)Projective ) ? Projective : 01260 (LhsMode == (int)Affine || RhsMode == (int)Affine ) ? Affine : 01261 (LhsMode == (int)AffineCompact || RhsMode == (int)AffineCompact ) ? AffineCompact : 01262 (LhsMode == (int)Isometry || RhsMode == (int)Isometry ) ? Isometry : Projective 01263 }; 01264 }; 01265 01266 template< typename TransformType, typename MatrixType > 01267 struct transform_right_product_impl< TransformType, MatrixType, 0 > 01268 { 01269 typedef typename MatrixType::PlainObject ResultType; 01270 01271 static EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other) 01272 { 01273 return T.matrix() * other; 01274 } 01275 }; 01276 01277 template< typename TransformType, typename MatrixType > 01278 struct transform_right_product_impl< TransformType, MatrixType, 1 > 01279 { 01280 enum { 01281 Dim = TransformType::Dim, 01282 HDim = TransformType::HDim, 01283 OtherRows = MatrixType::RowsAtCompileTime, 01284 OtherCols = MatrixType::ColsAtCompileTime 01285 }; 01286 01287 typedef typename MatrixType::PlainObject ResultType; 01288 01289 static EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other) 01290 { 01291 EIGEN_STATIC_ASSERT(OtherRows==HDim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES); 01292 01293 typedef Block<ResultType, Dim, OtherCols, int(MatrixType::RowsAtCompileTime)==Dim> TopLeftLhs; 01294 01295 ResultType res(other.rows(),other.cols()); 01296 TopLeftLhs(res, 0, 0, Dim, other.cols()).noalias() = T.affine() * other; 01297 res.row(OtherRows-1) = other.row(OtherRows-1); 01298 01299 return res; 01300 } 01301 }; 01302 01303 template< typename TransformType, typename MatrixType > 01304 struct transform_right_product_impl< TransformType, MatrixType, 2 > 01305 { 01306 enum { 01307 Dim = TransformType::Dim, 01308 HDim = TransformType::HDim, 01309 OtherRows = MatrixType::RowsAtCompileTime, 01310 OtherCols = MatrixType::ColsAtCompileTime 01311 }; 01312 01313 typedef typename MatrixType::PlainObject ResultType; 01314 01315 static EIGEN_STRONG_INLINE ResultType run(const TransformType& T, const MatrixType& other) 01316 { 01317 EIGEN_STATIC_ASSERT(OtherRows==Dim, YOU_MIXED_MATRICES_OF_DIFFERENT_SIZES); 01318 01319 typedef Block<ResultType, Dim, OtherCols, true> TopLeftLhs; 01320 ResultType res(Replicate<typename TransformType::ConstTranslationPart, 1, OtherCols>(T.translation(),1,other.cols())); 01321 TopLeftLhs(res, 0, 0, Dim, other.cols()).noalias() += T.linear() * other; 01322 01323 return res; 01324 } 01325 }; 01326 01327 /********************************************************** 01328 *** Specializations of operator* with lhs EigenBase *** 01329 **********************************************************/ 01330 01331 // generic HDim x HDim matrix * T => Projective 01332 template<typename Other,int Mode, int Options, int Dim, int HDim> 01333 struct transform_left_product_impl<Other,Mode,Options,Dim,HDim, HDim,HDim> 01334 { 01335 typedef Transform<typename Other::Scalar,Dim,Mode,Options> TransformType; 01336 typedef typename TransformType::MatrixType MatrixType; 01337 typedef Transform<typename Other::Scalar,Dim,Projective,Options> ResultType; 01338 static ResultType run(const Other& other,const TransformType& tr) 01339 { return ResultType(other * tr.matrix()); } 01340 }; 01341 01342 // generic HDim x HDim matrix * AffineCompact => Projective 01343 template<typename Other, int Options, int Dim, int HDim> 01344 struct transform_left_product_impl<Other,AffineCompact,Options,Dim,HDim, HDim,HDim> 01345 { 01346 typedef Transform<typename Other::Scalar,Dim,AffineCompact,Options> TransformType; 01347 typedef typename TransformType::MatrixType MatrixType; 01348 typedef Transform<typename Other::Scalar,Dim,Projective,Options> ResultType; 01349 static ResultType run(const Other& other,const TransformType& tr) 01350 { 01351 ResultType res; 01352 res.matrix().noalias() = other.template block<HDim,Dim>(0,0) * tr.matrix(); 01353 res.matrix().col(Dim) += other.col(Dim); 01354 return res; 01355 } 01356 }; 01357 01358 // affine matrix * T 01359 template<typename Other,int Mode, int Options, int Dim, int HDim> 01360 struct transform_left_product_impl<Other,Mode,Options,Dim,HDim, Dim,HDim> 01361 { 01362 typedef Transform<typename Other::Scalar,Dim,Mode,Options> TransformType; 01363 typedef typename TransformType::MatrixType MatrixType; 01364 typedef TransformType ResultType; 01365 static ResultType run(const Other& other,const TransformType& tr) 01366 { 01367 ResultType res; 01368 res.affine().noalias() = other * tr.matrix(); 01369 res.matrix().row(Dim) = tr.matrix().row(Dim); 01370 return res; 01371 } 01372 }; 01373 01374 // affine matrix * AffineCompact 01375 template<typename Other, int Options, int Dim, int HDim> 01376 struct transform_left_product_impl<Other,AffineCompact,Options,Dim,HDim, Dim,HDim> 01377 { 01378 typedef Transform<typename Other::Scalar,Dim,AffineCompact,Options> TransformType; 01379 typedef typename TransformType::MatrixType MatrixType; 01380 typedef TransformType ResultType; 01381 static ResultType run(const Other& other,const TransformType& tr) 01382 { 01383 ResultType res; 01384 res.matrix().noalias() = other.template block<Dim,Dim>(0,0) * tr.matrix(); 01385 res.translation() += other.col(Dim); 01386 return res; 01387 } 01388 }; 01389 01390 // linear matrix * T 01391 template<typename Other,int Mode, int Options, int Dim, int HDim> 01392 struct transform_left_product_impl<Other,Mode,Options,Dim,HDim, Dim,Dim> 01393 { 01394 typedef Transform<typename Other::Scalar,Dim,Mode,Options> TransformType; 01395 typedef typename TransformType::MatrixType MatrixType; 01396 typedef TransformType ResultType; 01397 static ResultType run(const Other& other, const TransformType& tr) 01398 { 01399 TransformType res; 01400 if(Mode!=int(AffineCompact)) 01401 res.matrix().row(Dim) = tr.matrix().row(Dim); 01402 res.matrix().template topRows<Dim>().noalias() 01403 = other * tr.matrix().template topRows<Dim>(); 01404 return res; 01405 } 01406 }; 01407 01408 /********************************************************** 01409 *** Specializations of operator* with another Transform *** 01410 **********************************************************/ 01411 01412 template<typename Scalar, int Dim, int LhsMode, int LhsOptions, int RhsMode, int RhsOptions> 01413 struct transform_transform_product_impl<Transform<Scalar,Dim,LhsMode,LhsOptions>,Transform<Scalar,Dim,RhsMode,RhsOptions>,false > 01414 { 01415 enum { ResultMode = transform_product_result<LhsMode,RhsMode>::Mode }; 01416 typedef Transform<Scalar,Dim,LhsMode,LhsOptions> Lhs; 01417 typedef Transform<Scalar,Dim,RhsMode,RhsOptions> Rhs; 01418 typedef Transform<Scalar,Dim,ResultMode,LhsOptions> ResultType; 01419 static ResultType run(const Lhs& lhs, const Rhs& rhs) 01420 { 01421 ResultType res; 01422 res.linear() = lhs.linear() * rhs.linear(); 01423 res.translation() = lhs.linear() * rhs.translation() + lhs.translation(); 01424 res.makeAffine(); 01425 return res; 01426 } 01427 }; 01428 01429 template<typename Scalar, int Dim, int LhsMode, int LhsOptions, int RhsMode, int RhsOptions> 01430 struct transform_transform_product_impl<Transform<Scalar,Dim,LhsMode,LhsOptions>,Transform<Scalar,Dim,RhsMode,RhsOptions>,true > 01431 { 01432 typedef Transform<Scalar,Dim,LhsMode,LhsOptions> Lhs; 01433 typedef Transform<Scalar,Dim,RhsMode,RhsOptions> Rhs; 01434 typedef Transform<Scalar,Dim,Projective> ResultType; 01435 static ResultType run(const Lhs& lhs, const Rhs& rhs) 01436 { 01437 return ResultType( lhs.matrix() * rhs.matrix() ); 01438 } 01439 }; 01440 01441 template<typename Scalar, int Dim, int LhsOptions, int RhsOptions> 01442 struct transform_transform_product_impl<Transform<Scalar,Dim,AffineCompact,LhsOptions>,Transform<Scalar,Dim,Projective,RhsOptions>,true > 01443 { 01444 typedef Transform<Scalar,Dim,AffineCompact,LhsOptions> Lhs; 01445 typedef Transform<Scalar,Dim,Projective,RhsOptions> Rhs; 01446 typedef Transform<Scalar,Dim,Projective> ResultType; 01447 static ResultType run(const Lhs& lhs, const Rhs& rhs) 01448 { 01449 ResultType res; 01450 res.matrix().template topRows<Dim>() = lhs.matrix() * rhs.matrix(); 01451 res.matrix().row(Dim) = rhs.matrix().row(Dim); 01452 return res; 01453 } 01454 }; 01455 01456 template<typename Scalar, int Dim, int LhsOptions, int RhsOptions> 01457 struct transform_transform_product_impl<Transform<Scalar,Dim,Projective,LhsOptions>,Transform<Scalar,Dim,AffineCompact,RhsOptions>,true > 01458 { 01459 typedef Transform<Scalar,Dim,Projective,LhsOptions> Lhs; 01460 typedef Transform<Scalar,Dim,AffineCompact,RhsOptions> Rhs; 01461 typedef Transform<Scalar,Dim,Projective> ResultType; 01462 static ResultType run(const Lhs& lhs, const Rhs& rhs) 01463 { 01464 ResultType res(lhs.matrix().template leftCols<Dim>() * rhs.matrix()); 01465 res.matrix().col(Dim) += lhs.matrix().col(Dim); 01466 return res; 01467 } 01468 }; 01469 01470 } // end namespace internal 01471 01472 } // end namespace Eigen 01473 01474 #endif // EIGEN_TRANSFORM_H
Generated on Tue Jul 12 2022 17:47:00 by 1.7.2