Eigen  3.2.1
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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  typedef RotationBase<Derived, 3> Base;
38 public:
39  using Base::operator*;
40  using Base::derived;
41 
42  typedef typename internal::traits<Derived>::Scalar Scalar;
43  typedef typename NumTraits<Scalar>::Real RealScalar;
44  typedef typename internal::traits<Derived>::Coefficients Coefficients;
45  enum {
46  Flags = Eigen::internal::traits<Derived>::Flags
47  };
48 
49  // typedef typename Matrix<Scalar,4,1> Coefficients;
56 
57 
58 
60  inline Scalar x() const { return this->derived().coeffs().coeff(0); }
62  inline Scalar y() const { return this->derived().coeffs().coeff(1); }
64  inline Scalar z() const { return this->derived().coeffs().coeff(2); }
66  inline Scalar w() const { return this->derived().coeffs().coeff(3); }
67 
69  inline Scalar& x() { return this->derived().coeffs().coeffRef(0); }
71  inline Scalar& y() { return this->derived().coeffs().coeffRef(1); }
73  inline Scalar& z() { return this->derived().coeffs().coeffRef(2); }
75  inline Scalar& w() { return this->derived().coeffs().coeffRef(3); }
76 
78  inline const VectorBlock<const Coefficients,3> vec() const { return coeffs().template head<3>(); }
79 
81  inline VectorBlock<Coefficients,3> vec() { return coeffs().template head<3>(); }
82 
84  inline const typename internal::traits<Derived>::Coefficients& coeffs() const { return derived().coeffs(); }
85 
87  inline typename internal::traits<Derived>::Coefficients& coeffs() { return derived().coeffs(); }
88 
89  EIGEN_STRONG_INLINE QuaternionBase<Derived>& operator=(const QuaternionBase<Derived>& other);
90  template<class OtherDerived> EIGEN_STRONG_INLINE Derived& operator=(const QuaternionBase<OtherDerived>& other);
91 
92 // disabled this copy operator as it is giving very strange compilation errors when compiling
93 // test_stdvector with GCC 4.4.2. This looks like a GCC bug though, so feel free to re-enable it if it's
94 // useful; however notice that we already have the templated operator= above and e.g. in MatrixBase
95 // we didn't have to add, in addition to templated operator=, such a non-templated copy operator.
96 // Derived& operator=(const QuaternionBase& other)
97 // { return operator=<Derived>(other); }
98 
99  Derived& operator=(const AngleAxisType& aa);
100  template<class OtherDerived> Derived& operator=(const MatrixBase<OtherDerived>& m);
101 
105  static inline Quaternion<Scalar> Identity() { return Quaternion<Scalar>(1, 0, 0, 0); }
106 
109  inline QuaternionBase& setIdentity() { coeffs() << 0, 0, 0, 1; return *this; }
110 
114  inline Scalar squaredNorm() const { return coeffs().squaredNorm(); }
115 
119  inline Scalar norm() const { return coeffs().norm(); }
120 
123  inline void normalize() { coeffs().normalize(); }
126  inline Quaternion<Scalar> normalized() const { return Quaternion<Scalar>(coeffs().normalized()); }
127 
133  template<class OtherDerived> inline Scalar dot(const QuaternionBase<OtherDerived>& other) const { return coeffs().dot(other.coeffs()); }
134 
135  template<class OtherDerived> Scalar angularDistance(const QuaternionBase<OtherDerived>& other) const;
136 
138  Matrix3 toRotationMatrix() const;
139 
141  template<typename Derived1, typename Derived2>
142  Derived& setFromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
143 
144  template<class OtherDerived> EIGEN_STRONG_INLINE Quaternion<Scalar> operator* (const QuaternionBase<OtherDerived>& q) const;
145  template<class OtherDerived> EIGEN_STRONG_INLINE Derived& operator*= (const QuaternionBase<OtherDerived>& q);
146 
148  Quaternion<Scalar> inverse() const;
149 
152 
153  template<class OtherDerived> Quaternion<Scalar> slerp(const Scalar& t, const QuaternionBase<OtherDerived>& other) const;
154 
159  template<class OtherDerived>
160  bool isApprox(const QuaternionBase<OtherDerived>& other, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const
161  { return coeffs().isApprox(other.coeffs(), prec); }
162 
164  EIGEN_STRONG_INLINE Vector3 _transformVector(Vector3 v) const;
165 
171  template<typename NewScalarType>
172  inline typename internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type cast() const
173  {
174  return typename internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type(derived());
175  }
176 
177 #ifdef EIGEN_QUATERNIONBASE_PLUGIN
178 # include EIGEN_QUATERNIONBASE_PLUGIN
179 #endif
180 };
181 
182 /***************************************************************************
183 * Definition/implementation of Quaternion<Scalar>
184 ***************************************************************************/
185 
209 namespace internal {
210 template<typename _Scalar,int _Options>
211 struct traits<Quaternion<_Scalar,_Options> >
212 {
213  typedef Quaternion<_Scalar,_Options> PlainObject;
214  typedef _Scalar Scalar;
215  typedef Matrix<_Scalar,4,1,_Options> Coefficients;
216  enum{
217  IsAligned = internal::traits<Coefficients>::Flags & AlignedBit,
218  Flags = IsAligned ? (AlignedBit | LvalueBit) : LvalueBit
219  };
220 };
221 }
222 
223 template<typename _Scalar, int _Options>
224 class Quaternion : public QuaternionBase<Quaternion<_Scalar,_Options> >
225 {
226  typedef QuaternionBase<Quaternion<_Scalar,_Options> > Base;
227  enum { IsAligned = internal::traits<Quaternion>::IsAligned };
228 
229 public:
230  typedef _Scalar Scalar;
231 
232  EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(Quaternion)
233  using Base::operator*=;
234 
235  typedef typename internal::traits<Quaternion>::Coefficients Coefficients;
236  typedef typename Base::AngleAxisType AngleAxisType;
237 
239  inline Quaternion() {}
240 
248  inline Quaternion(const Scalar& w, const Scalar& x, const Scalar& y, const Scalar& z) : m_coeffs(x, y, z, w){}
249 
251  inline Quaternion(const Scalar* data) : m_coeffs(data) {}
252 
254  template<class Derived> EIGEN_STRONG_INLINE Quaternion(const QuaternionBase<Derived>& other) { this->Base::operator=(other); }
255 
257  explicit inline Quaternion(const AngleAxisType& aa) { *this = aa; }
258 
263  template<typename Derived>
264  explicit inline Quaternion(const MatrixBase<Derived>& other) { *this = other; }
265 
267  template<typename OtherScalar, int OtherOptions>
268  explicit inline Quaternion(const Quaternion<OtherScalar, OtherOptions>& other)
269  { m_coeffs = other.coeffs().template cast<Scalar>(); }
270 
271  template<typename Derived1, typename Derived2>
272  static Quaternion FromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
273 
274  inline Coefficients& coeffs() { return m_coeffs;}
275  inline const Coefficients& coeffs() const { return m_coeffs;}
276 
277  EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(IsAligned)
278 
279 protected:
280  Coefficients m_coeffs;
281 
282 #ifndef EIGEN_PARSED_BY_DOXYGEN
283  static EIGEN_STRONG_INLINE void _check_template_params()
284  {
285  EIGEN_STATIC_ASSERT( (_Options & DontAlign) == _Options,
286  INVALID_MATRIX_TEMPLATE_PARAMETERS)
287  }
288 #endif
289 };
290 
297 
298 /***************************************************************************
299 * Specialization of Map<Quaternion<Scalar>>
300 ***************************************************************************/
301 
302 namespace internal {
303  template<typename _Scalar, int _Options>
304  struct traits<Map<Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
305  {
306  typedef Map<Matrix<_Scalar,4,1>, _Options> Coefficients;
307  };
308 }
309 
310 namespace internal {
311  template<typename _Scalar, int _Options>
312  struct traits<Map<const Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
313  {
314  typedef Map<const Matrix<_Scalar,4,1>, _Options> Coefficients;
315  typedef traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> > TraitsBase;
316  enum {
317  Flags = TraitsBase::Flags & ~LvalueBit
318  };
319  };
320 }
321 
333 template<typename _Scalar, int _Options>
334 class Map<const Quaternion<_Scalar>, _Options >
335  : public QuaternionBase<Map<const Quaternion<_Scalar>, _Options> >
336 {
338 
339  public:
340  typedef _Scalar Scalar;
341  typedef typename internal::traits<Map>::Coefficients Coefficients;
342  EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(Map)
343  using Base::operator*=;
344 
351  EIGEN_STRONG_INLINE Map(const Scalar* coeffs) : m_coeffs(coeffs) {}
352 
353  inline const Coefficients& coeffs() const { return m_coeffs;}
354 
355  protected:
356  const Coefficients m_coeffs;
357 };
358 
370 template<typename _Scalar, int _Options>
371 class Map<Quaternion<_Scalar>, _Options >
372  : public QuaternionBase<Map<Quaternion<_Scalar>, _Options> >
373 {
374  typedef QuaternionBase<Map<Quaternion<_Scalar>, _Options> > Base;
375 
376  public:
377  typedef _Scalar Scalar;
378  typedef typename internal::traits<Map>::Coefficients Coefficients;
379  EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(Map)
380  using Base::operator*=;
381 
388  EIGEN_STRONG_INLINE Map(Scalar* coeffs) : m_coeffs(coeffs) {}
389 
390  inline Coefficients& coeffs() { return m_coeffs; }
391  inline const Coefficients& coeffs() const { return m_coeffs; }
392 
393  protected:
394  Coefficients m_coeffs;
395 };
396 
409 
410 /***************************************************************************
411 * Implementation of QuaternionBase methods
412 ***************************************************************************/
413 
414 // Generic Quaternion * Quaternion product
415 // This product can be specialized for a given architecture via the Arch template argument.
416 namespace internal {
417 template<int Arch, class Derived1, class Derived2, typename Scalar, int _Options> struct quat_product
418 {
419  static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived1>& a, const QuaternionBase<Derived2>& b){
420  return Quaternion<Scalar>
421  (
422  a.w() * b.w() - a.x() * b.x() - a.y() * b.y() - a.z() * b.z(),
423  a.w() * b.x() + a.x() * b.w() + a.y() * b.z() - a.z() * b.y(),
424  a.w() * b.y() + a.y() * b.w() + a.z() * b.x() - a.x() * b.z(),
425  a.w() * b.z() + a.z() * b.w() + a.x() * b.y() - a.y() * b.x()
426  );
427  }
428 };
429 }
430 
432 template <class Derived>
433 template <class OtherDerived>
434 EIGEN_STRONG_INLINE Quaternion<typename internal::traits<Derived>::Scalar>
436 {
437  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename OtherDerived::Scalar>::value),
438  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
439  return internal::quat_product<Architecture::Target, Derived, OtherDerived,
440  typename internal::traits<Derived>::Scalar,
441  internal::traits<Derived>::IsAligned && internal::traits<OtherDerived>::IsAligned>::run(*this, other);
442 }
443 
445 template <class Derived>
446 template <class OtherDerived>
447 EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator*= (const QuaternionBase<OtherDerived>& other)
448 {
449  derived() = derived() * other.derived();
450  return derived();
451 }
452 
460 template <class Derived>
461 EIGEN_STRONG_INLINE typename QuaternionBase<Derived>::Vector3
463 {
464  // Note that this algorithm comes from the optimization by hand
465  // of the conversion to a Matrix followed by a Matrix/Vector product.
466  // It appears to be much faster than the common algorithm found
467  // in the litterature (30 versus 39 flops). It also requires two
468  // Vector3 as temporaries.
469  Vector3 uv = this->vec().cross(v);
470  uv += uv;
471  return v + this->w() * uv + this->vec().cross(uv);
472 }
473 
474 template<class Derived>
476 {
477  coeffs() = other.coeffs();
478  return derived();
479 }
480 
481 template<class Derived>
482 template<class OtherDerived>
483 EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const QuaternionBase<OtherDerived>& other)
484 {
485  coeffs() = other.coeffs();
486  return derived();
487 }
488 
491 template<class Derived>
492 EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const AngleAxisType& aa)
493 {
494  using std::cos;
495  using std::sin;
496  Scalar ha = Scalar(0.5)*aa.angle(); // Scalar(0.5) to suppress precision loss warnings
497  this->w() = cos(ha);
498  this->vec() = sin(ha) * aa.axis();
499  return derived();
500 }
501 
508 template<class Derived>
509 template<class MatrixDerived>
511 {
512  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename MatrixDerived::Scalar>::value),
513  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
514  internal::quaternionbase_assign_impl<MatrixDerived>::run(*this, xpr.derived());
515  return derived();
516 }
517 
521 template<class Derived>
522 inline typename QuaternionBase<Derived>::Matrix3
524 {
525  // NOTE if inlined, then gcc 4.2 and 4.4 get rid of the temporary (not gcc 4.3 !!)
526  // if not inlined then the cost of the return by value is huge ~ +35%,
527  // however, not inlining this function is an order of magnitude slower, so
528  // it has to be inlined, and so the return by value is not an issue
529  Matrix3 res;
530 
531  const Scalar tx = Scalar(2)*this->x();
532  const Scalar ty = Scalar(2)*this->y();
533  const Scalar tz = Scalar(2)*this->z();
534  const Scalar twx = tx*this->w();
535  const Scalar twy = ty*this->w();
536  const Scalar twz = tz*this->w();
537  const Scalar txx = tx*this->x();
538  const Scalar txy = ty*this->x();
539  const Scalar txz = tz*this->x();
540  const Scalar tyy = ty*this->y();
541  const Scalar tyz = tz*this->y();
542  const Scalar tzz = tz*this->z();
543 
544  res.coeffRef(0,0) = Scalar(1)-(tyy+tzz);
545  res.coeffRef(0,1) = txy-twz;
546  res.coeffRef(0,2) = txz+twy;
547  res.coeffRef(1,0) = txy+twz;
548  res.coeffRef(1,1) = Scalar(1)-(txx+tzz);
549  res.coeffRef(1,2) = tyz-twx;
550  res.coeffRef(2,0) = txz-twy;
551  res.coeffRef(2,1) = tyz+twx;
552  res.coeffRef(2,2) = Scalar(1)-(txx+tyy);
553 
554  return res;
555 }
556 
567 template<class Derived>
568 template<typename Derived1, typename Derived2>
570 {
571  using std::max;
572  using std::sqrt;
573  Vector3 v0 = a.normalized();
574  Vector3 v1 = b.normalized();
575  Scalar c = v1.dot(v0);
576 
577  // if dot == -1, vectors are nearly opposites
578  // => accurately compute the rotation axis by computing the
579  // intersection of the two planes. This is done by solving:
580  // x^T v0 = 0
581  // x^T v1 = 0
582  // under the constraint:
583  // ||x|| = 1
584  // which yields a singular value problem
585  if (c < Scalar(-1)+NumTraits<Scalar>::dummy_precision())
586  {
587  c = max<Scalar>(c,-1);
588  Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
590  Vector3 axis = svd.matrixV().col(2);
591 
592  Scalar w2 = (Scalar(1)+c)*Scalar(0.5);
593  this->w() = sqrt(w2);
594  this->vec() = axis * sqrt(Scalar(1) - w2);
595  return derived();
596  }
597  Vector3 axis = v0.cross(v1);
598  Scalar s = sqrt((Scalar(1)+c)*Scalar(2));
599  Scalar invs = Scalar(1)/s;
600  this->vec() = axis * invs;
601  this->w() = s * Scalar(0.5);
602 
603  return derived();
604 }
605 
606 
617 template<typename Scalar, int Options>
618 template<typename Derived1, typename Derived2>
620 {
621  Quaternion quat;
622  quat.setFromTwoVectors(a, b);
623  return quat;
624 }
625 
626 
633 template <class Derived>
635 {
636  // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite() ??
637  Scalar n2 = this->squaredNorm();
638  if (n2 > 0)
639  return Quaternion<Scalar>(conjugate().coeffs() / n2);
640  else
641  {
642  // return an invalid result to flag the error
643  return Quaternion<Scalar>(Coefficients::Zero());
644  }
645 }
646 
653 template <class Derived>
656 {
657  return Quaternion<Scalar>(this->w(),-this->x(),-this->y(),-this->z());
658 }
659 
663 template <class Derived>
664 template <class OtherDerived>
665 inline typename internal::traits<Derived>::Scalar
667 {
668  using std::acos;
669  using std::abs;
670  double d = abs(this->dot(other));
671  if (d>=1.0)
672  return Scalar(0);
673  return static_cast<Scalar>(2 * acos(d));
674 }
675 
676 
677 
684 template <class Derived>
685 template <class OtherDerived>
688 {
689  using std::acos;
690  using std::sin;
691  using std::abs;
692  static const Scalar one = Scalar(1) - NumTraits<Scalar>::epsilon();
693  Scalar d = this->dot(other);
694  Scalar absD = abs(d);
695 
696  Scalar scale0;
697  Scalar scale1;
698 
699  if(absD>=one)
700  {
701  scale0 = Scalar(1) - t;
702  scale1 = t;
703  }
704  else
705  {
706  // theta is the angle between the 2 quaternions
707  Scalar theta = acos(absD);
708  Scalar sinTheta = sin(theta);
709 
710  scale0 = sin( ( Scalar(1) - t ) * theta) / sinTheta;
711  scale1 = sin( ( t * theta) ) / sinTheta;
712  }
713  if(d<0) scale1 = -scale1;
714 
715  return Quaternion<Scalar>(scale0 * coeffs() + scale1 * other.coeffs());
716 }
717 
718 namespace internal {
719 
720 // set from a rotation matrix
721 template<typename Other>
722 struct quaternionbase_assign_impl<Other,3,3>
723 {
724  typedef typename Other::Scalar Scalar;
725  typedef DenseIndex Index;
726  template<class Derived> static inline void run(QuaternionBase<Derived>& q, const Other& mat)
727  {
728  using std::sqrt;
729  // This algorithm comes from "Quaternion Calculus and Fast Animation",
730  // Ken Shoemake, 1987 SIGGRAPH course notes
731  Scalar t = mat.trace();
732  if (t > Scalar(0))
733  {
734  t = sqrt(t + Scalar(1.0));
735  q.w() = Scalar(0.5)*t;
736  t = Scalar(0.5)/t;
737  q.x() = (mat.coeff(2,1) - mat.coeff(1,2)) * t;
738  q.y() = (mat.coeff(0,2) - mat.coeff(2,0)) * t;
739  q.z() = (mat.coeff(1,0) - mat.coeff(0,1)) * t;
740  }
741  else
742  {
743  DenseIndex i = 0;
744  if (mat.coeff(1,1) > mat.coeff(0,0))
745  i = 1;
746  if (mat.coeff(2,2) > mat.coeff(i,i))
747  i = 2;
748  DenseIndex j = (i+1)%3;
749  DenseIndex k = (j+1)%3;
750 
751  t = sqrt(mat.coeff(i,i)-mat.coeff(j,j)-mat.coeff(k,k) + Scalar(1.0));
752  q.coeffs().coeffRef(i) = Scalar(0.5) * t;
753  t = Scalar(0.5)/t;
754  q.w() = (mat.coeff(k,j)-mat.coeff(j,k))*t;
755  q.coeffs().coeffRef(j) = (mat.coeff(j,i)+mat.coeff(i,j))*t;
756  q.coeffs().coeffRef(k) = (mat.coeff(k,i)+mat.coeff(i,k))*t;
757  }
758  }
759 };
760 
761 // set from a vector of coefficients assumed to be a quaternion
762 template<typename Other>
763 struct quaternionbase_assign_impl<Other,4,1>
764 {
765  typedef typename Other::Scalar Scalar;
766  template<class Derived> static inline void run(QuaternionBase<Derived>& q, const Other& vec)
767  {
768  q.coeffs() = vec;
769  }
770 };
771 
772 } // end namespace internal
773 
774 } // end namespace Eigen
775 
776 #endif // EIGEN_QUATERNION_H
Map(const Scalar *coeffs)
Definition: Quaternion.h:351
Scalar x() const
Definition: Quaternion.h:60
Quaternion< float > Quaternionf
Definition: Quaternion.h:293
Scalar norm() const
Definition: Quaternion.h:119
Scalar dot(const QuaternionBase< OtherDerived > &other) const
Definition: Quaternion.h:133
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:104
Scalar squaredNorm() const
Definition: Quaternion.h:114
internal::traits< Quaternion< _Scalar, _Options > >::Scalar Scalar
Definition: RotationBase.h:34
const internal::traits< Derived >::Coefficients & coeffs() const
Definition: Quaternion.h:84
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
const PlainObject normalized() const
Definition: Dot.h:139
Scalar & w()
Definition: Quaternion.h:75
Definition: Constants.h:331
Quaternion< Scalar > inverse() const
Definition: Quaternion.h:634
VectorBlock< Coefficients, 3 > vec()
Definition: Quaternion.h:81
Map(Scalar *coeffs)
Definition: Quaternion.h:388
Expression of a fixed-size or dynamic-size sub-vector.
Definition: ForwardDeclarations.h:83
Quaternion(const AngleAxisType &aa)
Definition: Quaternion.h:257
Scalar y() const
Definition: Quaternion.h:62
Scalar & z()
Definition: Quaternion.h:73
internal::traits< Derived >::Coefficients & coeffs()
Definition: Quaternion.h:87
const unsigned int LvalueBit
Definition: Constants.h:131
Quaternion(const Scalar *data)
Definition: Quaternion.h:251
Quaternion< Scalar > conjugate() const
Definition: Quaternion.h:655
Definition: Constants.h:270
Quaternion< double > Quaterniond
Definition: Quaternion.h:296
Derived & setFromTwoVectors(const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
Definition: Quaternion.h:569
Quaternion< Scalar > normalized() const
Definition: Quaternion.h:126
internal::cast_return_type< Derived, Quaternion< NewScalarType > >::type cast() const
Definition: Quaternion.h:172
Quaternion(const Quaternion< OtherScalar, OtherOptions > &other)
Definition: Quaternion.h:268
Base class for quaternion expressions.
Definition: ForwardDeclarations.h:233
Matrix< Scalar, 3, 1 > Vector3
Definition: Quaternion.h:51
Quaternion(const Scalar &w, const Scalar &x, const Scalar &y, const Scalar &z)
Definition: Quaternion.h:248
AngleAxis< Scalar > AngleAxisType
Definition: Quaternion.h:55
Scalar & x()
Definition: Quaternion.h:69
Matrix3 toRotationMatrix() const
Definition: Quaternion.h:523
const MatrixVType & matrixV() const
Definition: JacobiSVD.h:621
Scalar z() const
Definition: Quaternion.h:64
Map< Quaternion< double >, Aligned > QuaternionMapAlignedd
Definition: Quaternion.h:408
bool isApprox(const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Quaternion.h:160
The quaternion class used to represent 3D orientations and rotations.
Definition: ForwardDeclarations.h:260
Vector3 _transformVector(Vector3 v) const
Definition: Quaternion.h:462
Definition: Constants.h:194
Quaternion(const MatrixBase< Derived > &other)
Definition: Quaternion.h:264
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: ForwardDeclarations.h:224
Map< Quaternion< float >, Aligned > QuaternionMapAlignedf
Definition: Quaternion.h:405
Scalar & y()
Definition: Quaternion.h:71
void normalize()
Definition: Quaternion.h:123
Scalar w() const
Definition: Quaternion.h:66
Matrix< Scalar, 3, 3 > Matrix3
Definition: Quaternion.h:53
Map< Quaternion< float >, 0 > QuaternionMapf
Definition: Quaternion.h:399
Map< Quaternion< double >, 0 > QuaternionMapd
Definition: Quaternion.h:402
const unsigned int AlignedBit
Definition: Constants.h:147
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Represents a 3D rotation as a rotation angle around an arbitrary 3D axis.
Definition: ForwardDeclarations.h:235
Quaternion(const QuaternionBase< Derived > &other)
Definition: Quaternion.h:254
Derived & operator*=(const QuaternionBase< OtherDerived > &q)
Definition: Quaternion.h:447
static Quaternion< Scalar > Identity()
Definition: Quaternion.h:105
QuaternionBase & setIdentity()
Definition: Quaternion.h:109
const VectorBlock< const Coefficients, 3 > vec() const
Definition: Quaternion.h:78