You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

84 lines
2.6 KiB

  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. template<typename Derived>
  29. inline Matrix<typename MatrixBase<Derived>::Scalar,3,1>
  30. MatrixBase<Derived>::eulerAngles(Index a0, Index a1, Index a2) const
  31. {
  32. /* Implemented from Graphics Gems IV */
  33. EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Derived,3,3)
  34. Matrix<Scalar,3,1> res;
  35. typedef Matrix<typename Derived::Scalar,2,1> Vector2;
  36. const Scalar epsilon = NumTraits<Scalar>::dummy_precision();
  37. const Index odd = ((a0+1)%3 == a1) ? 0 : 1;
  38. const Index i = a0;
  39. const Index j = (a0 + 1 + odd)%3;
  40. const Index k = (a0 + 2 - odd)%3;
  41. if (a0==a2)
  42. {
  43. Scalar s = Vector2(coeff(j,i) , coeff(k,i)).norm();
  44. res[1] = internal::atan2(s, coeff(i,i));
  45. if (s > epsilon)
  46. {
  47. res[0] = internal::atan2(coeff(j,i), coeff(k,i));
  48. res[2] = internal::atan2(coeff(i,j),-coeff(i,k));
  49. }
  50. else
  51. {
  52. res[0] = Scalar(0);
  53. res[2] = (coeff(i,i)>0?1:-1)*internal::atan2(-coeff(k,j), coeff(j,j));
  54. }
  55. }
  56. else
  57. {
  58. Scalar c = Vector2(coeff(i,i) , coeff(i,j)).norm();
  59. res[1] = internal::atan2(-coeff(i,k), c);
  60. if (c > epsilon)
  61. {
  62. res[0] = internal::atan2(coeff(j,k), coeff(k,k));
  63. res[2] = internal::atan2(coeff(i,j), coeff(i,i));
  64. }
  65. else
  66. {
  67. res[0] = Scalar(0);
  68. res[2] = (coeff(i,k)>0?1:-1)*internal::atan2(-coeff(k,j), coeff(j,j));
  69. }
  70. }
  71. if (!odd)
  72. res = -res;
  73. return res;
  74. }
  75. } // end namespace Eigen
  76. #endif // EIGEN_EULERANGLES_H