Quaternion.h 34 KB

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