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