Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
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 <[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// IWYU pragma: private
14#include "./InternalHeaderCheck.h"
15
16namespace Eigen {
17
18/***************************************************************************
19 * Definition of QuaternionBase<Derived>
20 * The implementation is at the end of the file
21 ***************************************************************************/
22
23namespace internal {
24template <typename Other, int OtherRows = Other::RowsAtCompileTime, int OtherCols = Other::ColsAtCompileTime>
25struct quaternionbase_assign_impl;
26}
27
34template <class Derived>
35class QuaternionBase : public RotationBase<Derived, 3> {
36 public:
38
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 typedef typename Coefficients::CoeffReturnType CoeffReturnType;
46 typedef std::conditional_t<bool(internal::traits<Derived>::Flags& LvalueBit), Scalar&, CoeffReturnType>
47 NonConstCoeffReturnType;
48
49 enum { Flags = Eigen::internal::traits<Derived>::Flags };
50
51 // typedef typename Matrix<Scalar,4,1> Coefficients;
58
60 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline CoeffReturnType x() const { return this->derived().coeffs().coeff(0); }
62 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline CoeffReturnType y() const { return this->derived().coeffs().coeff(1); }
64 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline CoeffReturnType z() const { return this->derived().coeffs().coeff(2); }
66 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline CoeffReturnType w() const { return this->derived().coeffs().coeff(3); }
67
69 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline NonConstCoeffReturnType x() { return this->derived().coeffs().x(); }
71 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline NonConstCoeffReturnType y() { return this->derived().coeffs().y(); }
73 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline NonConstCoeffReturnType z() { return this->derived().coeffs().z(); }
75 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline NonConstCoeffReturnType w() { return this->derived().coeffs().w(); }
76
78 EIGEN_DEVICE_FUNC inline const VectorBlock<const Coefficients, 3> vec() const { return coeffs().template head<3>(); }
79
81 EIGEN_DEVICE_FUNC inline VectorBlock<Coefficients, 3> vec() { return coeffs().template head<3>(); }
82
84 EIGEN_DEVICE_FUNC inline const typename internal::traits<Derived>::Coefficients& coeffs() const {
85 return derived().coeffs();
86 }
87
89 EIGEN_DEVICE_FUNC inline typename internal::traits<Derived>::Coefficients& coeffs() { return derived().coeffs(); }
90
91 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE QuaternionBase<Derived>& operator=(const QuaternionBase<Derived>& other);
92 template <class OtherDerived>
93 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& operator=(const QuaternionBase<OtherDerived>& other);
94
95 // disabled this copy operator as it is giving very strange compilation errors when compiling
96 // test_stdvector with GCC 4.4.2. This looks like a GCC bug though, so feel free to re-enable it if it's
97 // useful; however notice that we already have the templated operator= above and e.g. in MatrixBase
98 // we didn't have to add, in addition to templated operator=, such a non-templated copy operator.
99 // Derived& operator=(const QuaternionBase& other)
100 // { return operator=<Derived>(other); }
101
102 EIGEN_DEVICE_FUNC Derived& operator=(const AngleAxisType& aa);
103 template <class OtherDerived>
104 EIGEN_DEVICE_FUNC Derived& operator=(const MatrixBase<OtherDerived>& m);
105
109 EIGEN_DEVICE_FUNC static inline Quaternion<Scalar> Identity() {
110 return Quaternion<Scalar>(Scalar(1), Scalar(0), Scalar(0), Scalar(0));
111 }
112
115 EIGEN_DEVICE_FUNC inline QuaternionBase& setIdentity() {
116 coeffs() << Scalar(0), Scalar(0), Scalar(0), Scalar(1);
117 return *this;
118 }
119
123 EIGEN_DEVICE_FUNC inline Scalar squaredNorm() const { return coeffs().squaredNorm(); }
124
128 EIGEN_DEVICE_FUNC inline Scalar norm() const { return coeffs().norm(); }
129
132 EIGEN_DEVICE_FUNC inline void normalize() { coeffs().normalize(); }
135 EIGEN_DEVICE_FUNC inline Quaternion<Scalar> normalized() const { return Quaternion<Scalar>(coeffs().normalized()); }
136
142 template <class OtherDerived>
143 EIGEN_DEVICE_FUNC inline Scalar dot(const QuaternionBase<OtherDerived>& other) const {
144 return coeffs().dot(other.coeffs());
145 }
146
147 template <class OtherDerived>
148 EIGEN_DEVICE_FUNC Scalar angularDistance(const QuaternionBase<OtherDerived>& other) const;
149
151 EIGEN_DEVICE_FUNC inline Matrix3 toRotationMatrix() const;
152
154 template <typename Derived1, typename Derived2>
155 EIGEN_DEVICE_FUNC Derived& setFromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
156
157 template <class OtherDerived>
158 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion<Scalar> operator*(const QuaternionBase<OtherDerived>& q) const;
159 template <class OtherDerived>
160 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& operator*=(const QuaternionBase<OtherDerived>& q);
161
163 EIGEN_DEVICE_FUNC Quaternion<Scalar> inverse() const;
164
166 EIGEN_DEVICE_FUNC Quaternion<Scalar> conjugate() const;
167
168 template <class OtherDerived>
169 EIGEN_DEVICE_FUNC Quaternion<Scalar> slerp(const Scalar& t, const QuaternionBase<OtherDerived>& other) const;
170
175 template <class OtherDerived>
176 EIGEN_DEVICE_FUNC inline bool operator==(const QuaternionBase<OtherDerived>& other) const {
177 return coeffs() == other.coeffs();
178 }
179
184 template <class OtherDerived>
185 EIGEN_DEVICE_FUNC inline bool operator!=(const QuaternionBase<OtherDerived>& other) const {
186 return coeffs() != other.coeffs();
187 }
188
193 template <class OtherDerived>
194 EIGEN_DEVICE_FUNC bool isApprox(const QuaternionBase<OtherDerived>& other,
195 const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const {
196 return coeffs().isApprox(other.coeffs(), prec);
197 }
198
200 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Vector3 _transformVector(const Vector3& v) const;
201
202#ifdef EIGEN_PARSED_BY_DOXYGEN
208 template <typename NewScalarType>
209 EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<Derived, Quaternion<NewScalarType> >::type cast() const;
210
211#else
212
213 template <typename NewScalarType>
214 EIGEN_DEVICE_FUNC inline std::enable_if_t<internal::is_same<Scalar, NewScalarType>::value, const Derived&> cast()
215 const {
216 return derived();
217 }
218
219 template <typename NewScalarType>
220 EIGEN_DEVICE_FUNC inline std::enable_if_t<!internal::is_same<Scalar, NewScalarType>::value,
222 cast() const {
224 }
225#endif
226
227#ifndef EIGEN_NO_IO
228 friend std::ostream& operator<<(std::ostream& s, const QuaternionBase<Derived>& q) {
229 s << q.x() << "i + " << q.y() << "j + " << q.z() << "k"
230 << " + " << q.w();
231 return s;
232 }
233#endif
234
235#ifdef EIGEN_QUATERNIONBASE_PLUGIN
236#include EIGEN_QUATERNIONBASE_PLUGIN
237#endif
238 protected:
239 EIGEN_DEFAULT_COPY_CONSTRUCTOR(QuaternionBase)
240 EIGEN_DEFAULT_EMPTY_CONSTRUCTOR_AND_DESTRUCTOR(QuaternionBase)
241};
242
243/***************************************************************************
244 * Definition/implementation of Quaternion<Scalar>
245 ***************************************************************************/
246
274namespace internal {
275template <typename Scalar_, int Options_>
276struct traits<Quaternion<Scalar_, Options_> > {
277 typedef Quaternion<Scalar_, Options_> PlainObject;
278 typedef Scalar_ Scalar;
279 typedef Matrix<Scalar_, 4, 1, Options_> Coefficients;
280 enum { Alignment = internal::traits<Coefficients>::Alignment, Flags = LvalueBit };
281};
282} // namespace internal
283
284template <typename Scalar_, int Options_>
285class Quaternion : public QuaternionBase<Quaternion<Scalar_, Options_> > {
286 public:
288 enum { NeedsAlignment = internal::traits<Quaternion>::Alignment > 0 };
289
290 typedef Scalar_ Scalar;
291
292 EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Quaternion)
293 using Base::operator*=;
294
295 typedef typename internal::traits<Quaternion>::Coefficients Coefficients;
296 typedef typename Base::AngleAxisType AngleAxisType;
297
299 EIGEN_DEVICE_FUNC inline Quaternion() {}
300
308 EIGEN_DEVICE_FUNC inline Quaternion(const Scalar& w, const Scalar& x, const Scalar& y, const Scalar& z)
309 : m_coeffs(x, y, z, w) {}
310
314 template <typename Derived>
315 EIGEN_DEVICE_FUNC inline Quaternion(const Scalar& w, const Eigen::MatrixBase<Derived>& vec)
316 : m_coeffs(vec.x(), vec.y(), vec.z(), w) {
317 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived, 3);
318 }
319
321 EIGEN_DEVICE_FUNC explicit inline Quaternion(const Scalar* data) : m_coeffs(data) {}
322
324 template <class Derived>
325 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion(const QuaternionBase<Derived>& other) {
326 this->Base::operator=(other);
327 }
328
330 EIGEN_DEVICE_FUNC explicit inline Quaternion(const AngleAxisType& aa) { *this = aa; }
331
336 template <typename Derived>
337 EIGEN_DEVICE_FUNC explicit inline Quaternion(const MatrixBase<Derived>& other) {
338 *this = other;
339 }
340
342 template <typename OtherScalar, int OtherOptions>
343 EIGEN_DEVICE_FUNC explicit inline Quaternion(const Quaternion<OtherScalar, OtherOptions>& other) {
344 m_coeffs = other.coeffs().template cast<Scalar>();
345 }
346
347 // We define a copy constructor, which means we don't get an implicit move constructor or assignment operator.
349 EIGEN_DEVICE_FUNC inline Quaternion(Quaternion&& other)
350 EIGEN_NOEXCEPT_IF(std::is_nothrow_move_constructible<Scalar>::value)
351 : m_coeffs(std::move(other.coeffs())) {}
352
354 EIGEN_DEVICE_FUNC Quaternion& operator=(Quaternion&& other)
355 EIGEN_NOEXCEPT_IF(std::is_nothrow_move_assignable<Scalar>::value) {
356 m_coeffs = std::move(other.coeffs());
357 return *this;
358 }
359
360 EIGEN_DEVICE_FUNC static Quaternion UnitRandom();
361
362 template <typename Derived1, typename Derived2>
363 EIGEN_DEVICE_FUNC static Quaternion FromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
364
365 EIGEN_DEVICE_FUNC inline Coefficients& coeffs() { return m_coeffs; }
366 EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs; }
367
368 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(bool(NeedsAlignment))
369
370#ifdef EIGEN_QUATERNION_PLUGIN
371#include EIGEN_QUATERNION_PLUGIN
372#endif
373
374 protected:
375 Coefficients m_coeffs;
376
377#ifndef EIGEN_PARSED_BY_DOXYGEN
378 EIGEN_STATIC_ASSERT((Options_ & DontAlign) == Options_, INVALID_MATRIX_TEMPLATE_PARAMETERS)
379#endif
380};
381
388
389/***************************************************************************
390 * Specialization of Map<Quaternion<Scalar>>
391 ***************************************************************************/
392
393namespace internal {
394template <typename Scalar_, int Options_>
395struct traits<Map<Quaternion<Scalar_>, Options_> >
396 : traits<Quaternion<Scalar_, (int(Options_) & Aligned) == Aligned ? AutoAlign : DontAlign> > {
397 typedef Map<Matrix<Scalar_, 4, 1>, Options_> Coefficients;
398};
399} // namespace internal
400
401namespace internal {
402template <typename Scalar_, int Options_>
403struct traits<Map<const Quaternion<Scalar_>, Options_> >
404 : traits<Quaternion<Scalar_, (int(Options_) & Aligned) == Aligned ? AutoAlign : DontAlign> > {
405 typedef Map<const Matrix<Scalar_, 4, 1>, Options_> Coefficients;
406 typedef traits<Quaternion<Scalar_, (int(Options_) & Aligned) == Aligned ? AutoAlign : DontAlign> > TraitsBase;
407 enum { Flags = TraitsBase::Flags & ~LvalueBit };
408};
409} // namespace internal
410
422template <typename Scalar_, int Options_>
423class Map<const Quaternion<Scalar_>, Options_> : public QuaternionBase<Map<const Quaternion<Scalar_>, Options_> > {
424 public:
426
427 typedef Scalar_ Scalar;
428 typedef typename internal::traits<Map>::Coefficients Coefficients;
429 EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
430 using Base::operator*=;
431
438 EIGEN_DEVICE_FUNC explicit EIGEN_STRONG_INLINE Map(const Scalar* coeffs) : m_coeffs(coeffs) {}
439
440 EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs; }
441
442 protected:
443 const Coefficients m_coeffs;
444};
445
457template <typename Scalar_, int Options_>
458class Map<Quaternion<Scalar_>, Options_> : public QuaternionBase<Map<Quaternion<Scalar_>, Options_> > {
459 public:
460 typedef QuaternionBase<Map<Quaternion<Scalar_>, Options_> > Base;
461
462 typedef Scalar_ Scalar;
463 typedef typename internal::traits<Map>::Coefficients Coefficients;
464 EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
465 using Base::operator*=;
466
473 EIGEN_DEVICE_FUNC explicit EIGEN_STRONG_INLINE Map(Scalar* coeffs) : m_coeffs(coeffs) {}
474
475 EIGEN_DEVICE_FUNC inline Coefficients& coeffs() { return m_coeffs; }
476 EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs; }
477
478 protected:
479 Coefficients m_coeffs;
480};
481
494
495/***************************************************************************
496 * Implementation of QuaternionBase methods
497 ***************************************************************************/
498
499// Generic Quaternion * Quaternion product
500// This product can be specialized for a given architecture via the Arch template argument.
501namespace internal {
502template <int Arch, class Derived1, class Derived2, typename Scalar>
503struct quat_product {
504 EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived1>& a,
505 const QuaternionBase<Derived2>& b) {
506 return Quaternion<Scalar>(a.w() * b.w() - a.x() * b.x() - a.y() * b.y() - a.z() * b.z(),
507 a.w() * b.x() + a.x() * b.w() + a.y() * b.z() - a.z() * b.y(),
508 a.w() * b.y() + a.y() * b.w() + a.z() * b.x() - a.x() * b.z(),
509 a.w() * b.z() + a.z() * b.w() + a.x() * b.y() - a.y() * b.x());
510 }
511};
512} // namespace internal
513
515template <class Derived>
516template <class OtherDerived>
517EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion<typename internal::traits<Derived>::Scalar>
519 EIGEN_STATIC_ASSERT(
520 (internal::is_same<typename Derived::Scalar, typename OtherDerived::Scalar>::value),
521 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
522 return internal::quat_product<Architecture::Target, Derived, OtherDerived,
523 typename internal::traits<Derived>::Scalar>::run(*this, other);
524}
525
527template <class Derived>
528template <class OtherDerived>
529EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator*=(
530 const QuaternionBase<OtherDerived>& other) {
531 derived() = derived() * other.derived();
532 return derived();
533}
534
542template <class Derived>
543EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE typename QuaternionBase<Derived>::Vector3
545 // Note that this algorithm comes from the optimization by hand
546 // of the conversion to a Matrix followed by a Matrix/Vector product.
547 // It appears to be much faster than the common algorithm found
548 // in the literature (30 versus 39 flops). It also requires two
549 // Vector3 as temporaries.
550 Vector3 uv = this->vec().cross(v);
551 uv += uv;
552 return v + this->w() * uv + this->vec().cross(uv);
553}
554
555template <class Derived>
556EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE QuaternionBase<Derived>& QuaternionBase<Derived>::operator=(
557 const QuaternionBase<Derived>& other) {
558 coeffs() = other.coeffs();
559 return derived();
560}
561
562template <class Derived>
563template <class OtherDerived>
564EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(
565 const QuaternionBase<OtherDerived>& other) {
566 coeffs() = other.coeffs();
567 return derived();
568}
569
572template <class Derived>
573EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const AngleAxisType& aa) {
574 EIGEN_USING_STD(cos)
575 EIGEN_USING_STD(sin)
576 Scalar ha = Scalar(0.5) * aa.angle(); // Scalar(0.5) to suppress precision loss warnings
577 this->w() = cos(ha);
578 this->vec() = sin(ha) * aa.axis();
579 return derived();
580}
581
588template <class Derived>
589template <class MatrixDerived>
590EIGEN_DEVICE_FUNC inline Derived& QuaternionBase<Derived>::operator=(const MatrixBase<MatrixDerived>& xpr) {
591 EIGEN_STATIC_ASSERT(
592 (internal::is_same<typename Derived::Scalar, typename MatrixDerived::Scalar>::value),
593 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
594 internal::quaternionbase_assign_impl<MatrixDerived>::run(*this, xpr.derived());
595 return derived();
596}
597
601template <class Derived>
603 void) const {
604 // NOTE if inlined, then gcc 4.2 and 4.4 get rid of the temporary (not gcc 4.3 !!)
605 // if not inlined then the cost of the return by value is huge ~ +35%,
606 // however, not inlining this function is an order of magnitude slower, so
607 // it has to be inlined, and so the return by value is not an issue
608 Matrix3 res;
609
610 const Scalar tx = Scalar(2) * this->x();
611 const Scalar ty = Scalar(2) * this->y();
612 const Scalar tz = Scalar(2) * this->z();
613 const Scalar twx = tx * this->w();
614 const Scalar twy = ty * this->w();
615 const Scalar twz = tz * this->w();
616 const Scalar txx = tx * this->x();
617 const Scalar txy = ty * this->x();
618 const Scalar txz = tz * this->x();
619 const Scalar tyy = ty * this->y();
620 const Scalar tyz = tz * this->y();
621 const Scalar tzz = tz * this->z();
622
623 res.coeffRef(0, 0) = Scalar(1) - (tyy + tzz);
624 res.coeffRef(0, 1) = txy - twz;
625 res.coeffRef(0, 2) = txz + twy;
626 res.coeffRef(1, 0) = txy + twz;
627 res.coeffRef(1, 1) = Scalar(1) - (txx + tzz);
628 res.coeffRef(1, 2) = tyz - twx;
629 res.coeffRef(2, 0) = txz - twy;
630 res.coeffRef(2, 1) = tyz + twx;
631 res.coeffRef(2, 2) = Scalar(1) - (txx + tyy);
632
633 return res;
634}
635
646template <class Derived>
647template <typename Derived1, typename Derived2>
648EIGEN_DEVICE_FUNC inline Derived& QuaternionBase<Derived>::setFromTwoVectors(const MatrixBase<Derived1>& a,
649 const MatrixBase<Derived2>& b) {
650 EIGEN_USING_STD(sqrt)
651 Vector3 v0 = a.normalized();
652 Vector3 v1 = b.normalized();
653 Scalar c = v1.dot(v0);
654
655 // if dot == -1, vectors are nearly opposites
656 // => accurately compute the rotation axis by computing the
657 // intersection of the two planes. This is done by solving:
658 // x^T v0 = 0
659 // x^T v1 = 0
660 // under the constraint:
661 // ||x|| = 1
662 // which yields a singular value problem
663 if (c < Scalar(-1) + NumTraits<Scalar>::dummy_precision()) {
664 c = numext::maxi(c, Scalar(-1));
666 m << v0.transpose(), v1.transpose();
668 Vector3 axis = svd.matrixV().col(2);
669
670 Scalar w2 = (Scalar(1) + c) * Scalar(0.5);
671 this->w() = sqrt(w2);
672 this->vec() = axis * sqrt(Scalar(1) - w2);
673 return derived();
674 }
675 Vector3 axis = v0.cross(v1);
676 Scalar s = sqrt((Scalar(1) + c) * Scalar(2));
677 Scalar invs = Scalar(1) / s;
678 this->vec() = axis * invs;
679 this->w() = s * Scalar(0.5);
680
681 return derived();
682}
683
688template <typename Scalar, int Options>
690 EIGEN_USING_STD(sqrt)
691 EIGEN_USING_STD(sin)
692 EIGEN_USING_STD(cos)
693 const Scalar u1 = internal::random<Scalar>(0, 1), u2 = internal::random<Scalar>(0, 2 * EIGEN_PI),
694 u3 = internal::random<Scalar>(0, 2 * EIGEN_PI);
695 const Scalar a = sqrt(Scalar(1) - u1), b = sqrt(u1);
696 return Quaternion(a * sin(u2), a * cos(u2), b * sin(u3), b * cos(u3));
697}
698
709template <typename Scalar, int Options>
710template <typename Derived1, typename Derived2>
712 const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b) {
713 Quaternion quat;
714 quat.setFromTwoVectors(a, b);
715 return quat;
716}
717
724template <class Derived>
726 const {
727 // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite() ??
728 Scalar n2 = this->squaredNorm();
729 if (n2 > Scalar(0))
730 return Quaternion<Scalar>(conjugate().coeffs() / n2);
731 else {
732 // return an invalid result to flag the error
733 return Quaternion<Scalar>(Coefficients::Zero());
734 }
735}
736
737// Generic conjugate of a Quaternion
738namespace internal {
739template <int Arch, class Derived, typename Scalar>
740struct quat_conj {
741 EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived>& q) {
742 return Quaternion<Scalar>(q.w(), -q.x(), -q.y(), -q.z());
743 }
744};
745} // namespace internal
746
753template <class Derived>
755 const {
756 return internal::quat_conj<Architecture::Target, Derived, typename internal::traits<Derived>::Scalar>::run(*this);
757}
758
762template <class Derived>
763template <class OtherDerived>
764EIGEN_DEVICE_FUNC inline typename internal::traits<Derived>::Scalar QuaternionBase<Derived>::angularDistance(
765 const QuaternionBase<OtherDerived>& other) const {
766 EIGEN_USING_STD(atan2)
767 Quaternion<Scalar> d = (*this) * other.conjugate();
768 return Scalar(2) * atan2(d.vec().norm(), numext::abs(d.w()));
769}
770
777template <class Derived>
778template <class OtherDerived>
780 const Scalar& t, const QuaternionBase<OtherDerived>& other) const {
781 EIGEN_USING_STD(acos)
782 EIGEN_USING_STD(sin)
783 const Scalar one = Scalar(1) - NumTraits<Scalar>::epsilon();
784 Scalar d = this->dot(other);
785 Scalar absD = numext::abs(d);
786
787 Scalar scale0;
788 Scalar scale1;
789
790 if (absD >= one) {
791 scale0 = Scalar(1) - t;
792 scale1 = t;
793 } else {
794 // theta is the angle between the 2 quaternions
795 Scalar theta = acos(absD);
796 Scalar sinTheta = sin(theta);
797
798 scale0 = sin((Scalar(1) - t) * theta) / sinTheta;
799 scale1 = sin((t * theta)) / sinTheta;
800 }
801 if (d < Scalar(0)) scale1 = -scale1;
802
803 return Quaternion<Scalar>(scale0 * coeffs() + scale1 * other.coeffs());
804}
805
806namespace internal {
807
808// set from a rotation matrix
809template <typename Other>
810struct quaternionbase_assign_impl<Other, 3, 3> {
811 typedef typename Other::Scalar Scalar;
812 template <class Derived>
813 EIGEN_DEVICE_FUNC static inline void run(QuaternionBase<Derived>& q, const Other& a_mat) {
814 const typename internal::nested_eval<Other, 2>::type mat(a_mat);
815 EIGEN_USING_STD(sqrt)
816 // This algorithm comes from "Quaternion Calculus and Fast Animation",
817 // Ken Shoemake, 1987 SIGGRAPH course notes
818 Scalar t = mat.trace();
819 if (t > Scalar(0)) {
820 t = sqrt(t + Scalar(1.0));
821 q.w() = Scalar(0.5) * t;
822 t = Scalar(0.5) / t;
823 q.x() = (mat.coeff(2, 1) - mat.coeff(1, 2)) * t;
824 q.y() = (mat.coeff(0, 2) - mat.coeff(2, 0)) * t;
825 q.z() = (mat.coeff(1, 0) - mat.coeff(0, 1)) * t;
826 } else {
827 Index i = 0;
828 if (mat.coeff(1, 1) > mat.coeff(0, 0)) i = 1;
829 if (mat.coeff(2, 2) > mat.coeff(i, i)) i = 2;
830 Index j = (i + 1) % 3;
831 Index k = (j + 1) % 3;
832
833 t = sqrt(mat.coeff(i, i) - mat.coeff(j, j) - mat.coeff(k, k) + Scalar(1.0));
834 q.coeffs().coeffRef(i) = Scalar(0.5) * t;
835 t = Scalar(0.5) / t;
836 q.w() = (mat.coeff(k, j) - mat.coeff(j, k)) * t;
837 q.coeffs().coeffRef(j) = (mat.coeff(j, i) + mat.coeff(i, j)) * t;
838 q.coeffs().coeffRef(k) = (mat.coeff(k, i) + mat.coeff(i, k)) * t;
839 }
840 }
841};
842
843// set from a vector of coefficients assumed to be a quaternion
844template <typename Other>
845struct quaternionbase_assign_impl<Other, 4, 1> {
846 typedef typename Other::Scalar Scalar;
847 template <class Derived>
848 EIGEN_DEVICE_FUNC static inline void run(QuaternionBase<Derived>& q, const Other& vec) {
849 q.coeffs() = vec;
850 }
851};
852
853} // end namespace internal
854
855} // end namespace Eigen
856
857#endif // EIGEN_QUATERNION_H
Represents a 3D rotation as a rotation angle around an arbitrary 3D axis.
Definition AngleAxis.h:52
const Vector3 & axis() const
Definition AngleAxis.h:99
Scalar angle() const
Definition AngleAxis.h:94
Derived & derived()
Definition EigenBase.h:49
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition JacobiSVD.h:500
Map(Scalar *coeffs)
Definition Quaternion.h:473
Map(const Scalar *coeffs)
Definition Quaternion.h:438
A matrix or vector expression mapping an existing array of data.
Definition Map.h:96
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
const PlainObject normalized() const
Definition Dot.h:113
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
constexpr Scalar & coeffRef(Index rowId, Index colId)
Definition PlainObjectBase.h:217
Base class for quaternion expressions.
Definition Quaternion.h:35
Matrix< Scalar, 3, 3 > Matrix3
Definition Quaternion.h:55
Scalar squaredNorm() const
Definition Quaternion.h:123
QuaternionBase & setIdentity()
Definition Quaternion.h:115
Quaternion< Scalar > normalized() const
Definition Quaternion.h:135
EIGEN_CONSTEXPR CoeffReturnType x() const
Definition Quaternion.h:60
internal::traits< Derived >::Coefficients & coeffs()
Definition Quaternion.h:89
Matrix< Scalar, 3, 1 > Vector3
Definition Quaternion.h:53
internal::cast_return_type< Derived, Quaternion< NewScalarType > >::type cast() const
VectorBlock< Coefficients, 3 > vec()
Definition Quaternion.h:81
const VectorBlock< const Coefficients, 3 > vec() const
Definition Quaternion.h:78
static Quaternion< Scalar > Identity()
Definition Quaternion.h:109
EIGEN_CONSTEXPR NonConstCoeffReturnType z()
Definition Quaternion.h:73
EIGEN_CONSTEXPR NonConstCoeffReturnType y()
Definition Quaternion.h:71
bool operator!=(const QuaternionBase< OtherDerived > &other) const
Definition Quaternion.h:185
Quaternion< Scalar > conjugate() const
Definition Quaternion.h:754
bool isApprox(const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition Quaternion.h:194
EIGEN_CONSTEXPR NonConstCoeffReturnType w()
Definition Quaternion.h:75
EIGEN_CONSTEXPR NonConstCoeffReturnType x()
Definition Quaternion.h:69
void normalize()
Definition Quaternion.h:132
Derived & setFromTwoVectors(const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
Definition Quaternion.h:648
EIGEN_CONSTEXPR CoeffReturnType w() const
Definition Quaternion.h:66
Matrix3 toRotationMatrix() const
Definition Quaternion.h:602
EIGEN_CONSTEXPR CoeffReturnType y() const
Definition Quaternion.h:62
Scalar dot(const QuaternionBase< OtherDerived > &other) const
Definition Quaternion.h:143
Derived & operator=(const AngleAxisType &aa)
Definition Quaternion.h:573
Vector3 _transformVector(const Vector3 &v) const
Definition Quaternion.h:544
Scalar norm() const
Definition Quaternion.h:128
Quaternion< Scalar > inverse() const
Definition Quaternion.h:725
EIGEN_CONSTEXPR CoeffReturnType z() const
Definition Quaternion.h:64
const internal::traits< Derived >::Coefficients & coeffs() const
Definition Quaternion.h:84
bool operator==(const QuaternionBase< OtherDerived > &other) const
Definition Quaternion.h:176
AngleAxis< Scalar > AngleAxisType
Definition Quaternion.h:57
Derived & operator*=(const QuaternionBase< OtherDerived > &q)
Definition Quaternion.h:529
The quaternion class used to represent 3D orientations and rotations.
Definition Quaternion.h:285
Quaternion(const Quaternion< OtherScalar, OtherOptions > &other)
Definition Quaternion.h:343
Quaternion(const Scalar &w, const Scalar &x, const Scalar &y, const Scalar &z)
Definition Quaternion.h:308
Quaternion(const QuaternionBase< Derived > &other)
Definition Quaternion.h:325
Quaternion(const Scalar &w, const Eigen::MatrixBase< Derived > &vec)
Definition Quaternion.h:315
static Quaternion UnitRandom()
Definition Quaternion.h:689
Quaternion(const MatrixBase< Derived > &other)
Definition Quaternion.h:337
Quaternion(const AngleAxisType &aa)
Definition Quaternion.h:330
Quaternion()
Definition Quaternion.h:299
Quaternion & operator=(Quaternion &&other) EIGEN_NOEXCEPT_IF(std
Definition Quaternion.h:354
Quaternion(Quaternion &&other) EIGEN_NOEXCEPT_IF(std
Definition Quaternion.h:349
Quaternion(const Scalar *data)
Definition Quaternion.h:321
Common base class for compact rotation representations.
Definition RotationBase.h:32
friend RotationMatrixType operator*(const EigenBase< OtherDerived > &l, const Derived &r)
Definition RotationBase.h:83
const MatrixVType & matrixV() const
Definition SVDBase.h:189
Expression of a fixed-size or dynamic-size sub-vector.
Definition VectorBlock.h:58
@ Aligned
Definition Constants.h:242
@ DontAlign
Definition Constants.h:324
@ ComputeFullV
Definition Constants.h:393
const unsigned int LvalueBit
Definition Constants.h:148
Namespace containing all symbols from the Eigen library.
Definition Core:137
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523