Eigen-unsupported  5.0.1-dev+284dcc12
 
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) 2015 Tal Hadad <tal_hd@hotmail.com>
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_EULERANGLESCLASS_H // TODO: Fix previous "EIGEN_EULERANGLES_H" definition?
11#define EIGEN_EULERANGLESCLASS_H
12
13// IWYU pragma: private
14#include "./InternalHeaderCheck.h"
15
16namespace Eigen {
102template <typename Scalar_, class _System>
103class EulerAngles : public RotationBase<EulerAngles<Scalar_, _System>, 3> {
104 public:
106
108 typedef Scalar_ Scalar;
109 typedef typename NumTraits<Scalar>::Real RealScalar;
110
112 typedef _System System;
113
118
121 const Vector3& u = Vector3::Unit(System::AlphaAxisAbs - 1);
122 return System::IsAlphaOpposite ? -u : u;
123 }
124
127 const Vector3& u = Vector3::Unit(System::BetaAxisAbs - 1);
128 return System::IsBetaOpposite ? -u : u;
129 }
130
133 const Vector3& u = Vector3::Unit(System::GammaAxisAbs - 1);
134 return System::IsGammaOpposite ? -u : u;
135 }
136
137 private:
138 Vector3 m_angles;
139
140 public:
144 EulerAngles(const Scalar& alpha, const Scalar& beta, const Scalar& gamma) : m_angles(alpha, beta, gamma) {}
145
146 // TODO: Test this constructor
148 explicit EulerAngles(const Scalar* data) : m_angles(data) {}
149
162 template <typename Derived>
163 explicit EulerAngles(const MatrixBase<Derived>& other) {
164 *this = other;
165 }
166
178 template <typename Derived>
180 System::CalcEulerAngles(*this, rot.toRotationMatrix());
181 }
182
183 /*EulerAngles(const QuaternionType& q)
184 {
185 // TODO: Implement it in a faster way for quaternions
186 // According to http://www.euclideanspace.com/maths/geometry/rotations/conversions/quaternionToEuler/
187 // we can compute only the needed matrix cells and then convert to euler angles. (see ZYX example below)
188 // Currently we compute all matrix cells from quaternion.
189
190 // Special case only for ZYX
191 //Scalar y2 = q.y() * q.y();
192 //m_angles[0] = std::atan2(2*(q.w()*q.z() + q.x()*q.y()), (1 - 2*(y2 + q.z()*q.z())));
193 //m_angles[1] = std::asin( 2*(q.w()*q.y() - q.z()*q.x()));
194 //m_angles[2] = std::atan2(2*(q.w()*q.x() + q.y()*q.z()), (1 - 2*(q.x()*q.x() + y2)));
195 }*/
196
198 const Vector3& angles() const { return m_angles; }
200 Vector3& angles() { return m_angles; }
201
203 Scalar alpha() const { return m_angles[0]; }
205 Scalar& alpha() { return m_angles[0]; }
206
208 Scalar beta() const { return m_angles[1]; }
210 Scalar& beta() { return m_angles[1]; }
211
213 Scalar gamma() const { return m_angles[2]; }
215 Scalar& gamma() { return m_angles[2]; }
216
221 EulerAngles res;
222 res.m_angles = -m_angles;
223 return res;
224 }
225
229 EulerAngles operator-() const { return inverse(); }
230
238 template <class Derived>
240 EIGEN_STATIC_ASSERT(
241 (internal::is_same<Scalar, typename Derived::Scalar>::value),
242 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
243
244 internal::eulerangles_assign_impl<System, Derived>::run(*this, other.derived());
245 return *this;
246 }
247
248 // TODO: Assign and construct from another EulerAngles (with different system)
249
255 template <typename Derived>
257 System::CalcEulerAngles(*this, rot.toRotationMatrix());
258 return *this;
259 }
260
265 bool isApprox(const EulerAngles& other, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const {
266 return angles().isApprox(other.angles(), prec);
267 }
268
271 // TODO: Calc it faster
272 return static_cast<QuaternionType>(*this).toRotationMatrix();
273 }
274
280
281 friend std::ostream& operator<<(std::ostream& s, const EulerAngles<Scalar, System>& eulerAngles) {
282 s << eulerAngles.angles().transpose();
283 return s;
284 }
285
287 template <typename NewScalarType>
293};
294
295#define EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(AXES, SCALAR_TYPE, SCALAR_POSTFIX) \
296 \
297 typedef EulerAngles<SCALAR_TYPE, EulerSystem##AXES> EulerAngles##AXES##SCALAR_POSTFIX;
298
299#define EIGEN_EULER_ANGLES_TYPEDEFS(SCALAR_TYPE, SCALAR_POSTFIX) \
300 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XYZ, SCALAR_TYPE, SCALAR_POSTFIX) \
301 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XYX, SCALAR_TYPE, SCALAR_POSTFIX) \
302 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XZY, SCALAR_TYPE, SCALAR_POSTFIX) \
303 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XZX, SCALAR_TYPE, SCALAR_POSTFIX) \
304 \
305 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YZX, SCALAR_TYPE, SCALAR_POSTFIX) \
306 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YZY, SCALAR_TYPE, SCALAR_POSTFIX) \
307 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YXZ, SCALAR_TYPE, SCALAR_POSTFIX) \
308 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YXY, SCALAR_TYPE, SCALAR_POSTFIX) \
309 \
310 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZXY, SCALAR_TYPE, SCALAR_POSTFIX) \
311 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZXZ, SCALAR_TYPE, SCALAR_POSTFIX) \
312 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZYX, SCALAR_TYPE, SCALAR_POSTFIX) \
313 EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZYZ, SCALAR_TYPE, SCALAR_POSTFIX)
314
315EIGEN_EULER_ANGLES_TYPEDEFS(float, f)
316EIGEN_EULER_ANGLES_TYPEDEFS(double, d)
317
318// Specifically-referenced instantiations.
327
328namespace internal {
329template <typename Scalar_, class _System>
330struct traits<EulerAngles<Scalar_, _System> > {
331 typedef Scalar_ Scalar;
332};
333
334// set from a rotation matrix
335template <class System, class Other>
336struct eulerangles_assign_impl<System, Other, 3, 3> {
337 typedef typename Other::Scalar Scalar;
338 static void run(EulerAngles<Scalar, System>& e, const Other& m) { System::CalcEulerAngles(e, m); }
339};
340
341// set from a vector of Euler angles
342template <class System, class Other>
343struct eulerangles_assign_impl<System, Other, 3, 1> {
344 typedef typename Other::Scalar Scalar;
345 static void run(EulerAngles<Scalar, System>& e, const Other& vec) { e.angles() = vec; }
346};
347} // namespace internal
348} // namespace Eigen
349
350#endif // EIGEN_EULERANGLESCLASS_H
Represents a rotation in a 3 dimensional space as three Euler angles.
Definition EulerAngles.h:103
Scalar & gamma()
Definition EulerAngles.h:215
_System System
Definition EulerAngles.h:112
Scalar & alpha()
Definition EulerAngles.h:205
static Vector3 GammaAxisVector()
Definition EulerAngles.h:132
Scalar alpha() const
Definition EulerAngles.h:203
static Vector3 AlphaAxisVector()
Definition EulerAngles.h:120
Vector3 & angles()
Definition EulerAngles.h:200
Quaternion< Scalar > QuaternionType
Definition EulerAngles.h:116
const Vector3 & angles() const
Definition EulerAngles.h:198
EulerAngles(const Scalar &alpha, const Scalar &beta, const Scalar &gamma)
Definition EulerAngles.h:144
EulerAngles(const Scalar *data)
Definition EulerAngles.h:148
EulerAngles< NewScalarType, System > cast() const
Definition EulerAngles.h:288
AngleAxis< Scalar > AngleAxisType
Definition EulerAngles.h:117
Matrix3 toRotationMatrix() const
Definition EulerAngles.h:270
EulerAngles & operator=(const MatrixBase< Derived > &other)
Definition EulerAngles.h:239
Scalar gamma() const
Definition EulerAngles.h:213
Scalar_ Scalar
Definition EulerAngles.h:108
Scalar & beta()
Definition EulerAngles.h:210
static Vector3 BetaAxisVector()
Definition EulerAngles.h:126
Scalar beta() const
Definition EulerAngles.h:208
EulerAngles inverse() const
Definition EulerAngles.h:220
EulerAngles(const MatrixBase< Derived > &other)
Definition EulerAngles.h:163
EulerAngles operator-() const
Definition EulerAngles.h:229
Matrix< Scalar, 3, 3 > Matrix3
Definition EulerAngles.h:114
EulerAngles(const RotationBase< Derived, 3 > &rot)
Definition EulerAngles.h:179
EulerAngles & operator=(const RotationBase< Derived, 3 > &rot)
Definition EulerAngles.h:256
bool isApprox(const EulerAngles &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition EulerAngles.h:265
EulerAngles()
Definition EulerAngles.h:142
Matrix< Scalar, 3, 1 > Vector3
Definition EulerAngles.h:115
@ GammaAxisAbs
Definition EulerSystem.h:138
@ IsGammaOpposite
Definition EulerSystem.h:142
@ IsAlphaOpposite
Definition EulerSystem.h:140
@ AlphaAxisAbs
Definition EulerSystem.h:136
@ BetaAxisAbs
Definition EulerSystem.h:137
@ IsBetaOpposite
Definition EulerSystem.h:141
RotationMatrixType toRotationMatrix() const
Matrix< Scalar, 3, 1 > eulerAngles(Index a0, Index a1, Index a2) const
Matrix< Type, 3, 1 > Vector3
Namespace containing all symbols from the Eigen library.