Eigen  3.3.7
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:
38  typedef RotationBase<Derived, 3> Base;
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 
177  template<typename NewScalarType>
178  EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type cast() const
179  {
180  return typename internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type(derived());
181  }
182 
183 #ifdef EIGEN_QUATERNIONBASE_PLUGIN
184 # include EIGEN_QUATERNIONBASE_PLUGIN
185 #endif
186 protected:
187  EIGEN_DEFAULT_COPY_CONSTRUCTOR(QuaternionBase)
188  EIGEN_DEFAULT_EMPTY_CONSTRUCTOR_AND_DESTRUCTOR(QuaternionBase)
189 };
190 
191 /***************************************************************************
192 * Definition/implementation of Quaternion<Scalar>
193 ***************************************************************************/
194 
220 namespace internal {
221 template<typename _Scalar,int _Options>
222 struct traits<Quaternion<_Scalar,_Options> >
223 {
224  typedef Quaternion<_Scalar,_Options> PlainObject;
225  typedef _Scalar Scalar;
226  typedef Matrix<_Scalar,4,1,_Options> Coefficients;
227  enum{
228  Alignment = internal::traits<Coefficients>::Alignment,
229  Flags = LvalueBit
230  };
231 };
232 }
233 
234 template<typename _Scalar, int _Options>
235 class Quaternion : public QuaternionBase<Quaternion<_Scalar,_Options> >
236 {
237 public:
238  typedef QuaternionBase<Quaternion<_Scalar,_Options> > Base;
239  enum { NeedsAlignment = internal::traits<Quaternion>::Alignment>0 };
240 
241  typedef _Scalar Scalar;
242 
243  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Quaternion)
244  using Base::operator*=;
245 
246  typedef typename internal::traits<Quaternion>::Coefficients Coefficients;
247  typedef typename Base::AngleAxisType AngleAxisType;
248 
250  EIGEN_DEVICE_FUNC inline Quaternion() {}
251 
259  EIGEN_DEVICE_FUNC inline Quaternion(const Scalar& w, const Scalar& x, const Scalar& y, const Scalar& z) : m_coeffs(x, y, z, w){}
260 
262  EIGEN_DEVICE_FUNC explicit inline Quaternion(const Scalar* data) : m_coeffs(data) {}
263 
265  template<class Derived> EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion(const QuaternionBase<Derived>& other) { this->Base::operator=(other); }
266 
268  EIGEN_DEVICE_FUNC explicit inline Quaternion(const AngleAxisType& aa) { *this = aa; }
269 
274  template<typename Derived>
275  EIGEN_DEVICE_FUNC explicit inline Quaternion(const MatrixBase<Derived>& other) { *this = other; }
276 
278  template<typename OtherScalar, int OtherOptions>
279  EIGEN_DEVICE_FUNC explicit inline Quaternion(const Quaternion<OtherScalar, OtherOptions>& other)
280  { m_coeffs = other.coeffs().template cast<Scalar>(); }
281 
282  EIGEN_DEVICE_FUNC static Quaternion UnitRandom();
283 
284  template<typename Derived1, typename Derived2>
285  EIGEN_DEVICE_FUNC static Quaternion FromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
286 
287  EIGEN_DEVICE_FUNC inline Coefficients& coeffs() { return m_coeffs;}
288  EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs;}
289 
290  EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(bool(NeedsAlignment))
291 
292 #ifdef EIGEN_QUATERNION_PLUGIN
293 # include EIGEN_QUATERNION_PLUGIN
294 #endif
295 
296 protected:
297  Coefficients m_coeffs;
298 
299 #ifndef EIGEN_PARSED_BY_DOXYGEN
300  static EIGEN_STRONG_INLINE void _check_template_params()
301  {
302  EIGEN_STATIC_ASSERT( (_Options & DontAlign) == _Options,
303  INVALID_MATRIX_TEMPLATE_PARAMETERS)
304  }
305 #endif
306 };
307 
314 
315 /***************************************************************************
316 * Specialization of Map<Quaternion<Scalar>>
317 ***************************************************************************/
318 
319 namespace internal {
320  template<typename _Scalar, int _Options>
321  struct traits<Map<Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
322  {
323  typedef Map<Matrix<_Scalar,4,1>, _Options> Coefficients;
324  };
325 }
326 
327 namespace internal {
328  template<typename _Scalar, int _Options>
329  struct traits<Map<const Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
330  {
331  typedef Map<const Matrix<_Scalar,4,1>, _Options> Coefficients;
332  typedef traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> > TraitsBase;
333  enum {
334  Flags = TraitsBase::Flags & ~LvalueBit
335  };
336  };
337 }
338 
350 template<typename _Scalar, int _Options>
351 class Map<const Quaternion<_Scalar>, _Options >
352  : public QuaternionBase<Map<const Quaternion<_Scalar>, _Options> >
353 {
354  public:
356 
357  typedef _Scalar Scalar;
358  typedef typename internal::traits<Map>::Coefficients Coefficients;
359  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
360  using Base::operator*=;
361 
368  EIGEN_DEVICE_FUNC explicit EIGEN_STRONG_INLINE Map(const Scalar* coeffs) : m_coeffs(coeffs) {}
369 
370  EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs;}
371 
372  protected:
373  const Coefficients m_coeffs;
374 };
375 
387 template<typename _Scalar, int _Options>
388 class Map<Quaternion<_Scalar>, _Options >
389  : public QuaternionBase<Map<Quaternion<_Scalar>, _Options> >
390 {
391  public:
392  typedef QuaternionBase<Map<Quaternion<_Scalar>, _Options> > Base;
393 
394  typedef _Scalar Scalar;
395  typedef typename internal::traits<Map>::Coefficients Coefficients;
396  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
397  using Base::operator*=;
398 
405  EIGEN_DEVICE_FUNC explicit EIGEN_STRONG_INLINE Map(Scalar* coeffs) : m_coeffs(coeffs) {}
406 
407  EIGEN_DEVICE_FUNC inline Coefficients& coeffs() { return m_coeffs; }
408  EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs; }
409 
410  protected:
411  Coefficients m_coeffs;
412 };
413 
426 
427 /***************************************************************************
428 * Implementation of QuaternionBase methods
429 ***************************************************************************/
430 
431 // Generic Quaternion * Quaternion product
432 // This product can be specialized for a given architecture via the Arch template argument.
433 namespace internal {
434 template<int Arch, class Derived1, class Derived2, typename Scalar> struct quat_product
435 {
436  EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived1>& a, const QuaternionBase<Derived2>& b){
437  return Quaternion<Scalar>
438  (
439  a.w() * b.w() - a.x() * b.x() - a.y() * b.y() - a.z() * b.z(),
440  a.w() * b.x() + a.x() * b.w() + a.y() * b.z() - a.z() * b.y(),
441  a.w() * b.y() + a.y() * b.w() + a.z() * b.x() - a.x() * b.z(),
442  a.w() * b.z() + a.z() * b.w() + a.x() * b.y() - a.y() * b.x()
443  );
444  }
445 };
446 }
447 
449 template <class Derived>
450 template <class OtherDerived>
451 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion<typename internal::traits<Derived>::Scalar>
453 {
454  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename OtherDerived::Scalar>::value),
455  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
456  return internal::quat_product<Architecture::Target, Derived, OtherDerived,
457  typename internal::traits<Derived>::Scalar>::run(*this, other);
458 }
459 
461 template <class Derived>
462 template <class OtherDerived>
463 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator*= (const QuaternionBase<OtherDerived>& other)
464 {
465  derived() = derived() * other.derived();
466  return derived();
467 }
468 
476 template <class Derived>
477 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE typename QuaternionBase<Derived>::Vector3
479 {
480  // Note that this algorithm comes from the optimization by hand
481  // of the conversion to a Matrix followed by a Matrix/Vector product.
482  // It appears to be much faster than the common algorithm found
483  // in the literature (30 versus 39 flops). It also requires two
484  // Vector3 as temporaries.
485  Vector3 uv = this->vec().cross(v);
486  uv += uv;
487  return v + this->w() * uv + this->vec().cross(uv);
488 }
489 
490 template<class Derived>
491 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE QuaternionBase<Derived>& QuaternionBase<Derived>::operator=(const QuaternionBase<Derived>& other)
492 {
493  coeffs() = other.coeffs();
494  return derived();
495 }
496 
497 template<class Derived>
498 template<class OtherDerived>
499 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const QuaternionBase<OtherDerived>& other)
500 {
501  coeffs() = other.coeffs();
502  return derived();
503 }
504 
507 template<class Derived>
508 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const AngleAxisType& aa)
509 {
510  EIGEN_USING_STD_MATH(cos)
511  EIGEN_USING_STD_MATH(sin)
512  Scalar ha = Scalar(0.5)*aa.angle(); // Scalar(0.5) to suppress precision loss warnings
513  this->w() = cos(ha);
514  this->vec() = sin(ha) * aa.axis();
515  return derived();
516 }
517 
524 template<class Derived>
525 template<class MatrixDerived>
526 EIGEN_DEVICE_FUNC inline Derived& QuaternionBase<Derived>::operator=(const MatrixBase<MatrixDerived>& xpr)
527 {
528  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename MatrixDerived::Scalar>::value),
529  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
530  internal::quaternionbase_assign_impl<MatrixDerived>::run(*this, xpr.derived());
531  return derived();
532 }
533 
537 template<class Derived>
538 EIGEN_DEVICE_FUNC inline typename QuaternionBase<Derived>::Matrix3
540 {
541  // NOTE if inlined, then gcc 4.2 and 4.4 get rid of the temporary (not gcc 4.3 !!)
542  // if not inlined then the cost of the return by value is huge ~ +35%,
543  // however, not inlining this function is an order of magnitude slower, so
544  // it has to be inlined, and so the return by value is not an issue
545  Matrix3 res;
546 
547  const Scalar tx = Scalar(2)*this->x();
548  const Scalar ty = Scalar(2)*this->y();
549  const Scalar tz = Scalar(2)*this->z();
550  const Scalar twx = tx*this->w();
551  const Scalar twy = ty*this->w();
552  const Scalar twz = tz*this->w();
553  const Scalar txx = tx*this->x();
554  const Scalar txy = ty*this->x();
555  const Scalar txz = tz*this->x();
556  const Scalar tyy = ty*this->y();
557  const Scalar tyz = tz*this->y();
558  const Scalar tzz = tz*this->z();
559 
560  res.coeffRef(0,0) = Scalar(1)-(tyy+tzz);
561  res.coeffRef(0,1) = txy-twz;
562  res.coeffRef(0,2) = txz+twy;
563  res.coeffRef(1,0) = txy+twz;
564  res.coeffRef(1,1) = Scalar(1)-(txx+tzz);
565  res.coeffRef(1,2) = tyz-twx;
566  res.coeffRef(2,0) = txz-twy;
567  res.coeffRef(2,1) = tyz+twx;
568  res.coeffRef(2,2) = Scalar(1)-(txx+tyy);
569 
570  return res;
571 }
572 
583 template<class Derived>
584 template<typename Derived1, typename Derived2>
585 EIGEN_DEVICE_FUNC inline Derived& QuaternionBase<Derived>::setFromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b)
586 {
587  EIGEN_USING_STD_MATH(sqrt)
588  Vector3 v0 = a.normalized();
589  Vector3 v1 = b.normalized();
590  Scalar c = v1.dot(v0);
591 
592  // if dot == -1, vectors are nearly opposites
593  // => accurately compute the rotation axis by computing the
594  // intersection of the two planes. This is done by solving:
595  // x^T v0 = 0
596  // x^T v1 = 0
597  // under the constraint:
598  // ||x|| = 1
599  // which yields a singular value problem
600  if (c < Scalar(-1)+NumTraits<Scalar>::dummy_precision())
601  {
602  c = numext::maxi(c,Scalar(-1));
603  Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
605  Vector3 axis = svd.matrixV().col(2);
606 
607  Scalar w2 = (Scalar(1)+c)*Scalar(0.5);
608  this->w() = sqrt(w2);
609  this->vec() = axis * sqrt(Scalar(1) - w2);
610  return derived();
611  }
612  Vector3 axis = v0.cross(v1);
613  Scalar s = sqrt((Scalar(1)+c)*Scalar(2));
614  Scalar invs = Scalar(1)/s;
615  this->vec() = axis * invs;
616  this->w() = s * Scalar(0.5);
617 
618  return derived();
619 }
620 
625 template<typename Scalar, int Options>
627 {
628  EIGEN_USING_STD_MATH(sqrt)
629  EIGEN_USING_STD_MATH(sin)
630  EIGEN_USING_STD_MATH(cos)
631  const Scalar u1 = internal::random<Scalar>(0, 1),
632  u2 = internal::random<Scalar>(0, 2*EIGEN_PI),
633  u3 = internal::random<Scalar>(0, 2*EIGEN_PI);
634  const Scalar a = sqrt(1 - u1),
635  b = sqrt(u1);
636  return Quaternion (a * sin(u2), a * cos(u2), b * sin(u3), b * cos(u3));
637 }
638 
639 
650 template<typename Scalar, int Options>
651 template<typename Derived1, typename Derived2>
653 {
654  Quaternion quat;
655  quat.setFromTwoVectors(a, b);
656  return quat;
657 }
658 
659 
666 template <class Derived>
668 {
669  // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite() ??
670  Scalar n2 = this->squaredNorm();
671  if (n2 > Scalar(0))
672  return Quaternion<Scalar>(conjugate().coeffs() / n2);
673  else
674  {
675  // return an invalid result to flag the error
676  return Quaternion<Scalar>(Coefficients::Zero());
677  }
678 }
679 
680 // Generic conjugate of a Quaternion
681 namespace internal {
682 template<int Arch, class Derived, typename Scalar> struct quat_conj
683 {
684  EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived>& q){
685  return Quaternion<Scalar>(q.w(),-q.x(),-q.y(),-q.z());
686  }
687 };
688 }
689 
696 template <class Derived>
697 EIGEN_DEVICE_FUNC inline Quaternion<typename internal::traits<Derived>::Scalar>
699 {
700  return internal::quat_conj<Architecture::Target, Derived,
701  typename internal::traits<Derived>::Scalar>::run(*this);
702 
703 }
704 
708 template <class Derived>
709 template <class OtherDerived>
710 EIGEN_DEVICE_FUNC inline typename internal::traits<Derived>::Scalar
712 {
713  EIGEN_USING_STD_MATH(atan2)
714  Quaternion<Scalar> d = (*this) * other.conjugate();
715  return Scalar(2) * atan2( d.vec().norm(), numext::abs(d.w()) );
716 }
717 
718 
719 
726 template <class Derived>
727 template <class OtherDerived>
730 {
731  EIGEN_USING_STD_MATH(acos)
732  EIGEN_USING_STD_MATH(sin)
733  const Scalar one = Scalar(1) - NumTraits<Scalar>::epsilon();
734  Scalar d = this->dot(other);
735  Scalar absD = numext::abs(d);
736 
737  Scalar scale0;
738  Scalar scale1;
739 
740  if(absD>=one)
741  {
742  scale0 = Scalar(1) - t;
743  scale1 = t;
744  }
745  else
746  {
747  // theta is the angle between the 2 quaternions
748  Scalar theta = acos(absD);
749  Scalar sinTheta = sin(theta);
750 
751  scale0 = sin( ( Scalar(1) - t ) * theta) / sinTheta;
752  scale1 = sin( ( t * theta) ) / sinTheta;
753  }
754  if(d<Scalar(0)) scale1 = -scale1;
755 
756  return Quaternion<Scalar>(scale0 * coeffs() + scale1 * other.coeffs());
757 }
758 
759 namespace internal {
760 
761 // set from a rotation matrix
762 template<typename Other>
763 struct quaternionbase_assign_impl<Other,3,3>
764 {
765  typedef typename Other::Scalar Scalar;
766  template<class Derived> EIGEN_DEVICE_FUNC static inline void run(QuaternionBase<Derived>& q, const Other& a_mat)
767  {
768  const typename internal::nested_eval<Other,2>::type mat(a_mat);
769  EIGEN_USING_STD_MATH(sqrt)
770  // This algorithm comes from "Quaternion Calculus and Fast Animation",
771  // Ken Shoemake, 1987 SIGGRAPH course notes
772  Scalar t = mat.trace();
773  if (t > Scalar(0))
774  {
775  t = sqrt(t + Scalar(1.0));
776  q.w() = Scalar(0.5)*t;
777  t = Scalar(0.5)/t;
778  q.x() = (mat.coeff(2,1) - mat.coeff(1,2)) * t;
779  q.y() = (mat.coeff(0,2) - mat.coeff(2,0)) * t;
780  q.z() = (mat.coeff(1,0) - mat.coeff(0,1)) * t;
781  }
782  else
783  {
784  Index i = 0;
785  if (mat.coeff(1,1) > mat.coeff(0,0))
786  i = 1;
787  if (mat.coeff(2,2) > mat.coeff(i,i))
788  i = 2;
789  Index j = (i+1)%3;
790  Index k = (j+1)%3;
791 
792  t = sqrt(mat.coeff(i,i)-mat.coeff(j,j)-mat.coeff(k,k) + Scalar(1.0));
793  q.coeffs().coeffRef(i) = Scalar(0.5) * t;
794  t = Scalar(0.5)/t;
795  q.w() = (mat.coeff(k,j)-mat.coeff(j,k))*t;
796  q.coeffs().coeffRef(j) = (mat.coeff(j,i)+mat.coeff(i,j))*t;
797  q.coeffs().coeffRef(k) = (mat.coeff(k,i)+mat.coeff(i,k))*t;
798  }
799  }
800 };
801 
802 // set from a vector of coefficients assumed to be a quaternion
803 template<typename Other>
804 struct quaternionbase_assign_impl<Other,4,1>
805 {
806  typedef typename Other::Scalar Scalar;
807  template<class Derived> EIGEN_DEVICE_FUNC static inline void run(QuaternionBase<Derived>& q, const Other& vec)
808  {
809  q.coeffs() = vec;
810  }
811 };
812 
813 } // end namespace internal
814 
815 } // end namespace Eigen
816 
817 #endif // EIGEN_QUATERNION_H
Eigen::QuaternionBase::y
NonConstCoeffReturnType y()
Definition: Quaternion.h:77
Eigen
Namespace containing all symbols from the Eigen library.
Definition: Core:306
Eigen::QuaternionBase::inverse
Quaternion< Scalar > inverse() const
Definition: Quaternion.h:667
Eigen::AngleAxis::axis
const Vector3 & axis() const
Definition: AngleAxis.h:96
Eigen::ComputeFullV
@ ComputeFullV
Definition: Constants.h:387
Eigen::QuaternionBase::normalize
void normalize()
Definition: Quaternion.h:129
Eigen::sqrt
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Eigen::QuaternionBase::x
CoeffReturnType x() const
Definition: Quaternion.h:66
Eigen::Quaternion::Quaternion
Quaternion(const MatrixBase< Derived > &other)
Definition: Quaternion.h:275
Eigen::QuaternionMapAlignedf
Map< Quaternion< float >, Aligned > QuaternionMapAlignedf
Definition: Quaternion.h:422
Eigen::QuaternionBase::z
NonConstCoeffReturnType z()
Definition: Quaternion.h:79
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: ForwardDeclarations.h:270
Eigen::QuaternionBase::dot
Scalar dot(const QuaternionBase< OtherDerived > &other) const
Definition: Quaternion.h:139
Eigen::QuaternionMapd
Map< Quaternion< double >, 0 > QuaternionMapd
Definition: Quaternion.h:419
Eigen::Quaternion::Quaternion
Quaternion(const AngleAxisType &aa)
Definition: Quaternion.h:268
Eigen::QuaternionBase::z
CoeffReturnType z() const
Definition: Quaternion.h:70
Eigen::DontAlign
@ DontAlign
Definition: Constants.h:326
Eigen::QuaternionBase::AngleAxisType
AngleAxis< Scalar > AngleAxisType
Definition: Quaternion.h:61
Eigen::QuaternionBase::_transformVector
Vector3 _transformVector(const Vector3 &v) const
Definition: Quaternion.h:478
Eigen::Quaternion::Quaternion
Quaternion(const QuaternionBase< Derived > &other)
Definition: Quaternion.h:265
Eigen::QuaternionBase::setFromTwoVectors
Derived & setFromTwoVectors(const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
Definition: Quaternion.h:585
Eigen::QuaternionBase::setIdentity
QuaternionBase & setIdentity()
Definition: Quaternion.h:115
Eigen::Quaternion::Quaternion
Quaternion(const Scalar *data)
Definition: Quaternion.h:262
Eigen::QuaternionBase::vec
const VectorBlock< const Coefficients, 3 > vec() const
Definition: Quaternion.h:84
Eigen::QuaternionBase::toRotationMatrix
Matrix3 toRotationMatrix() const
Definition: Quaternion.h:539
Eigen::QuaternionBase::normalized
Quaternion< Scalar > normalized() const
Definition: Quaternion.h:132
Eigen::QuaternionBase::conjugate
Quaternion< Scalar > conjugate() const
Definition: Quaternion.h:698
Eigen::LvalueBit
const unsigned int LvalueBit
Definition: Constants.h:139
Eigen::Quaternion::UnitRandom
static Quaternion UnitRandom()
Definition: Quaternion.h:626
Eigen::QuaternionBase::isApprox
bool isApprox(const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Quaternion.h:166
Eigen::acos
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_acos_op< typename Derived::Scalar >, const Derived > acos(const Eigen::ArrayBase< Derived > &x)
Eigen::QuaternionBase::y
CoeffReturnType y() const
Definition: Quaternion.h:68
Eigen::Quaternion::Quaternion
Quaternion(const Scalar &w, const Scalar &x, const Scalar &y, const Scalar &z)
Definition: Quaternion.h:259
Eigen::QuaternionBase::x
NonConstCoeffReturnType x()
Definition: Quaternion.h:75
Eigen::QuaternionBase::w
NonConstCoeffReturnType w()
Definition: Quaternion.h:81
Eigen::Map< Quaternion< _Scalar >, _Options >::Map
Map(Scalar *coeffs)
Definition: Quaternion.h:405
Eigen::Quaternionf
Quaternion< float > Quaternionf
Definition: Quaternion.h:310
Eigen::PlainObjectBase::coeffRef
Scalar & coeffRef(Index rowId, Index colId)
Definition: PlainObjectBase.h:183
Eigen::QuaternionBase::cast
internal::cast_return_type< Derived, Quaternion< NewScalarType > >::type cast() const
Definition: Quaternion.h:178
Eigen::Quaterniond
Quaternion< double > Quaterniond
Definition: Quaternion.h:313
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::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: ForwardDeclarations.h:258
Eigen::QuaternionBase::norm
Scalar norm() const
Definition: Quaternion.h:125
Eigen::Quaternion
The quaternion class used to represent 3D orientations and rotations.
Definition: ForwardDeclarations.h:273
Eigen::AngleAxis::angle
Scalar angle() const
Definition: AngleAxis.h:91
Eigen::QuaternionBase::coeffs
const internal::traits< Derived >::Coefficients & coeffs() const
Definition: Quaternion.h:90
Eigen::Quaternion::Quaternion
Quaternion()
Definition: Quaternion.h:250
Eigen::QuaternionBase::Vector3
Matrix< Scalar, 3, 1 > Vector3
Definition: Quaternion.h:57
Eigen::QuaternionBase::Matrix3
Matrix< Scalar, 3, 3 > Matrix3
Definition: Quaternion.h:59
Eigen::QuaternionBase::Identity
static Quaternion< Scalar > Identity()
Definition: Quaternion.h:111
Eigen::QuaternionBase::operator*=
Derived & operator*=(const QuaternionBase< OtherDerived > &q)
Definition: Quaternion.h:463
Eigen::QuaternionMapf
Map< Quaternion< float >, 0 > QuaternionMapf
Definition: Quaternion.h:416
Eigen::VectorBlock
Expression of a fixed-size or dynamic-size sub-vector.
Definition: ForwardDeclarations.h:87
Eigen::Quaternion::Quaternion
Quaternion(const Quaternion< OtherScalar, OtherOptions > &other)
Definition: Quaternion.h:279
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::QuaternionBase::coeffs
internal::traits< Derived >::Coefficients & coeffs()
Definition: Quaternion.h:93
Eigen::QuaternionBase::w
CoeffReturnType w() const
Definition: Quaternion.h:72
Eigen::QuaternionBase::squaredNorm
Scalar squaredNorm() const
Definition: Quaternion.h:120
Eigen::Map< const Quaternion< _Scalar >, _Options >::Map
Map(const Scalar *coeffs)
Definition: Quaternion.h:368
Eigen::QuaternionMapAlignedd
Map< Quaternion< double >, Aligned > QuaternionMapAlignedd
Definition: Quaternion.h:425
Eigen::QuaternionBase
Base class for quaternion expressions.
Definition: ForwardDeclarations.h:268
Eigen::NumTraits
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:151
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::MatrixBase::normalized
const PlainObject normalized() const
Definition: Dot.h:124