Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/eigen3/Eigen/src/Geometry/EulerAngles.h was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
0005 //
0006 // This Source Code Form is subject to the terms of the Mozilla
0007 // Public License v. 2.0. If a copy of the MPL was not distributed
0008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0009 
0010 #ifndef EIGEN_EULERANGLES_H
0011 #define EIGEN_EULERANGLES_H
0012 
0013 namespace Eigen { 
0014 
0015 /** \geometry_module \ingroup Geometry_Module
0016   *
0017   *
0018   * \returns the Euler-angles of the rotation matrix \c *this using the convention defined by the triplet (\a a0,\a a1,\a a2)
0019   *
0020   * Each of the three parameters \a a0,\a a1,\a a2 represents the respective rotation axis as an integer in {0,1,2}.
0021   * For instance, in:
0022   * \code Vector3f ea = mat.eulerAngles(2, 0, 2); \endcode
0023   * "2" represents the z axis and "0" the x axis, etc. The returned angles are such that
0024   * we have the following equality:
0025   * \code
0026   * mat == AngleAxisf(ea[0], Vector3f::UnitZ())
0027   *      * AngleAxisf(ea[1], Vector3f::UnitX())
0028   *      * AngleAxisf(ea[2], Vector3f::UnitZ()); \endcode
0029   * This corresponds to the right-multiply conventions (with right hand side frames).
0030   * 
0031   * The returned angles are in the ranges [0:pi]x[-pi:pi]x[-pi:pi].
0032   * 
0033   * \sa class AngleAxis
0034   */
0035 template<typename Derived>
0036 EIGEN_DEVICE_FUNC inline Matrix<typename MatrixBase<Derived>::Scalar,3,1>
0037 MatrixBase<Derived>::eulerAngles(Index a0, Index a1, Index a2) const
0038 {
0039   EIGEN_USING_STD(atan2)
0040   EIGEN_USING_STD(sin)
0041   EIGEN_USING_STD(cos)
0042   /* Implemented from Graphics Gems IV */
0043   EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Derived,3,3)
0044 
0045   Matrix<Scalar,3,1> res;
0046   typedef Matrix<typename Derived::Scalar,2,1> Vector2;
0047 
0048   const Index odd = ((a0+1)%3 == a1) ? 0 : 1;
0049   const Index i = a0;
0050   const Index j = (a0 + 1 + odd)%3;
0051   const Index k = (a0 + 2 - odd)%3;
0052   
0053   if (a0==a2)
0054   {
0055     res[0] = atan2(coeff(j,i), coeff(k,i));
0056     if((odd && res[0]<Scalar(0)) || ((!odd) && res[0]>Scalar(0)))
0057     {
0058       if(res[0] > Scalar(0)) {
0059         res[0] -= Scalar(EIGEN_PI);
0060       }
0061       else {
0062         res[0] += Scalar(EIGEN_PI);
0063       }
0064       Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm();
0065       res[1] = -atan2(s2, coeff(i,i));
0066     }
0067     else
0068     {
0069       Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm();
0070       res[1] = atan2(s2, coeff(i,i));
0071     }
0072     
0073     // With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles,
0074     // we can compute their respective rotation, and apply its inverse to M. Since the result must
0075     // be a rotation around x, we have:
0076     //
0077     //  c2  s1.s2 c1.s2                   1  0   0 
0078     //  0   c1    -s1       *    M    =   0  c3  s3
0079     //  -s2 s1.c2 c1.c2                   0 -s3  c3
0080     //
0081     //  Thus:  m11.c1 - m21.s1 = c3  &   m12.c1 - m22.s1 = s3
0082     
0083     Scalar s1 = sin(res[0]);
0084     Scalar c1 = cos(res[0]);
0085     res[2] = atan2(c1*coeff(j,k)-s1*coeff(k,k), c1*coeff(j,j) - s1 * coeff(k,j));
0086   } 
0087   else
0088   {
0089     res[0] = atan2(coeff(j,k), coeff(k,k));
0090     Scalar c2 = Vector2(coeff(i,i), coeff(i,j)).norm();
0091     if((odd && res[0]<Scalar(0)) || ((!odd) && res[0]>Scalar(0))) {
0092       if(res[0] > Scalar(0)) {
0093         res[0] -= Scalar(EIGEN_PI);
0094       }
0095       else {
0096         res[0] += Scalar(EIGEN_PI);
0097       }
0098       res[1] = atan2(-coeff(i,k), -c2);
0099     }
0100     else
0101       res[1] = atan2(-coeff(i,k), c2);
0102     Scalar s1 = sin(res[0]);
0103     Scalar c1 = cos(res[0]);
0104     res[2] = atan2(s1*coeff(k,i)-c1*coeff(j,i), c1*coeff(j,j) - s1 * coeff(k,j));
0105   }
0106   if (!odd)
0107     res = -res;
0108   
0109   return res;
0110 }
0111 
0112 } // end namespace Eigen
0113 
0114 #endif // EIGEN_EULERANGLES_H