Eigen/src/Geometry/EulerAngles.h
Go to the documentation of this file.
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 // Copyright (C) 2023 Juraj Oršulić, University of Zagreb <juraj.orsulic@fer.hr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_EULERANGLES_H
12 #define EIGEN_EULERANGLES_H
13 
14 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
44 template <typename Derived>
46  Index a0, Index a1, Index a2) const {
47  /* Implemented from Graphics Gems IV */
49 
51 
52  const Index odd = ((a0 + 1) % 3 == a1) ? 0 : 1;
53  const Index i = a0;
54  const Index j = (a0 + 1 + odd) % 3;
55  const Index k = (a0 + 2 - odd) % 3;
56 
57  if (a0 == a2) {
58  // Proper Euler angles (same first and last axis).
59  // The i, j, k indices enable addressing the input matrix as the XYX archetype matrix (see Graphics Gems IV),
60  // where e.g. coeff(k, i) means third column, first row in the XYX archetype matrix:
61  // c2 s2s1 s2c1
62  // s2s3 -c2s1s3 + c1c3 -c2c1s3 - s1c3
63  // -s2c3 c2s1c3 + c1s3 c2c1c3 - s1s3
64 
65  // Note: s2 is always positive.
66  Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
67  if (odd) {
68  res[0] = numext::atan2(coeff(j, i), coeff(k, i));
69  // s2 is always positive, so res[1] will be within the canonical [0, pi] range
70  res[1] = numext::atan2(s2, coeff(i, i));
71  } else {
72  // In the !odd case, signs of all three angles are flipped at the very end. To keep the solution within the
73  // canonical range, we flip the solution and make res[1] always negative here (since s2 is always positive,
74  // -atan2(s2, c2) will always be negative). The final flip at the end due to !odd will thus make res[1] positive
75  // and canonical. NB: in the general case, there are two correct solutions, but only one is canonical. For proper
76  // Euler angles, flipping from one solution to the other involves flipping the sign of the second angle res[1] and
77  // adding/subtracting pi to the first and third angles. The addition/subtraction of pi to the first angle res[0]
78  // is handled here by flipping the signs of arguments to atan2, while the calculation of the third angle does not
79  // need special adjustment since it uses the adjusted res[0] as the input and produces a correct result.
80  res[0] = numext::atan2(-coeff(j, i), -coeff(k, i));
81  res[1] = -numext::atan2(s2, coeff(i, i));
82  }
83 
84  // With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles,
85  // we can compute their respective rotation, and apply its inverse to M. Since the result must
86  // be a rotation around x, we have:
87  //
88  // c2 s1.s2 c1.s2 1 0 0
89  // 0 c1 -s1 * M = 0 c3 s3
90  // -s2 s1.c2 c1.c2 0 -s3 c3
91  //
92  // Thus: m11.c1 - m21.s1 = c3 & m12.c1 - m22.s1 = s3
93 
94  Scalar s1 = numext::sin(res[0]);
95  Scalar c1 = numext::cos(res[0]);
96  res[2] = numext::atan2(c1 * coeff(j, k) - s1 * coeff(k, k), c1 * coeff(j, j) - s1 * coeff(k, j));
97  } else {
98  // Tait-Bryan angles (all three axes are different; typically used for yaw-pitch-roll calculations).
99  // The i, j, k indices enable addressing the input matrix as the XYZ archetype matrix (see Graphics Gems IV),
100  // where e.g. coeff(k, i) means third column, first row in the XYZ archetype matrix:
101  // c2c3 s2s1c3 - c1s3 s2c1c3 + s1s3
102  // c2s3 s2s1s3 + c1c3 s2c1s3 - s1c3
103  // -s2 c2s1 c2c1
104 
105  res[0] = numext::atan2(coeff(j, k), coeff(k, k));
106 
107  Scalar c2 = numext::hypot(coeff(i, i), coeff(i, j));
108  // c2 is always positive, so the following atan2 will always return a result in the correct canonical middle angle
109  // range [-pi/2, pi/2]
110  res[1] = numext::atan2(-coeff(i, k), c2);
111 
112  Scalar s1 = numext::sin(res[0]);
113  Scalar c1 = numext::cos(res[0]);
114  res[2] = numext::atan2(s1 * coeff(k, i) - c1 * coeff(j, i), c1 * coeff(j, j) - s1 * coeff(k, j));
115  }
116  if (!odd) {
117  res = -res;
118  }
119 
120  return res;
121 }
122 
135 template <typename Derived>
138  /* Implemented from Graphics Gems IV */
140 
142 
143  const Index odd = ((a0 + 1) % 3 == a1) ? 0 : 1;
144  const Index i = a0;
145  const Index j = (a0 + 1 + odd) % 3;
146  const Index k = (a0 + 2 - odd) % 3;
147 
148  if (a0 == a2) {
149  res[0] = numext::atan2(coeff(j, i), coeff(k, i));
150  if ((odd && res[0] < Scalar(0)) || ((!odd) && res[0] > Scalar(0))) {
151  if (res[0] > Scalar(0)) {
152  res[0] -= Scalar(EIGEN_PI);
153  } else {
154  res[0] += Scalar(EIGEN_PI);
155  }
156 
157  Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
158  res[1] = -numext::atan2(s2, coeff(i, i));
159  } else {
160  Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
161  res[1] = numext::atan2(s2, coeff(i, i));
162  }
163 
164  // With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles,
165  // we can compute their respective rotation, and apply its inverse to M. Since the result must
166  // be a rotation around x, we have:
167  //
168  // c2 s1.s2 c1.s2 1 0 0
169  // 0 c1 -s1 * M = 0 c3 s3
170  // -s2 s1.c2 c1.c2 0 -s3 c3
171  //
172  // Thus: m11.c1 - m21.s1 = c3 & m12.c1 - m22.s1 = s3
173 
174  Scalar s1 = numext::sin(res[0]);
175  Scalar c1 = numext::cos(res[0]);
176  res[2] = numext::atan2(c1 * coeff(j, k) - s1 * coeff(k, k), c1 * coeff(j, j) - s1 * coeff(k, j));
177  } else {
178  res[0] = numext::atan2(coeff(j, k), coeff(k, k));
179  Scalar c2 = numext::hypot(coeff(i, i), coeff(i, j));
180  if ((odd && res[0] < Scalar(0)) || ((!odd) && res[0] > Scalar(0))) {
181  if (res[0] > Scalar(0)) {
182  res[0] -= Scalar(EIGEN_PI);
183  } else {
184  res[0] += Scalar(EIGEN_PI);
185  }
186  res[1] = numext::atan2(-coeff(i, k), -c2);
187  } else {
188  res[1] = numext::atan2(-coeff(i, k), c2);
189  }
190  Scalar s1 = numext::sin(res[0]);
191  Scalar c1 = numext::cos(res[0]);
192  res[2] = numext::atan2(s1 * coeff(k, i) - c1 * coeff(j, i), c1 * coeff(j, j) - s1 * coeff(k, j));
193  }
194  if (!odd) {
195  res = -res;
196  }
197 
198  return res;
199 }
200 
201 } // end namespace Eigen
202 
203 #endif // EIGEN_EULERANGLES_H
int i
Definition: BiCGSTAB_step_by_step.cpp:9
#define EIGEN_DEPRECATED
Definition: Macros.h:931
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:892
#define EIGEN_PI
Definition: MathFunctions.h:16
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition: PartialRedux_count.cpp:3
#define EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(TYPE, ROWS, COLS)
Definition: StaticAssert.h:55
SCALAR Scalar
Definition: bench_gemm.cpp:45
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:62
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
EIGEN_DEVICE_FUNC Matrix< Scalar, 3, 1 > canonicalEulerAngles(Index a0, Index a1, Index a2) const
Definition: Eigen/src/Geometry/EulerAngles.h:45
EIGEN_DEPRECATED EIGEN_DEVICE_FUNC Matrix< Scalar, 3, 1 > eulerAngles(Index a0, Index a1, Index a2) const
Definition: Eigen/src/Geometry/EulerAngles.h:137
char char char int int * k
Definition: level2_impl.h:374
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T cos(const T &x)
Definition: MathFunctions.h:1559
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T sin(const T &x)
Definition: MathFunctions.h:1581
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T atan2(const T &y, const T &x)
Definition: MathFunctions.h:1689
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2