EulerAngles.h 3.6 KB

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