Eigen  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
EulerAngles.h
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 
10 #ifndef EIGEN_EULERANGLES_H
11 #define EIGEN_EULERANGLES_H
12 
13 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
37 template<typename Derived>
38 EIGEN_DEVICE_FUNC inline Matrix<typename MatrixBase<Derived>::Scalar,3,1>
40 {
41  EIGEN_USING_STD(atan2)
42  EIGEN_USING_STD(sin)
43  EIGEN_USING_STD(cos)
44  /* Implemented from Graphics Gems IV */
45  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Derived,3,3)
46 
49 
50  const Index odd = ((a0+1)%3 == a1) ? 0 : 1;
51  const Index i = a0;
52  const Index j = (a0 + 1 + odd)%3;
53  const Index k = (a0 + 2 - odd)%3;
54 
55  if (a0==a2)
56  {
57  res[0] = atan2(coeff(j,i), coeff(k,i));
58  if((odd && res[0]<Scalar(0)) || ((!odd) && res[0]>Scalar(0)))
59  {
60  if(res[0] > Scalar(0)) {
61  res[0] -= Scalar(EIGEN_PI);
62  }
63  else {
64  res[0] += Scalar(EIGEN_PI);
65  }
66  Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm();
67  res[1] = -atan2(s2, coeff(i,i));
68  }
69  else
70  {
71  Scalar s2 = Vector2(coeff(j,i), coeff(k,i)).norm();
72  res[1] = atan2(s2, coeff(i,i));
73  }
74 
75  // With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles,
76  // we can compute their respective rotation, and apply its inverse to M. Since the result must
77  // be a rotation around x, we have:
78  //
79  // c2 s1.s2 c1.s2 1 0 0
80  // 0 c1 -s1 * M = 0 c3 s3
81  // -s2 s1.c2 c1.c2 0 -s3 c3
82  //
83  // Thus: m11.c1 - m21.s1 = c3 & m12.c1 - m22.s1 = s3
84 
85  Scalar s1 = sin(res[0]);
86  Scalar c1 = cos(res[0]);
87  res[2] = atan2(c1*coeff(j,k)-s1*coeff(k,k), c1*coeff(j,j) - s1 * coeff(k,j));
88  }
89  else
90  {
91  res[0] = atan2(coeff(j,k), coeff(k,k));
92  Scalar c2 = Vector2(coeff(i,i), coeff(i,j)).norm();
93  if((odd && res[0]<Scalar(0)) || ((!odd) && res[0]>Scalar(0))) {
94  if(res[0] > Scalar(0)) {
95  res[0] -= Scalar(EIGEN_PI);
96  }
97  else {
98  res[0] += Scalar(EIGEN_PI);
99  }
100  res[1] = atan2(-coeff(i,k), -c2);
101  }
102  else
103  res[1] = atan2(-coeff(i,k), c2);
104  Scalar s1 = sin(res[0]);
105  Scalar c1 = cos(res[0]);
106  res[2] = atan2(s1*coeff(k,i)-c1*coeff(j,i), c1*coeff(j,j) - s1 * coeff(k,j));
107  }
108  if (!odd)
109  res = -res;
110 
111  return res;
112 }
113 
114 } // end namespace Eigen
115 
116 #endif // EIGEN_EULERANGLES_H
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:61
Matrix< Scalar, 3, 1 > eulerAngles(Index a0, Index a1, Index a2) const
Definition: EulerAngles.h:39
Matrix< Type, 2, 1 > Vector2
[c++11] 2×1 vector of type Type.
Definition: Matrix.h:533
Namespace containing all symbols from the Eigen library.
Definition: Core:139
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_cos_op< typename Derived::Scalar >, const Derived > cos(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sin_op< typename Derived::Scalar >, const Derived > sin(const Eigen::ArrayBase< Derived > &x)
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41