From 611d2b0b1d64880ff207368a901c7faba73a78c5 Mon Sep 17 00:00:00 2001 From: Mathieu Gautier Date: Tue, 27 Oct 2009 13:19:16 +0000 Subject: Quaternion could now map an array of 4 scalars : new classes : * QuaternionBase * Map --- Eigen/src/Core/util/ForwardDeclarations.h | 1 + Eigen/src/Geometry/Quaternion.h | 545 ++++++++++++++++++++++++++++++ Eigen/src/Geometry/arch/Geometry_SSE.h | 21 ++ 3 files changed, 567 insertions(+) diff --git a/Eigen/src/Core/util/ForwardDeclarations.h b/Eigen/src/Core/util/ForwardDeclarations.h index 65e5ce687..025904734 100644 --- a/Eigen/src/Core/util/ForwardDeclarations.h +++ b/Eigen/src/Core/util/ForwardDeclarations.h @@ -129,6 +129,7 @@ template class PlanarRotation; // Geometry module: template class RotationBase; template class Cross; +template class QuaternionBase; template class Quaternion; template class Rotation2D; template class AngleAxis; diff --git a/Eigen/src/Geometry/Quaternion.h b/Eigen/src/Geometry/Quaternion.h index 2f9f97807..f3f60a08a 100644 --- a/Eigen/src/Geometry/Quaternion.h +++ b/Eigen/src/Geometry/Quaternion.h @@ -507,4 +507,549 @@ struct ei_quaternion_assign_impl } }; +/*################################################################### + QuaternionBase and Map and Quat + ###################################################################*/ + +template +struct ei_quaternionbase_assign_impl; + +template class Quat; // [XXX] => remove when Quat becomes Quaternion + +template +struct ei_traits > +{ + typedef typename ei_traits::Scalar Scalar; + enum { + PacketAccess = ei_traits::PacketAccess + }; +}; + +template +class QuaternionBase : public RotationBase { + typedef RotationBase Base; +public: + using Base::operator*; + + typedef typename ei_traits >::Scalar Scalar; + typedef typename NumTraits::Real RealScalar; + + // typedef typename Matrix Coefficients; + /** the type of a 3D vector */ + typedef Matrix Vector3; + /** the equivalent rotation matrix type */ + typedef Matrix Matrix3; + /** the equivalent angle-axis type */ + typedef AngleAxis AngleAxisType; + + /** \returns the \c x coefficient */ + inline Scalar x() const { return this->derived().coeffs().coeff(0); } + /** \returns the \c y coefficient */ + inline Scalar y() const { return this->derived().coeffs().coeff(1); } + /** \returns the \c z coefficient */ + inline Scalar z() const { return this->derived().coeffs().coeff(2); } + /** \returns the \c w coefficient */ + inline Scalar w() const { return this->derived().coeffs().coeff(3); } + + /** \returns a reference to the \c x coefficient */ + inline Scalar& x() { return this->derived().coeffs().coeffRef(0); } + /** \returns a reference to the \c y coefficient */ + inline Scalar& y() { return this->derived().coeffs().coeffRef(1); } + /** \returns a reference to the \c z coefficient */ + inline Scalar& z() { return this->derived().coeffs().coeffRef(2); } + /** \returns a reference to the \c w coefficient */ + inline Scalar& w() { return this->derived().coeffs().coeffRef(3); } + + /** \returns a read-only vector expression of the imaginary part (x,y,z) */ + inline const VectorBlock::Coefficients,3,1> vec() const { return this->derived().coeffs().template start<3>(); } + + /** \returns a vector expression of the imaginary part (x,y,z) */ + inline VectorBlock::Coefficients,3,1> vec() { return this->derived().coeffs().template start<3>(); } + + /** \returns a read-only vector expression of the coefficients (x,y,z,w) */ + inline const typename ei_traits::Coefficients& coeffs() const { return this->derived().coeffs(); } + + /** \returns a vector expression of the coefficients (x,y,z,w) */ + inline typename ei_traits::Coefficients& coeffs() { return this->derived().coeffs(); } + + template QuaternionBase& operator=(const QuaternionBase& other); + QuaternionBase& operator=(const AngleAxisType& aa); + template + QuaternionBase& operator=(const MatrixBase& m); + + /** \returns a quaternion representing an identity rotation + * \sa MatrixBase::Identity() + */ + inline static Quat Identity() { return Quat(1, 0, 0, 0); } + + /** \sa Quaternion2::Identity(), MatrixBase::setIdentity() + */ + inline QuaternionBase& setIdentity() { coeffs() << 0, 0, 0, 1; return *this; } + + /** \returns the squared norm of the quaternion's coefficients + * \sa Quaternion2::norm(), MatrixBase::squaredNorm() + */ + inline Scalar squaredNorm() const { return coeffs().squaredNorm(); } + + /** \returns the norm of the quaternion's coefficients + * \sa Quaternion2::squaredNorm(), MatrixBase::norm() + */ + inline Scalar norm() const { return coeffs().norm(); } + + /** Normalizes the quaternion \c *this + * \sa normalized(), MatrixBase::normalize() */ + inline void normalize() { coeffs().normalize(); } + /** \returns a normalized version of \c *this + * \sa normalize(), MatrixBase::normalized() */ + inline Quat normalized() const { return Quat(coeffs().normalized()); } + + /** \returns the dot product of \c *this and \a other + * Geometrically speaking, the dot product of two unit quaternions + * corresponds to the cosine of half the angle between the two rotations. + * \sa angularDistance() + */ + template inline Scalar dot(const QuaternionBase& other) const { return coeffs().dot(other.coeffs()); } + + template inline Scalar angularDistance(const QuaternionBase& other) const; + + Matrix3 toRotationMatrix(void) const; + + template + QuaternionBase& setFromTwoVectors(const MatrixBase& a, const MatrixBase& b); + + template inline Quat operator* (const QuaternionBase& q) const; + template inline QuaternionBase& operator*= (const QuaternionBase& q); + + Quat inverse(void) const; + Quat conjugate(void) const; + + template Quat slerp(Scalar t, const QuaternionBase& other) const; + + /** \returns \c true if \c *this is approximately equal to \a other, within the precision + * determined by \a prec. + * + * \sa MatrixBase::isApprox() */ + bool isApprox(const QuaternionBase& other, typename RealScalar prec = precision()) const + { return coeffs().isApprox(other.coeffs(), prec); } + + Vector3 _transformVector(Vector3 v) const; + +}; + +/* ########### Quat -> Quaternion */ + +template +struct ei_traits > +{ + typedef _Scalar Scalar; + typedef Matrix<_Scalar,4,1> Coefficients; + enum{ + PacketAccess = Aligned + }; +}; + +template +class Quat : public QuaternionBase >{ + typedef QuaternionBase > Base; +public: + using Base::operator=; + + typedef _Scalar Scalar; + + typedef typename ei_traits >::Coefficients Coefficients; + typedef typename Base::AngleAxisType AngleAxisType; + + /** Default constructor leaving the quaternion uninitialized. */ + inline Quat() {} + + /** Constructs and initializes the quaternion \f$ w+xi+yj+zk \f$ from + * its four coefficients \a w, \a x, \a y and \a z. + * + * \warning Note the order of the arguments: the real \a w coefficient first, + * while internally the coefficients are stored in the following order: + * [\c x, \c y, \c z, \c w] + */ + inline Quat(Scalar w, Scalar x, Scalar y, Scalar z) + { coeffs() << x, y, z, w; } + + /** Constructs and initialize a quaternion from the array data + * This constructor is also used to map an array */ + inline Quat(const Scalar* data) : m_coeffs(data) {} + + /** Copy constructor */ +// template inline Quat(const QuaternionBase& other) { m_coeffs = other.coeffs(); } [XXX] redundant with 703 + + /** Constructs and initializes a quaternion from the angle-axis \a aa */ + explicit inline Quat(const AngleAxisType& aa) { *this = aa; } + + /** Constructs and initializes a quaternion from either: + * - a rotation matrix expression, + * - a 4D vector expression representing quaternion coefficients. + */ + template + explicit inline Quat(const MatrixBase& other) { *this = other; } + + /** \returns \c *this with scalar type casted to \a NewScalarType + * + * Note that if \a NewScalarType is equal to the current scalar type of \c *this + * then this function smartly returns a const reference to \c *this. + */ + template + inline typename ei_cast_return_type >::type cast() const + { return typename ei_cast_return_type >::type(*this); } + + /** Copy constructor with scalar type conversion */ + template + inline explicit Quat(const QuaternionBase& other) + { m_coeffs = other.coeffs().template cast(); } + + inline Coefficients& coeffs() { return m_coeffs;} + inline const Coefficients& coeffs() const { return m_coeffs;} + +protected: + Coefficients m_coeffs; +}; + +/* ########### Map */ + +/** \class Map + * \nonstableyet + * + * \brief Expression of a quaternion + * + * \param Scalar the type of the vector of diagonal coefficients + * + * \sa class Quaternion, class QuaternionBase + */ +template +struct ei_traits, _PacketAccess> >: +ei_traits > +{ + typedef _Scalar Scalar; + typedef Map > Coefficients; + enum { + PacketAccess = _PacketAccess + }; +}; + +template +class Map, PacketAccess > : public QuaternionBase, PacketAccess> >, ei_no_assignment_operator { + public: + + typedef _Scalar Scalar; + + typedef typename ei_traits, PacketAccess> >::Coefficients Coefficients; + + inline Map, PacketAccess >(const Scalar* coeffs) : m_coeffs(coeffs) {} + + inline Coefficients& coeffs() { return m_coeffs;} + inline const Coefficients& coeffs() const { return m_coeffs;} + + protected: + Coefficients m_coeffs; +}; + +typedef Map > QuaternionMapd; +typedef Map > QuaternionMapf; +typedef Map, Aligned> QuaternionMapAlignedd; +typedef Map, Aligned> QuaternionMapAlignedf; + +// Generic Quaternion * Quaternion product +template struct ei_quat_product +{ + inline static Quat run(const QuaternionBase& a, const QuaternionBase& b){ + return Quat + ( + a.w() * b.w() - a.x() * b.x() - a.y() * b.y() - a.z() * b.z(), + a.w() * b.x() + a.x() * b.w() + a.y() * b.z() - a.z() * b.y(), + a.w() * b.y() + a.y() * b.w() + a.z() * b.x() - a.x() * b.z(), + a.w() * b.z() + a.z() * b.w() + a.x() * b.y() - a.y() * b.x() + ); + } +}; + +/** \returns the concatenation of two rotations as a quaternion-quaternion product */ +template +template +inline Quat >::Scalar> QuaternionBase::operator* (const QuaternionBase& other) const +{ + EIGEN_STATIC_ASSERT((ei_is_same_type::ret), + YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY) + return ei_quat_product::Scalar, + ei_traits::PacketAccess && ei_traits::PacketAccess>::run(*this, other); +} + +/** \sa operator*(Quaternion) */ +template +template +inline QuaternionBase& QuaternionBase::operator*= (const QuaternionBase& other) +{ + return (*this = *this * other); +} + +/** Rotation of a vector by a quaternion. + * \remarks If the quaternion is used to rotate several points (>1) + * then it is much more efficient to first convert it to a 3x3 Matrix. + * Comparison of the operation cost for n transformations: + * - Quaternion2: 30n + * - Via a Matrix3: 24 + 15n + */ +template +inline typename QuaternionBase::Vector3 +QuaternionBase::_transformVector(Vector3 v) const +{ + // Note that this algorithm comes from the optimization by hand + // of the conversion to a Matrix followed by a Matrix/Vector product. + // It appears to be much faster than the common algorithm found + // in the litterature (30 versus 39 flops). It also requires two + // Vector3 as temporaries. + Vector3 uv = Scalar(2) * this->vec().cross(v); + return v + this->w() * uv + this->vec().cross(uv); +} + +template +template +inline QuaternionBase& QuaternionBase::operator=(const QuaternionBase& other) +{ + coeffs() = other.coeffs(); + return *this; +} + +/** Set \c *this from an angle-axis \a aa and returns a reference to \c *this + */ +template +inline QuaternionBase& QuaternionBase::operator=(const AngleAxisType& aa) +{ + Scalar ha = Scalar(0.5)*aa.angle(); // Scalar(0.5) to suppress precision loss warnings + this->w() = ei_cos(ha); + this->vec() = ei_sin(ha) * aa.axis(); + return *this; +} + +/** Set \c *this from the expression \a xpr: + * - if \a xpr is a 4x1 vector, then \a xpr is assumed to be a quaternion + * - if \a xpr is a 3x3 matrix, then \a xpr is assumed to be rotation matrix + * and \a xpr is converted to a quaternion + */ + +template +template +inline QuaternionBase& QuaternionBase::operator=(const MatrixBase& xpr) +{ + EIGEN_STATIC_ASSERT((ei_is_same_type::ret), + YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY) + ei_quaternionbase_assign_impl::run(*this, xpr.derived()); + return *this; +} + +/** Convert the quaternion to a 3x3 rotation matrix. The quaternion is required to + * be normalized, otherwise the result is undefined. + */ +template +inline typename QuaternionBase::Matrix3 +QuaternionBase::toRotationMatrix(void) const +{ + // NOTE if inlined, then gcc 4.2 and 4.4 get rid of the temporary (not gcc 4.3 !!) + // if not inlined then the cost of the return by value is huge ~ +35%, + // however, not inlining this function is an order of magnitude slower, so + // it has to be inlined, and so the return by value is not an issue + Matrix3 res; + + const Scalar tx = 2*this->x(); + const Scalar ty = 2*this->y(); + const Scalar tz = 2*this->z(); + const Scalar twx = tx*this->w(); + const Scalar twy = ty*this->w(); + const Scalar twz = tz*this->w(); + const Scalar txx = tx*this->x(); + const Scalar txy = ty*this->x(); + const Scalar txz = tz*this->x(); + const Scalar tyy = ty*this->y(); + const Scalar tyz = tz*this->y(); + const Scalar tzz = tz*this->z(); + + res.coeffRef(0,0) = 1-(tyy+tzz); + res.coeffRef(0,1) = txy-twz; + res.coeffRef(0,2) = txz+twy; + res.coeffRef(1,0) = txy+twz; + res.coeffRef(1,1) = 1-(txx+tzz); + res.coeffRef(1,2) = tyz-twx; + res.coeffRef(2,0) = txz-twy; + res.coeffRef(2,1) = tyz+twx; + res.coeffRef(2,2) = 1-(txx+tyy); + + return res; +} + +/** Sets \c *this to be a quaternion representing a rotation between + * the two arbitrary vectors \a a and \a b. In other words, the built + * rotation represent a rotation sending the line of direction \a a + * to the line of direction \a b, both lines passing through the origin. + * + * \returns a reference to \c *this. + * + * Note that the two input vectors do \b not have to be normalized, and + * do not need to have the same norm. + */ +template +template +inline QuaternionBase& QuaternionBase::setFromTwoVectors(const MatrixBase& a, const MatrixBase& b) +{ + Vector3 v0 = a.normalized(); + Vector3 v1 = b.normalized(); + Scalar c = v1.dot(v0); + + // if dot == -1, vectors are nearly opposites + // => accuraletly compute the rotation axis by computing the + // intersection of the two planes. This is done by solving: + // x^T v0 = 0 + // x^T v1 = 0 + // under the constraint: + // ||x|| = 1 + // which yields a singular value problem + if (c < Scalar(-1)+precision()) + { + c = std::max(c,-1); + Matrix m; m << v0.transpose(), v1.transpose(); + SVD > svd(m); + Vector3 axis = svd.matrixV().col(2); + + Scalar w2 = (Scalar(1)+c)*Scalar(0.5); + this->w() = ei_sqrt(w2); + this->vec() = axis * ei_sqrt(Scalar(1) - w2); + return *this; + } + Vector3 axis = v0.cross(v1); + Scalar s = ei_sqrt((Scalar(1)+c)*Scalar(2)); + Scalar invs = Scalar(1)/s; + this->vec() = axis * invs; + this->w() = s * Scalar(0.5); + + return *this; +} + +/** \returns the multiplicative inverse of \c *this + * Note that in most cases, i.e., if you simply want the opposite rotation, + * and/or the quaternion is normalized, then it is enough to use the conjugate. + * + * \sa Quaternion2::conjugate() + */ +template +inline Quat >::Scalar> QuaternionBase::inverse() const +{ + // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite() ?? + Scalar n2 = this->squaredNorm(); + if (n2 > 0) + return Quat(conjugate().coeffs() / n2); + else + { + // return an invalid result to flag the error + return Quat(ei_traits::Coefficients::Zero()); + } +} + +/** \returns the conjugate of the \c *this which is equal to the multiplicative inverse + * if the quaternion is normalized. + * The conjugate of a quaternion represents the opposite rotation. + * + * \sa Quaternion2::inverse() + */ +template +inline Quat >::Scalar> QuaternionBase::conjugate() const +{ + return Quat(this->w(),-this->x(),-this->y(),-this->z()); +} + +/** \returns the angle (in radian) between two rotations + * \sa dot() + */ +template +template +inline typename ei_traits >::Scalar QuaternionBase::angularDistance(const QuaternionBase& other) const +{ + double d = ei_abs(this->dot(other)); + if (d>=1.0) + return 0; + return Scalar(2) * std::acos(d); +} + +/** \returns the spherical linear interpolation between the two quaternions + * \c *this and \a other at the parameter \a t + */ +template +template +Quat >::Scalar> QuaternionBase::slerp(Scalar t, const QuaternionBase& other) const +{ + static const Scalar one = Scalar(1) - precision(); + Scalar d = this->dot(other); + Scalar absD = ei_abs(d); + if (absD>=one) + return Quat(*this); + + // theta is the angle between the 2 quaternions + Scalar theta = std::acos(absD); + Scalar sinTheta = ei_sin(theta); + + Scalar scale0 = ei_sin( ( Scalar(1) - t ) * theta) / sinTheta; + Scalar scale1 = ei_sin( ( t * theta) ) / sinTheta; + if (d<0) + scale1 = -scale1; + + return Quat(scale0 * coeffs() + scale1 * other.coeffs()); +} + +// set from a rotation matrix +template +struct ei_quaternionbase_assign_impl +{ + typedef typename Other::Scalar Scalar; + template inline static void run(QuaternionBase& q, const Other& mat) + { + // This algorithm comes from "Quaternion Calculus and Fast Animation", + // Ken Shoemake, 1987 SIGGRAPH course notes + Scalar t = mat.trace(); + if (t > 0) + { + t = ei_sqrt(t + Scalar(1.0)); + q.w() = Scalar(0.5)*t; + t = Scalar(0.5)/t; + q.x() = (mat.coeff(2,1) - mat.coeff(1,2)) * t; + q.y() = (mat.coeff(0,2) - mat.coeff(2,0)) * t; + q.z() = (mat.coeff(1,0) - mat.coeff(0,1)) * t; + } + else + { + int i = 0; + if (mat.coeff(1,1) > mat.coeff(0,0)) + i = 1; + if (mat.coeff(2,2) > mat.coeff(i,i)) + i = 2; + int j = (i+1)%3; + int k = (j+1)%3; + + t = ei_sqrt(mat.coeff(i,i)-mat.coeff(j,j)-mat.coeff(k,k) + Scalar(1.0)); + q.coeffs().coeffRef(i) = Scalar(0.5) * t; + t = Scalar(0.5)/t; + q.w() = (mat.coeff(k,j)-mat.coeff(j,k))*t; + q.coeffs().coeffRef(j) = (mat.coeff(j,i)+mat.coeff(i,j))*t; + q.coeffs().coeffRef(k) = (mat.coeff(k,i)+mat.coeff(i,k))*t; + } + } +}; + +// set from a vector of coefficients assumed to be a quaternion +template +struct ei_quaternionbase_assign_impl +{ + typedef typename Other::Scalar Scalar; + template inline static void run(QuaternionBase& q, const Other& vec) + { + q.coeffs() = vec; + } +}; + + #endif // EIGEN_QUATERNION_H diff --git a/Eigen/src/Geometry/arch/Geometry_SSE.h b/Eigen/src/Geometry/arch/Geometry_SSE.h index d0342febc..c608e4843 100644 --- a/Eigen/src/Geometry/arch/Geometry_SSE.h +++ b/Eigen/src/Geometry/arch/Geometry_SSE.h @@ -45,6 +45,27 @@ ei_quaternion_product(const Quaternion& _a, const Quate return res; } +template struct ei_quat_product +{ + inline static Quat run(const QuaternionBase& _a, const QuaternionBase& _b) + { + const __m128 mask = _mm_castsi128_ps(_mm_setr_epi32(0,0,0,0x80000000)); + Quat res; + __m128 a = _a.coeffs().packet(0); + __m128 b = _b.coeffs().packet(0); + __m128 flip1 = _mm_xor_ps(_mm_mul_ps(ei_vec4f_swizzle1(a,1,2,0,2), + ei_vec4f_swizzle1(b,2,0,1,2)),mask); + __m128 flip2 = _mm_xor_ps(_mm_mul_ps(ei_vec4f_swizzle1(a,3,3,3,1), + ei_vec4f_swizzle1(b,0,1,2,1)),mask); + ei_pstore(&res.x(), + _mm_add_ps(_mm_sub_ps(_mm_mul_ps(a,ei_vec4f_swizzle1(b,3,3,3,3)), + _mm_mul_ps(ei_vec4f_swizzle1(a,2,0,1,0), + ei_vec4f_swizzle1(b,1,2,0,0))), + _mm_add_ps(flip1,flip2))); + return res; + } +}; + template struct ei_cross3_impl { inline static typename ei_plain_matrix_type::type -- cgit v1.2.3