MOAB  4.9.3pre
Eigen::QuaternionBase< Derived > Class Template Reference

Base class for quaternion expressions. More...

#include <Quaternion.h>

Inheritance diagram for Eigen::QuaternionBase< Derived >:
Collaboration diagram for Eigen::QuaternionBase< Derived >:

List of all members.

Public Types

enum  { Flags = Eigen::internal::traits<Derived>::Flags }
typedef RotationBase< Derived, 3 > Base
typedef internal::traits
< Derived >::Scalar 
Scalar
typedef NumTraits< Scalar >::Real RealScalar
typedef internal::traits
< Derived >::Coefficients 
Coefficients
typedef Matrix< Scalar, 3, 1 > Vector3
typedef Matrix< Scalar, 3, 3 > Matrix3
typedef AngleAxis< ScalarAngleAxisType

Public Member Functions

Scalar x () const
Scalar y () const
Scalar z () const
Scalar w () const
Scalarx ()
Scalary ()
Scalarz ()
Scalarw ()
const VectorBlock< const
Coefficients, 3 > 
vec () const
VectorBlock< Coefficients, 3 > vec ()
const internal::traits
< Derived >::Coefficients
coeffs () const
internal::traits< Derived >
::Coefficients
coeffs ()
EIGEN_STRONG_INLINE
QuaternionBase< Derived > & 
operator= (const QuaternionBase< Derived > &other)
template<class OtherDerived >
EIGEN_STRONG_INLINE Derived & operator= (const QuaternionBase< OtherDerived > &other)
Derived & operator= (const AngleAxisType &aa)
template<class OtherDerived >
Derived & operator= (const MatrixBase< OtherDerived > &m)
QuaternionBasesetIdentity ()
Scalar squaredNorm () const
Scalar norm () const
void normalize ()
Quaternion< Scalarnormalized () const
template<class OtherDerived >
Scalar dot (const QuaternionBase< OtherDerived > &other) const
template<class OtherDerived >
Scalar angularDistance (const QuaternionBase< OtherDerived > &other) const
Matrix3 toRotationMatrix () const
template<typename Derived1 , typename Derived2 >
Derived & setFromTwoVectors (const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)
template<class OtherDerived >
EIGEN_STRONG_INLINE Quaternion
< Scalar
operator* (const QuaternionBase< OtherDerived > &q) const
template<class OtherDerived >
EIGEN_STRONG_INLINE Derived & operator*= (const QuaternionBase< OtherDerived > &q)
Quaternion< Scalarinverse () const
Quaternion< Scalarconjugate () const
template<class OtherDerived >
Quaternion< Scalarslerp (const Scalar &t, const QuaternionBase< OtherDerived > &other) const
template<class OtherDerived >
bool isApprox (const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
EIGEN_STRONG_INLINE Vector3 _transformVector (const Vector3 &v) const
template<typename NewScalarType >
internal::cast_return_type
< Derived, Quaternion
< NewScalarType > >::type 
cast () const
template<class MatrixDerived >
Derived & operator= (const MatrixBase< MatrixDerived > &xpr)

Static Public Member Functions

static Quaternion< ScalarIdentity ()

Detailed Description

template<class Derived>
class Eigen::QuaternionBase< Derived >

Base class for quaternion expressions.

Template Parameters:
Derivedderived type (CRTP)
See also:
class Quaternion

Definition at line 35 of file Quaternion.h.


Member Typedef Documentation

template<class Derived>
typedef AngleAxis<Scalar> Eigen::QuaternionBase< Derived >::AngleAxisType

the equivalent angle-axis type

Reimplemented in Eigen::Quaternion< _Scalar, _Options >.

Definition at line 56 of file Quaternion.h.

template<class Derived>
typedef internal::traits<Derived>::Coefficients Eigen::QuaternionBase< Derived >::Coefficients
template<class Derived>
typedef Matrix<Scalar,3,3> Eigen::QuaternionBase< Derived >::Matrix3

the equivalent rotation matrix type

Definition at line 54 of file Quaternion.h.

template<class Derived>
typedef NumTraits<Scalar>::Real Eigen::QuaternionBase< Derived >::RealScalar

Definition at line 44 of file Quaternion.h.

template<class Derived>
typedef internal::traits<Derived>::Scalar Eigen::QuaternionBase< Derived >::Scalar
template<class Derived>
typedef Matrix<Scalar,3,1> Eigen::QuaternionBase< Derived >::Vector3

the type of a 3D vector

Definition at line 52 of file Quaternion.h.


Member Enumeration Documentation

template<class Derived>
anonymous enum
Enumerator:
Flags 

Definition at line 46 of file Quaternion.h.

       {
    Flags = Eigen::internal::traits<Derived>::Flags
  };

Member Function Documentation

template<class Derived >
EIGEN_STRONG_INLINE QuaternionBase< Derived >::Vector3 Eigen::QuaternionBase< Derived >::_transformVector ( const Vector3 v) const

return the result vector of v through the rotation

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

Definition at line 469 of file Quaternion.h.

{
    // 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 literature (30 versus 39 flops). It also requires two
    // Vector3 as temporaries.
    Vector3 uv = this->vec().cross(v);
    uv += uv;
    return v + this->w() * uv + this->vec().cross(uv);
}
template<class Derived >
template<class OtherDerived >
internal::traits< Derived >::Scalar Eigen::QuaternionBase< Derived >::angularDistance ( const QuaternionBase< OtherDerived > &  other) const [inline]
Returns:
the angle (in radian) between two rotations
See also:
dot()

Definition at line 685 of file Quaternion.h.

{
  using std::atan2;
  using std::abs;
  Quaternion<Scalar> d = (*this) * other.conjugate();
  return Scalar(2) * atan2( d.vec().norm(), abs(d.w()) );
}
template<class Derived>
template<typename NewScalarType >
internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type Eigen::QuaternionBase< Derived >::cast ( ) const [inline]
Returns:
*this with scalar type casted to NewScalarType

Note that if NewScalarType is equal to the current scalar type of *this then this function smartly returns a const reference to *this.

Definition at line 173 of file Quaternion.h.

  {
    return typename internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type(derived());
  }
template<class Derived>
const internal::traits<Derived>::Coefficients& Eigen::QuaternionBase< Derived >::coeffs ( ) const [inline]
Returns:
a read-only vector expression of the coefficients (x,y,z,w)

Reimplemented in Eigen::Map< Quaternion< _Scalar >, _Options >, Eigen::Map< const Quaternion< _Scalar >, _Options >, and Eigen::Quaternion< _Scalar, _Options >.

Definition at line 85 of file Quaternion.h.

{ return derived().coeffs(); }
template<class Derived>
internal::traits<Derived>::Coefficients& Eigen::QuaternionBase< Derived >::coeffs ( ) [inline]
Returns:
a vector expression of the coefficients (x,y,z,w)

Reimplemented in Eigen::Map< Quaternion< _Scalar >, _Options >, and Eigen::Quaternion< _Scalar, _Options >.

Definition at line 88 of file Quaternion.h.

{ return derived().coeffs(); }
template<class Derived >
Quaternion< typename internal::traits< Derived >::Scalar > Eigen::QuaternionBase< Derived >::conjugate ( ) const [inline]
Returns:
the conjugated quaternion
the conjugate of the *this which is equal to the multiplicative inverse if the quaternion is normalized. The conjugate of a quaternion represents the opposite rotation.
See also:
Quaternion2::inverse()

Definition at line 671 of file Quaternion.h.

{
  return internal::quat_conj<Architecture::Target, Derived,
                         typename internal::traits<Derived>::Scalar,
                         internal::traits<Derived>::Alignment>::run(*this);
                         
}
template<class Derived>
template<class OtherDerived >
Scalar Eigen::QuaternionBase< Derived >::dot ( const QuaternionBase< OtherDerived > &  other) const [inline]
Returns:
the dot product of *this and other Geometrically speaking, the dot product of two unit quaternions corresponds to the cosine of half the angle between the two rotations.
See also:
angularDistance()

Definition at line 134 of file Quaternion.h.

{ return coeffs().dot(other.coeffs()); }
template<class Derived>
static Quaternion<Scalar> Eigen::QuaternionBase< Derived >::Identity ( ) [inline, static]
Returns:
a quaternion representing an identity rotation
See also:
MatrixBase::Identity()

Definition at line 106 of file Quaternion.h.

{ return Quaternion<Scalar>(Scalar(1), Scalar(0), Scalar(0), Scalar(0)); }
template<class Derived >
Quaternion< typename internal::traits< Derived >::Scalar > Eigen::QuaternionBase< Derived >::inverse ( ) const [inline]
Returns:
the quaternion describing the inverse rotation
the multiplicative inverse of *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.
See also:
QuaternionBase::conjugate()

Reimplemented from Eigen::RotationBase< Derived, 3 >.

Definition at line 640 of file Quaternion.h.

{
  // FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite()  ??
  Scalar n2 = this->squaredNorm();
  if (n2 > Scalar(0))
    return Quaternion<Scalar>(conjugate().coeffs() / n2);
  else
  {
    // return an invalid result to flag the error
    return Quaternion<Scalar>(Coefficients::Zero());
  }
}
template<class Derived>
template<class OtherDerived >
bool Eigen::QuaternionBase< Derived >::isApprox ( const QuaternionBase< OtherDerived > &  other,
const RealScalar prec = NumTraits<Scalar>::dummy_precision() 
) const [inline]
Returns:
true if *this is approximately equal to other, within the precision determined by prec.
See also:
MatrixBase::isApprox()

Definition at line 161 of file Quaternion.h.

  { return coeffs().isApprox(other.coeffs(), prec); }
template<class Derived>
Scalar Eigen::QuaternionBase< Derived >::norm ( ) const [inline]
Returns:
the norm of the quaternion's coefficients
See also:
QuaternionBase::squaredNorm(), MatrixBase::norm()

Definition at line 120 of file Quaternion.h.

{ return coeffs().norm(); }
template<class Derived>
void Eigen::QuaternionBase< Derived >::normalize ( void  ) [inline]

Normalizes the quaternion *this

See also:
normalized(), MatrixBase::normalize()

Definition at line 124 of file Quaternion.h.

{ coeffs().normalize(); }
template<class Derived>
Quaternion<Scalar> Eigen::QuaternionBase< Derived >::normalized ( ) const [inline]
Returns:
a normalized copy of *this
See also:
normalize(), MatrixBase::normalized()

Definition at line 127 of file Quaternion.h.

{ return Quaternion<Scalar>(coeffs().normalized()); }
template<class Derived >
template<class OtherDerived >
EIGEN_STRONG_INLINE Quaternion< typename internal::traits< Derived >::Scalar > Eigen::QuaternionBase< Derived >::operator* ( const QuaternionBase< OtherDerived > &  other) const
Returns:
the concatenation of two rotations as a quaternion-quaternion product

Definition at line 442 of file Quaternion.h.

{
  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename OtherDerived::Scalar>::value),
   YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
  return internal::quat_product<Architecture::Target, Derived, OtherDerived,
                         typename internal::traits<Derived>::Scalar,
                         EIGEN_PLAIN_ENUM_MIN(internal::traits<Derived>::Alignment, internal::traits<OtherDerived>::Alignment)>::run(*this, other);
}
template<class Derived >
template<class OtherDerived >
EIGEN_STRONG_INLINE Derived & Eigen::QuaternionBase< Derived >::operator*= ( const QuaternionBase< OtherDerived > &  other)
See also:
operator*(Quaternion)

Definition at line 454 of file Quaternion.h.

{
  derived() = derived() * other.derived();
  return derived();
}
template<class Derived>
EIGEN_STRONG_INLINE QuaternionBase< Derived > & Eigen::QuaternionBase< Derived >::operator= ( const QuaternionBase< Derived > &  other)

Definition at line 482 of file Quaternion.h.

{
  coeffs() = other.coeffs();
  return derived();
}
template<class Derived >
template<class OtherDerived >
EIGEN_STRONG_INLINE Derived & Eigen::QuaternionBase< Derived >::operator= ( const QuaternionBase< OtherDerived > &  other)

Definition at line 490 of file Quaternion.h.

{
  coeffs() = other.coeffs();
  return derived();
}
template<class Derived>
EIGEN_STRONG_INLINE Derived & Eigen::QuaternionBase< Derived >::operator= ( const AngleAxisType aa)

Set *this from an angle-axis aa and returns a reference to *this

Definition at line 499 of file Quaternion.h.

{
  using std::cos;
  using std::sin;
  Scalar ha = Scalar(0.5)*aa.angle(); // Scalar(0.5) to suppress precision loss warnings
  this->w() = cos(ha);
  this->vec() = sin(ha) * aa.axis();
  return derived();
}
template<class Derived>
template<class OtherDerived >
Derived& Eigen::QuaternionBase< Derived >::operator= ( const MatrixBase< OtherDerived > &  m)
template<class Derived>
template<class MatrixDerived >
Derived& Eigen::QuaternionBase< Derived >::operator= ( const MatrixBase< MatrixDerived > &  xpr) [inline]

Set *this from the expression xpr:

  • if xpr is a 4x1 vector, then xpr is assumed to be a quaternion
  • if xpr is a 3x3 matrix, then xpr is assumed to be rotation matrix and xpr is converted to a quaternion

Definition at line 517 of file Quaternion.h.

{
  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename MatrixDerived::Scalar>::value),
   YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
  internal::quaternionbase_assign_impl<MatrixDerived>::run(*this, xpr.derived());
  return derived();
}
template<class Derived >
template<typename Derived1 , typename Derived2 >
Derived & Eigen::QuaternionBase< Derived >::setFromTwoVectors ( const MatrixBase< Derived1 > &  a,
const MatrixBase< Derived2 > &  b 
) [inline]
Returns:
the quaternion which transform a into b through a rotation

Sets *this to be a quaternion representing a rotation between the two arbitrary vectors a and b. In other words, the built rotation represent a rotation sending the line of direction a to the line of direction b, both lines passing through the origin.

Returns:
a reference to *this.

Note that the two input vectors do not have to be normalized, and do not need to have the same norm.

Definition at line 576 of file Quaternion.h.

{
  using std::sqrt;
  Vector3 v0 = a.normalized();
  Vector3 v1 = b.normalized();
  Scalar c = v1.dot(v0);

  // if dot == -1, vectors are nearly opposites
  // => accurately 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)+NumTraits<Scalar>::dummy_precision())
  {
    c = numext::maxi(c,Scalar(-1));
    Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
    JacobiSVD<Matrix<Scalar,2,3> > svd(m, ComputeFullV);
    Vector3 axis = svd.matrixV().col(2);

    Scalar w2 = (Scalar(1)+c)*Scalar(0.5);
    this->w() = sqrt(w2);
    this->vec() = axis * sqrt(Scalar(1) - w2);
    return derived();
  }
  Vector3 axis = v0.cross(v1);
  Scalar s = sqrt((Scalar(1)+c)*Scalar(2));
  Scalar invs = Scalar(1)/s;
  this->vec() = axis * invs;
  this->w() = s * Scalar(0.5);

  return derived();
}
template<class Derived>
QuaternionBase& Eigen::QuaternionBase< Derived >::setIdentity ( ) [inline]
See also:
QuaternionBase::Identity(), MatrixBase::setIdentity()

Definition at line 110 of file Quaternion.h.

{ coeffs() << Scalar(0), Scalar(0), Scalar(0), Scalar(1); return *this; }
template<class Derived >
template<class OtherDerived >
Quaternion< typename internal::traits< Derived >::Scalar > Eigen::QuaternionBase< Derived >::slerp ( const Scalar t,
const QuaternionBase< OtherDerived > &  other 
) const
Returns:
the spherical linear interpolation between the two quaternions *this and other at the parameter t in [0;1].

This represents an interpolation for a constant motion between *this and other, see also http://en.wikipedia.org/wiki/Slerp.

Definition at line 704 of file Quaternion.h.

{
  using std::acos;
  using std::sin;
  using std::abs;
  static const Scalar one = Scalar(1) - NumTraits<Scalar>::epsilon();
  Scalar d = this->dot(other);
  Scalar absD = abs(d);

  Scalar scale0;
  Scalar scale1;

  if(absD>=one)
  {
    scale0 = Scalar(1) - t;
    scale1 = t;
  }
  else
  {
    // theta is the angle between the 2 quaternions
    Scalar theta = acos(absD);
    Scalar sinTheta = sin(theta);

    scale0 = sin( ( Scalar(1) - t ) * theta) / sinTheta;
    scale1 = sin( ( t * theta) ) / sinTheta;
  }
  if(d<Scalar(0)) scale1 = -scale1;

  return Quaternion<Scalar>(scale0 * coeffs() + scale1 * other.coeffs());
}
template<class Derived>
Scalar Eigen::QuaternionBase< Derived >::squaredNorm ( ) const [inline]
Returns:
the squared norm of the quaternion's coefficients
See also:
QuaternionBase::norm(), MatrixBase::squaredNorm()

Definition at line 115 of file Quaternion.h.

{ return coeffs().squaredNorm(); }
template<class Derived >
QuaternionBase< Derived >::Matrix3 Eigen::QuaternionBase< Derived >::toRotationMatrix ( void  ) const [inline]
Returns:
an equivalent 3x3 rotation matrix

Convert the quaternion to a 3x3 rotation matrix. The quaternion is required to be normalized, otherwise the result is undefined.

Reimplemented from Eigen::RotationBase< Derived, 3 >.

Definition at line 530 of file Quaternion.h.

{
  // 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  = Scalar(2)*this->x();
  const Scalar ty  = Scalar(2)*this->y();
  const Scalar tz  = Scalar(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) = Scalar(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) = Scalar(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) = Scalar(1)-(txx+tyy);

  return res;
}
template<class Derived>
const VectorBlock<const Coefficients,3> Eigen::QuaternionBase< Derived >::vec ( ) const [inline]
Returns:
a read-only vector expression of the imaginary part (x,y,z)

Definition at line 79 of file Quaternion.h.

{ return coeffs().template head<3>(); }
template<class Derived>
VectorBlock<Coefficients,3> Eigen::QuaternionBase< Derived >::vec ( ) [inline]
Returns:
a vector expression of the imaginary part (x,y,z)

Definition at line 82 of file Quaternion.h.

{ return coeffs().template head<3>(); }
template<class Derived>
Scalar Eigen::QuaternionBase< Derived >::w ( ) const [inline]
Returns:
the w coefficient

Definition at line 67 of file Quaternion.h.

{ return this->derived().coeffs().coeff(3); }
template<class Derived>
Scalar& Eigen::QuaternionBase< Derived >::w ( ) [inline]
Returns:
a reference to the w coefficient

Definition at line 76 of file Quaternion.h.

{ return this->derived().coeffs().coeffRef(3); }
template<class Derived>
Scalar Eigen::QuaternionBase< Derived >::x ( ) const [inline]
Returns:
the x coefficient

Definition at line 61 of file Quaternion.h.

{ return this->derived().coeffs().coeff(0); }
template<class Derived>
Scalar& Eigen::QuaternionBase< Derived >::x ( ) [inline]
Returns:
a reference to the x coefficient

Definition at line 70 of file Quaternion.h.

{ return this->derived().coeffs().coeffRef(0); }
template<class Derived>
Scalar Eigen::QuaternionBase< Derived >::y ( ) const [inline]
Returns:
the y coefficient

Definition at line 63 of file Quaternion.h.

{ return this->derived().coeffs().coeff(1); }
template<class Derived>
Scalar& Eigen::QuaternionBase< Derived >::y ( ) [inline]
Returns:
a reference to the y coefficient

Definition at line 72 of file Quaternion.h.

{ return this->derived().coeffs().coeffRef(1); }
template<class Derived>
Scalar Eigen::QuaternionBase< Derived >::z ( ) const [inline]
Returns:
the z coefficient

Definition at line 65 of file Quaternion.h.

{ return this->derived().coeffs().coeff(2); }
template<class Derived>
Scalar& Eigen::QuaternionBase< Derived >::z ( ) [inline]
Returns:
a reference to the z coefficient

Definition at line 74 of file Quaternion.h.

{ return this->derived().coeffs().coeffRef(2); }

The documentation for this class was generated from the following file:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines