Important changes to repositories hosted on mbed.com
Mbed hosted mercurial repositories are deprecated and are due to be permanently deleted in July 2026.
To keep a copy of this software download the repository Zip archive or clone locally using Mercurial.
It is also possible to export all your personal repositories from the account settings page.
EulerAngles.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 // 00006 // This Source Code Form is subject to the terms of the Mozilla 00007 // Public License v. 2.0. If a copy of the MPL was not distributed 00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 00009 00010 #ifndef EIGEN_EULERANGLES_H 00011 #define EIGEN_EULERANGLES_H 00012 00013 namespace Eigen { 00014 00015 /** \geometry_module \ingroup Geometry_Module 00016 * 00017 * 00018 * \returns the Euler-angles of the rotation matrix \c *this using the convention defined by the triplet (\a a0,\a a1,\a a2) 00019 * 00020 * Each of the three parameters \a a0,\a a1,\a a2 represents the respective rotation axis as an integer in {0,1,2}. 00021 * For instance, in: 00022 * \code Vector3f ea = mat.eulerAngles(2, 0, 2); \endcode 00023 * "2" represents the z axis and "0" the x axis, etc. The returned angles are such that 00024 * we have the following equality: 00025 * \code 00026 * mat == AngleAxisf(ea[0], Vector3f::UnitZ()) 00027 * * AngleAxisf(ea[1], Vector3f::UnitX()) 00028 * * AngleAxisf(ea[2], Vector3f::UnitZ()); \endcode 00029 * This corresponds to the right-multiply conventions (with right hand side frames). 00030 * 00031 * The returned angles are in the ranges [0:pi]x[-pi:pi]x[-pi:pi]. 00032 * 00033 * \sa class AngleAxis 00034 */ 00035 template<typename Derived> 00036 inline Matrix<typename MatrixBase<Derived>::Scalar,3,1> 00037 MatrixBase<Derived>::eulerAngles (Index a0, Index a1, Index a2) const 00038 { 00039 using std::atan2; 00040 using std::sin; 00041 using std::cos; 00042 /* Implemented from Graphics Gems IV */ 00043 EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Derived,3,3) 00044 00045 Matrix<Scalar,3,1> res; 00046 typedef Matrix<typename Derived::Scalar,2,1> Vector2; 00047 00048 const Index odd = ((a0+1)%3 == a1) ? 0 : 1; 00049 const Index i = a0; 00050 const Index j = (a0 + 1 + odd)%3; 00051 const Index k = (a0 + 2 - odd)%3; 00052 00053 if (a0==a2) 00054 { 00055 res[0] = atan2(coeff(j,i), coeff(k,i)); 00056 if((odd && res[0]<Scalar(0)) || ((!odd) && res[0]>Scalar(0))) 00057 { 00058 res[0] = (res[0] > Scalar(0)) ? res[0] - Scalar(M_PI) : res[0] + Scalar(M_PI); 00059 Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm(); 00060 res[1] = -atan2(s2, coeff(i,i)); 00061 } 00062 else 00063 { 00064 Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm(); 00065 res[1] = atan2(s2, coeff(i,i)); 00066 } 00067 00068 // With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles, 00069 // we can compute their respective rotation, and apply its inverse to M. Since the result must 00070 // be a rotation around x, we have: 00071 // 00072 // c2 s1.s2 c1.s2 1 0 0 00073 // 0 c1 -s1 * M = 0 c3 s3 00074 // -s2 s1.c2 c1.c2 0 -s3 c3 00075 // 00076 // Thus: m11.c1 - m21.s1 = c3 & m12.c1 - m22.s1 = s3 00077 00078 Scalar s1 = sin(res[0]); 00079 Scalar c1 = cos(res[0]); 00080 res[2] = atan2(c1*coeff(j,k)-s1*coeff(k,k), c1*coeff(j,j) - s1 * coeff(k,j)); 00081 } 00082 else 00083 { 00084 res[0] = atan2(coeff(j,k), coeff(k,k)); 00085 Scalar c2 = Vector2(coeff(i,i), coeff(i,j)).norm(); 00086 if((odd && res[0]<Scalar(0)) || ((!odd) && res[0]>Scalar(0))) { 00087 res[0] = (res[0] > Scalar(0)) ? res[0] - Scalar(M_PI) : res[0] + Scalar(M_PI); 00088 res[1] = atan2(-coeff(i,k), -c2); 00089 } 00090 else 00091 res[1] = atan2(-coeff(i,k), c2); 00092 Scalar s1 = sin(res[0]); 00093 Scalar c1 = cos(res[0]); 00094 res[2] = atan2(s1*coeff(k,i)-c1*coeff(j,i), c1*coeff(j,j) - s1 * coeff(k,j)); 00095 } 00096 if (!odd) 00097 res = -res; 00098 00099 return res; 00100 } 00101 00102 } // end namespace Eigen 00103 00104 #endif // EIGEN_EULERANGLES_H
Generated on Thu Nov 17 2022 22:01:28 by
1.7.2