Home | History | Annotate | Download | only in Geometry
      1 // This file is part of Eigen, a lightweight C++ template library
      2 // for linear algebra.
      3 //
      4 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      5 //
      6 // This Source Code Form is subject to the terms of the Mozilla
      7 // Public License v. 2.0. If a copy of the MPL was not distributed
      8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
      9 
     10 #ifndef EIGEN_EULERANGLES_H
     11 #define EIGEN_EULERANGLES_H
     12 
     13 namespace Eigen {
     14 
     15 /** \geometry_module \ingroup Geometry_Module
     16   *
     17   *
     18   * \returns the Euler-angles of the rotation matrix \c *this using the convention defined by the triplet (\a a0,\a a1,\a a2)
     19   *
     20   * Each of the three parameters \a a0,\a a1,\a a2 represents the respective rotation axis as an integer in {0,1,2}.
     21   * For instance, in:
     22   * \code Vector3f ea = mat.eulerAngles(2, 0, 2); \endcode
     23   * "2" represents the z axis and "0" the x axis, etc. The returned angles are such that
     24   * we have the following equality:
     25   * \code
     26   * mat == AngleAxisf(ea[0], Vector3f::UnitZ())
     27   *      * AngleAxisf(ea[1], Vector3f::UnitX())
     28   *      * AngleAxisf(ea[2], Vector3f::UnitZ()); \endcode
     29   * This corresponds to the right-multiply conventions (with right hand side frames).
     30   *
     31   * The returned angles are in the ranges [0:pi]x[-pi:pi]x[-pi:pi].
     32   *
     33   * \sa class AngleAxis
     34   */
     35 template<typename Derived>
     36 EIGEN_DEVICE_FUNC inline Matrix<typename MatrixBase<Derived>::Scalar,3,1>
     37 MatrixBase<Derived>::eulerAngles(Index a0, Index a1, Index a2) const
     38 {
     39   EIGEN_USING_STD_MATH(atan2)
     40   EIGEN_USING_STD_MATH(sin)
     41   EIGEN_USING_STD_MATH(cos)
     42   /* Implemented from Graphics Gems IV */
     43   EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Derived,3,3)
     44 
     45   Matrix<Scalar,3,1> res;
     46   typedef Matrix<typename Derived::Scalar,2,1> Vector2;
     47 
     48   const Index odd = ((a0+1)%3 == a1) ? 0 : 1;
     49   const Index i = a0;
     50   const Index j = (a0 + 1 + odd)%3;
     51   const Index k = (a0 + 2 - odd)%3;
     52 
     53   if (a0==a2)
     54   {
     55     res[0] = atan2(coeff(j,i), coeff(k,i));
     56     if((odd && res[0]<Scalar(0)) || ((!odd) && res[0]>Scalar(0)))
     57     {
     58       if(res[0] > Scalar(0)) {
     59         res[0] -= Scalar(EIGEN_PI);
     60       }
     61       else {
     62         res[0] += Scalar(EIGEN_PI);
     63       }
     64       Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm();
     65       res[1] = -atan2(s2, coeff(i,i));
     66     }
     67     else
     68     {
     69       Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm();
     70       res[1] = atan2(s2, coeff(i,i));
     71     }
     72 
     73     // With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles,
     74     // we can compute their respective rotation, and apply its inverse to M. Since the result must
     75     // be a rotation around x, we have:
     76     //
     77     //  c2  s1.s2 c1.s2                   1  0   0
     78     //  0   c1    -s1       *    M    =   0  c3  s3
     79     //  -s2 s1.c2 c1.c2                   0 -s3  c3
     80     //
     81     //  Thus:  m11.c1 - m21.s1 = c3  &   m12.c1 - m22.s1 = s3
     82 
     83     Scalar s1 = sin(res[0]);
     84     Scalar c1 = cos(res[0]);
     85     res[2] = atan2(c1*coeff(j,k)-s1*coeff(k,k), c1*coeff(j,j) - s1 * coeff(k,j));
     86   }
     87   else
     88   {
     89     res[0] = atan2(coeff(j,k), coeff(k,k));
     90     Scalar c2 = Vector2(coeff(i,i), coeff(i,j)).norm();
     91     if((odd && res[0]<Scalar(0)) || ((!odd) && res[0]>Scalar(0))) {
     92       if(res[0] > Scalar(0)) {
     93         res[0] -= Scalar(EIGEN_PI);
     94       }
     95       else {
     96         res[0] += Scalar(EIGEN_PI);
     97       }
     98       res[1] = atan2(-coeff(i,k), -c2);
     99     }
    100     else
    101       res[1] = atan2(-coeff(i,k), c2);
    102     Scalar s1 = sin(res[0]);
    103     Scalar c1 = cos(res[0]);
    104     res[2] = atan2(s1*coeff(k,i)-c1*coeff(j,i), c1*coeff(j,j) - s1 * coeff(k,j));
    105   }
    106   if (!odd)
    107     res = -res;
    108 
    109   return res;
    110 }
    111 
    112 } // end namespace Eigen
    113 
    114 #endif // EIGEN_EULERANGLES_H
    115