Eigen  3.2.7
 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 <[email protected]>
5 // Copyright (C) 2009 Mathieu Gautier <[email protected]>
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>(Scalar(1), Scalar(0), Scalar(0), Scalar(0)); }
106 
109  inline QuaternionBase& setIdentity() { coeffs() << Scalar(0), Scalar(0), Scalar(0), Scalar(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(const 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 
211 namespace internal {
212 template<typename _Scalar,int _Options>
213 struct traits<Quaternion<_Scalar,_Options> >
214 {
215  typedef Quaternion<_Scalar,_Options> PlainObject;
216  typedef _Scalar Scalar;
217  typedef Matrix<_Scalar,4,1,_Options> Coefficients;
218  enum{
219  IsAligned = internal::traits<Coefficients>::Flags & AlignedBit,
220  Flags = IsAligned ? (AlignedBit | LvalueBit) : LvalueBit
221  };
222 };
223 }
224 
225 template<typename _Scalar, int _Options>
226 class Quaternion : public QuaternionBase<Quaternion<_Scalar,_Options> >
227 {
228  typedef QuaternionBase<Quaternion<_Scalar,_Options> > Base;
229  enum { IsAligned = internal::traits<Quaternion>::IsAligned };
230 
231 public:
232  typedef _Scalar Scalar;
233 
234  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Quaternion)
235  using Base::operator*=;
236 
237  typedef typename internal::traits<Quaternion>::Coefficients Coefficients;
238  typedef typename Base::AngleAxisType AngleAxisType;
239 
241  inline Quaternion() {}
242 
250  inline Quaternion(const Scalar& w, const Scalar& x, const Scalar& y, const Scalar& z) : m_coeffs(x, y, z, w){}
251 
253  inline Quaternion(const Scalar* data) : m_coeffs(data) {}
254 
256  template<class Derived> EIGEN_STRONG_INLINE Quaternion(const QuaternionBase<Derived>& other) { this->Base::operator=(other); }
257 
259  explicit inline Quaternion(const AngleAxisType& aa) { *this = aa; }
260 
265  template<typename Derived>
266  explicit inline Quaternion(const MatrixBase<Derived>& other) { *this = other; }
267 
269  template<typename OtherScalar, int OtherOptions>
270  explicit inline Quaternion(const Quaternion<OtherScalar, OtherOptions>& other)
271  { m_coeffs = other.coeffs().template cast<Scalar>(); }
272 
273  template<typename Derived1, typename Derived2>
274  static Quaternion FromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
275 
276  inline Coefficients& coeffs() { return m_coeffs;}
277  inline const Coefficients& coeffs() const { return m_coeffs;}
278 
279  EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(IsAligned)
280 
281 protected:
282  Coefficients m_coeffs;
283 
284 #ifndef EIGEN_PARSED_BY_DOXYGEN
285  static EIGEN_STRONG_INLINE void _check_template_params()
286  {
287  EIGEN_STATIC_ASSERT( (_Options & DontAlign) == _Options,
288  INVALID_MATRIX_TEMPLATE_PARAMETERS)
289  }
290 #endif
291 };
292 
299 
300 /***************************************************************************
301 * Specialization of Map<Quaternion<Scalar>>
302 ***************************************************************************/
303 
304 namespace internal {
305  template<typename _Scalar, int _Options>
306  struct traits<Map<Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
307  {
308  typedef Map<Matrix<_Scalar,4,1>, _Options> Coefficients;
309  };
310 }
311 
312 namespace internal {
313  template<typename _Scalar, int _Options>
314  struct traits<Map<const Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
315  {
316  typedef Map<const Matrix<_Scalar,4,1>, _Options> Coefficients;
317  typedef traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> > TraitsBase;
318  enum {
319  Flags = TraitsBase::Flags & ~LvalueBit
320  };
321  };
322 }
323 
335 template<typename _Scalar, int _Options>
336 class Map<const Quaternion<_Scalar>, _Options >
337  : public QuaternionBase<Map<const Quaternion<_Scalar>, _Options> >
338 {
340 
341  public:
342  typedef _Scalar Scalar;
343  typedef typename internal::traits<Map>::Coefficients Coefficients;
344  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
345  using Base::operator*=;
346 
353  EIGEN_STRONG_INLINE Map(const Scalar* coeffs) : m_coeffs(coeffs) {}
354 
355  inline const Coefficients& coeffs() const { return m_coeffs;}
356 
357  protected:
358  const Coefficients m_coeffs;
359 };
360 
372 template<typename _Scalar, int _Options>
373 class Map<Quaternion<_Scalar>, _Options >
374  : public QuaternionBase<Map<Quaternion<_Scalar>, _Options> >
375 {
376  typedef QuaternionBase<Map<Quaternion<_Scalar>, _Options> > Base;
377 
378  public:
379  typedef _Scalar Scalar;
380  typedef typename internal::traits<Map>::Coefficients Coefficients;
381  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
382  using Base::operator*=;
383 
390  EIGEN_STRONG_INLINE Map(Scalar* coeffs) : m_coeffs(coeffs) {}
391 
392  inline Coefficients& coeffs() { return m_coeffs; }
393  inline const Coefficients& coeffs() const { return m_coeffs; }
394 
395  protected:
396  Coefficients m_coeffs;
397 };
398 
411 
412 /***************************************************************************
413 * Implementation of QuaternionBase methods
414 ***************************************************************************/
415 
416 // Generic Quaternion * Quaternion product
417 // This product can be specialized for a given architecture via the Arch template argument.
418 namespace internal {
419 template<int Arch, class Derived1, class Derived2, typename Scalar, int _Options> struct quat_product
420 {
421  static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived1>& a, const QuaternionBase<Derived2>& b){
422  return Quaternion<Scalar>
423  (
424  a.w() * b.w() - a.x() * b.x() - a.y() * b.y() - a.z() * b.z(),
425  a.w() * b.x() + a.x() * b.w() + a.y() * b.z() - a.z() * b.y(),
426  a.w() * b.y() + a.y() * b.w() + a.z() * b.x() - a.x() * b.z(),
427  a.w() * b.z() + a.z() * b.w() + a.x() * b.y() - a.y() * b.x()
428  );
429  }
430 };
431 }
432 
434 template <class Derived>
435 template <class OtherDerived>
436 EIGEN_STRONG_INLINE Quaternion<typename internal::traits<Derived>::Scalar>
438 {
439  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename OtherDerived::Scalar>::value),
440  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
441  return internal::quat_product<Architecture::Target, Derived, OtherDerived,
442  typename internal::traits<Derived>::Scalar,
443  internal::traits<Derived>::IsAligned && internal::traits<OtherDerived>::IsAligned>::run(*this, other);
444 }
445 
447 template <class Derived>
448 template <class OtherDerived>
449 EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator*= (const QuaternionBase<OtherDerived>& other)
450 {
451  derived() = derived() * other.derived();
452  return derived();
453 }
454 
462 template <class Derived>
463 EIGEN_STRONG_INLINE typename QuaternionBase<Derived>::Vector3
465 {
466  // Note that this algorithm comes from the optimization by hand
467  // of the conversion to a Matrix followed by a Matrix/Vector product.
468  // It appears to be much faster than the common algorithm found
469  // in the literature (30 versus 39 flops). It also requires two
470  // Vector3 as temporaries.
471  Vector3 uv = this->vec().cross(v);
472  uv += uv;
473  return v + this->w() * uv + this->vec().cross(uv);
474 }
475 
476 template<class Derived>
478 {
479  coeffs() = other.coeffs();
480  return derived();
481 }
482 
483 template<class Derived>
484 template<class OtherDerived>
485 EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const QuaternionBase<OtherDerived>& other)
486 {
487  coeffs() = other.coeffs();
488  return derived();
489 }
490 
493 template<class Derived>
494 EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const AngleAxisType& aa)
495 {
496  using std::cos;
497  using std::sin;
498  Scalar ha = Scalar(0.5)*aa.angle(); // Scalar(0.5) to suppress precision loss warnings
499  this->w() = cos(ha);
500  this->vec() = sin(ha) * aa.axis();
501  return derived();
502 }
503 
510 template<class Derived>
511 template<class MatrixDerived>
513 {
514  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename MatrixDerived::Scalar>::value),
515  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
516  internal::quaternionbase_assign_impl<MatrixDerived>::run(*this, xpr.derived());
517  return derived();
518 }
519 
523 template<class Derived>
524 inline typename QuaternionBase<Derived>::Matrix3
526 {
527  // NOTE if inlined, then gcc 4.2 and 4.4 get rid of the temporary (not gcc 4.3 !!)
528  // if not inlined then the cost of the return by value is huge ~ +35%,
529  // however, not inlining this function is an order of magnitude slower, so
530  // it has to be inlined, and so the return by value is not an issue
531  Matrix3 res;
532 
533  const Scalar tx = Scalar(2)*this->x();
534  const Scalar ty = Scalar(2)*this->y();
535  const Scalar tz = Scalar(2)*this->z();
536  const Scalar twx = tx*this->w();
537  const Scalar twy = ty*this->w();
538  const Scalar twz = tz*this->w();
539  const Scalar txx = tx*this->x();
540  const Scalar txy = ty*this->x();
541  const Scalar txz = tz*this->x();
542  const Scalar tyy = ty*this->y();
543  const Scalar tyz = tz*this->y();
544  const Scalar tzz = tz*this->z();
545 
546  res.coeffRef(0,0) = Scalar(1)-(tyy+tzz);
547  res.coeffRef(0,1) = txy-twz;
548  res.coeffRef(0,2) = txz+twy;
549  res.coeffRef(1,0) = txy+twz;
550  res.coeffRef(1,1) = Scalar(1)-(txx+tzz);
551  res.coeffRef(1,2) = tyz-twx;
552  res.coeffRef(2,0) = txz-twy;
553  res.coeffRef(2,1) = tyz+twx;
554  res.coeffRef(2,2) = Scalar(1)-(txx+tyy);
555 
556  return res;
557 }
558 
569 template<class Derived>
570 template<typename Derived1, typename Derived2>
572 {
573  using std::max;
574  using std::sqrt;
575  Vector3 v0 = a.normalized();
576  Vector3 v1 = b.normalized();
577  Scalar c = v1.dot(v0);
578 
579  // if dot == -1, vectors are nearly opposites
580  // => accurately compute the rotation axis by computing the
581  // intersection of the two planes. This is done by solving:
582  // x^T v0 = 0
583  // x^T v1 = 0
584  // under the constraint:
585  // ||x|| = 1
586  // which yields a singular value problem
587  if (c < Scalar(-1)+NumTraits<Scalar>::dummy_precision())
588  {
589  c = (max)(c,Scalar(-1));
590  Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
592  Vector3 axis = svd.matrixV().col(2);
593 
594  Scalar w2 = (Scalar(1)+c)*Scalar(0.5);
595  this->w() = sqrt(w2);
596  this->vec() = axis * sqrt(Scalar(1) - w2);
597  return derived();
598  }
599  Vector3 axis = v0.cross(v1);
600  Scalar s = sqrt((Scalar(1)+c)*Scalar(2));
601  Scalar invs = Scalar(1)/s;
602  this->vec() = axis * invs;
603  this->w() = s * Scalar(0.5);
604 
605  return derived();
606 }
607 
608 
619 template<typename Scalar, int Options>
620 template<typename Derived1, typename Derived2>
622 {
623  Quaternion quat;
624  quat.setFromTwoVectors(a, b);
625  return quat;
626 }
627 
628 
635 template <class Derived>
637 {
638  // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite() ??
639  Scalar n2 = this->squaredNorm();
640  if (n2 > Scalar(0))
641  return Quaternion<Scalar>(conjugate().coeffs() / n2);
642  else
643  {
644  // return an invalid result to flag the error
645  return Quaternion<Scalar>(Coefficients::Zero());
646  }
647 }
648 
655 template <class Derived>
658 {
659  return Quaternion<Scalar>(this->w(),-this->x(),-this->y(),-this->z());
660 }
661 
665 template <class Derived>
666 template <class OtherDerived>
667 inline typename internal::traits<Derived>::Scalar
669 {
670  using std::atan2;
671  using std::abs;
672  Quaternion<Scalar> d = (*this) * other.conjugate();
673  return Scalar(2) * atan2( d.vec().norm(), abs(d.w()) );
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<Scalar(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
Definition: Constants.h:270
Scalar angle() const
Definition: AngleAxis.h:87
Map(const Scalar *coeffs)
Definition: Quaternion.h:353
Scalar x() const
Definition: Quaternion.h:60
Quaternion< float > Quaternionf
Definition: Quaternion.h:295
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:636
Definition: Constants.h:194
VectorBlock< Coefficients, 3 > vec()
Definition: Quaternion.h:81
Map(Scalar *coeffs)
Definition: Quaternion.h:390
Expression of a fixed-size or dynamic-size sub-vector.
Definition: ForwardDeclarations.h:83
Vector3 _transformVector(const Vector3 &v) const
Definition: Quaternion.h:464
Quaternion(const AngleAxisType &aa)
Definition: Quaternion.h:259
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:253
Quaternion< Scalar > conjugate() const
Definition: Quaternion.h:657
Quaternion< double > Quaterniond
Definition: Quaternion.h:298
Derived & setFromTwoVectors(const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
Definition: Quaternion.h:571
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:270
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:250
AngleAxis< Scalar > AngleAxisType
Definition: Quaternion.h:55
Scalar & x()
Definition: Quaternion.h:69
Matrix3 toRotationMatrix() const
Definition: Quaternion.h:525
const MatrixVType & matrixV() const
Definition: JacobiSVD.h:629
Scalar z() const
Definition: Quaternion.h:64
Map< Quaternion< double >, Aligned > QuaternionMapAlignedd
Definition: Quaternion.h:410
bool isApprox(const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Quaternion.h:160
const Vector3 & axis() const
Definition: AngleAxis.h:92
The quaternion class used to represent 3D orientations and rotations.
Definition: ForwardDeclarations.h:263
Quaternion(const MatrixBase< Derived > &other)
Definition: Quaternion.h:266
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: ForwardDeclarations.h:224
Map< Quaternion< float >, Aligned > QuaternionMapAlignedf
Definition: Quaternion.h:407
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:401
Map< Quaternion< double >, 0 > QuaternionMapd
Definition: Quaternion.h:404
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:256
Derived & operator*=(const QuaternionBase< OtherDerived > &q)
Definition: Quaternion.h:449
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