Eigen  3.4.90 (git rev 21cd3fe20990a5ac1d683806f605110962aac3f1)
 
Loading...
Searching...
No Matches
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// 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#include "./InternalHeaderCheck.h"
15
16namespace Eigen {
17
41template<typename Derived>
42EIGEN_DEVICE_FUNC inline Matrix<typename MatrixBase<Derived>::Scalar,3,1>
44{
45 /* Implemented from Graphics Gems IV */
46 EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Derived, 3, 3)
47
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 // Proper Euler angles (same first and last axis).
58 // The i, j, k indices enable addressing the input matrix as the XYX archetype matrix (see Graphics Gems IV),
59 // where e.g. coeff(k, i) means third column, first row in the XYX archetype matrix:
60 // c2 s2s1 s2c1
61 // s2s3 -c2s1s3 + c1c3 -c2c1s3 - s1c3
62 // -s2c3 c2s1c3 + c1s3 c2c1c3 - s1s3
63
64 // Note: s2 is always positive.
65 Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
66 if (odd)
67 {
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 }
72 else
73 {
74 // In the !odd case, signs of all three angles are flipped at the very end. To keep the solution within the canonical range,
75 // we flip the solution and make res[1] always negative here (since s2 is always positive, -atan2(s2, c2) will always be negative).
76 // The final flip at the end due to !odd will thus make res[1] positive and canonical.
77 // NB: in the general case, there are two correct solutions, but only one is canonical. For proper Euler angles,
78 // flipping from one solution to the other involves flipping the sign of the second angle res[1] and adding/subtracting pi
79 // to the first and third angles. The addition/subtraction of pi to the first angle res[0] is handled here by flipping
80 // the signs of arguments to atan2, while the calculation of the third angle does not need special adjustment since
81 // it uses the adjusted res[0] as the input and produces a correct result.
82 res[0] = numext::atan2(-coeff(j, i), -coeff(k, i));
83 res[1] = -numext::atan2(s2, coeff(i, i));
84 }
85
86 // With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles,
87 // we can compute their respective rotation, and apply its inverse to M. Since the result must
88 // be a rotation around x, we have:
89 //
90 // c2 s1.s2 c1.s2 1 0 0
91 // 0 c1 -s1 * M = 0 c3 s3
92 // -s2 s1.c2 c1.c2 0 -s3 c3
93 //
94 // Thus: m11.c1 - m21.s1 = c3 & m12.c1 - m22.s1 = s3
95
96 Scalar s1 = numext::sin(res[0]);
97 Scalar c1 = numext::cos(res[0]);
98 res[2] = numext::atan2(c1 * coeff(j, k) - s1 * coeff(k, k), c1 * coeff(j, j) - s1 * coeff(k, j));
99 }
100 else
101 {
102 // Tait-Bryan angles (all three axes are different; typically used for yaw-pitch-roll calculations).
103 // The i, j, k indices enable addressing the input matrix as the XYZ archetype matrix (see Graphics Gems IV),
104 // where e.g. coeff(k, i) means third column, first row in the XYZ archetype matrix:
105 // c2c3 s2s1c3 - c1s3 s2c1c3 + s1s3
106 // c2s3 s2s1s3 + c1c3 s2c1s3 - s1c3
107 // -s2 c2s1 c2c1
108
109 res[0] = numext::atan2(coeff(j, k), coeff(k, k));
110
111 Scalar c2 = numext::hypot(coeff(i, i), coeff(i, j));
112 // c2 is always positive, so the following atan2 will always return a result in the correct canonical middle angle range [-pi/2, pi/2]
113 res[1] = numext::atan2(-coeff(i, k), c2);
114
115 Scalar s1 = numext::sin(res[0]);
116 Scalar c1 = numext::cos(res[0]);
117 res[2] = numext::atan2(s1 * coeff(k, i) - c1 * coeff(j, i), c1 * coeff(j, j) - s1 * coeff(k, j));
118 }
119 if (!odd)
120 {
121 res = -res;
122 }
123
124 return res;
125}
126
137template<typename Derived>
138EIGEN_DEPRECATED EIGEN_DEVICE_FUNC inline Matrix<typename MatrixBase<Derived>::Scalar,3,1>
140{
141 /* Implemented from Graphics Gems IV */
142 EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Derived, 3, 3)
143
145
146 const Index odd = ((a0 + 1) % 3 == a1) ? 0 : 1;
147 const Index i = a0;
148 const Index j = (a0 + 1 + odd) % 3;
149 const Index k = (a0 + 2 - odd) % 3;
150
151 if (a0 == a2)
152 {
153 res[0] = numext::atan2(coeff(j, i), coeff(k, i));
154 if ((odd && res[0] < Scalar(0)) || ((!odd) && res[0] > Scalar(0)))
155 {
156 if (res[0] > Scalar(0))
157 {
158 res[0] -= Scalar(EIGEN_PI);
159 }
160 else
161 {
162 res[0] += Scalar(EIGEN_PI);
163 }
164
165 Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
166 res[1] = -numext::atan2(s2, coeff(i, i));
167 }
168 else
169 {
170 Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
171 res[1] = numext::atan2(s2, coeff(i, i));
172 }
173
174 // With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles,
175 // we can compute their respective rotation, and apply its inverse to M. Since the result must
176 // be a rotation around x, we have:
177 //
178 // c2 s1.s2 c1.s2 1 0 0
179 // 0 c1 -s1 * M = 0 c3 s3
180 // -s2 s1.c2 c1.c2 0 -s3 c3
181 //
182 // Thus: m11.c1 - m21.s1 = c3 & m12.c1 - m22.s1 = s3
183
184 Scalar s1 = numext::sin(res[0]);
185 Scalar c1 = numext::cos(res[0]);
186 res[2] = numext::atan2(c1 * coeff(j, k) - s1 * coeff(k, k), c1 * coeff(j, j) - s1 * coeff(k, j));
187 }
188 else
189 {
190 res[0] = numext::atan2(coeff(j, k), coeff(k, k));
191 Scalar c2 = numext::hypot(coeff(i, i), coeff(i, j));
192 if ((odd && res[0] < Scalar(0)) || ((!odd) && res[0] > Scalar(0)))
193 {
194 if (res[0] > Scalar(0))
195 {
196 res[0] -= Scalar(EIGEN_PI);
197 }
198 else
199 {
200 res[0] += Scalar(EIGEN_PI);
201 }
202 res[1] = numext::atan2(-coeff(i, k), -c2);
203 }
204 else
205 {
206 res[1] = numext::atan2(-coeff(i, k), c2);
207 }
208 Scalar s1 = numext::sin(res[0]);
209 Scalar c1 = numext::cos(res[0]);
210 res[2] = numext::atan2(s1 * coeff(k, i) - c1 * coeff(j, i), c1 * coeff(j, j) - s1 * coeff(k, j));
211 }
212 if (!odd)
213 {
214 res = -res;
215 }
216
217 return res;
218}
219
220} // end namespace Eigen
221
222#endif // EIGEN_EULERANGLES_H
internal::traits< Derived >::Scalar Scalar
Definition DenseBase.h:61
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:182
Namespace containing all symbols from the Eigen library.
Definition Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82