![]() |
Eigen
3.3.3
|
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 00035 template<typename Derived> 00036 EIGEN_DEVICE_FUNC inline Matrix<typename MatrixBase<Derived>::Scalar,3,1> 00037 MatrixBase<Derived>::eulerAngles(Index a0, Index a1, Index a2) const 00038 { 00039 EIGEN_USING_STD_MATH(atan2) 00040 EIGEN_USING_STD_MATH(sin) 00041 EIGEN_USING_STD_MATH(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 if(res[0] > Scalar(0)) { 00059 res[0] -= Scalar(EIGEN_PI); 00060 } 00061 else { 00062 res[0] += Scalar(EIGEN_PI); 00063 } 00064 Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm(); 00065 res[1] = -atan2(s2, coeff(i,i)); 00066 } 00067 else 00068 { 00069 Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm(); 00070 res[1] = atan2(s2, coeff(i,i)); 00071 } 00072 00073 // With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles, 00074 // we can compute their respective rotation, and apply its inverse to M. Since the result must 00075 // be a rotation around x, we have: 00076 // 00077 // c2 s1.s2 c1.s2 1 0 0 00078 // 0 c1 -s1 * M = 0 c3 s3 00079 // -s2 s1.c2 c1.c2 0 -s3 c3 00080 // 00081 // Thus: m11.c1 - m21.s1 = c3 & m12.c1 - m22.s1 = s3 00082 00083 Scalar s1 = sin(res[0]); 00084 Scalar c1 = cos(res[0]); 00085 res[2] = atan2(c1*coeff(j,k)-s1*coeff(k,k), c1*coeff(j,j) - s1 * coeff(k,j)); 00086 } 00087 else 00088 { 00089 res[0] = atan2(coeff(j,k), coeff(k,k)); 00090 Scalar c2 = Vector2(coeff(i,i), coeff(i,j)).norm(); 00091 if((odd && res[0]<Scalar(0)) || ((!odd) && res[0]>Scalar(0))) { 00092 if(res[0] > Scalar(0)) { 00093 res[0] -= Scalar(EIGEN_PI); 00094 } 00095 else { 00096 res[0] += Scalar(EIGEN_PI); 00097 } 00098 res[1] = atan2(-coeff(i,k), -c2); 00099 } 00100 else 00101 res[1] = atan2(-coeff(i,k), c2); 00102 Scalar s1 = sin(res[0]); 00103 Scalar c1 = cos(res[0]); 00104 res[2] = atan2(s1*coeff(k,i)-c1*coeff(j,i), c1*coeff(j,j) - s1 * coeff(k,j)); 00105 } 00106 if (!odd) 00107 res = -res; 00108 00109 return res; 00110 } 00111 00112 } // end namespace Eigen 00113 00114 #endif // EIGEN_EULERANGLES_H