xref: /aosp_15_r20/external/eigen/Eigen/src/Geometry/Quaternion.h (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2008-2010 Gael Guennebaud <[email protected]>
5*bf2c3715SXin Li // Copyright (C) 2009 Mathieu Gautier <[email protected]>
6*bf2c3715SXin Li //
7*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
8*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
9*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10*bf2c3715SXin Li 
11*bf2c3715SXin Li #ifndef EIGEN_QUATERNION_H
12*bf2c3715SXin Li #define EIGEN_QUATERNION_H
13*bf2c3715SXin Li namespace Eigen {
14*bf2c3715SXin Li 
15*bf2c3715SXin Li 
16*bf2c3715SXin Li /***************************************************************************
17*bf2c3715SXin Li * Definition of QuaternionBase<Derived>
18*bf2c3715SXin Li * The implementation is at the end of the file
19*bf2c3715SXin Li ***************************************************************************/
20*bf2c3715SXin Li 
21*bf2c3715SXin Li namespace internal {
22*bf2c3715SXin Li template<typename Other,
23*bf2c3715SXin Li          int OtherRows=Other::RowsAtCompileTime,
24*bf2c3715SXin Li          int OtherCols=Other::ColsAtCompileTime>
25*bf2c3715SXin Li struct quaternionbase_assign_impl;
26*bf2c3715SXin Li }
27*bf2c3715SXin Li 
28*bf2c3715SXin Li /** \geometry_module \ingroup Geometry_Module
29*bf2c3715SXin Li   * \class QuaternionBase
30*bf2c3715SXin Li   * \brief Base class for quaternion expressions
31*bf2c3715SXin Li   * \tparam Derived derived type (CRTP)
32*bf2c3715SXin Li   * \sa class Quaternion
33*bf2c3715SXin Li   */
34*bf2c3715SXin Li template<class Derived>
35*bf2c3715SXin Li class QuaternionBase : public RotationBase<Derived, 3>
36*bf2c3715SXin Li {
37*bf2c3715SXin Li  public:
38*bf2c3715SXin Li   typedef RotationBase<Derived, 3> Base;
39*bf2c3715SXin Li 
40*bf2c3715SXin Li   using Base::operator*;
41*bf2c3715SXin Li   using Base::derived;
42*bf2c3715SXin Li 
43*bf2c3715SXin Li   typedef typename internal::traits<Derived>::Scalar Scalar;
44*bf2c3715SXin Li   typedef typename NumTraits<Scalar>::Real RealScalar;
45*bf2c3715SXin Li   typedef typename internal::traits<Derived>::Coefficients Coefficients;
46*bf2c3715SXin Li   typedef typename Coefficients::CoeffReturnType CoeffReturnType;
47*bf2c3715SXin Li   typedef typename internal::conditional<bool(internal::traits<Derived>::Flags&LvalueBit),
48*bf2c3715SXin Li                                         Scalar&, CoeffReturnType>::type NonConstCoeffReturnType;
49*bf2c3715SXin Li 
50*bf2c3715SXin Li 
51*bf2c3715SXin Li   enum {
52*bf2c3715SXin Li     Flags = Eigen::internal::traits<Derived>::Flags
53*bf2c3715SXin Li   };
54*bf2c3715SXin Li 
55*bf2c3715SXin Li  // typedef typename Matrix<Scalar,4,1> Coefficients;
56*bf2c3715SXin Li   /** the type of a 3D vector */
57*bf2c3715SXin Li   typedef Matrix<Scalar,3,1> Vector3;
58*bf2c3715SXin Li   /** the equivalent rotation matrix type */
59*bf2c3715SXin Li   typedef Matrix<Scalar,3,3> Matrix3;
60*bf2c3715SXin Li   /** the equivalent angle-axis type */
61*bf2c3715SXin Li   typedef AngleAxis<Scalar> AngleAxisType;
62*bf2c3715SXin Li 
63*bf2c3715SXin Li 
64*bf2c3715SXin Li 
65*bf2c3715SXin Li   /** \returns the \c x coefficient */
x()66*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline CoeffReturnType x() const { return this->derived().coeffs().coeff(0); }
67*bf2c3715SXin Li   /** \returns the \c y coefficient */
y()68*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline CoeffReturnType y() const { return this->derived().coeffs().coeff(1); }
69*bf2c3715SXin Li   /** \returns the \c z coefficient */
z()70*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline CoeffReturnType z() const { return this->derived().coeffs().coeff(2); }
71*bf2c3715SXin Li   /** \returns the \c w coefficient */
w()72*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline CoeffReturnType w() const { return this->derived().coeffs().coeff(3); }
73*bf2c3715SXin Li 
74*bf2c3715SXin Li   /** \returns a reference to the \c x coefficient (if Derived is a non-const lvalue) */
x()75*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline NonConstCoeffReturnType x() { return this->derived().coeffs().x(); }
76*bf2c3715SXin Li   /** \returns a reference to the \c y coefficient (if Derived is a non-const lvalue) */
y()77*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline NonConstCoeffReturnType y() { return this->derived().coeffs().y(); }
78*bf2c3715SXin Li   /** \returns a reference to the \c z coefficient (if Derived is a non-const lvalue) */
z()79*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline NonConstCoeffReturnType z() { return this->derived().coeffs().z(); }
80*bf2c3715SXin Li   /** \returns a reference to the \c w coefficient (if Derived is a non-const lvalue) */
w()81*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline NonConstCoeffReturnType w() { return this->derived().coeffs().w(); }
82*bf2c3715SXin Li 
83*bf2c3715SXin Li   /** \returns a read-only vector expression of the imaginary part (x,y,z) */
vec()84*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline const VectorBlock<const Coefficients,3> vec() const { return coeffs().template head<3>(); }
85*bf2c3715SXin Li 
86*bf2c3715SXin Li   /** \returns a vector expression of the imaginary part (x,y,z) */
vec()87*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline VectorBlock<Coefficients,3> vec() { return coeffs().template head<3>(); }
88*bf2c3715SXin Li 
89*bf2c3715SXin Li   /** \returns a read-only vector expression of the coefficients (x,y,z,w) */
coeffs()90*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline const typename internal::traits<Derived>::Coefficients& coeffs() const { return derived().coeffs(); }
91*bf2c3715SXin Li 
92*bf2c3715SXin Li   /** \returns a vector expression of the coefficients (x,y,z,w) */
coeffs()93*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline typename internal::traits<Derived>::Coefficients& coeffs() { return derived().coeffs(); }
94*bf2c3715SXin Li 
95*bf2c3715SXin Li   EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE QuaternionBase<Derived>& operator=(const QuaternionBase<Derived>& other);
96*bf2c3715SXin Li   template<class OtherDerived> EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& operator=(const QuaternionBase<OtherDerived>& other);
97*bf2c3715SXin Li 
98*bf2c3715SXin Li // disabled this copy operator as it is giving very strange compilation errors when compiling
99*bf2c3715SXin Li // test_stdvector with GCC 4.4.2. This looks like a GCC bug though, so feel free to re-enable it if it's
100*bf2c3715SXin Li // useful; however notice that we already have the templated operator= above and e.g. in MatrixBase
101*bf2c3715SXin Li // we didn't have to add, in addition to templated operator=, such a non-templated copy operator.
102*bf2c3715SXin Li //  Derived& operator=(const QuaternionBase& other)
103*bf2c3715SXin Li //  { return operator=<Derived>(other); }
104*bf2c3715SXin Li 
105*bf2c3715SXin Li   EIGEN_DEVICE_FUNC Derived& operator=(const AngleAxisType& aa);
106*bf2c3715SXin Li   template<class OtherDerived> EIGEN_DEVICE_FUNC Derived& operator=(const MatrixBase<OtherDerived>& m);
107*bf2c3715SXin Li 
108*bf2c3715SXin Li   /** \returns a quaternion representing an identity rotation
109*bf2c3715SXin Li     * \sa MatrixBase::Identity()
110*bf2c3715SXin Li     */
Identity()111*bf2c3715SXin Li   EIGEN_DEVICE_FUNC static inline Quaternion<Scalar> Identity() { return Quaternion<Scalar>(Scalar(1), Scalar(0), Scalar(0), Scalar(0)); }
112*bf2c3715SXin Li 
113*bf2c3715SXin Li   /** \sa QuaternionBase::Identity(), MatrixBase::setIdentity()
114*bf2c3715SXin Li     */
setIdentity()115*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline QuaternionBase& setIdentity() { coeffs() << Scalar(0), Scalar(0), Scalar(0), Scalar(1); return *this; }
116*bf2c3715SXin Li 
117*bf2c3715SXin Li   /** \returns the squared norm of the quaternion's coefficients
118*bf2c3715SXin Li     * \sa QuaternionBase::norm(), MatrixBase::squaredNorm()
119*bf2c3715SXin Li     */
squaredNorm()120*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline Scalar squaredNorm() const { return coeffs().squaredNorm(); }
121*bf2c3715SXin Li 
122*bf2c3715SXin Li   /** \returns the norm of the quaternion's coefficients
123*bf2c3715SXin Li     * \sa QuaternionBase::squaredNorm(), MatrixBase::norm()
124*bf2c3715SXin Li     */
norm()125*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline Scalar norm() const { return coeffs().norm(); }
126*bf2c3715SXin Li 
127*bf2c3715SXin Li   /** Normalizes the quaternion \c *this
128*bf2c3715SXin Li     * \sa normalized(), MatrixBase::normalize() */
normalize()129*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline void normalize() { coeffs().normalize(); }
130*bf2c3715SXin Li   /** \returns a normalized copy of \c *this
131*bf2c3715SXin Li     * \sa normalize(), MatrixBase::normalized() */
normalized()132*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline Quaternion<Scalar> normalized() const { return Quaternion<Scalar>(coeffs().normalized()); }
133*bf2c3715SXin Li 
134*bf2c3715SXin Li     /** \returns the dot product of \c *this and \a other
135*bf2c3715SXin Li     * Geometrically speaking, the dot product of two unit quaternions
136*bf2c3715SXin Li     * corresponds to the cosine of half the angle between the two rotations.
137*bf2c3715SXin Li     * \sa angularDistance()
138*bf2c3715SXin Li     */
dot(const QuaternionBase<OtherDerived> & other)139*bf2c3715SXin Li   template<class OtherDerived> EIGEN_DEVICE_FUNC inline Scalar dot(const QuaternionBase<OtherDerived>& other) const { return coeffs().dot(other.coeffs()); }
140*bf2c3715SXin Li 
141*bf2c3715SXin Li   template<class OtherDerived> EIGEN_DEVICE_FUNC Scalar angularDistance(const QuaternionBase<OtherDerived>& other) const;
142*bf2c3715SXin Li 
143*bf2c3715SXin Li   /** \returns an equivalent 3x3 rotation matrix */
144*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline Matrix3 toRotationMatrix() const;
145*bf2c3715SXin Li 
146*bf2c3715SXin Li   /** \returns the quaternion which transform \a a into \a b through a rotation */
147*bf2c3715SXin Li   template<typename Derived1, typename Derived2>
148*bf2c3715SXin Li   EIGEN_DEVICE_FUNC Derived& setFromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
149*bf2c3715SXin Li 
150*bf2c3715SXin Li   template<class OtherDerived> EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion<Scalar> operator* (const QuaternionBase<OtherDerived>& q) const;
151*bf2c3715SXin Li   template<class OtherDerived> EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& operator*= (const QuaternionBase<OtherDerived>& q);
152*bf2c3715SXin Li 
153*bf2c3715SXin Li   /** \returns the quaternion describing the inverse rotation */
154*bf2c3715SXin Li   EIGEN_DEVICE_FUNC Quaternion<Scalar> inverse() const;
155*bf2c3715SXin Li 
156*bf2c3715SXin Li   /** \returns the conjugated quaternion */
157*bf2c3715SXin Li   EIGEN_DEVICE_FUNC Quaternion<Scalar> conjugate() const;
158*bf2c3715SXin Li 
159*bf2c3715SXin Li   template<class OtherDerived> EIGEN_DEVICE_FUNC Quaternion<Scalar> slerp(const Scalar& t, const QuaternionBase<OtherDerived>& other) const;
160*bf2c3715SXin Li 
161*bf2c3715SXin Li   /** \returns true if each coefficients of \c *this and \a other are all exactly equal.
162*bf2c3715SXin Li     * \warning When using floating point scalar values you probably should rather use a
163*bf2c3715SXin Li     *          fuzzy comparison such as isApprox()
164*bf2c3715SXin Li     * \sa isApprox(), operator!= */
165*bf2c3715SXin Li   template<class OtherDerived>
166*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline bool operator==(const QuaternionBase<OtherDerived>& other) const
167*bf2c3715SXin Li   { return coeffs() == other.coeffs(); }
168*bf2c3715SXin Li 
169*bf2c3715SXin Li   /** \returns true if at least one pair of coefficients of \c *this and \a other are not exactly equal to each other.
170*bf2c3715SXin Li     * \warning When using floating point scalar values you probably should rather use a
171*bf2c3715SXin Li     *          fuzzy comparison such as isApprox()
172*bf2c3715SXin Li     * \sa isApprox(), operator== */
173*bf2c3715SXin Li   template<class OtherDerived>
174*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline bool operator!=(const QuaternionBase<OtherDerived>& other) const
175*bf2c3715SXin Li   { return coeffs() != other.coeffs(); }
176*bf2c3715SXin Li 
177*bf2c3715SXin Li   /** \returns \c true if \c *this is approximately equal to \a other, within the precision
178*bf2c3715SXin Li     * determined by \a prec.
179*bf2c3715SXin Li     *
180*bf2c3715SXin Li     * \sa MatrixBase::isApprox() */
181*bf2c3715SXin Li   template<class OtherDerived>
182*bf2c3715SXin Li   EIGEN_DEVICE_FUNC bool isApprox(const QuaternionBase<OtherDerived>& other, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const
183*bf2c3715SXin Li   { return coeffs().isApprox(other.coeffs(), prec); }
184*bf2c3715SXin Li 
185*bf2c3715SXin Li   /** return the result vector of \a v through the rotation*/
186*bf2c3715SXin Li   EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Vector3 _transformVector(const Vector3& v) const;
187*bf2c3715SXin Li 
188*bf2c3715SXin Li   #ifdef EIGEN_PARSED_BY_DOXYGEN
189*bf2c3715SXin Li   /** \returns \c *this with scalar type casted to \a NewScalarType
190*bf2c3715SXin Li     *
191*bf2c3715SXin Li     * Note that if \a NewScalarType is equal to the current scalar type of \c *this
192*bf2c3715SXin Li     * then this function smartly returns a const reference to \c *this.
193*bf2c3715SXin Li     */
194*bf2c3715SXin Li   template<typename NewScalarType>
195*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline typename internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type cast() const;
196*bf2c3715SXin Li 
197*bf2c3715SXin Li   #else
198*bf2c3715SXin Li 
199*bf2c3715SXin Li   template<typename NewScalarType>
200*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline
cast()201*bf2c3715SXin Li   typename internal::enable_if<internal::is_same<Scalar,NewScalarType>::value,const Derived&>::type cast() const
202*bf2c3715SXin Li   {
203*bf2c3715SXin Li     return derived();
204*bf2c3715SXin Li   }
205*bf2c3715SXin Li 
206*bf2c3715SXin Li   template<typename NewScalarType>
207*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline
cast()208*bf2c3715SXin Li   typename internal::enable_if<!internal::is_same<Scalar,NewScalarType>::value,Quaternion<NewScalarType> >::type cast() const
209*bf2c3715SXin Li   {
210*bf2c3715SXin Li     return Quaternion<NewScalarType>(coeffs().template cast<NewScalarType>());
211*bf2c3715SXin Li   }
212*bf2c3715SXin Li   #endif
213*bf2c3715SXin Li 
214*bf2c3715SXin Li #ifndef EIGEN_NO_IO
215*bf2c3715SXin Li   friend std::ostream& operator<<(std::ostream& s, const QuaternionBase<Derived>& q) {
216*bf2c3715SXin Li     s << q.x() << "i + " << q.y() << "j + " << q.z() << "k" << " + " << q.w();
217*bf2c3715SXin Li     return s;
218*bf2c3715SXin Li   }
219*bf2c3715SXin Li #endif
220*bf2c3715SXin Li 
221*bf2c3715SXin Li #ifdef EIGEN_QUATERNIONBASE_PLUGIN
222*bf2c3715SXin Li # include EIGEN_QUATERNIONBASE_PLUGIN
223*bf2c3715SXin Li #endif
224*bf2c3715SXin Li protected:
225*bf2c3715SXin Li   EIGEN_DEFAULT_COPY_CONSTRUCTOR(QuaternionBase)
226*bf2c3715SXin Li   EIGEN_DEFAULT_EMPTY_CONSTRUCTOR_AND_DESTRUCTOR(QuaternionBase)
227*bf2c3715SXin Li };
228*bf2c3715SXin Li 
229*bf2c3715SXin Li /***************************************************************************
230*bf2c3715SXin Li * Definition/implementation of Quaternion<Scalar>
231*bf2c3715SXin Li ***************************************************************************/
232*bf2c3715SXin Li 
233*bf2c3715SXin Li /** \geometry_module \ingroup Geometry_Module
234*bf2c3715SXin Li   *
235*bf2c3715SXin Li   * \class Quaternion
236*bf2c3715SXin Li   *
237*bf2c3715SXin Li   * \brief The quaternion class used to represent 3D orientations and rotations
238*bf2c3715SXin Li   *
239*bf2c3715SXin Li   * \tparam _Scalar the scalar type, i.e., the type of the coefficients
240*bf2c3715SXin Li   * \tparam _Options controls the memory alignment of the coefficients. Can be \# AutoAlign or \# DontAlign. Default is AutoAlign.
241*bf2c3715SXin Li   *
242*bf2c3715SXin Li   * This class represents a quaternion \f$ w+xi+yj+zk \f$ that is a convenient representation of
243*bf2c3715SXin Li   * orientations and rotations of objects in three dimensions. Compared to other representations
244*bf2c3715SXin Li   * like Euler angles or 3x3 matrices, quaternions offer the following advantages:
245*bf2c3715SXin Li   * \li \b compact storage (4 scalars)
246*bf2c3715SXin Li   * \li \b efficient to compose (28 flops),
247*bf2c3715SXin Li   * \li \b stable spherical interpolation
248*bf2c3715SXin Li   *
249*bf2c3715SXin Li   * The following two typedefs are provided for convenience:
250*bf2c3715SXin Li   * \li \c Quaternionf for \c float
251*bf2c3715SXin Li   * \li \c Quaterniond for \c double
252*bf2c3715SXin Li   *
253*bf2c3715SXin Li   * \warning Operations interpreting the quaternion as rotation have undefined behavior if the quaternion is not normalized.
254*bf2c3715SXin Li   *
255*bf2c3715SXin Li   * \sa  class AngleAxis, class Transform
256*bf2c3715SXin Li   */
257*bf2c3715SXin Li 
258*bf2c3715SXin Li namespace internal {
259*bf2c3715SXin Li template<typename _Scalar,int _Options>
260*bf2c3715SXin Li struct traits<Quaternion<_Scalar,_Options> >
261*bf2c3715SXin Li {
262*bf2c3715SXin Li   typedef Quaternion<_Scalar,_Options> PlainObject;
263*bf2c3715SXin Li   typedef _Scalar Scalar;
264*bf2c3715SXin Li   typedef Matrix<_Scalar,4,1,_Options> Coefficients;
265*bf2c3715SXin Li   enum{
266*bf2c3715SXin Li     Alignment = internal::traits<Coefficients>::Alignment,
267*bf2c3715SXin Li     Flags = LvalueBit
268*bf2c3715SXin Li   };
269*bf2c3715SXin Li };
270*bf2c3715SXin Li }
271*bf2c3715SXin Li 
272*bf2c3715SXin Li template<typename _Scalar, int _Options>
273*bf2c3715SXin Li class Quaternion : public QuaternionBase<Quaternion<_Scalar,_Options> >
274*bf2c3715SXin Li {
275*bf2c3715SXin Li public:
276*bf2c3715SXin Li   typedef QuaternionBase<Quaternion<_Scalar,_Options> > Base;
277*bf2c3715SXin Li   enum { NeedsAlignment = internal::traits<Quaternion>::Alignment>0 };
278*bf2c3715SXin Li 
279*bf2c3715SXin Li   typedef _Scalar Scalar;
280*bf2c3715SXin Li 
281*bf2c3715SXin Li   EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Quaternion)
282*bf2c3715SXin Li   using Base::operator*=;
283*bf2c3715SXin Li 
284*bf2c3715SXin Li   typedef typename internal::traits<Quaternion>::Coefficients Coefficients;
285*bf2c3715SXin Li   typedef typename Base::AngleAxisType AngleAxisType;
286*bf2c3715SXin Li 
287*bf2c3715SXin Li   /** Default constructor leaving the quaternion uninitialized. */
288*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline Quaternion() {}
289*bf2c3715SXin Li 
290*bf2c3715SXin Li   /** Constructs and initializes the quaternion \f$ w+xi+yj+zk \f$ from
291*bf2c3715SXin Li     * its four coefficients \a w, \a x, \a y and \a z.
292*bf2c3715SXin Li     *
293*bf2c3715SXin Li     * \warning Note the order of the arguments: the real \a w coefficient first,
294*bf2c3715SXin Li     * while internally the coefficients are stored in the following order:
295*bf2c3715SXin Li     * [\c x, \c y, \c z, \c w]
296*bf2c3715SXin Li     */
297*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline Quaternion(const Scalar& w, const Scalar& x, const Scalar& y, const Scalar& z) : m_coeffs(x, y, z, w){}
298*bf2c3715SXin Li 
299*bf2c3715SXin Li   /** Constructs and initialize a quaternion from the array data */
300*bf2c3715SXin Li   EIGEN_DEVICE_FUNC explicit inline Quaternion(const Scalar* data) : m_coeffs(data) {}
301*bf2c3715SXin Li 
302*bf2c3715SXin Li   /** Copy constructor */
303*bf2c3715SXin Li   template<class Derived> EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion(const QuaternionBase<Derived>& other) { this->Base::operator=(other); }
304*bf2c3715SXin Li 
305*bf2c3715SXin Li   /** Constructs and initializes a quaternion from the angle-axis \a aa */
306*bf2c3715SXin Li   EIGEN_DEVICE_FUNC explicit inline Quaternion(const AngleAxisType& aa) { *this = aa; }
307*bf2c3715SXin Li 
308*bf2c3715SXin Li   /** Constructs and initializes a quaternion from either:
309*bf2c3715SXin Li     *  - a rotation matrix expression,
310*bf2c3715SXin Li     *  - a 4D vector expression representing quaternion coefficients.
311*bf2c3715SXin Li     */
312*bf2c3715SXin Li   template<typename Derived>
313*bf2c3715SXin Li   EIGEN_DEVICE_FUNC explicit inline Quaternion(const MatrixBase<Derived>& other) { *this = other; }
314*bf2c3715SXin Li 
315*bf2c3715SXin Li   /** Explicit copy constructor with scalar conversion */
316*bf2c3715SXin Li   template<typename OtherScalar, int OtherOptions>
317*bf2c3715SXin Li   EIGEN_DEVICE_FUNC explicit inline Quaternion(const Quaternion<OtherScalar, OtherOptions>& other)
318*bf2c3715SXin Li   { m_coeffs = other.coeffs().template cast<Scalar>(); }
319*bf2c3715SXin Li 
320*bf2c3715SXin Li #if EIGEN_HAS_RVALUE_REFERENCES
321*bf2c3715SXin Li   // We define a copy constructor, which means we don't get an implicit move constructor or assignment operator.
322*bf2c3715SXin Li   /** Default move constructor */
323*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline Quaternion(Quaternion&& other) EIGEN_NOEXCEPT_IF(std::is_nothrow_move_constructible<Scalar>::value)
324*bf2c3715SXin Li     : m_coeffs(std::move(other.coeffs()))
325*bf2c3715SXin Li   {}
326*bf2c3715SXin Li 
327*bf2c3715SXin Li   /** Default move assignment operator */
328*bf2c3715SXin Li   EIGEN_DEVICE_FUNC Quaternion& operator=(Quaternion&& other) EIGEN_NOEXCEPT_IF(std::is_nothrow_move_assignable<Scalar>::value)
329*bf2c3715SXin Li   {
330*bf2c3715SXin Li     m_coeffs = std::move(other.coeffs());
331*bf2c3715SXin Li     return *this;
332*bf2c3715SXin Li   }
333*bf2c3715SXin Li #endif
334*bf2c3715SXin Li 
335*bf2c3715SXin Li   EIGEN_DEVICE_FUNC static Quaternion UnitRandom();
336*bf2c3715SXin Li 
337*bf2c3715SXin Li   template<typename Derived1, typename Derived2>
338*bf2c3715SXin Li   EIGEN_DEVICE_FUNC static Quaternion FromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b);
339*bf2c3715SXin Li 
340*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline Coefficients& coeffs() { return m_coeffs;}
341*bf2c3715SXin Li   EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs;}
342*bf2c3715SXin Li 
343*bf2c3715SXin Li   EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(bool(NeedsAlignment))
344*bf2c3715SXin Li 
345*bf2c3715SXin Li #ifdef EIGEN_QUATERNION_PLUGIN
346*bf2c3715SXin Li # include EIGEN_QUATERNION_PLUGIN
347*bf2c3715SXin Li #endif
348*bf2c3715SXin Li 
349*bf2c3715SXin Li protected:
350*bf2c3715SXin Li   Coefficients m_coeffs;
351*bf2c3715SXin Li 
352*bf2c3715SXin Li #ifndef EIGEN_PARSED_BY_DOXYGEN
353*bf2c3715SXin Li     static EIGEN_STRONG_INLINE void _check_template_params()
354*bf2c3715SXin Li     {
355*bf2c3715SXin Li       EIGEN_STATIC_ASSERT( (_Options & DontAlign) == _Options,
356*bf2c3715SXin Li         INVALID_MATRIX_TEMPLATE_PARAMETERS)
357*bf2c3715SXin Li     }
358*bf2c3715SXin Li #endif
359*bf2c3715SXin Li };
360*bf2c3715SXin Li 
361*bf2c3715SXin Li /** \ingroup Geometry_Module
362*bf2c3715SXin Li   * single precision quaternion type */
363*bf2c3715SXin Li typedef Quaternion<float> Quaternionf;
364*bf2c3715SXin Li /** \ingroup Geometry_Module
365*bf2c3715SXin Li   * double precision quaternion type */
366*bf2c3715SXin Li typedef Quaternion<double> Quaterniond;
367*bf2c3715SXin Li 
368*bf2c3715SXin Li /***************************************************************************
369*bf2c3715SXin Li * Specialization of Map<Quaternion<Scalar>>
370*bf2c3715SXin Li ***************************************************************************/
371*bf2c3715SXin Li 
372*bf2c3715SXin Li namespace internal {
373*bf2c3715SXin Li   template<typename _Scalar, int _Options>
374*bf2c3715SXin Li   struct traits<Map<Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
375*bf2c3715SXin Li   {
376*bf2c3715SXin Li     typedef Map<Matrix<_Scalar,4,1>, _Options> Coefficients;
377*bf2c3715SXin Li   };
378*bf2c3715SXin Li }
379*bf2c3715SXin Li 
380*bf2c3715SXin Li namespace internal {
381*bf2c3715SXin Li   template<typename _Scalar, int _Options>
382*bf2c3715SXin Li   struct traits<Map<const Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> >
383*bf2c3715SXin Li   {
384*bf2c3715SXin Li     typedef Map<const Matrix<_Scalar,4,1>, _Options> Coefficients;
385*bf2c3715SXin Li     typedef traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> > TraitsBase;
386*bf2c3715SXin Li     enum {
387*bf2c3715SXin Li       Flags = TraitsBase::Flags & ~LvalueBit
388*bf2c3715SXin Li     };
389*bf2c3715SXin Li   };
390*bf2c3715SXin Li }
391*bf2c3715SXin Li 
392*bf2c3715SXin Li /** \ingroup Geometry_Module
393*bf2c3715SXin Li   * \brief Quaternion expression mapping a constant memory buffer
394*bf2c3715SXin Li   *
395*bf2c3715SXin Li   * \tparam _Scalar the type of the Quaternion coefficients
396*bf2c3715SXin Li   * \tparam _Options see class Map
397*bf2c3715SXin Li   *
398*bf2c3715SXin Li   * This is a specialization of class Map for Quaternion. This class allows to view
399*bf2c3715SXin Li   * a 4 scalar memory buffer as an Eigen's Quaternion object.
400*bf2c3715SXin Li   *
401*bf2c3715SXin Li   * \sa class Map, class Quaternion, class QuaternionBase
402*bf2c3715SXin Li   */
403*bf2c3715SXin Li template<typename _Scalar, int _Options>
404*bf2c3715SXin Li class Map<const Quaternion<_Scalar>, _Options >
405*bf2c3715SXin Li   : public QuaternionBase<Map<const Quaternion<_Scalar>, _Options> >
406*bf2c3715SXin Li {
407*bf2c3715SXin Li   public:
408*bf2c3715SXin Li     typedef QuaternionBase<Map<const Quaternion<_Scalar>, _Options> > Base;
409*bf2c3715SXin Li 
410*bf2c3715SXin Li     typedef _Scalar Scalar;
411*bf2c3715SXin Li     typedef typename internal::traits<Map>::Coefficients Coefficients;
412*bf2c3715SXin Li     EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
413*bf2c3715SXin Li     using Base::operator*=;
414*bf2c3715SXin Li 
415*bf2c3715SXin Li     /** Constructs a Mapped Quaternion object from the pointer \a coeffs
416*bf2c3715SXin Li       *
417*bf2c3715SXin Li       * The pointer \a coeffs must reference the four coefficients of Quaternion in the following order:
418*bf2c3715SXin Li       * \code *coeffs == {x, y, z, w} \endcode
419*bf2c3715SXin Li       *
420*bf2c3715SXin Li       * If the template parameter _Options is set to #Aligned, then the pointer coeffs must be aligned. */
421*bf2c3715SXin Li     EIGEN_DEVICE_FUNC explicit EIGEN_STRONG_INLINE Map(const Scalar* coeffs) : m_coeffs(coeffs) {}
422*bf2c3715SXin Li 
423*bf2c3715SXin Li     EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs;}
424*bf2c3715SXin Li 
425*bf2c3715SXin Li   protected:
426*bf2c3715SXin Li     const Coefficients m_coeffs;
427*bf2c3715SXin Li };
428*bf2c3715SXin Li 
429*bf2c3715SXin Li /** \ingroup Geometry_Module
430*bf2c3715SXin Li   * \brief Expression of a quaternion from a memory buffer
431*bf2c3715SXin Li   *
432*bf2c3715SXin Li   * \tparam _Scalar the type of the Quaternion coefficients
433*bf2c3715SXin Li   * \tparam _Options see class Map
434*bf2c3715SXin Li   *
435*bf2c3715SXin Li   * This is a specialization of class Map for Quaternion. This class allows to view
436*bf2c3715SXin Li   * a 4 scalar memory buffer as an Eigen's  Quaternion object.
437*bf2c3715SXin Li   *
438*bf2c3715SXin Li   * \sa class Map, class Quaternion, class QuaternionBase
439*bf2c3715SXin Li   */
440*bf2c3715SXin Li template<typename _Scalar, int _Options>
441*bf2c3715SXin Li class Map<Quaternion<_Scalar>, _Options >
442*bf2c3715SXin Li   : public QuaternionBase<Map<Quaternion<_Scalar>, _Options> >
443*bf2c3715SXin Li {
444*bf2c3715SXin Li   public:
445*bf2c3715SXin Li     typedef QuaternionBase<Map<Quaternion<_Scalar>, _Options> > Base;
446*bf2c3715SXin Li 
447*bf2c3715SXin Li     typedef _Scalar Scalar;
448*bf2c3715SXin Li     typedef typename internal::traits<Map>::Coefficients Coefficients;
449*bf2c3715SXin Li     EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
450*bf2c3715SXin Li     using Base::operator*=;
451*bf2c3715SXin Li 
452*bf2c3715SXin Li     /** Constructs a Mapped Quaternion object from the pointer \a coeffs
453*bf2c3715SXin Li       *
454*bf2c3715SXin Li       * The pointer \a coeffs must reference the four coefficients of Quaternion in the following order:
455*bf2c3715SXin Li       * \code *coeffs == {x, y, z, w} \endcode
456*bf2c3715SXin Li       *
457*bf2c3715SXin Li       * If the template parameter _Options is set to #Aligned, then the pointer coeffs must be aligned. */
458*bf2c3715SXin Li     EIGEN_DEVICE_FUNC explicit EIGEN_STRONG_INLINE Map(Scalar* coeffs) : m_coeffs(coeffs) {}
459*bf2c3715SXin Li 
460*bf2c3715SXin Li     EIGEN_DEVICE_FUNC inline Coefficients& coeffs() { return m_coeffs; }
461*bf2c3715SXin Li     EIGEN_DEVICE_FUNC inline const Coefficients& coeffs() const { return m_coeffs; }
462*bf2c3715SXin Li 
463*bf2c3715SXin Li   protected:
464*bf2c3715SXin Li     Coefficients m_coeffs;
465*bf2c3715SXin Li };
466*bf2c3715SXin Li 
467*bf2c3715SXin Li /** \ingroup Geometry_Module
468*bf2c3715SXin Li   * Map an unaligned array of single precision scalars as a quaternion */
469*bf2c3715SXin Li typedef Map<Quaternion<float>, 0>         QuaternionMapf;
470*bf2c3715SXin Li /** \ingroup Geometry_Module
471*bf2c3715SXin Li   * Map an unaligned array of double precision scalars as a quaternion */
472*bf2c3715SXin Li typedef Map<Quaternion<double>, 0>        QuaternionMapd;
473*bf2c3715SXin Li /** \ingroup Geometry_Module
474*bf2c3715SXin Li   * Map a 16-byte aligned array of single precision scalars as a quaternion */
475*bf2c3715SXin Li typedef Map<Quaternion<float>, Aligned>   QuaternionMapAlignedf;
476*bf2c3715SXin Li /** \ingroup Geometry_Module
477*bf2c3715SXin Li   * Map a 16-byte aligned array of double precision scalars as a quaternion */
478*bf2c3715SXin Li typedef Map<Quaternion<double>, Aligned>  QuaternionMapAlignedd;
479*bf2c3715SXin Li 
480*bf2c3715SXin Li /***************************************************************************
481*bf2c3715SXin Li * Implementation of QuaternionBase methods
482*bf2c3715SXin Li ***************************************************************************/
483*bf2c3715SXin Li 
484*bf2c3715SXin Li // Generic Quaternion * Quaternion product
485*bf2c3715SXin Li // This product can be specialized for a given architecture via the Arch template argument.
486*bf2c3715SXin Li namespace internal {
487*bf2c3715SXin Li template<int Arch, class Derived1, class Derived2, typename Scalar> struct quat_product
488*bf2c3715SXin Li {
489*bf2c3715SXin Li   EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived1>& a, const QuaternionBase<Derived2>& b){
490*bf2c3715SXin Li     return Quaternion<Scalar>
491*bf2c3715SXin Li     (
492*bf2c3715SXin Li       a.w() * b.w() - a.x() * b.x() - a.y() * b.y() - a.z() * b.z(),
493*bf2c3715SXin Li       a.w() * b.x() + a.x() * b.w() + a.y() * b.z() - a.z() * b.y(),
494*bf2c3715SXin Li       a.w() * b.y() + a.y() * b.w() + a.z() * b.x() - a.x() * b.z(),
495*bf2c3715SXin Li       a.w() * b.z() + a.z() * b.w() + a.x() * b.y() - a.y() * b.x()
496*bf2c3715SXin Li     );
497*bf2c3715SXin Li   }
498*bf2c3715SXin Li };
499*bf2c3715SXin Li }
500*bf2c3715SXin Li 
501*bf2c3715SXin Li /** \returns the concatenation of two rotations as a quaternion-quaternion product */
502*bf2c3715SXin Li template <class Derived>
503*bf2c3715SXin Li template <class OtherDerived>
504*bf2c3715SXin Li EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Quaternion<typename internal::traits<Derived>::Scalar>
505*bf2c3715SXin Li QuaternionBase<Derived>::operator* (const QuaternionBase<OtherDerived>& other) const
506*bf2c3715SXin Li {
507*bf2c3715SXin Li   EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename OtherDerived::Scalar>::value),
508*bf2c3715SXin Li    YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
509*bf2c3715SXin Li   return internal::quat_product<Architecture::Target, Derived, OtherDerived,
510*bf2c3715SXin Li                          typename internal::traits<Derived>::Scalar>::run(*this, other);
511*bf2c3715SXin Li }
512*bf2c3715SXin Li 
513*bf2c3715SXin Li /** \sa operator*(Quaternion) */
514*bf2c3715SXin Li template <class Derived>
515*bf2c3715SXin Li template <class OtherDerived>
516*bf2c3715SXin Li EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator*= (const QuaternionBase<OtherDerived>& other)
517*bf2c3715SXin Li {
518*bf2c3715SXin Li   derived() = derived() * other.derived();
519*bf2c3715SXin Li   return derived();
520*bf2c3715SXin Li }
521*bf2c3715SXin Li 
522*bf2c3715SXin Li /** Rotation of a vector by a quaternion.
523*bf2c3715SXin Li   * \remarks If the quaternion is used to rotate several points (>1)
524*bf2c3715SXin Li   * then it is much more efficient to first convert it to a 3x3 Matrix.
525*bf2c3715SXin Li   * Comparison of the operation cost for n transformations:
526*bf2c3715SXin Li   *   - Quaternion2:    30n
527*bf2c3715SXin Li   *   - Via a Matrix3: 24 + 15n
528*bf2c3715SXin Li   */
529*bf2c3715SXin Li template <class Derived>
530*bf2c3715SXin Li EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE typename QuaternionBase<Derived>::Vector3
531*bf2c3715SXin Li QuaternionBase<Derived>::_transformVector(const Vector3& v) const
532*bf2c3715SXin Li {
533*bf2c3715SXin Li     // Note that this algorithm comes from the optimization by hand
534*bf2c3715SXin Li     // of the conversion to a Matrix followed by a Matrix/Vector product.
535*bf2c3715SXin Li     // It appears to be much faster than the common algorithm found
536*bf2c3715SXin Li     // in the literature (30 versus 39 flops). It also requires two
537*bf2c3715SXin Li     // Vector3 as temporaries.
538*bf2c3715SXin Li     Vector3 uv = this->vec().cross(v);
539*bf2c3715SXin Li     uv += uv;
540*bf2c3715SXin Li     return v + this->w() * uv + this->vec().cross(uv);
541*bf2c3715SXin Li }
542*bf2c3715SXin Li 
543*bf2c3715SXin Li template<class Derived>
544*bf2c3715SXin Li EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE QuaternionBase<Derived>& QuaternionBase<Derived>::operator=(const QuaternionBase<Derived>& other)
545*bf2c3715SXin Li {
546*bf2c3715SXin Li   coeffs() = other.coeffs();
547*bf2c3715SXin Li   return derived();
548*bf2c3715SXin Li }
549*bf2c3715SXin Li 
550*bf2c3715SXin Li template<class Derived>
551*bf2c3715SXin Li template<class OtherDerived>
552*bf2c3715SXin Li EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const QuaternionBase<OtherDerived>& other)
553*bf2c3715SXin Li {
554*bf2c3715SXin Li   coeffs() = other.coeffs();
555*bf2c3715SXin Li   return derived();
556*bf2c3715SXin Li }
557*bf2c3715SXin Li 
558*bf2c3715SXin Li /** Set \c *this from an angle-axis \a aa and returns a reference to \c *this
559*bf2c3715SXin Li   */
560*bf2c3715SXin Li template<class Derived>
561*bf2c3715SXin Li EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(const AngleAxisType& aa)
562*bf2c3715SXin Li {
563*bf2c3715SXin Li   EIGEN_USING_STD(cos)
564*bf2c3715SXin Li   EIGEN_USING_STD(sin)
565*bf2c3715SXin Li   Scalar ha = Scalar(0.5)*aa.angle(); // Scalar(0.5) to suppress precision loss warnings
566*bf2c3715SXin Li   this->w() = cos(ha);
567*bf2c3715SXin Li   this->vec() = sin(ha) * aa.axis();
568*bf2c3715SXin Li   return derived();
569*bf2c3715SXin Li }
570*bf2c3715SXin Li 
571*bf2c3715SXin Li /** Set \c *this from the expression \a xpr:
572*bf2c3715SXin Li   *   - if \a xpr is a 4x1 vector, then \a xpr is assumed to be a quaternion
573*bf2c3715SXin Li   *   - if \a xpr is a 3x3 matrix, then \a xpr is assumed to be rotation matrix
574*bf2c3715SXin Li   *     and \a xpr is converted to a quaternion
575*bf2c3715SXin Li   */
576*bf2c3715SXin Li 
577*bf2c3715SXin Li template<class Derived>
578*bf2c3715SXin Li template<class MatrixDerived>
579*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline Derived& QuaternionBase<Derived>::operator=(const MatrixBase<MatrixDerived>& xpr)
580*bf2c3715SXin Li {
581*bf2c3715SXin Li   EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename MatrixDerived::Scalar>::value),
582*bf2c3715SXin Li    YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
583*bf2c3715SXin Li   internal::quaternionbase_assign_impl<MatrixDerived>::run(*this, xpr.derived());
584*bf2c3715SXin Li   return derived();
585*bf2c3715SXin Li }
586*bf2c3715SXin Li 
587*bf2c3715SXin Li /** Convert the quaternion to a 3x3 rotation matrix. The quaternion is required to
588*bf2c3715SXin Li   * be normalized, otherwise the result is undefined.
589*bf2c3715SXin Li   */
590*bf2c3715SXin Li template<class Derived>
591*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline typename QuaternionBase<Derived>::Matrix3
592*bf2c3715SXin Li QuaternionBase<Derived>::toRotationMatrix(void) const
593*bf2c3715SXin Li {
594*bf2c3715SXin Li   // NOTE if inlined, then gcc 4.2 and 4.4 get rid of the temporary (not gcc 4.3 !!)
595*bf2c3715SXin Li   // if not inlined then the cost of the return by value is huge ~ +35%,
596*bf2c3715SXin Li   // however, not inlining this function is an order of magnitude slower, so
597*bf2c3715SXin Li   // it has to be inlined, and so the return by value is not an issue
598*bf2c3715SXin Li   Matrix3 res;
599*bf2c3715SXin Li 
600*bf2c3715SXin Li   const Scalar tx  = Scalar(2)*this->x();
601*bf2c3715SXin Li   const Scalar ty  = Scalar(2)*this->y();
602*bf2c3715SXin Li   const Scalar tz  = Scalar(2)*this->z();
603*bf2c3715SXin Li   const Scalar twx = tx*this->w();
604*bf2c3715SXin Li   const Scalar twy = ty*this->w();
605*bf2c3715SXin Li   const Scalar twz = tz*this->w();
606*bf2c3715SXin Li   const Scalar txx = tx*this->x();
607*bf2c3715SXin Li   const Scalar txy = ty*this->x();
608*bf2c3715SXin Li   const Scalar txz = tz*this->x();
609*bf2c3715SXin Li   const Scalar tyy = ty*this->y();
610*bf2c3715SXin Li   const Scalar tyz = tz*this->y();
611*bf2c3715SXin Li   const Scalar tzz = tz*this->z();
612*bf2c3715SXin Li 
613*bf2c3715SXin Li   res.coeffRef(0,0) = Scalar(1)-(tyy+tzz);
614*bf2c3715SXin Li   res.coeffRef(0,1) = txy-twz;
615*bf2c3715SXin Li   res.coeffRef(0,2) = txz+twy;
616*bf2c3715SXin Li   res.coeffRef(1,0) = txy+twz;
617*bf2c3715SXin Li   res.coeffRef(1,1) = Scalar(1)-(txx+tzz);
618*bf2c3715SXin Li   res.coeffRef(1,2) = tyz-twx;
619*bf2c3715SXin Li   res.coeffRef(2,0) = txz-twy;
620*bf2c3715SXin Li   res.coeffRef(2,1) = tyz+twx;
621*bf2c3715SXin Li   res.coeffRef(2,2) = Scalar(1)-(txx+tyy);
622*bf2c3715SXin Li 
623*bf2c3715SXin Li   return res;
624*bf2c3715SXin Li }
625*bf2c3715SXin Li 
626*bf2c3715SXin Li /** Sets \c *this to be a quaternion representing a rotation between
627*bf2c3715SXin Li   * the two arbitrary vectors \a a and \a b. In other words, the built
628*bf2c3715SXin Li   * rotation represent a rotation sending the line of direction \a a
629*bf2c3715SXin Li   * to the line of direction \a b, both lines passing through the origin.
630*bf2c3715SXin Li   *
631*bf2c3715SXin Li   * \returns a reference to \c *this.
632*bf2c3715SXin Li   *
633*bf2c3715SXin Li   * Note that the two input vectors do \b not have to be normalized, and
634*bf2c3715SXin Li   * do not need to have the same norm.
635*bf2c3715SXin Li   */
636*bf2c3715SXin Li template<class Derived>
637*bf2c3715SXin Li template<typename Derived1, typename Derived2>
638*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline Derived& QuaternionBase<Derived>::setFromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b)
639*bf2c3715SXin Li {
640*bf2c3715SXin Li   EIGEN_USING_STD(sqrt)
641*bf2c3715SXin Li   Vector3 v0 = a.normalized();
642*bf2c3715SXin Li   Vector3 v1 = b.normalized();
643*bf2c3715SXin Li   Scalar c = v1.dot(v0);
644*bf2c3715SXin Li 
645*bf2c3715SXin Li   // if dot == -1, vectors are nearly opposites
646*bf2c3715SXin Li   // => accurately compute the rotation axis by computing the
647*bf2c3715SXin Li   //    intersection of the two planes. This is done by solving:
648*bf2c3715SXin Li   //       x^T v0 = 0
649*bf2c3715SXin Li   //       x^T v1 = 0
650*bf2c3715SXin Li   //    under the constraint:
651*bf2c3715SXin Li   //       ||x|| = 1
652*bf2c3715SXin Li   //    which yields a singular value problem
653*bf2c3715SXin Li   if (c < Scalar(-1)+NumTraits<Scalar>::dummy_precision())
654*bf2c3715SXin Li   {
655*bf2c3715SXin Li     c = numext::maxi(c,Scalar(-1));
656*bf2c3715SXin Li     Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
657*bf2c3715SXin Li     JacobiSVD<Matrix<Scalar,2,3> > svd(m, ComputeFullV);
658*bf2c3715SXin Li     Vector3 axis = svd.matrixV().col(2);
659*bf2c3715SXin Li 
660*bf2c3715SXin Li     Scalar w2 = (Scalar(1)+c)*Scalar(0.5);
661*bf2c3715SXin Li     this->w() = sqrt(w2);
662*bf2c3715SXin Li     this->vec() = axis * sqrt(Scalar(1) - w2);
663*bf2c3715SXin Li     return derived();
664*bf2c3715SXin Li   }
665*bf2c3715SXin Li   Vector3 axis = v0.cross(v1);
666*bf2c3715SXin Li   Scalar s = sqrt((Scalar(1)+c)*Scalar(2));
667*bf2c3715SXin Li   Scalar invs = Scalar(1)/s;
668*bf2c3715SXin Li   this->vec() = axis * invs;
669*bf2c3715SXin Li   this->w() = s * Scalar(0.5);
670*bf2c3715SXin Li 
671*bf2c3715SXin Li   return derived();
672*bf2c3715SXin Li }
673*bf2c3715SXin Li 
674*bf2c3715SXin Li /** \returns a random unit quaternion following a uniform distribution law on SO(3)
675*bf2c3715SXin Li   *
676*bf2c3715SXin Li   * \note The implementation is based on http://planning.cs.uiuc.edu/node198.html
677*bf2c3715SXin Li   */
678*bf2c3715SXin Li template<typename Scalar, int Options>
679*bf2c3715SXin Li EIGEN_DEVICE_FUNC Quaternion<Scalar,Options> Quaternion<Scalar,Options>::UnitRandom()
680*bf2c3715SXin Li {
681*bf2c3715SXin Li   EIGEN_USING_STD(sqrt)
682*bf2c3715SXin Li   EIGEN_USING_STD(sin)
683*bf2c3715SXin Li   EIGEN_USING_STD(cos)
684*bf2c3715SXin Li   const Scalar u1 = internal::random<Scalar>(0, 1),
685*bf2c3715SXin Li                u2 = internal::random<Scalar>(0, 2*EIGEN_PI),
686*bf2c3715SXin Li                u3 = internal::random<Scalar>(0, 2*EIGEN_PI);
687*bf2c3715SXin Li   const Scalar a = sqrt(Scalar(1) - u1),
688*bf2c3715SXin Li                b = sqrt(u1);
689*bf2c3715SXin Li   return Quaternion (a * sin(u2), a * cos(u2), b * sin(u3), b * cos(u3));
690*bf2c3715SXin Li }
691*bf2c3715SXin Li 
692*bf2c3715SXin Li 
693*bf2c3715SXin Li /** Returns a quaternion representing a rotation between
694*bf2c3715SXin Li   * the two arbitrary vectors \a a and \a b. In other words, the built
695*bf2c3715SXin Li   * rotation represent a rotation sending the line of direction \a a
696*bf2c3715SXin Li   * to the line of direction \a b, both lines passing through the origin.
697*bf2c3715SXin Li   *
698*bf2c3715SXin Li   * \returns resulting quaternion
699*bf2c3715SXin Li   *
700*bf2c3715SXin Li   * Note that the two input vectors do \b not have to be normalized, and
701*bf2c3715SXin Li   * do not need to have the same norm.
702*bf2c3715SXin Li   */
703*bf2c3715SXin Li template<typename Scalar, int Options>
704*bf2c3715SXin Li template<typename Derived1, typename Derived2>
705*bf2c3715SXin Li EIGEN_DEVICE_FUNC Quaternion<Scalar,Options> Quaternion<Scalar,Options>::FromTwoVectors(const MatrixBase<Derived1>& a, const MatrixBase<Derived2>& b)
706*bf2c3715SXin Li {
707*bf2c3715SXin Li     Quaternion quat;
708*bf2c3715SXin Li     quat.setFromTwoVectors(a, b);
709*bf2c3715SXin Li     return quat;
710*bf2c3715SXin Li }
711*bf2c3715SXin Li 
712*bf2c3715SXin Li 
713*bf2c3715SXin Li /** \returns the multiplicative inverse of \c *this
714*bf2c3715SXin Li   * Note that in most cases, i.e., if you simply want the opposite rotation,
715*bf2c3715SXin Li   * and/or the quaternion is normalized, then it is enough to use the conjugate.
716*bf2c3715SXin Li   *
717*bf2c3715SXin Li   * \sa QuaternionBase::conjugate()
718*bf2c3715SXin Li   */
719*bf2c3715SXin Li template <class Derived>
720*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline Quaternion<typename internal::traits<Derived>::Scalar> QuaternionBase<Derived>::inverse() const
721*bf2c3715SXin Li {
722*bf2c3715SXin Li   // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite()  ??
723*bf2c3715SXin Li   Scalar n2 = this->squaredNorm();
724*bf2c3715SXin Li   if (n2 > Scalar(0))
725*bf2c3715SXin Li     return Quaternion<Scalar>(conjugate().coeffs() / n2);
726*bf2c3715SXin Li   else
727*bf2c3715SXin Li   {
728*bf2c3715SXin Li     // return an invalid result to flag the error
729*bf2c3715SXin Li     return Quaternion<Scalar>(Coefficients::Zero());
730*bf2c3715SXin Li   }
731*bf2c3715SXin Li }
732*bf2c3715SXin Li 
733*bf2c3715SXin Li // Generic conjugate of a Quaternion
734*bf2c3715SXin Li namespace internal {
735*bf2c3715SXin Li template<int Arch, class Derived, typename Scalar> struct quat_conj
736*bf2c3715SXin Li {
737*bf2c3715SXin Li   EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE Quaternion<Scalar> run(const QuaternionBase<Derived>& q){
738*bf2c3715SXin Li     return Quaternion<Scalar>(q.w(),-q.x(),-q.y(),-q.z());
739*bf2c3715SXin Li   }
740*bf2c3715SXin Li };
741*bf2c3715SXin Li }
742*bf2c3715SXin Li 
743*bf2c3715SXin Li /** \returns the conjugate of the \c *this which is equal to the multiplicative inverse
744*bf2c3715SXin Li   * if the quaternion is normalized.
745*bf2c3715SXin Li   * The conjugate of a quaternion represents the opposite rotation.
746*bf2c3715SXin Li   *
747*bf2c3715SXin Li   * \sa Quaternion2::inverse()
748*bf2c3715SXin Li   */
749*bf2c3715SXin Li template <class Derived>
750*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline Quaternion<typename internal::traits<Derived>::Scalar>
751*bf2c3715SXin Li QuaternionBase<Derived>::conjugate() const
752*bf2c3715SXin Li {
753*bf2c3715SXin Li   return internal::quat_conj<Architecture::Target, Derived,
754*bf2c3715SXin Li                          typename internal::traits<Derived>::Scalar>::run(*this);
755*bf2c3715SXin Li 
756*bf2c3715SXin Li }
757*bf2c3715SXin Li 
758*bf2c3715SXin Li /** \returns the angle (in radian) between two rotations
759*bf2c3715SXin Li   * \sa dot()
760*bf2c3715SXin Li   */
761*bf2c3715SXin Li template <class Derived>
762*bf2c3715SXin Li template <class OtherDerived>
763*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline typename internal::traits<Derived>::Scalar
764*bf2c3715SXin Li QuaternionBase<Derived>::angularDistance(const QuaternionBase<OtherDerived>& other) const
765*bf2c3715SXin Li {
766*bf2c3715SXin Li   EIGEN_USING_STD(atan2)
767*bf2c3715SXin Li   Quaternion<Scalar> d = (*this) * other.conjugate();
768*bf2c3715SXin Li   return Scalar(2) * atan2( d.vec().norm(), numext::abs(d.w()) );
769*bf2c3715SXin Li }
770*bf2c3715SXin Li 
771*bf2c3715SXin Li 
772*bf2c3715SXin Li 
773*bf2c3715SXin Li /** \returns the spherical linear interpolation between the two quaternions
774*bf2c3715SXin Li   * \c *this and \a other at the parameter \a t in [0;1].
775*bf2c3715SXin Li   *
776*bf2c3715SXin Li   * This represents an interpolation for a constant motion between \c *this and \a other,
777*bf2c3715SXin Li   * see also http://en.wikipedia.org/wiki/Slerp.
778*bf2c3715SXin Li   */
779*bf2c3715SXin Li template <class Derived>
780*bf2c3715SXin Li template <class OtherDerived>
781*bf2c3715SXin Li EIGEN_DEVICE_FUNC Quaternion<typename internal::traits<Derived>::Scalar>
782*bf2c3715SXin Li QuaternionBase<Derived>::slerp(const Scalar& t, const QuaternionBase<OtherDerived>& other) const
783*bf2c3715SXin Li {
784*bf2c3715SXin Li   EIGEN_USING_STD(acos)
785*bf2c3715SXin Li   EIGEN_USING_STD(sin)
786*bf2c3715SXin Li   const Scalar one = Scalar(1) - NumTraits<Scalar>::epsilon();
787*bf2c3715SXin Li   Scalar d = this->dot(other);
788*bf2c3715SXin Li   Scalar absD = numext::abs(d);
789*bf2c3715SXin Li 
790*bf2c3715SXin Li   Scalar scale0;
791*bf2c3715SXin Li   Scalar scale1;
792*bf2c3715SXin Li 
793*bf2c3715SXin Li   if(absD>=one)
794*bf2c3715SXin Li   {
795*bf2c3715SXin Li     scale0 = Scalar(1) - t;
796*bf2c3715SXin Li     scale1 = t;
797*bf2c3715SXin Li   }
798*bf2c3715SXin Li   else
799*bf2c3715SXin Li   {
800*bf2c3715SXin Li     // theta is the angle between the 2 quaternions
801*bf2c3715SXin Li     Scalar theta = acos(absD);
802*bf2c3715SXin Li     Scalar sinTheta = sin(theta);
803*bf2c3715SXin Li 
804*bf2c3715SXin Li     scale0 = sin( ( Scalar(1) - t ) * theta) / sinTheta;
805*bf2c3715SXin Li     scale1 = sin( ( t * theta) ) / sinTheta;
806*bf2c3715SXin Li   }
807*bf2c3715SXin Li   if(d<Scalar(0)) scale1 = -scale1;
808*bf2c3715SXin Li 
809*bf2c3715SXin Li   return Quaternion<Scalar>(scale0 * coeffs() + scale1 * other.coeffs());
810*bf2c3715SXin Li }
811*bf2c3715SXin Li 
812*bf2c3715SXin Li namespace internal {
813*bf2c3715SXin Li 
814*bf2c3715SXin Li // set from a rotation matrix
815*bf2c3715SXin Li template<typename Other>
816*bf2c3715SXin Li struct quaternionbase_assign_impl<Other,3,3>
817*bf2c3715SXin Li {
818*bf2c3715SXin Li   typedef typename Other::Scalar Scalar;
819*bf2c3715SXin Li   template<class Derived> EIGEN_DEVICE_FUNC static inline void run(QuaternionBase<Derived>& q, const Other& a_mat)
820*bf2c3715SXin Li   {
821*bf2c3715SXin Li     const typename internal::nested_eval<Other,2>::type mat(a_mat);
822*bf2c3715SXin Li     EIGEN_USING_STD(sqrt)
823*bf2c3715SXin Li     // This algorithm comes from  "Quaternion Calculus and Fast Animation",
824*bf2c3715SXin Li     // Ken Shoemake, 1987 SIGGRAPH course notes
825*bf2c3715SXin Li     Scalar t = mat.trace();
826*bf2c3715SXin Li     if (t > Scalar(0))
827*bf2c3715SXin Li     {
828*bf2c3715SXin Li       t = sqrt(t + Scalar(1.0));
829*bf2c3715SXin Li       q.w() = Scalar(0.5)*t;
830*bf2c3715SXin Li       t = Scalar(0.5)/t;
831*bf2c3715SXin Li       q.x() = (mat.coeff(2,1) - mat.coeff(1,2)) * t;
832*bf2c3715SXin Li       q.y() = (mat.coeff(0,2) - mat.coeff(2,0)) * t;
833*bf2c3715SXin Li       q.z() = (mat.coeff(1,0) - mat.coeff(0,1)) * t;
834*bf2c3715SXin Li     }
835*bf2c3715SXin Li     else
836*bf2c3715SXin Li     {
837*bf2c3715SXin Li       Index i = 0;
838*bf2c3715SXin Li       if (mat.coeff(1,1) > mat.coeff(0,0))
839*bf2c3715SXin Li         i = 1;
840*bf2c3715SXin Li       if (mat.coeff(2,2) > mat.coeff(i,i))
841*bf2c3715SXin Li         i = 2;
842*bf2c3715SXin Li       Index j = (i+1)%3;
843*bf2c3715SXin Li       Index k = (j+1)%3;
844*bf2c3715SXin Li 
845*bf2c3715SXin Li       t = sqrt(mat.coeff(i,i)-mat.coeff(j,j)-mat.coeff(k,k) + Scalar(1.0));
846*bf2c3715SXin Li       q.coeffs().coeffRef(i) = Scalar(0.5) * t;
847*bf2c3715SXin Li       t = Scalar(0.5)/t;
848*bf2c3715SXin Li       q.w() = (mat.coeff(k,j)-mat.coeff(j,k))*t;
849*bf2c3715SXin Li       q.coeffs().coeffRef(j) = (mat.coeff(j,i)+mat.coeff(i,j))*t;
850*bf2c3715SXin Li       q.coeffs().coeffRef(k) = (mat.coeff(k,i)+mat.coeff(i,k))*t;
851*bf2c3715SXin Li     }
852*bf2c3715SXin Li   }
853*bf2c3715SXin Li };
854*bf2c3715SXin Li 
855*bf2c3715SXin Li // set from a vector of coefficients assumed to be a quaternion
856*bf2c3715SXin Li template<typename Other>
857*bf2c3715SXin Li struct quaternionbase_assign_impl<Other,4,1>
858*bf2c3715SXin Li {
859*bf2c3715SXin Li   typedef typename Other::Scalar Scalar;
860*bf2c3715SXin Li   template<class Derived> EIGEN_DEVICE_FUNC static inline void run(QuaternionBase<Derived>& q, const Other& vec)
861*bf2c3715SXin Li   {
862*bf2c3715SXin Li     q.coeffs() = vec;
863*bf2c3715SXin Li   }
864*bf2c3715SXin Li };
865*bf2c3715SXin Li 
866*bf2c3715SXin Li } // end namespace internal
867*bf2c3715SXin Li 
868*bf2c3715SXin Li } // end namespace Eigen
869*bf2c3715SXin Li 
870*bf2c3715SXin Li #endif // EIGEN_QUATERNION_H
871