Medial Code Documentation
Loading...
Searching...
No Matches
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
13namespace Eigen {
14
15
16/***************************************************************************
17* Definition of QuaternionBase<Derived>
18* The implementation is at the end of the file
19***************************************************************************/
20
21namespace internal {
22template<typename Other,
23 int OtherRows=Other::RowsAtCompileTime,
24 int OtherCols=Other::ColsAtCompileTime>
26}
27
34template<class Derived>
35class QuaternionBase : public RotationBase<Derived, 3>
36{
37 public:
39
40 using Base::operator*;
41 using Base::derived;
42
43 typedef typename internal::traits<Derived>::Scalar Scalar;
44 typedef typename NumTraits<Scalar>::Real RealScalar;
45 typedef typename internal::traits<Derived>::Coefficients Coefficients;
46 enum {
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(); }
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
140
142 template<typename Derived1, typename Derived2>
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
150
153
154 template<class OtherDerived> Quaternion<Scalar> slerp(const Scalar& t, const QuaternionBase<OtherDerived>& other) const;
155
160 template<class OtherDerived>
162 { return coeffs().isApprox(other.coeffs(), prec); }
163
165 EIGEN_STRONG_INLINE Vector3 _transformVector(const Vector3& v) const;
166
172 template<typename NewScalarType>
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
212namespace internal {
213template<typename _Scalar,int _Options>
214struct traits<Quaternion<_Scalar,_Options> >
215{
217 typedef _Scalar Scalar;
219 enum{
221 Flags = LvalueBit
222 };
223};
224}
225
226template<typename _Scalar, int _Options>
227class Quaternion : public QuaternionBase<Quaternion<_Scalar,_Options> >
228{
229public:
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>
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
286protected:
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
309namespace 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
317namespace 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 {
322 typedef traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> > TraitsBase;
323 enum {
324 Flags = TraitsBase::Flags & ~LvalueBit
325 };
326 };
327}
328
340template<typename _Scalar, int _Options>
341class Map<const Quaternion<_Scalar>, _Options >
342 : public QuaternionBase<Map<const Quaternion<_Scalar>, _Options> >
343{
344 public:
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
377template<typename _Scalar, int _Options>
378class 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
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.
423namespace internal {
424template<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
439template <class Derived>
440template <class OtherDerived>
443{
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,
449}
450
452template <class Derived>
453template <class OtherDerived>
455{
456 derived() = derived() * other.derived();
457 return derived();
458}
459
467template <class Derived>
468EIGEN_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
481template<class Derived>
483{
484 coeffs() = other.coeffs();
485 return derived();
486}
487
488template<class Derived>
489template<class OtherDerived>
490EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const QuaternionBase<OtherDerived>& other)
491{
492 coeffs() = other.coeffs();
493 return derived();
494}
495
498template<class Derived>
499EIGEN_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
515template<class Derived>
516template<class MatrixDerived>
518{
520 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
522 return derived();
523}
524
528template<class Derived>
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
574template<class Derived>
575template<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
623template<typename Scalar, int Options>
624template<typename Derived1, typename Derived2>
631
632
639template <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
654namespace internal {
655template<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
669template <class Derived>
672{
673 return internal::quat_conj<Architecture::Target, Derived,
676
677}
678
682template <class Derived>
683template <class OtherDerived>
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
701template <class Derived>
702template <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
735namespace internal {
736
737// set from a rotation matrix
738template<typename Other>
740{
741 typedef typename Other::Scalar Scalar;
742 template<class Derived> static inline void run(QuaternionBase<Derived>& q, const Other& a_mat)
743 {
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
779template<typename Other>
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
\geometry_module
Definition AngleAxis.h:50
EIGEN_STRONG_INLINE Map(Scalar *coeffs)
Constructs a Mapped Quaternion object from the pointer coeffs.
Definition Quaternion.h:395
EIGEN_STRONG_INLINE Map(const Scalar *coeffs)
Constructs a Mapped Quaternion object from the pointer coeffs.
Definition Quaternion.h:358
A matrix or vector expression mapping an existing array of data.
Definition Map.h:91
\geometry_module
Definition Quaternion.h:36
Scalar y() const
Definition Quaternion.h:63
Scalar squaredNorm() const
Definition Quaternion.h:115
Scalar & y()
Definition Quaternion.h:72
QuaternionBase & setIdentity()
Definition Quaternion.h:110
Quaternion< Scalar > normalized() const
Definition Quaternion.h:127
internal::traits< Derived >::Coefficients & coeffs()
Definition Quaternion.h:88
internal::cast_return_type< Derived, Quaternion< NewScalarType > >::type cast() const
Definition Quaternion.h:173
VectorBlock< Coefficients, 3 > vec()
Definition Quaternion.h:82
const VectorBlock< const Coefficients, 3 > vec() const
Definition Quaternion.h:79
static Quaternion< Scalar > Identity()
Definition Quaternion.h:106
Scalar w() const
Definition Quaternion.h:67
Quaternion< Scalar > conjugate() const
Definition Quaternion.h:671
bool isApprox(const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition Quaternion.h:161
Scalar & z()
Definition Quaternion.h:74
void normalize()
Normalizes the quaternion *this.
Definition Quaternion.h:124
Derived & setFromTwoVectors(const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
Sets *this to be a quaternion representing a rotation between the two arbitrary vectors a and b.
Definition Quaternion.h:576
Derived & operator=(const AngleAxisType &aa)
Set *this from an angle-axis aa and returns a reference to *this.
Definition Quaternion.h:499
Scalar & x()
Definition Quaternion.h:70
Matrix3 toRotationMatrix() const
Convert the quaternion to a 3x3 rotation matrix.
Definition Quaternion.h:530
Matrix< Scalar, 3, 1 > Vector3
the type of a 3D vector
Definition Quaternion.h:52
Scalar dot(const QuaternionBase< OtherDerived > &other) const
Definition Quaternion.h:134
Scalar norm() const
Definition Quaternion.h:120
Quaternion< Scalar > inverse() const
Definition Quaternion.h:640
Scalar z() const
Definition Quaternion.h:65
Matrix< Scalar, 3, 3 > Matrix3
the equivalent rotation matrix type
Definition Quaternion.h:54
const internal::traits< Derived >::Coefficients & coeffs() const
Definition Quaternion.h:85
EIGEN_STRONG_INLINE Derived & operator*=(const QuaternionBase< OtherDerived > &q)
Definition Quaternion.h:454
EIGEN_STRONG_INLINE Vector3 _transformVector(const Vector3 &v) const
return the result vector of v through the rotation
Definition Quaternion.h:469
AngleAxis< Scalar > AngleAxisType
the equivalent angle-axis type
Definition Quaternion.h:56
Scalar x() const
Definition Quaternion.h:61
Scalar & w()
Definition Quaternion.h:76
\geometry_module
Definition Quaternion.h:228
EIGEN_STRONG_INLINE Quaternion(const QuaternionBase< Derived > &other)
Copy constructor.
Definition Quaternion.h:257
Quaternion(const AngleAxisType &aa)
Constructs and initializes a quaternion from the angle-axis aa.
Definition Quaternion.h:260
Quaternion(const Quaternion< OtherScalar, OtherOptions > &other)
Explicit copy constructor with scalar conversion.
Definition Quaternion.h:271
Quaternion(const MatrixBase< Derived > &other)
Constructs and initializes a quaternion from either:
Definition Quaternion.h:267
Quaternion()
Default constructor leaving the quaternion uninitialized.
Definition Quaternion.h:242
Quaternion(const Scalar &w, const Scalar &x, const Scalar &y, const Scalar &z)
Constructs and initializes the quaternion from its four coefficients w, x, y and z.
Definition Quaternion.h:251
Quaternion(const Scalar *data)
Constructs and initialize a quaternion from the array data.
Definition Quaternion.h:254
Common base class for compact rotation representations.
Definition RotationBase.h:30
friend RotationMatrixType operator*(const EigenBase< OtherDerived > &l, const Derived &r)
Definition RotationBase.h:76
Pseudo expression representing a solving operation.
Definition Solve.h:63
@ Aligned
Definition Constants.h:235
@ DontAlign
Don't require alignment for the matrix itself (the array of coefficients, if dynamically allocated,...
Definition Constants.h:326
@ AutoAlign
Align the matrix itself if it is vectorizable fixed-size.
Definition Constants.h:324
@ ComputeFullV
Used in JacobiSVD to indicate that the square matrix V is to be computed.
Definition Constants.h:387
const unsigned int LvalueBit
Means the expression has a coeffRef() method, i.e.
Definition Constants.h:138
Holds information about the various numeric (i.e.
Definition NumTraits.h:108
Definition XprHelper.h:500
Definition Meta.h:39
Definition Quaternion.h:656
Definition Quaternion.h:425
Definition ForwardDeclarations.h:17
Definition Meta.h:30