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) 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 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen
16 {
101  template <typename Scalar_, class _System>
102  class EulerAngles : public RotationBase<EulerAngles<Scalar_, _System>, 3>
103  {
104  public:
106 
108  typedef Scalar_ Scalar;
110 
112  typedef _System System;
113 
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:
139 
140  public:
144  EulerAngles(const Scalar& alpha, const Scalar& beta, const Scalar& gamma) :
145  m_angles(alpha, beta, gamma) {}
146 
147  // TODO: Test this constructor
149  explicit EulerAngles(const Scalar* data) : m_angles(data) {}
150 
163  template<typename Derived>
164  explicit EulerAngles(const MatrixBase<Derived>& other) { *this = other; }
165 
177  template<typename Derived>
178  EulerAngles(const RotationBase<Derived, 3>& rot) { System::CalcEulerAngles(*this, rot.toRotationMatrix()); }
179 
180  /*EulerAngles(const QuaternionType& q)
181  {
182  // TODO: Implement it in a faster way for quaternions
183  // According to http://www.euclideanspace.com/maths/geometry/rotations/conversions/quaternionToEuler/
184  // we can compute only the needed matrix cells and then convert to euler angles. (see ZYX example below)
185  // Currently we compute all matrix cells from quaternion.
186 
187  // Special case only for ZYX
188  //Scalar y2 = q.y() * q.y();
189  //m_angles[0] = std::atan2(2*(q.w()*q.z() + q.x()*q.y()), (1 - 2*(y2 + q.z()*q.z())));
190  //m_angles[1] = std::asin( 2*(q.w()*q.y() - q.z()*q.x()));
191  //m_angles[2] = std::atan2(2*(q.w()*q.x() + q.y()*q.z()), (1 - 2*(q.x()*q.x() + y2)));
192  }*/
193 
195  const Vector3& angles() const { return m_angles; }
197  Vector3& angles() { return m_angles; }
198 
200  Scalar alpha() const { return m_angles[0]; }
202  Scalar& alpha() { return m_angles[0]; }
203 
205  Scalar beta() const { return m_angles[1]; }
207  Scalar& beta() { return m_angles[1]; }
208 
210  Scalar gamma() const { return m_angles[2]; }
212  Scalar& gamma() { return m_angles[2]; }
213 
218  {
220  res.m_angles = -m_angles;
221  return res;
222  }
223 
228  {
229  return inverse();
230  }
231 
239  template<class Derived>
241  {
242  EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename Derived::Scalar>::value),
243  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
244 
245  internal::eulerangles_assign_impl<System, Derived>::run(*this, other.derived());
246  return *this;
247  }
248 
249  // TODO: Assign and construct from another EulerAngles (with different system)
250 
256  template<typename Derived>
258  System::CalcEulerAngles(*this, rot.toRotationMatrix());
259  return *this;
260  }
261 
266  bool isApprox(const EulerAngles& other,
267  const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const
268  { return angles().isApprox(other.angles(), prec); }
269 
272  {
273  // TODO: Calc it faster
274  return static_cast<QuaternionType>(*this).toRotationMatrix();
275  }
276 
278  operator QuaternionType() const
279  {
280  return
284  }
285 
286  friend std::ostream& operator<<(std::ostream& s, const EulerAngles<Scalar, System>& eulerAngles)
287  {
288  s << eulerAngles.angles().transpose();
289  return s;
290  }
291 
293  template <typename NewScalarType>
295  {
297  e.angles() = angles().template cast<NewScalarType>();
298  return e;
299  }
300  };
301 
302 #define EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(AXES, SCALAR_TYPE, SCALAR_POSTFIX) \
303  \
304  typedef EulerAngles<SCALAR_TYPE, EulerSystem##AXES> EulerAngles##AXES##SCALAR_POSTFIX;
305 
306 #define EIGEN_EULER_ANGLES_TYPEDEFS(SCALAR_TYPE, SCALAR_POSTFIX) \
307  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XYZ, SCALAR_TYPE, SCALAR_POSTFIX) \
308  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XYX, SCALAR_TYPE, SCALAR_POSTFIX) \
309  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XZY, SCALAR_TYPE, SCALAR_POSTFIX) \
310  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(XZX, SCALAR_TYPE, SCALAR_POSTFIX) \
311  \
312  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YZX, SCALAR_TYPE, SCALAR_POSTFIX) \
313  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YZY, SCALAR_TYPE, SCALAR_POSTFIX) \
314  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YXZ, SCALAR_TYPE, SCALAR_POSTFIX) \
315  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(YXY, SCALAR_TYPE, SCALAR_POSTFIX) \
316  \
317  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZXY, SCALAR_TYPE, SCALAR_POSTFIX) \
318  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZXZ, SCALAR_TYPE, SCALAR_POSTFIX) \
319  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZYX, SCALAR_TYPE, SCALAR_POSTFIX) \
320  EIGEN_EULER_ANGLES_SINGLE_TYPEDEF(ZYZ, SCALAR_TYPE, SCALAR_POSTFIX)
321 
324 
325  namespace internal
326  {
327  template<typename Scalar_, class _System>
328  struct traits<EulerAngles<Scalar_, _System> >
329  {
330  typedef Scalar_ Scalar;
331  };
332 
333  // set from a rotation matrix
334  template<class System, class Other>
335  struct eulerangles_assign_impl<System,Other,3,3>
336  {
337  typedef typename Other::Scalar Scalar;
338  static void run(EulerAngles<Scalar, System>& e, const Other& m)
339  {
340  System::CalcEulerAngles(e, m);
341  }
342  };
343 
344  // set from a vector of Euler angles
345  template<class System, class Other>
346  struct eulerangles_assign_impl<System,Other,3,1>
347  {
348  typedef typename Other::Scalar Scalar;
349  static void run(EulerAngles<Scalar, System>& e, const Other& vec)
350  {
351  e.angles() = vec;
352  }
353  };
354  }
355 }
356 
357 #endif // EIGEN_EULERANGLESCLASS_H
Matrix3f m
Array< double, 1, 3 > e(1./3., 0.5, 2.)
#define EIGEN_EULER_ANGLES_TYPEDEFS(SCALAR_TYPE, SCALAR_POSTFIX)
Definition: EulerAngles.h:306
int data[]
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
#define EIGEN_STATIC_ASSERT(X, MSG)
Represents a rotation in a 3 dimensional space as three Euler angles.
Definition: EulerAngles.h:103
EulerAngles & operator=(const RotationBase< Derived, 3 > &rot)
Definition: EulerAngles.h:257
EulerAngles & operator=(const MatrixBase< Derived > &other)
Definition: EulerAngles.h:240
static Vector3 GammaAxisVector()
Definition: EulerAngles.h:132
Scalar alpha() const
Definition: EulerAngles.h:200
static Vector3 AlphaAxisVector()
Definition: EulerAngles.h:120
EulerAngles< NewScalarType, System > cast() const
Definition: EulerAngles.h:294
Quaternion< Scalar > QuaternionType
Definition: EulerAngles.h:116
EulerAngles(const Scalar &alpha, const Scalar &beta, const Scalar &gamma)
Definition: EulerAngles.h:144
const Vector3 & angles() const
Definition: EulerAngles.h:195
EulerAngles(const Scalar *data)
Definition: EulerAngles.h:149
Scalar & gamma()
Definition: EulerAngles.h:212
AngleAxis< Scalar > AngleAxisType
Definition: EulerAngles.h:117
Matrix3 toRotationMatrix() const
Definition: EulerAngles.h:271
friend std::ostream & operator<<(std::ostream &s, const EulerAngles< Scalar, System > &eulerAngles)
Definition: EulerAngles.h:286
NumTraits< Scalar >::Real RealScalar
Definition: EulerAngles.h:109
Scalar gamma() const
Definition: EulerAngles.h:210
static Vector3 BetaAxisVector()
Definition: EulerAngles.h:126
Vector3 & angles()
Definition: EulerAngles.h:197
Scalar beta() const
Definition: EulerAngles.h:205
EulerAngles inverse() const
Definition: EulerAngles.h:217
EulerAngles(const MatrixBase< Derived > &other)
Definition: EulerAngles.h:164
EulerAngles operator-() const
Definition: EulerAngles.h:227
Matrix< Scalar, 3, 3 > Matrix3
Definition: EulerAngles.h:114
EulerAngles(const RotationBase< Derived, 3 > &rot)
Definition: EulerAngles.h:178
bool isApprox(const EulerAngles &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: EulerAngles.h:266
Scalar & alpha()
Definition: EulerAngles.h:202
RotationBase< EulerAngles< Scalar_, _System >, 3 > Base
Definition: EulerAngles.h:105
Matrix< Scalar, 3, 1 > Vector3
Definition: EulerAngles.h:115
static const BasisReturnType Unit(Index i)
RotationMatrixType toRotationMatrix() const
EIGEN_DEPRECATED Matrix< Scalar, 3, 1 > eulerAngles(Index a0, Index a1, Index a2) const
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend