Please, help us to better know about our user community by answering the following short survey: https://forms.gle/wpyrxWi18ox9Z5ae9
Eigen  3.3.8
Quaternion.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2009 Mathieu Gautier <mathieu.gautier@cea.fr>
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_QUATERNION_H
12 #define EIGEN_QUATERNION_H
13 namespace Eigen {
14 
15 
16 /***************************************************************************
17 * Definition of QuaternionBase<Derived>
18 * The implementation is at the end of the file
19 ***************************************************************************/
20 
21 namespace internal {
22 template<typename Other,
23  int OtherRows=Other::RowsAtCompileTime,
24  int OtherCols=Other::ColsAtCompileTime>
25 struct quaternionbase_assign_impl;
26 }
27 
34 template<class Derived>
35 class QuaternionBase : public RotationBase<Derived, 3>
36 {
37  public:
39 
40  using Base::operator*;
41  using Base::derived;
42 
43  typedef typename internal::traits<Derived>::Scalar Scalar;
44  typedef typename NumTraits<Scalar>::Real RealScalar;
45  typedef typename internal::traits<Derived>::Coefficients Coefficients;
46  typedef typename Coefficients::CoeffReturnType CoeffReturnType;
47  typedef typename internal::conditional<bool(internal::traits<Derived>::Flags&LvalueBit),
48  Scalar&, CoeffReturnType>::type NonConstCoeffReturnType;
49 
50 
51  enum {
52  Flags = Eigen::internal::traits<Derived>::Flags
53  };
54 
55  // typedef typename Matrix<Scalar,4,1> Coefficients;
62 
63 
64 
66  EIGEN_DEVICE_FUNC inline CoeffReturnType x() const { return this->derived().coeffs().coeff(0); }
68  EIGEN_DEVICE_FUNC inline CoeffReturnType y() const { return this->derived().coeffs().coeff(1); }
70  EIGEN_DEVICE_FUNC inline CoeffReturnType z() const { return this->derived().coeffs().coeff(2); }
72  EIGEN_DEVICE_FUNC inline CoeffReturnType w() const { return this->derived().coeffs().coeff(3); }
73 
75  EIGEN_DEVICE_FUNC inline NonConstCoeffReturnType x() { return this->derived().coeffs().x(); }
77  EIGEN_DEVICE_FUNC inline NonConstCoeffReturnType y() { return this->derived().coeffs().y(); }
79  EIGEN_DEVICE_FUNC inline NonConstCoeffReturnType z() { return this->derived().coeffs().z(); }
81  EIGEN_DEVICE_FUNC inline NonConstCoeffReturnType w() { return this->derived().coeffs().w(); }
82 
84  EIGEN_DEVICE_FUNC inline const VectorBlock<const Coefficients,3> vec() const { return coeffs().template head<3>(); }
85 
87  EIGEN_DEVICE_FUNC inline VectorBlock<Coefficients,3> vec() { return coeffs().template head<3>(); }
88 
90  EIGEN_DEVICE_FUNC inline const typename internal::traits<Derived>::Coefficients& coeffs() const { return derived().coeffs(); }
91 
93  EIGEN_DEVICE_FUNC inline typename internal::traits<Derived>::Coefficients& coeffs() { return derived().coeffs(); }
94 
95  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE QuaternionBase<Derived>& operator=(const QuaternionBase<Derived>& other);
96  template<class OtherDerived> EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& operator=(const QuaternionBase<OtherDerived>& other);
97 
98 // disabled this copy operator as it is giving very strange compilation errors when compiling
99 // test_stdvector with GCC 4.4.2. This looks like a GCC bug though, so feel free to re-enable it if it's
100 // useful; however notice that we already have the templated operator= above and e.g. in MatrixBase
101 // we didn't have to add, in addition to templated operator=, such a non-templated copy operator.
102 // Derived& operator=(const QuaternionBase& other)
103 // { return operator=<Derived>(other); }
104 
105  EIGEN_DEVICE_FUNC Derived& operator=(const AngleAxisType& aa);
106  template<class OtherDerived> EIGEN_DEVICE_FUNC Derived& operator=(const MatrixBase<OtherDerived>& m);
107 
111  EIGEN_DEVICE_FUNC static inline Quaternion<Scalar> Identity() { return Quaternion<Scalar>(Scalar(1), Scalar(0), Scalar(0), Scalar(0)); }
112 
115  EIGEN_DEVICE_FUNC inline QuaternionBase& setIdentity() { coeffs() << Scalar(0), Scalar(0), Scalar(0), Scalar(1); return *this; }
116 
120  EIGEN_DEVICE_FUNC inline Scalar squaredNorm() const { return coeffs().squaredNorm(); }
121 
125  EIGEN_DEVICE_FUNC inline Scalar norm() const { return coeffs().norm(); }
126 
129  EIGEN_DEVICE_FUNC inline void normalize() { coeffs().normalize(); }
132  EIGEN_DEVICE_FUNC inline Quaternion<Scalar> normalized() const { return Quaternion<Scalar>(coeffs().normalized()); }
133 
139  template<class OtherDerived> EIGEN_DEVICE_FUNC inline Scalar dot(const QuaternionBase<OtherDerived>& other) const { return coeffs().dot(other.coeffs()); }
140 
141  template<class OtherDerived> EIGEN_DEVICE_FUNC Scalar angularDistance(const QuaternionBase<OtherDerived>& other) const;
142 
144  EIGEN_DEVICE_FUNC Matrix3 toRotationMatrix() const;
145 
147  template<typename Derived1, typename Derived2>
148  EIGEN_DEVICE_FUNC Derived& setFromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
149 
150  template<class OtherDerived> EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion<Scalar> operator* (const QuaternionBase<OtherDerived>& q) const;
151  template<class OtherDerived> EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& operator*= (const QuaternionBase<OtherDerived>& q);
152 
154  EIGEN_DEVICE_FUNC Quaternion<Scalar> inverse() const;
155 
157  EIGEN_DEVICE_FUNC Quaternion<Scalar> conjugate() const;
158 
159  template<class OtherDerived> EIGEN_DEVICE_FUNC Quaternion<Scalar> slerp(const Scalar& t, const QuaternionBase<OtherDerived>& other) const;
160 
165  template<class OtherDerived>
166  EIGEN_DEVICE_FUNC bool isApprox(const QuaternionBase<OtherDerived>& other, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const
167  { return coeffs().isApprox(other.coeffs(), prec); }
168 
170  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Vector3 _transformVector(const Vector3& v) const;
171 
172  #ifdef EIGEN_PARSED_BY_DOXYGEN
173 
178  template<typename NewScalarType>
179  EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type cast() const;
180 
181  #else
182 
183  template<typename NewScalarType>
184  EIGEN_DEVICE_FUNC inline
185  typename internal::enable_if<internal::is_same<Scalar,NewScalarType>::value,const Derived&>::type cast() const
186  {
187  return derived();
188  }
189 
190  template<typename NewScalarType>
191  EIGEN_DEVICE_FUNC inline
192  typename internal::enable_if<!internal::is_same<Scalar,NewScalarType>::value,Quaternion<NewScalarType> >::type cast() const
193  {
194  return Quaternion<NewScalarType>(coeffs().template cast<NewScalarType>());
195  }
196  #endif
197 
198 #ifdef EIGEN_QUATERNIONBASE_PLUGIN
199 # include EIGEN_QUATERNIONBASE_PLUGIN
200 #endif
201 protected:
202  EIGEN_DEFAULT_COPY_CONSTRUCTOR(QuaternionBase)
203  EIGEN_DEFAULT_EMPTY_CONSTRUCTOR_AND_DESTRUCTOR(QuaternionBase)
204 };
205 
206 /***************************************************************************
207 * Definition/implementation of Quaternion<Scalar>
208 ***************************************************************************/
209 
235 namespace internal {
236 template<typename _Scalar,int _Options>
237 struct traits<Quaternion<_Scalar,_Options> >
238 {
239  typedef Quaternion<_Scalar,_Options> PlainObject;
240  typedef _Scalar Scalar;
241  typedef Matrix<_Scalar,4,1,_Options> Coefficients;
242  enum{
243  Alignment = internal::traits<Coefficients>::Alignment,
244  Flags = LvalueBit
245  };
246 };
247 }
248 
249 template<typename _Scalar, int _Options>
250 class Quaternion : public QuaternionBase<Quaternion<_Scalar,_Options> >
251 {
252 public:
254  enum { NeedsAlignment = internal::traits<Quaternion>::Alignment>0 };
255 
256  typedef _Scalar Scalar;
257 
258  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Quaternion)
259  using Base::operator*=;
260 
261  typedef typename internal::traits<Quaternion>::Coefficients Coefficients;
262  typedef typename Base::AngleAxisType AngleAxisType;
263 
265  EIGEN_DEVICE_FUNC inline Quaternion() {}
266 
274  EIGEN_DEVICE_FUNC inline Quaternion(const Scalar& w, const Scalar& x, const Scalar& y, const Scalar& z) : m_coeffs(x, y, z, w){}
275 
277  EIGEN_DEVICE_FUNC explicit inline Quaternion(const Scalar* data) : m_coeffs(data) {}
278 
280  template<class Derived> EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion(const QuaternionBase<Derived>& other) { this->Base::operator=(other); }
281 
283  EIGEN_DEVICE_FUNC explicit inline Quaternion(const AngleAxisType& aa) { *this = aa; }
284 
289  template<typename Derived>
290  EIGEN_DEVICE_FUNC explicit inline Quaternion(const MatrixBase<Derived>& other) { *this = other; }
291 
293  template<typename OtherScalar, int OtherOptions>
294  EIGEN_DEVICE_FUNC explicit inline Quaternion(const Quaternion<OtherScalar, OtherOptions>& other)
295  { m_coeffs = other.coeffs().template cast<Scalar>(); }
296 
297  EIGEN_DEVICE_FUNC static Quaternion UnitRandom();
298 
299  template<typename Derived1, typename Derived2>
300  EIGEN_DEVICE_FUNC static Quaternion FromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
301 
302  EIGEN_DEVICE_FUNC inline Coefficients& coeffs() { return m_coeffs;}
303  EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs;}
304 
305  EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(bool(NeedsAlignment))
306 
307 #ifdef EIGEN_QUATERNION_PLUGIN
308 # include EIGEN_QUATERNION_PLUGIN
309 #endif
310 
311 protected:
312  Coefficients m_coeffs;
313 
314 #ifndef EIGEN_PARSED_BY_DOXYGEN
315  static EIGEN_STRONG_INLINE void _check_template_params()
316  {
317  EIGEN_STATIC_ASSERT( (_Options & DontAlign) == _Options,
318  INVALID_MATRIX_TEMPLATE_PARAMETERS)
319  }
320 #endif
321 };
322 
329 
330 /***************************************************************************
331 * Specialization of Map<Quaternion<Scalar>>
332 ***************************************************************************/
333 
334 namespace internal {
335  template<typename _Scalar, int _Options>
336  struct traits<Map<Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
337  {
338  typedef Map<Matrix<_Scalar,4,1>, _Options> Coefficients;
339  };
340 }
341 
342 namespace internal {
343  template<typename _Scalar, int _Options>
344  struct traits<Map<const Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
345  {
346  typedef Map<const Matrix<_Scalar,4,1>, _Options> Coefficients;
347  typedef traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> > TraitsBase;
348  enum {
349  Flags = TraitsBase::Flags & ~LvalueBit
350  };
351  };
352 }
353 
365 template<typename _Scalar, int _Options>
366 class Map<const Quaternion<_Scalar>, _Options >
367  : public QuaternionBase<Map<const Quaternion<_Scalar>, _Options> >
368 {
369  public:
371 
372  typedef _Scalar Scalar;
373  typedef typename internal::traits<Map>::Coefficients Coefficients;
374  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
375  using Base::operator*=;
376 
383  EIGEN_DEVICE_FUNC explicit EIGEN_STRONG_INLINE Map(const Scalar* coeffs) : m_coeffs(coeffs) {}
384 
385  EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs;}
386 
387  protected:
388  const Coefficients m_coeffs;
389 };
390 
402 template<typename _Scalar, int _Options>
403 class Map<Quaternion<_Scalar>, _Options >
404  : public QuaternionBase<Map<Quaternion<_Scalar>, _Options> >
405 {
406  public:
407  typedef QuaternionBase<Map<Quaternion<_Scalar>, _Options> > Base;
408 
409  typedef _Scalar Scalar;
410  typedef typename internal::traits<Map>::Coefficients Coefficients;
411  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
412  using Base::operator*=;
413 
420  EIGEN_DEVICE_FUNC explicit EIGEN_STRONG_INLINE Map(Scalar* coeffs) : m_coeffs(coeffs) {}
421 
422  EIGEN_DEVICE_FUNC inline Coefficients& coeffs() { return m_coeffs; }
423  EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs; }
424 
425  protected:
426  Coefficients m_coeffs;
427 };
428 
441 
442 /***************************************************************************
443 * Implementation of QuaternionBase methods
444 ***************************************************************************/
445 
446 // Generic Quaternion * Quaternion product
447 // This product can be specialized for a given architecture via the Arch template argument.
448 namespace internal {
449 template<int Arch, class Derived1, class Derived2, typename Scalar> struct quat_product
450 {
451  EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived1>& a, const QuaternionBase<Derived2>& b){
452  return Quaternion<Scalar>
453  (
454  a.w() * b.w() - a.x() * b.x() - a.y() * b.y() - a.z() * b.z(),
455  a.w() * b.x() + a.x() * b.w() + a.y() * b.z() - a.z() * b.y(),
456  a.w() * b.y() + a.y() * b.w() + a.z() * b.x() - a.x() * b.z(),
457  a.w() * b.z() + a.z() * b.w() + a.x() * b.y() - a.y() * b.x()
458  );
459  }
460 };
461 }
462 
464 template <class Derived>
465 template <class OtherDerived>
466 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion<typename internal::traits<Derived>::Scalar>
468 {
469  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename OtherDerived::Scalar>::value),
470  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
471  return internal::quat_product<Architecture::Target, Derived, OtherDerived,
472  typename internal::traits<Derived>::Scalar>::run(*this, other);
473 }
474 
476 template <class Derived>
477 template <class OtherDerived>
478 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator*= (const QuaternionBase<OtherDerived>& other)
479 {
480  derived() = derived() * other.derived();
481  return derived();
482 }
483 
491 template <class Derived>
492 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE typename QuaternionBase<Derived>::Vector3
494 {
495  // Note that this algorithm comes from the optimization by hand
496  // of the conversion to a Matrix followed by a Matrix/Vector product.
497  // It appears to be much faster than the common algorithm found
498  // in the literature (30 versus 39 flops). It also requires two
499  // Vector3 as temporaries.
500  Vector3 uv = this->vec().cross(v);
501  uv += uv;
502  return v + this->w() * uv + this->vec().cross(uv);
503 }
504 
505 template<class Derived>
506 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE QuaternionBase<Derived>& QuaternionBase<Derived>::operator=(const QuaternionBase<Derived>& other)
507 {
508  coeffs() = other.coeffs();
509  return derived();
510 }
511 
512 template<class Derived>
513 template<class OtherDerived>
514 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const QuaternionBase<OtherDerived>& other)
515 {
516  coeffs() = other.coeffs();
517  return derived();
518 }
519 
522 template<class Derived>
523 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const AngleAxisType& aa)
524 {
525  EIGEN_USING_STD_MATH(cos)
526  EIGEN_USING_STD_MATH(sin)
527  Scalar ha = Scalar(0.5)*aa.angle(); // Scalar(0.5) to suppress precision loss warnings
528  this->w() = cos(ha);
529  this->vec() = sin(ha) * aa.axis();
530  return derived();
531 }
532 
539 template<class Derived>
540 template<class MatrixDerived>
541 EIGEN_DEVICE_FUNC inline Derived& QuaternionBase<Derived>::operator=(const MatrixBase<MatrixDerived>& xpr)
542 {
543  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename MatrixDerived::Scalar>::value),
544  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
545  internal::quaternionbase_assign_impl<MatrixDerived>::run(*this, xpr.derived());
546  return derived();
547 }
548 
552 template<class Derived>
553 EIGEN_DEVICE_FUNC inline typename QuaternionBase<Derived>::Matrix3
555 {
556  // NOTE if inlined, then gcc 4.2 and 4.4 get rid of the temporary (not gcc 4.3 !!)
557  // if not inlined then the cost of the return by value is huge ~ +35%,
558  // however, not inlining this function is an order of magnitude slower, so
559  // it has to be inlined, and so the return by value is not an issue
560  Matrix3 res;
561 
562  const Scalar tx = Scalar(2)*this->x();
563  const Scalar ty = Scalar(2)*this->y();
564  const Scalar tz = Scalar(2)*this->z();
565  const Scalar twx = tx*this->w();
566  const Scalar twy = ty*this->w();
567  const Scalar twz = tz*this->w();
568  const Scalar txx = tx*this->x();
569  const Scalar txy = ty*this->x();
570  const Scalar txz = tz*this->x();
571  const Scalar tyy = ty*this->y();
572  const Scalar tyz = tz*this->y();
573  const Scalar tzz = tz*this->z();
574 
575  res.coeffRef(0,0) = Scalar(1)-(tyy+tzz);
576  res.coeffRef(0,1) = txy-twz;
577  res.coeffRef(0,2) = txz+twy;
578  res.coeffRef(1,0) = txy+twz;
579  res.coeffRef(1,1) = Scalar(1)-(txx+tzz);
580  res.coeffRef(1,2) = tyz-twx;
581  res.coeffRef(2,0) = txz-twy;
582  res.coeffRef(2,1) = tyz+twx;
583  res.coeffRef(2,2) = Scalar(1)-(txx+tyy);
584 
585  return res;
586 }
587 
598 template<class Derived>
599 template<typename Derived1, typename Derived2>
600 EIGEN_DEVICE_FUNC inline Derived& QuaternionBase<Derived>::setFromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b)
601 {
602  EIGEN_USING_STD_MATH(sqrt)
603  Vector3 v0 = a.normalized();
604  Vector3 v1 = b.normalized();
605  Scalar c = v1.dot(v0);
606 
607  // if dot == -1, vectors are nearly opposites
608  // => accurately compute the rotation axis by computing the
609  // intersection of the two planes. This is done by solving:
610  // x^T v0 = 0
611  // x^T v1 = 0
612  // under the constraint:
613  // ||x|| = 1
614  // which yields a singular value problem
615  if (c < Scalar(-1)+NumTraits<Scalar>::dummy_precision())
616  {
617  c = numext::maxi(c,Scalar(-1));
618  Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
620  Vector3 axis = svd.matrixV().col(2);
621 
622  Scalar w2 = (Scalar(1)+c)*Scalar(0.5);
623  this->w() = sqrt(w2);
624  this->vec() = axis * sqrt(Scalar(1) - w2);
625  return derived();
626  }
627  Vector3 axis = v0.cross(v1);
628  Scalar s = sqrt((Scalar(1)+c)*Scalar(2));
629  Scalar invs = Scalar(1)/s;
630  this->vec() = axis * invs;
631  this->w() = s * Scalar(0.5);
632 
633  return derived();
634 }
635 
640 template<typename Scalar, int Options>
642 {
643  EIGEN_USING_STD_MATH(sqrt)
644  EIGEN_USING_STD_MATH(sin)
645  EIGEN_USING_STD_MATH(cos)
646  const Scalar u1 = internal::random<Scalar>(0, 1),
647  u2 = internal::random<Scalar>(0, 2*EIGEN_PI),
648  u3 = internal::random<Scalar>(0, 2*EIGEN_PI);
649  const Scalar a = sqrt(1 - u1),
650  b = sqrt(u1);
651  return Quaternion (a * sin(u2), a * cos(u2), b * sin(u3), b * cos(u3));
652 }
653 
654 
665 template<typename Scalar, int Options>
666 template<typename Derived1, typename Derived2>
668 {
669  Quaternion quat;
670  quat.setFromTwoVectors(a, b);
671  return quat;
672 }
673 
674 
681 template <class Derived>
683 {
684  // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite() ??
685  Scalar n2 = this->squaredNorm();
686  if (n2 > Scalar(0))
687  return Quaternion<Scalar>(conjugate().coeffs() / n2);
688  else
689  {
690  // return an invalid result to flag the error
691  return Quaternion<Scalar>(Coefficients::Zero());
692  }
693 }
694 
695 // Generic conjugate of a Quaternion
696 namespace internal {
697 template<int Arch, class Derived, typename Scalar> struct quat_conj
698 {
699  EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived>& q){
700  return Quaternion<Scalar>(q.w(),-q.x(),-q.y(),-q.z());
701  }
702 };
703 }
704 
711 template <class Derived>
712 EIGEN_DEVICE_FUNC inline Quaternion<typename internal::traits<Derived>::Scalar>
714 {
715  return internal::quat_conj<Architecture::Target, Derived,
716  typename internal::traits<Derived>::Scalar>::run(*this);
717 
718 }
719 
723 template <class Derived>
724 template <class OtherDerived>
725 EIGEN_DEVICE_FUNC inline typename internal::traits<Derived>::Scalar
727 {
728  EIGEN_USING_STD_MATH(atan2)
729  Quaternion<Scalar> d = (*this) * other.conjugate();
730  return Scalar(2) * atan2( d.vec().norm(), numext::abs(d.w()) );
731 }
732 
733 
734 
741 template <class Derived>
742 template <class OtherDerived>
745 {
746  EIGEN_USING_STD_MATH(acos)
747  EIGEN_USING_STD_MATH(sin)
748  const Scalar one = Scalar(1) - NumTraits<Scalar>::epsilon();
749  Scalar d = this->dot(other);
750  Scalar absD = numext::abs(d);
751 
752  Scalar scale0;
753  Scalar scale1;
754 
755  if(absD>=one)
756  {
757  scale0 = Scalar(1) - t;
758  scale1 = t;
759  }
760  else
761  {
762  // theta is the angle between the 2 quaternions
763  Scalar theta = acos(absD);
764  Scalar sinTheta = sin(theta);
765 
766  scale0 = sin( ( Scalar(1) - t ) * theta) / sinTheta;
767  scale1 = sin( ( t * theta) ) / sinTheta;
768  }
769  if(d<Scalar(0)) scale1 = -scale1;
770 
771  return Quaternion<Scalar>(scale0 * coeffs() + scale1 * other.coeffs());
772 }
773 
774 namespace internal {
775 
776 // set from a rotation matrix
777 template<typename Other>
778 struct quaternionbase_assign_impl<Other,3,3>
779 {
780  typedef typename Other::Scalar Scalar;
781  template<class Derived> EIGEN_DEVICE_FUNC static inline void run(QuaternionBase<Derived>& q, const Other& a_mat)
782  {
783  const typename internal::nested_eval<Other,2>::type mat(a_mat);
784  EIGEN_USING_STD_MATH(sqrt)
785  // This algorithm comes from "Quaternion Calculus and Fast Animation",
786  // Ken Shoemake, 1987 SIGGRAPH course notes
787  Scalar t = mat.trace();
788  if (t > Scalar(0))
789  {
790  t = sqrt(t + Scalar(1.0));
791  q.w() = Scalar(0.5)*t;
792  t = Scalar(0.5)/t;
793  q.x() = (mat.coeff(2,1) - mat.coeff(1,2)) * t;
794  q.y() = (mat.coeff(0,2) - mat.coeff(2,0)) * t;
795  q.z() = (mat.coeff(1,0) - mat.coeff(0,1)) * t;
796  }
797  else
798  {
799  Index i = 0;
800  if (mat.coeff(1,1) > mat.coeff(0,0))
801  i = 1;
802  if (mat.coeff(2,2) > mat.coeff(i,i))
803  i = 2;
804  Index j = (i+1)%3;
805  Index k = (j+1)%3;
806 
807  t = sqrt(mat.coeff(i,i)-mat.coeff(j,j)-mat.coeff(k,k) + Scalar(1.0));
808  q.coeffs().coeffRef(i) = Scalar(0.5) * t;
809  t = Scalar(0.5)/t;
810  q.w() = (mat.coeff(k,j)-mat.coeff(j,k))*t;
811  q.coeffs().coeffRef(j) = (mat.coeff(j,i)+mat.coeff(i,j))*t;
812  q.coeffs().coeffRef(k) = (mat.coeff(k,i)+mat.coeff(i,k))*t;
813  }
814  }
815 };
816 
817 // set from a vector of coefficients assumed to be a quaternion
818 template<typename Other>
819 struct quaternionbase_assign_impl<Other,4,1>
820 {
821  typedef typename Other::Scalar Scalar;
822  template<class Derived> EIGEN_DEVICE_FUNC static inline void run(QuaternionBase<Derived>& q, const Other& vec)
823  {
824  q.coeffs() = vec;
825  }
826 };
827 
828 } // end namespace internal
829 
830 } // end namespace Eigen
831 
832 #endif // EIGEN_QUATERNION_H
Eigen::QuaternionBase::Vector3
Matrix< Scalar, 3, 1 > Vector3
Definition: Quaternion.h:57
Eigen::QuaternionBase::x
NonConstCoeffReturnType x()
Definition: Quaternion.h:75
Eigen
Namespace containing all symbols from the Eigen library.
Definition: Core:309
Eigen::QuaternionBase::dot
Scalar dot(const QuaternionBase< OtherDerived > &other) const
Definition: Quaternion.h:139
Eigen::QuaternionBase::w
CoeffReturnType w() const
Definition: Quaternion.h:72
Eigen::ComputeFullV
@ ComputeFullV
Definition: Constants.h:387
Eigen::sqrt
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Eigen::QuaternionMapAlignedf
Map< Quaternion< float >, Aligned > QuaternionMapAlignedf
Definition: Quaternion.h:437
Eigen::QuaternionBase::operator=
Derived & operator=(const AngleAxisType &aa)
Definition: Quaternion.h:523
Eigen::DenseCoeffsBase< Derived, DirectWriteAccessors >::derived
Derived & derived()
Definition: EigenBase.h:45
Eigen::sin
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sin_op< typename Derived::Scalar >, const Derived > sin(const Eigen::ArrayBase< Derived > &x)
Eigen::AngleAxis
Represents a 3D rotation as a rotation angle around an arbitrary 3D axis.
Definition: AngleAxis.h:50
Eigen::QuaternionMapd
Map< Quaternion< double >, 0 > QuaternionMapd
Definition: Quaternion.h:434
Eigen::QuaternionBase::norm
Scalar norm() const
Definition: Quaternion.h:125
Eigen::QuaternionBase::cast
internal::cast_return_type< Derived, Quaternion< NewScalarType > >::type cast() const
Eigen::Quaternion::Quaternion
Quaternion(const AngleAxisType &aa)
Definition: Quaternion.h:283
Eigen::QuaternionBase::Identity
static Quaternion< Scalar > Identity()
Definition: Quaternion.h:111
Eigen::QuaternionBase::z
NonConstCoeffReturnType z()
Definition: Quaternion.h:79
Eigen::Quaternion::Quaternion
Quaternion(const MatrixBase< Derived > &other)
Definition: Quaternion.h:290
Eigen::AngleAxis::angle
Scalar angle() const
Definition: AngleAxis.h:91
Eigen::QuaternionBase::y
NonConstCoeffReturnType y()
Definition: Quaternion.h:77
Eigen::DontAlign
@ DontAlign
Definition: Constants.h:326
Eigen::Quaternion::Quaternion
Quaternion(const Quaternion< OtherScalar, OtherOptions > &other)
Definition: Quaternion.h:294
Eigen::QuaternionBase::isApprox
bool isApprox(const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Quaternion.h:166
Eigen::Quaternion::UnitRandom
static Quaternion UnitRandom()
Definition: Quaternion.h:641
Eigen::QuaternionBase::z
CoeffReturnType z() const
Definition: Quaternion.h:70
Eigen::QuaternionBase::operator*=
Derived & operator*=(const QuaternionBase< OtherDerived > &q)
Definition: Quaternion.h:478
Eigen::QuaternionBase::squaredNorm
Scalar squaredNorm() const
Definition: Quaternion.h:120
Eigen::QuaternionBase::vec
const VectorBlock< const Coefficients, 3 > vec() const
Definition: Quaternion.h:84
Eigen::Quaternion::Quaternion
Quaternion(const Scalar &w, const Scalar &x, const Scalar &y, const Scalar &z)
Definition: Quaternion.h:274
Eigen::QuaternionBase::normalized
Quaternion< Scalar > normalized() const
Definition: Quaternion.h:132
Eigen::QuaternionBase::AngleAxisType
AngleAxis< Scalar > AngleAxisType
Definition: Quaternion.h:61
Eigen::LvalueBit
const unsigned int LvalueBit
Definition: Constants.h:139
Eigen::QuaternionBase::w
NonConstCoeffReturnType w()
Definition: Quaternion.h:81
Eigen::acos
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_acos_op< typename Derived::Scalar >, const Derived > acos(const Eigen::ArrayBase< Derived > &x)
Eigen::RotationBase
Common base class for compact rotation representations.
Definition: RotationBase.h:30
Eigen::QuaternionBase::setFromTwoVectors
Derived & setFromTwoVectors(const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
Definition: Quaternion.h:600
Eigen::Map< Quaternion< _Scalar >, _Options >::Map
Map(Scalar *coeffs)
Definition: Quaternion.h:420
Eigen::QuaternionBase::normalize
void normalize()
Definition: Quaternion.h:129
Eigen::Quaternionf
Quaternion< float > Quaternionf
Definition: Quaternion.h:325
Eigen::PlainObjectBase::coeffRef
Scalar & coeffRef(Index rowId, Index colId)
Definition: PlainObjectBase.h:183
Eigen::Quaterniond
Quaternion< double > Quaterniond
Definition: Quaternion.h:328
Eigen::QuaternionBase::inverse
Quaternion< Scalar > inverse() const
Definition: Quaternion.h:682
Eigen::AutoAlign
@ AutoAlign
Definition: Constants.h:324
Eigen::Map
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
Eigen::QuaternionBase::coeffs
internal::traits< Derived >::Coefficients & coeffs()
Definition: Quaternion.h:93
Eigen::QuaternionBase::operator*
Transform< Scalar, Dim, Isometry > operator*(const Translation< Scalar, Dim > &t) const
Definition: RotationBase.h:56
Eigen::cos
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_cos_op< typename Derived::Scalar >, const Derived > cos(const Eigen::ArrayBase< Derived > &x)
Eigen::JacobiSVD
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: JacobiSVD.h:489
Eigen::Quaternion
The quaternion class used to represent 3D orientations and rotations.
Definition: Quaternion.h:251
Eigen::QuaternionBase::y
CoeffReturnType y() const
Definition: Quaternion.h:68
Eigen::QuaternionMapf
Map< Quaternion< float >, 0 > QuaternionMapf
Definition: Quaternion.h:431
Eigen::QuaternionBase::_transformVector
Vector3 _transformVector(const Vector3 &v) const
Definition: Quaternion.h:493
Eigen::QuaternionBase::x
CoeffReturnType x() const
Definition: Quaternion.h:66
Eigen::QuaternionBase::setIdentity
QuaternionBase & setIdentity()
Definition: Quaternion.h:115
Eigen::QuaternionBase::conjugate
Quaternion< Scalar > conjugate() const
Definition: Quaternion.h:713
Eigen::QuaternionBase::toRotationMatrix
Matrix3 toRotationMatrix() const
Definition: Quaternion.h:554
Eigen::VectorBlock
Expression of a fixed-size or dynamic-size sub-vector.
Definition: VectorBlock.h:60
Eigen::QuaternionBase::Matrix3
Matrix< Scalar, 3, 3 > Matrix3
Definition: Quaternion.h:59
Eigen::Quaternion::Quaternion
Quaternion()
Definition: Quaternion.h:265
Eigen::QuaternionBase::vec
VectorBlock< Coefficients, 3 > vec()
Definition: Quaternion.h:87
Eigen::Matrix< Scalar, 3, 1 >
Eigen::MatrixBase
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:50
Eigen::Quaternion::Quaternion
Quaternion(const QuaternionBase< Derived > &other)
Definition: Quaternion.h:280
Eigen::QuaternionBase::coeffs
const internal::traits< Derived >::Coefficients & coeffs() const
Definition: Quaternion.h:90
Eigen::Map< const Quaternion< _Scalar >, _Options >::Map
Map(const Scalar *coeffs)
Definition: Quaternion.h:383
Eigen::QuaternionMapAlignedd
Map< Quaternion< double >, Aligned > QuaternionMapAlignedd
Definition: Quaternion.h:440
Eigen::QuaternionBase
Base class for quaternion expressions.
Definition: Quaternion.h:36
Eigen::NumTraits
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:151
Eigen::AngleAxis::axis
const Vector3 & axis() const
Definition: AngleAxis.h:96
Eigen::Index
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
Eigen::Aligned
@ Aligned
Definition: Constants.h:235
Eigen::Quaternion::Quaternion
Quaternion(const Scalar *data)
Definition: Quaternion.h:277
Eigen::MatrixBase::normalized
const PlainObject normalized() const
Definition: Dot.h:124