Eigen  3.2.92
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  enum {
47  Flags = Eigen::internal::traits<Derived>::Flags
48  };
49 
50  // typedef typename Matrix<Scalar,4,1> Coefficients;
57 
58 
59 
61  inline Scalar x() const { return this->derived().coeffs().coeff(0); }
63  inline Scalar y() const { return this->derived().coeffs().coeff(1); }
65  inline Scalar z() const { return this->derived().coeffs().coeff(2); }
67  inline Scalar w() const { return this->derived().coeffs().coeff(3); }
68 
70  inline Scalar& x() { return this->derived().coeffs().coeffRef(0); }
72  inline Scalar& y() { return this->derived().coeffs().coeffRef(1); }
74  inline Scalar& z() { return this->derived().coeffs().coeffRef(2); }
76  inline Scalar& w() { return this->derived().coeffs().coeffRef(3); }
77 
79  inline const VectorBlock<const Coefficients,3> vec() const { return coeffs().template head<3>(); }
80 
82  inline VectorBlock<Coefficients,3> vec() { return coeffs().template head<3>(); }
83 
85  inline const typename internal::traits<Derived>::Coefficients& coeffs() const { return derived().coeffs(); }
86 
88  inline typename internal::traits<Derived>::Coefficients& coeffs() { return derived().coeffs(); }
89 
90  EIGEN_STRONG_INLINE QuaternionBase<Derived>& operator=(const QuaternionBase<Derived>& other);
91  template<class OtherDerived> EIGEN_STRONG_INLINE Derived& operator=(const QuaternionBase<OtherDerived>& other);
92 
93 // disabled this copy operator as it is giving very strange compilation errors when compiling
94 // test_stdvector with GCC 4.4.2. This looks like a GCC bug though, so feel free to re-enable it if it's
95 // useful; however notice that we already have the templated operator= above and e.g. in MatrixBase
96 // we didn't have to add, in addition to templated operator=, such a non-templated copy operator.
97 // Derived& operator=(const QuaternionBase& other)
98 // { return operator=<Derived>(other); }
99 
100  Derived& operator=(const AngleAxisType& aa);
101  template<class OtherDerived> Derived& operator=(const MatrixBase<OtherDerived>& m);
102 
106  static inline Quaternion<Scalar> Identity() { return Quaternion<Scalar>(Scalar(1), Scalar(0), Scalar(0), Scalar(0)); }
107 
110  inline QuaternionBase& setIdentity() { coeffs() << Scalar(0), Scalar(0), Scalar(0), Scalar(1); return *this; }
111 
115  inline Scalar squaredNorm() const { return coeffs().squaredNorm(); }
116 
120  inline Scalar norm() const { return coeffs().norm(); }
121 
124  inline void normalize() { coeffs().normalize(); }
127  inline Quaternion<Scalar> normalized() const { return Quaternion<Scalar>(coeffs().normalized()); }
128 
134  template<class OtherDerived> inline Scalar dot(const QuaternionBase<OtherDerived>& other) const { return coeffs().dot(other.coeffs()); }
135 
136  template<class OtherDerived> Scalar angularDistance(const QuaternionBase<OtherDerived>& other) const;
137 
139  Matrix3 toRotationMatrix() const;
140 
142  template<typename Derived1, typename Derived2>
143  Derived& setFromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
144 
145  template<class OtherDerived> EIGEN_STRONG_INLINE Quaternion<Scalar> operator* (const QuaternionBase<OtherDerived>& q) const;
146  template<class OtherDerived> EIGEN_STRONG_INLINE Derived& operator*= (const QuaternionBase<OtherDerived>& q);
147 
149  Quaternion<Scalar> inverse() const;
150 
152  Quaternion<Scalar> conjugate() const;
153 
154  template<class OtherDerived> Quaternion<Scalar> slerp(const Scalar& t, const QuaternionBase<OtherDerived>& other) const;
155 
160  template<class OtherDerived>
161  bool isApprox(const QuaternionBase<OtherDerived>& other, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const
162  { return coeffs().isApprox(other.coeffs(), prec); }
163 
165  EIGEN_STRONG_INLINE Vector3 _transformVector(const Vector3& v) const;
166 
172  template<typename NewScalarType>
173  inline typename internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type cast() const
174  {
175  return typename internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type(derived());
176  }
177 
178 #ifdef EIGEN_QUATERNIONBASE_PLUGIN
179 # include EIGEN_QUATERNIONBASE_PLUGIN
180 #endif
181 };
182 
183 /***************************************************************************
184 * Definition/implementation of Quaternion<Scalar>
185 ***************************************************************************/
186 
212 namespace internal {
213 template<typename _Scalar,int _Options>
214 struct traits<Quaternion<_Scalar,_Options> >
215 {
216  typedef Quaternion<_Scalar,_Options> PlainObject;
217  typedef _Scalar Scalar;
218  typedef Matrix<_Scalar,4,1,_Options> Coefficients;
219  enum{
220  Alignment = internal::traits<Coefficients>::Alignment,
221  Flags = LvalueBit
222  };
223 };
224 }
225 
226 template<typename _Scalar, int _Options>
227 class Quaternion : public QuaternionBase<Quaternion<_Scalar,_Options> >
228 {
229 public:
231  enum { NeedsAlignment = internal::traits<Quaternion>::Alignment>0 };
232 
233  typedef _Scalar Scalar;
234 
235  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Quaternion)
236  using Base::operator*=;
237 
238  typedef typename internal::traits<Quaternion>::Coefficients Coefficients;
239  typedef typename Base::AngleAxisType AngleAxisType;
240 
242  inline Quaternion() {}
243 
251  inline Quaternion(const Scalar& w, const Scalar& x, const Scalar& y, const Scalar& z) : m_coeffs(x, y, z, w){}
252 
254  explicit inline Quaternion(const Scalar* data) : m_coeffs(data) {}
255 
257  template<class Derived> EIGEN_STRONG_INLINE Quaternion(const QuaternionBase<Derived>& other) { this->Base::operator=(other); }
258 
260  explicit inline Quaternion(const AngleAxisType& aa) { *this = aa; }
261 
266  template<typename Derived>
267  explicit inline Quaternion(const MatrixBase<Derived>& other) { *this = other; }
268 
270  template<typename OtherScalar, int OtherOptions>
271  explicit inline Quaternion(const Quaternion<OtherScalar, OtherOptions>& other)
272  { m_coeffs = other.coeffs().template cast<Scalar>(); }
273 
274  template<typename Derived1, typename Derived2>
275  static Quaternion FromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
276 
277  inline Coefficients& coeffs() { return m_coeffs;}
278  inline const Coefficients& coeffs() const { return m_coeffs;}
279 
280  EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(NeedsAlignment)
281 
282 #ifdef EIGEN_QUATERNION_PLUGIN
283 # include EIGEN_QUATERNION_PLUGIN
284 #endif
285 
286 protected:
287  Coefficients m_coeffs;
288 
289 #ifndef EIGEN_PARSED_BY_DOXYGEN
290  static EIGEN_STRONG_INLINE void _check_template_params()
291  {
292  EIGEN_STATIC_ASSERT( (_Options & DontAlign) == _Options,
293  INVALID_MATRIX_TEMPLATE_PARAMETERS)
294  }
295 #endif
296 };
297 
304 
305 /***************************************************************************
306 * Specialization of Map<Quaternion<Scalar>>
307 ***************************************************************************/
308 
309 namespace internal {
310  template<typename _Scalar, int _Options>
311  struct traits<Map<Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
312  {
313  typedef Map<Matrix<_Scalar,4,1>, _Options> Coefficients;
314  };
315 }
316 
317 namespace internal {
318  template<typename _Scalar, int _Options>
319  struct traits<Map<const Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
320  {
321  typedef Map<const Matrix<_Scalar,4,1>, _Options> Coefficients;
322  typedef traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> > TraitsBase;
323  enum {
324  Flags = TraitsBase::Flags & ~LvalueBit
325  };
326  };
327 }
328 
340 template<typename _Scalar, int _Options>
341 class Map<const Quaternion<_Scalar>, _Options >
342  : public QuaternionBase<Map<const Quaternion<_Scalar>, _Options> >
343 {
344  public:
345  typedef QuaternionBase<Map<const Quaternion<_Scalar>, _Options> > Base;
346 
347  typedef _Scalar Scalar;
348  typedef typename internal::traits<Map>::Coefficients Coefficients;
349  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
350  using Base::operator*=;
351 
358  explicit EIGEN_STRONG_INLINE Map(const Scalar* coeffs) : m_coeffs(coeffs) {}
359 
360  inline const Coefficients& coeffs() const { return m_coeffs;}
361 
362  protected:
363  const Coefficients m_coeffs;
364 };
365 
377 template<typename _Scalar, int _Options>
378 class Map<Quaternion<_Scalar>, _Options >
379  : public QuaternionBase<Map<Quaternion<_Scalar>, _Options> >
380 {
381  public:
382  typedef QuaternionBase<Map<Quaternion<_Scalar>, _Options> > Base;
383 
384  typedef _Scalar Scalar;
385  typedef typename internal::traits<Map>::Coefficients Coefficients;
386  EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
387  using Base::operator*=;
388 
395  explicit EIGEN_STRONG_INLINE Map(Scalar* coeffs) : m_coeffs(coeffs) {}
396 
397  inline Coefficients& coeffs() { return m_coeffs; }
398  inline const Coefficients& coeffs() const { return m_coeffs; }
399 
400  protected:
401  Coefficients m_coeffs;
402 };
403 
412 typedef Map<Quaternion<float>, Aligned> QuaternionMapAlignedf;
415 typedef Map<Quaternion<double>, Aligned> QuaternionMapAlignedd;
416 
417 /***************************************************************************
418 * Implementation of QuaternionBase methods
419 ***************************************************************************/
420 
421 // Generic Quaternion * Quaternion product
422 // This product can be specialized for a given architecture via the Arch template argument.
423 namespace internal {
424 template<int Arch, class Derived1, class Derived2, typename Scalar, int _Options> struct quat_product
425 {
426  static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived1>& a, const QuaternionBase<Derived2>& b){
427  return Quaternion<Scalar>
428  (
429  a.w() * b.w() - a.x() * b.x() - a.y() * b.y() - a.z() * b.z(),
430  a.w() * b.x() + a.x() * b.w() + a.y() * b.z() - a.z() * b.y(),
431  a.w() * b.y() + a.y() * b.w() + a.z() * b.x() - a.x() * b.z(),
432  a.w() * b.z() + a.z() * b.w() + a.x() * b.y() - a.y() * b.x()
433  );
434  }
435 };
436 }
437 
439 template <class Derived>
440 template <class OtherDerived>
443 {
444  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename OtherDerived::Scalar>::value),
445  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
446  return internal::quat_product<Architecture::Target, Derived, OtherDerived,
447  typename internal::traits<Derived>::Scalar,
448  EIGEN_PLAIN_ENUM_MIN(internal::traits<Derived>::Alignment, internal::traits<OtherDerived>::Alignment)>::run(*this, other);
449 }
450 
452 template <class Derived>
453 template <class OtherDerived>
454 EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator*= (const QuaternionBase<OtherDerived>& other)
455 {
456  derived() = derived() * other.derived();
457  return derived();
458 }
459 
467 template <class Derived>
468 EIGEN_STRONG_INLINE typename QuaternionBase<Derived>::Vector3
470 {
471  // Note that this algorithm comes from the optimization by hand
472  // of the conversion to a Matrix followed by a Matrix/Vector product.
473  // It appears to be much faster than the common algorithm found
474  // in the literature (30 versus 39 flops). It also requires two
475  // Vector3 as temporaries.
476  Vector3 uv = this->vec().cross(v);
477  uv += uv;
478  return v + this->w() * uv + this->vec().cross(uv);
479 }
480 
481 template<class Derived>
483 {
484  coeffs() = other.coeffs();
485  return derived();
486 }
487 
488 template<class Derived>
489 template<class OtherDerived>
490 EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const QuaternionBase<OtherDerived>& other)
491 {
492  coeffs() = other.coeffs();
493  return derived();
494 }
495 
498 template<class Derived>
499 EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const AngleAxisType& aa)
500 {
501  using std::cos;
502  using std::sin;
503  Scalar ha = Scalar(0.5)*aa.angle(); // Scalar(0.5) to suppress precision loss warnings
504  this->w() = cos(ha);
505  this->vec() = sin(ha) * aa.axis();
506  return derived();
507 }
508 
515 template<class Derived>
516 template<class MatrixDerived>
518 {
519  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename MatrixDerived::Scalar>::value),
520  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
521  internal::quaternionbase_assign_impl<MatrixDerived>::run(*this, xpr.derived());
522  return derived();
523 }
524 
528 template<class Derived>
529 inline typename QuaternionBase<Derived>::Matrix3
531 {
532  // NOTE if inlined, then gcc 4.2 and 4.4 get rid of the temporary (not gcc 4.3 !!)
533  // if not inlined then the cost of the return by value is huge ~ +35%,
534  // however, not inlining this function is an order of magnitude slower, so
535  // it has to be inlined, and so the return by value is not an issue
536  Matrix3 res;
537 
538  const Scalar tx = Scalar(2)*this->x();
539  const Scalar ty = Scalar(2)*this->y();
540  const Scalar tz = Scalar(2)*this->z();
541  const Scalar twx = tx*this->w();
542  const Scalar twy = ty*this->w();
543  const Scalar twz = tz*this->w();
544  const Scalar txx = tx*this->x();
545  const Scalar txy = ty*this->x();
546  const Scalar txz = tz*this->x();
547  const Scalar tyy = ty*this->y();
548  const Scalar tyz = tz*this->y();
549  const Scalar tzz = tz*this->z();
550 
551  res.coeffRef(0,0) = Scalar(1)-(tyy+tzz);
552  res.coeffRef(0,1) = txy-twz;
553  res.coeffRef(0,2) = txz+twy;
554  res.coeffRef(1,0) = txy+twz;
555  res.coeffRef(1,1) = Scalar(1)-(txx+tzz);
556  res.coeffRef(1,2) = tyz-twx;
557  res.coeffRef(2,0) = txz-twy;
558  res.coeffRef(2,1) = tyz+twx;
559  res.coeffRef(2,2) = Scalar(1)-(txx+tyy);
560 
561  return res;
562 }
563 
574 template<class Derived>
575 template<typename Derived1, typename Derived2>
577 {
578  using std::sqrt;
579  Vector3 v0 = a.normalized();
580  Vector3 v1 = b.normalized();
581  Scalar c = v1.dot(v0);
582 
583  // if dot == -1, vectors are nearly opposites
584  // => accurately compute the rotation axis by computing the
585  // intersection of the two planes. This is done by solving:
586  // x^T v0 = 0
587  // x^T v1 = 0
588  // under the constraint:
589  // ||x|| = 1
590  // which yields a singular value problem
591  if (c < Scalar(-1)+NumTraits<Scalar>::dummy_precision())
592  {
593  c = numext::maxi(c,Scalar(-1));
594  Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
596  Vector3 axis = svd.matrixV().col(2);
597 
598  Scalar w2 = (Scalar(1)+c)*Scalar(0.5);
599  this->w() = sqrt(w2);
600  this->vec() = axis * sqrt(Scalar(1) - w2);
601  return derived();
602  }
603  Vector3 axis = v0.cross(v1);
604  Scalar s = sqrt((Scalar(1)+c)*Scalar(2));
605  Scalar invs = Scalar(1)/s;
606  this->vec() = axis * invs;
607  this->w() = s * Scalar(0.5);
608 
609  return derived();
610 }
611 
612 
623 template<typename Scalar, int Options>
624 template<typename Derived1, typename Derived2>
626 {
627  Quaternion quat;
628  quat.setFromTwoVectors(a, b);
629  return quat;
630 }
631 
632 
639 template <class Derived>
641 {
642  // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite() ??
643  Scalar n2 = this->squaredNorm();
644  if (n2 > Scalar(0))
645  return Quaternion<Scalar>(conjugate().coeffs() / n2);
646  else
647  {
648  // return an invalid result to flag the error
649  return Quaternion<Scalar>(Coefficients::Zero());
650  }
651 }
652 
653 // Generic conjugate of a Quaternion
654 namespace internal {
655 template<int Arch, class Derived, typename Scalar, int _Options> struct quat_conj
656 {
657  static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived>& q){
658  return Quaternion<Scalar>(q.w(),-q.x(),-q.y(),-q.z());
659  }
660 };
661 }
662 
669 template <class Derived>
672 {
673  return internal::quat_conj<Architecture::Target, Derived,
674  typename internal::traits<Derived>::Scalar,
675  internal::traits<Derived>::Alignment>::run(*this);
676 
677 }
678 
682 template <class Derived>
683 template <class OtherDerived>
684 inline typename internal::traits<Derived>::Scalar
686 {
687  using std::atan2;
688  using std::abs;
689  Quaternion<Scalar> d = (*this) * other.conjugate();
690  return Scalar(2) * atan2( d.vec().norm(), abs(d.w()) );
691 }
692 
693 
694 
701 template <class Derived>
702 template <class OtherDerived>
705 {
706  using std::acos;
707  using std::sin;
708  using std::abs;
709  static const Scalar one = Scalar(1) - NumTraits<Scalar>::epsilon();
710  Scalar d = this->dot(other);
711  Scalar absD = abs(d);
712 
713  Scalar scale0;
714  Scalar scale1;
715 
716  if(absD>=one)
717  {
718  scale0 = Scalar(1) - t;
719  scale1 = t;
720  }
721  else
722  {
723  // theta is the angle between the 2 quaternions
724  Scalar theta = acos(absD);
725  Scalar sinTheta = sin(theta);
726 
727  scale0 = sin( ( Scalar(1) - t ) * theta) / sinTheta;
728  scale1 = sin( ( t * theta) ) / sinTheta;
729  }
730  if(d<Scalar(0)) scale1 = -scale1;
731 
732  return Quaternion<Scalar>(scale0 * coeffs() + scale1 * other.coeffs());
733 }
734 
735 namespace internal {
736 
737 // set from a rotation matrix
738 template<typename Other>
739 struct quaternionbase_assign_impl<Other,3,3>
740 {
741  typedef typename Other::Scalar Scalar;
742  template<class Derived> static inline void run(QuaternionBase<Derived>& q, const Other& a_mat)
743  {
744  const typename internal::nested_eval<Other,2>::type mat(a_mat);
745  using std::sqrt;
746  // This algorithm comes from "Quaternion Calculus and Fast Animation",
747  // Ken Shoemake, 1987 SIGGRAPH course notes
748  Scalar t = mat.trace();
749  if (t > Scalar(0))
750  {
751  t = sqrt(t + Scalar(1.0));
752  q.w() = Scalar(0.5)*t;
753  t = Scalar(0.5)/t;
754  q.x() = (mat.coeff(2,1) - mat.coeff(1,2)) * t;
755  q.y() = (mat.coeff(0,2) - mat.coeff(2,0)) * t;
756  q.z() = (mat.coeff(1,0) - mat.coeff(0,1)) * t;
757  }
758  else
759  {
760  Index i = 0;
761  if (mat.coeff(1,1) > mat.coeff(0,0))
762  i = 1;
763  if (mat.coeff(2,2) > mat.coeff(i,i))
764  i = 2;
765  Index j = (i+1)%3;
766  Index k = (j+1)%3;
767 
768  t = sqrt(mat.coeff(i,i)-mat.coeff(j,j)-mat.coeff(k,k) + Scalar(1.0));
769  q.coeffs().coeffRef(i) = Scalar(0.5) * t;
770  t = Scalar(0.5)/t;
771  q.w() = (mat.coeff(k,j)-mat.coeff(j,k))*t;
772  q.coeffs().coeffRef(j) = (mat.coeff(j,i)+mat.coeff(i,j))*t;
773  q.coeffs().coeffRef(k) = (mat.coeff(k,i)+mat.coeff(i,k))*t;
774  }
775  }
776 };
777 
778 // set from a vector of coefficients assumed to be a quaternion
779 template<typename Other>
780 struct quaternionbase_assign_impl<Other,4,1>
781 {
782  typedef typename Other::Scalar Scalar;
783  template<class Derived> static inline void run(QuaternionBase<Derived>& q, const Other& vec)
784  {
785  q.coeffs() = vec;
786  }
787 };
788 
789 } // end namespace internal
790 
791 } // end namespace Eigen
792 
793 #endif // EIGEN_QUATERNION_H
Scalar x() const
Definition: Quaternion.h:61
Scalar dot(const QuaternionBase< OtherDerived > &other) const
Definition: Quaternion.h:134
Vector3 _transformVector(const Vector3 &v) const
Definition: Quaternion.h:469
Scalar y() const
Definition: Quaternion.h:63
Scalar norm() const
Definition: Quaternion.h:120
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:89
Matrix< Scalar, 3, 1 > Vector3
Definition: Quaternion.h:52
Scalar angle() const
Definition: AngleAxis.h:89
Quaternion< double > Quaterniond
Definition: Quaternion.h:303
const unsigned int LvalueBit
Definition: Constants.h:138
Definition: LDLT.h:16
const VectorBlock< const Coefficients, 3 > vec() const
Definition: Quaternion.h:79
internal::cast_return_type< Derived, Quaternion< NewScalarType > >::type cast() const
Definition: Quaternion.h:173
const internal::traits< Derived >::Coefficients & coeffs() const
Definition: Quaternion.h:85
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:107
AngleAxis< Scalar > AngleAxisType
Definition: Quaternion.h:56
Scalar z() const
Definition: Quaternion.h:65
Quaternion()
Definition: Quaternion.h:242
Derived & operator*=(const QuaternionBase< OtherDerived > &q)
Definition: Quaternion.h:454
static Quaternion< Scalar > Identity()
Definition: Quaternion.h:106
Scalar squaredNorm() const
Definition: Quaternion.h:115
Map(const Scalar *coeffs)
Definition: Quaternion.h:358
const Vector3 & axis() const
Definition: AngleAxis.h:94
void normalize()
Definition: Quaternion.h:124
Derived & setFromTwoVectors(const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
Definition: Quaternion.h:576
Expression of a fixed-size or dynamic-size sub-vector.
Definition: ForwardDeclarations.h:87
internal::traits< Derived >::Coefficients & coeffs()
Definition: Quaternion.h:88
Scalar & x()
Definition: Quaternion.h:70
Map(Scalar *coeffs)
Definition: Quaternion.h:395
bool isApprox(const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Quaternion.h:161
const PlainObject normalized() const
Definition: Dot.h:114
const MatrixVType & matrixV() const
Definition: SVDBase.h:99
Quaternion(const Quaternion< OtherScalar, OtherOptions > &other)
Definition: Quaternion.h:271
Quaternion< Scalar > conjugate() const
Definition: Quaternion.h:671
Scalar & z()
Definition: Quaternion.h:74
Base class for quaternion expressions.
Definition: ForwardDeclarations.h:265
Quaternion< Scalar > inverse() const
Definition: Quaternion.h:640
Map< Quaternion< float >, 0 > QuaternionMapf
Definition: Quaternion.h:406
Definition: Eigen_Colamd.h:54
Definition: Constants.h:326
Quaternion(const AngleAxisType &aa)
Definition: Quaternion.h:260
The quaternion class used to represent 3D orientations and rotations.
Definition: ForwardDeclarations.h:270
Quaternion(const MatrixBase< Derived > &other)
Definition: Quaternion.h:267
Map< Quaternion< double >, 0 > QuaternionMapd
Definition: Quaternion.h:409
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: ForwardDeclarations.h:255
Scalar w() const
Definition: Quaternion.h:67
Matrix3 toRotationMatrix() const
Definition: Quaternion.h:530
Map< Quaternion< float >, Aligned > QuaternionMapAlignedf
Definition: Quaternion.h:412
QuaternionBase & setIdentity()
Definition: Quaternion.h:110
Scalar & w()
Definition: Quaternion.h:76
Definition: Constants.h:235
Quaternion< Scalar > normalized() const
Definition: Quaternion.h:127
Scalar & y()
Definition: Quaternion.h:72
Matrix< Scalar, 3, 3 > Matrix3
Definition: Quaternion.h:54
VectorBlock< Coefficients, 3 > vec()
Definition: Quaternion.h:82
Definition: Constants.h:387
Quaternion(const Scalar &w, const Scalar &x, const Scalar &y, const Scalar &z)
Definition: Quaternion.h:251
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:267
Quaternion(const Scalar *data)
Definition: Quaternion.h:254
Quaternion< float > Quaternionf
Definition: Quaternion.h:300
Map< Quaternion< double >, Aligned > QuaternionMapAlignedd
Definition: Quaternion.h:415
Quaternion(const QuaternionBase< Derived > &other)
Definition: Quaternion.h:257