Rotation Class Reference

The Rotation class is a Mat33 that guarantees that the matrix is a legitimate 3x3 array associated with the relative orientation of two right-handed, orthogonal, unit vector bases. More...

#include <Rotation.h>

Inheritance diagram for Rotation:
Mat< M, N, ELT, CS, RS >

List of all members.

Public Types

typedef UnitVec
< Mat33::RowSpacing > 
ColType
typedef UnitRow
< Mat33::ColSpacing > 
RowType

Public Member Functions

 Rotation ()
RotationsetRotationToIdentityMatrix ()
RotationsetRotationToNaN ()
 Rotation (const Rotation &R)
Rotationoperator= (const Rotation &R)
 Rotation (BodyOrSpaceType bodyOrSpace, Real angle1, const CoordinateAxis &axis1, Real angle2, const CoordinateAxis &axis2)
 Constructor for two-angle, two-axes, Body-fixed or Space-fixed rotation sequences (angles are in radians).
 Rotation (BodyOrSpaceType bodyOrSpace, Real angle1, const CoordinateAxis &axis1, Real angle2, const CoordinateAxis &axis2, Real angle3, const CoordinateAxis &axis3)
 Constructor for three-angle Body-fixed or Space-fixed rotation sequences (angles are in radians).
RotationsetRotationFromTwoAnglesTwoAxes (BodyOrSpaceType bodyOrSpace, Real angle1, const CoordinateAxis &axis1, Real angle2, const CoordinateAxis &axis2)
 Set this Rotation object to a two-angle, two-axes, Body-fixed or Space-fixed rotation sequences (angles are in radians).
RotationsetRotationFromThreeAnglesThreeAxes (BodyOrSpaceType bodyOrSpace, Real angle1, const CoordinateAxis &axis1, Real angle2, const CoordinateAxis &axis2, Real angle3, const CoordinateAxis &axis3)
 Set this Rotation object to a three-angle Body-fixed or Space-fixed rotation sequences (angles are in radians).
void setRotationToBodyFixedXY (const Vec2 &v)
 Set this Rotation to represent a rotation characterized by subsequent rotations of: +v[0] about the body frame's X axis, followed by a rotation of +v[1] about the body frame's NEW Y axis.
void setRotationToBodyFixedXYZ (const Vec3 &v)
 Set this Rotation to represent a rotation characterized by subsequent rotations of: +v[0] about the body frame's X axis, followed by a rotation of +v[1] about the body frame's NEW Y axis, followed by a rotation of +v[2] about the body frame's NEW Z axis.
 Rotation (const Quaternion &q)
 Constructor for relating a rotation matrix to a quaternion.
RotationsetRotationFromQuaternion (const Quaternion &q)
 Method for relating a rotation matrix to a quaternion.
 Rotation (const Mat33 &m, bool)
 Construct a Rotation directly from a Mat33 (we trust that m is a valid Rotation!).
 Rotation (const Mat33 &m)
 Constructs an (hopefully nearby) orthogonal rotation matrix from a generic Mat33.
RotationsetRotationFromApproximateMat33 (const Mat33 &m)
 Set this Rotation object to an (hopefully nearby) orthogonal rotation matrix from a generic Mat33.
Real convertOneAxisRotationToOneAngle (const CoordinateAxis &axis1) const
 Converts rotation matrix to a orientation angle.
Vec2 convertTwoAxesRotationToTwoAngles (BodyOrSpaceType bodyOrSpace, const CoordinateAxis &axis1, const CoordinateAxis &axis2) const
 Converts rotation matrix to two orientation angles.
Vec3 convertThreeAxesRotationToThreeAngles (BodyOrSpaceType bodyOrSpace, const CoordinateAxis &axis1, const CoordinateAxis &axis2, const CoordinateAxis &axis3) const
 Converts rotation matrix to three orientation angles.
Quaternion convertRotationToQuaternion () const
 Converts rotation matrix to a quaternion.
Vec4 convertRotationToAngleAxis () const
 Converts rotation matrix to angle axis form.
Vec2 convertRotationToBodyFixedXY () const
 Special case of convertTwoAxesRotationToTwoAngles().
Vec3 convertRotationToBodyFixedXYZ () const
 Special case of convertThreeAxesRotationToThreeAngles().
Real getMaxAbsDifferenceInRotationElements (const Rotation &R) const
bool areAllRotationElementsSameToEpsilon (const Rotation &R, Real epsilon) const
bool areAllRotationElementsSameToMachinePrecision (const Rotation &R) const
 Rotation (const InverseRotation &)
 Like copy constructor but for inverse rotation. This allows implicit conversion from InverseRotation to Rotation.
Rotationoperator= (const InverseRotation &)
 Like copy assign but for inverse rotation.
const InverseRotationinvert () const
 Convert from Rotation to InverseRotation (no cost).
InverseRotationupdInvert ()
 Convert from Rotation to InverseRotation (no cost).
const InverseRotationtranspose () const
const InverseRotationoperator~ () const
InverseRotationupdTranspose ()
InverseRotationoperator~ ()
Rotationoperator*= (const Rotation &R)
Rotationoperator/= (const Rotation &R)
Rotationoperator*= (const InverseRotation &)
Rotationoperator/= (const InverseRotation &)
const RowTyperow (int i) const
const ColTypecol (int j) const
const ColTypex () const
const ColTypey () const
const ColTypez () const
const RowTypeoperator[] (int i) const
const ColTypeoperator() (int j) const



 Rotation (Real angle, const CoordinateAxis &axis)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
 Rotation (Real angle, const CoordinateAxis::XCoordinateAxis)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
 Rotation (Real angle, const CoordinateAxis::YCoordinateAxis)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
 Rotation (Real angle, const CoordinateAxis::ZCoordinateAxis)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
RotationsetRotationFromAngleAboutAxis (Real angle, const CoordinateAxis &axis)
 Set this Rotation object to a right-handed rotation of an angle (in radians) about a coordinate axis.
RotationsetRotationFromAngleAboutX (Real angle)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
RotationsetRotationFromAngleAboutY (Real angle)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
RotationsetRotationFromAngleAboutZ (Real angle)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
RotationsetRotationFromAngleAboutX (Real cosAngle, Real sinAngle)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
RotationsetRotationFromAngleAboutY (Real cosAngle, Real sinAngle)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
RotationsetRotationFromAngleAboutZ (Real cosAngle, Real sinAngle)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
 Rotation (Real angle, const UnitVec3 &unitVector)
 Constructor for right-handed rotation of an angle (in radians) about an arbitrary vector.
 Rotation (Real angle, const Vec3 &nonUnitVector)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
RotationsetRotationFromAngleAboutNonUnitVector (Real angle, const Vec3 &nonUnitVector)
 Set this Rotation object to a right-handed rotation of an angle (in radians) about an arbitrary vector.
RotationsetRotationFromAngleAboutUnitVector (Real angle, const UnitVec3 &unitVector)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
 Rotation (const UnitVec3 &uvec, const CoordinateAxis axis)
 Calculate A.RotationMatrix.B by knowing one of B's unit vector expressed in A.
RotationsetRotationFromOneAxis (const UnitVec3 &uvec, const CoordinateAxis axis)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
 Rotation (const UnitVec3 &uveci, const CoordinateAxis &axisi, const Vec3 &vecjApprox, const CoordinateAxis &axisjApprox)
 Calculate A.RotationMatrix.B by knowing one of B's unit vectors expressed in A and another vector that may be perpendicular If the 2nd vector is not perpendicular, no worries - we'll make it so it is perpendicular.
RotationsetRotationFromTwoAxes (const UnitVec3 &uveci, const CoordinateAxis &axisi, const Vec3 &vecjApprox, const CoordinateAxis &axisjApprox)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
bool isSameRotationToWithinAngle (const Rotation &R, Real okPointingAngleErrorRads) const
 Return true if "this" Rotation is nearly identical to "R" within a specified pointing angle error.
bool isSameRotationToWithinAngleOfMachinePrecision (const Rotation &R) const
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
const Mat33asMat33 () const
 Conversion from Rotation to Mat33.
Mat33 toMat33 () const
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
RotationsetRotationFromMat33TrustMe (const Mat33 &m)
 Set the Rotation matrix directly - but you had better know what you are doing!
RotationsetRotationColFromUnitVecTrustMe (int coli, const UnitVec3 &uveci)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.
RotationsetRotationFromUnitVecsTrustMe (const UnitVec3 &colA, const UnitVec3 &colB, const UnitVec3 &colC)
 Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

Static Public Member Functions

static Vec3 convertAngVelToBodyFixed321Dot (const Vec3 &q, const Vec3 &w_PB_B)
 Given Euler angles forming a body-fixed 3-2-1 sequence, and the relative angular velocity vector of B in the parent frame, *BUT EXPRESSED IN THE BODY FRAME*, return the Euler angle derivatives.
static Vec3 convertBodyFixed321DotToAngVel (const Vec3 &q, const Vec3 &qd)
 Inverse of the above routine.
static Vec3 convertAngVelDotToBodyFixed321DotDot (const Vec3 &q, const Vec3 &w_PB_B, const Vec3 &wdot)
static Mat33 calcQBlockForBodyXYZInBodyFrame (const Vec3 &q)
 Given Euler angles forming a body-fixed X-Y-Z sequence q return the block of the Q matrix such that qdot=Q(q)*w where w is the angular velocity of B in P EXPRESSED IN *B*!!! This matrix will be singular if Y (q[1]) gets near 90 degrees! See Kane's Spacecraft Dynamics, page 427, body-three: 1-2-3.
static Mat33 calcQInvBlockForBodyXYZInBodyFrame (const Vec3 &q)
 Inverse of the above routine.
static Vec3 convertAngVelToBodyFixed123Dot (const Vec3 &q, const Vec3 &w_PB_B)
 Given Euler angles forming a body-fixed 1-2-3 sequence, and the relative angular velocity vector of B in the parent frame, *BUT EXPRESSED IN THE BODY FRAME*, return the Euler angle derivatives.
static Vec3 convertBodyFixed123DotToAngVel (const Vec3 &q, const Vec3 &qd)
 Inverse of the above routine.
static Vec3 convertAngVelDotToBodyFixed123DotDot (const Vec3 &q, const Vec3 &w_PB_B, const Vec3 &wdot)
static Mat43 calcUnnormalizedQBlockForQuaternion (const Vec4 &q)
 Given a possibly unnormalized quaternion q, calculate the 4x3 matrix Q which maps angular velocity w to quaternion derivatives qdot.
static Mat34 calcUnnormalizedQInvBlockForQuaternion (const Vec4 &q)
 Given a (possibly unnormalized) quaternion q, calculate the 3x4 matrix QInv (= Q^-1) which maps quaternion derivatives qdot to angular velocity w, where the angular velocity is in the parent frame, i.e.
static Vec4 convertAngVelToQuaternionDot (const Vec4 &q, const Vec3 &w_PB_P)
 Given a possibly unnormalized quaternion (0th element is the scalar) and the relative angular velocity vector of B in its parent, expressed in the *PARENT*, return the quaternion derivatives.
static Vec3 convertQuaternionDotToAngVel (const Vec4 &q, const Vec4 &qd)
 Inverse of the above routine.
static Vec4 convertAngVelDotToQuaternionDotDot (const Vec4 &q, const Vec3 &w_PB_P, const Vec3 &wdot)
 Everything is measured and expressed in the parent.

Detailed Description

The Rotation class is a Mat33 that guarantees that the matrix is a legitimate 3x3 array associated with the relative orientation of two right-handed, orthogonal, unit vector bases.

The Rotation class takes advantage of knowledge-specific information and properties of orthogonal matrices. For example, multiplication by a rotation matrix preserves a vector's length so unit vectors are still unit vectors afterwards and don't need to be re-normalized.

A rotation is an orthogonal matrix whose columns and rows are directions (that is, unit vectors) which are mutually orthogonal. Furthermore, if the columns (or rows) are labeled x,y,z it always holds that z = x X y (rather than -(x X y)) ensuring that this is a right-handed rotation matrix and not a reflection.

Suppose there is a vector vF expressed in terms of the right-handed, orthogonal unit vectors Fx, Fy, Fz and one would like to express vF in terms of the right-handed, orthogonal unit vectors Gx, Gy, Gz. To calculate it, form vG = G_rotationMatrix_F * vF. Because a rotation is orthogonal, its transpose is its inverse. Hence F_rotationMatrix_G = ~G_rotationMatrix_F (where ~ is the SimTK "transpose"). This transpose matrix can be used to expressed vG in terms of Fx, Fy, Fz as vF = F_rotationMatrix_G * vG or vF = ~(G_rotationMatrix_F) * vG.


Member Typedef Documentation

typedef UnitVec<Mat33::RowSpacing> ColType
typedef UnitRow<Mat33::ColSpacing> RowType

Constructor & Destructor Documentation

Rotation (  )  [inline]
Rotation ( const Rotation R  )  [inline]
Rotation ( Real  angle,
const CoordinateAxis axis 
) [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

References Rotation::setRotationFromAngleAboutAxis().

Rotation ( Real  angle,
const CoordinateAxis::XCoordinateAxis   
) [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

References Rotation::setRotationFromAngleAboutX().

Rotation ( Real  angle,
const CoordinateAxis::YCoordinateAxis   
) [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

References Rotation::setRotationFromAngleAboutY().

Rotation ( Real  angle,
const CoordinateAxis::ZCoordinateAxis   
) [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

References Rotation::setRotationFromAngleAboutZ().

Rotation ( Real  angle,
const UnitVec3 unitVector 
) [inline]

Constructor for right-handed rotation of an angle (in radians) about an arbitrary vector.

References Rotation::setRotationFromAngleAboutUnitVector().

Rotation ( Real  angle,
const Vec3 nonUnitVector 
) [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

References Rotation::setRotationFromAngleAboutNonUnitVector().

Rotation ( BodyOrSpaceType  bodyOrSpace,
Real  angle1,
const CoordinateAxis axis1,
Real  angle2,
const CoordinateAxis axis2 
) [inline]

Constructor for two-angle, two-axes, Body-fixed or Space-fixed rotation sequences (angles are in radians).

References Rotation::setRotationFromTwoAnglesTwoAxes().

Rotation ( BodyOrSpaceType  bodyOrSpace,
Real  angle1,
const CoordinateAxis axis1,
Real  angle2,
const CoordinateAxis axis2,
Real  angle3,
const CoordinateAxis axis3 
) [inline]

Constructor for three-angle Body-fixed or Space-fixed rotation sequences (angles are in radians).

References Rotation::setRotationFromThreeAnglesThreeAxes().

Rotation ( const Quaternion q  )  [inline, explicit]

Constructor for relating a rotation matrix to a quaternion.

References Rotation::setRotationFromQuaternion().

Rotation ( const Mat33 m,
bool   
) [inline]

Construct a Rotation directly from a Mat33 (we trust that m is a valid Rotation!).

Rotation ( const Mat33 m  )  [inline, explicit]

Constructs an (hopefully nearby) orthogonal rotation matrix from a generic Mat33.

References Rotation::setRotationFromApproximateMat33().

Rotation ( const UnitVec3 uvec,
const CoordinateAxis  axis 
) [inline]

Calculate A.RotationMatrix.B by knowing one of B's unit vector expressed in A.

Note: The other vectors are perpendicular (but somewhat arbitrarily so).

References Rotation::setRotationFromOneAxis().

Rotation ( const UnitVec3 uveci,
const CoordinateAxis axisi,
const Vec3 vecjApprox,
const CoordinateAxis axisjApprox 
) [inline]

Calculate A.RotationMatrix.B by knowing one of B's unit vectors expressed in A and another vector that may be perpendicular If the 2nd vector is not perpendicular, no worries - we'll make it so it is perpendicular.

References Rotation::setRotationFromTwoAxes().

Rotation ( const InverseRotation R  )  [inline]

Like copy constructor but for inverse rotation. This allows implicit conversion from InverseRotation to Rotation.


Member Function Documentation

bool areAllRotationElementsSameToEpsilon ( const Rotation R,
Real  epsilon 
) const [inline]
bool areAllRotationElementsSameToMachinePrecision ( const Rotation R  )  const [inline]
const Mat33& asMat33 (  )  const [inline]

Conversion from Rotation to Mat33.

Note: asMat33 is more efficient than toMat33() (no copy), but you have to know the internal layout.

Referenced by Rotation::col(), Rotation::getMaxAbsDifferenceInRotationElements(), SimTK::operator*(), Rotation::operator*=(), Rotation::operator/=(), Rotation::operator=(), Rotation::row(), and Rotation::toMat33().

static Mat33 calcQBlockForBodyXYZInBodyFrame ( const Vec3 q  )  [inline, static]

Given Euler angles forming a body-fixed X-Y-Z sequence q return the block of the Q matrix such that qdot=Q(q)*w where w is the angular velocity of B in P EXPRESSED IN *B*!!! This matrix will be singular if Y (q[1]) gets near 90 degrees! See Kane's Spacecraft Dynamics, page 427, body-three: 1-2-3.

Referenced by Rotation::convertAngVelToBodyFixed123Dot().

static Mat33 calcQInvBlockForBodyXYZInBodyFrame ( const Vec3 q  )  [inline, static]

Inverse of the above routine.

Return the inverse QInv of the Q block computed above, such that w=QInv(q)*qdot where w is the angular velocity of B in P EXPRESSED IN *B*!!! This matrix is never singular.

Referenced by Rotation::convertBodyFixed123DotToAngVel().

static Mat43 calcUnnormalizedQBlockForQuaternion ( const Vec4 q  )  [inline, static]

Given a possibly unnormalized quaternion q, calculate the 4x3 matrix Q which maps angular velocity w to quaternion derivatives qdot.

We expect the angular velocity in the parent frame, i.e. w==w_PB_P. We don't normalize, so Q=|q|Q' where Q' is the normalized version.

Referenced by Rotation::convertAngVelDotToQuaternionDotDot(), and Rotation::convertAngVelToQuaternionDot().

static Mat34 calcUnnormalizedQInvBlockForQuaternion ( const Vec4 q  )  [inline, static]

Given a (possibly unnormalized) quaternion q, calculate the 3x4 matrix QInv (= Q^-1) which maps quaternion derivatives qdot to angular velocity w, where the angular velocity is in the parent frame, i.e.

w==w_PB_P. Note: when the quaternion is not normalized, this is not precisely the (pseudo)inverse of Q. inv(Q)=inv(Q')/|q| but we're returning |q|*inv(Q')=|q|^2*inv(Q). That is, QInv*Q =|q|^2*I, which is I if the original q was normalized. (Note: Q*QInv != I, not even close.)

Referenced by Rotation::convertQuaternionDotToAngVel().

const ColType& col ( int  j  )  const [inline]
static Vec3 convertAngVelDotToBodyFixed123DotDot ( const Vec3 q,
const Vec3 w_PB_B,
const Vec3 wdot 
) [inline, static]
static Vec3 convertAngVelDotToBodyFixed321DotDot ( const Vec3 q,
const Vec3 w_PB_B,
const Vec3 wdot 
) [inline, static]
static Vec4 convertAngVelDotToQuaternionDotDot ( const Vec4 q,
const Vec3 w_PB_P,
const Vec3 wdot 
) [inline, static]

Everything is measured and expressed in the parent.

References Rotation::calcUnnormalizedQBlockForQuaternion().

static Vec3 convertAngVelToBodyFixed123Dot ( const Vec3 q,
const Vec3 w_PB_B 
) [inline, static]

Given Euler angles forming a body-fixed 1-2-3 sequence, and the relative angular velocity vector of B in the parent frame, *BUT EXPRESSED IN THE BODY FRAME*, return the Euler angle derivatives.

You are dead if q[1] gets near 90 degrees! See Kane's Spacecraft Dynamics, page 427, body-three: 1-2-3.

References Rotation::calcQBlockForBodyXYZInBodyFrame().

static Vec3 convertAngVelToBodyFixed321Dot ( const Vec3 q,
const Vec3 w_PB_B 
) [inline, static]

Given Euler angles forming a body-fixed 3-2-1 sequence, and the relative angular velocity vector of B in the parent frame, *BUT EXPRESSED IN THE BODY FRAME*, return the Euler angle derivatives.

You are dead if q[1] gets near 90 degrees! See Kane's Spacecraft Dynamics, page 428, body-three: 3-2-1.

static Vec4 convertAngVelToQuaternionDot ( const Vec4 q,
const Vec3 w_PB_P 
) [inline, static]

Given a possibly unnormalized quaternion (0th element is the scalar) and the relative angular velocity vector of B in its parent, expressed in the *PARENT*, return the quaternion derivatives.

This is never singular.

References Rotation::calcUnnormalizedQBlockForQuaternion().

static Vec3 convertBodyFixed123DotToAngVel ( const Vec3 q,
const Vec3 qd 
) [inline, static]

Inverse of the above routine.

Returned angular velocity is B in P, expressed in *B*: w_PB_B.

References Rotation::calcQInvBlockForBodyXYZInBodyFrame().

static Vec3 convertBodyFixed321DotToAngVel ( const Vec3 q,
const Vec3 qd 
) [inline, static]

Inverse of the above routine.

Returned angular velocity is B in P, expressed in *B*: w_PB_B.

Real convertOneAxisRotationToOneAngle ( const CoordinateAxis axis1  )  const

Converts rotation matrix to a orientation angle.

Note: The result is most meaningful if the Rotation matrix is one that can be produced by such a sequence.

static Vec3 convertQuaternionDotToAngVel ( const Vec4 q,
const Vec4 qd 
) [inline, static]

Inverse of the above routine.

Returned AngVel is expressed in the *PARENT* frame: w_PB_P.

References Rotation::calcUnnormalizedQInvBlockForQuaternion().

Vec4 convertRotationToAngleAxis (  )  const [inline]

Converts rotation matrix to angle axis form.

References Quaternion::convertQuaternionToAngleAxis(), and Rotation::convertRotationToQuaternion().

Vec2 convertRotationToBodyFixedXY (  )  const [inline]
Vec3 convertRotationToBodyFixedXYZ (  )  const [inline]
Quaternion convertRotationToQuaternion (  )  const
Vec3 convertThreeAxesRotationToThreeAngles ( BodyOrSpaceType  bodyOrSpace,
const CoordinateAxis axis1,
const CoordinateAxis axis2,
const CoordinateAxis axis3 
) const

Converts rotation matrix to three orientation angles.

Note: The result is most meaningful if the Rotation matrix is one that can be produced by such a sequence.

Referenced by Rotation::convertRotationToBodyFixedXYZ().

Vec2 convertTwoAxesRotationToTwoAngles ( BodyOrSpaceType  bodyOrSpace,
const CoordinateAxis axis1,
const CoordinateAxis axis2 
) const

Converts rotation matrix to two orientation angles.

Note: The result is most meaningful if the Rotation matrix is one that can be produced by such a sequence.

Referenced by Rotation::convertRotationToBodyFixedXY().

Real getMaxAbsDifferenceInRotationElements ( const Rotation R  )  const [inline]
const InverseRotation& invert (  )  const [inline]

Convert from Rotation to InverseRotation (no cost).

Reimplemented from Mat< M, N, ELT, CS, RS >.

Referenced by Rotation::operator~(), and Rotation::transpose().

bool isSameRotationToWithinAngle ( const Rotation R,
Real  okPointingAngleErrorRads 
) const

Return true if "this" Rotation is nearly identical to "R" within a specified pointing angle error.

Referenced by Rotation::isSameRotationToWithinAngleOfMachinePrecision(), and Test::numericallyEqual().

bool isSameRotationToWithinAngleOfMachinePrecision ( const Rotation R  )  const [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

References Rotation::isSameRotationToWithinAngle().

const ColType& operator() ( int  j  )  const [inline]

Reimplemented from Mat< M, N, ELT, CS, RS >.

References Rotation::col().

Rotation & operator*= ( const InverseRotation R  )  [inline]
Rotation & operator*= ( const Rotation R  )  [inline]

Reimplemented from Mat< M, N, ELT, CS, RS >.

References Rotation::asMat33().

Rotation & operator/= ( const InverseRotation R  )  [inline]

References Rotation::asMat33().

Rotation & operator/= ( const Rotation R  )  [inline]

References Rotation::asMat33().

Rotation & operator= ( const InverseRotation R  )  [inline]

Like copy assign but for inverse rotation.

References InverseRotation::asMat33().

Rotation& operator= ( const Rotation R  )  [inline]
const RowType& operator[] ( int  i  )  const [inline]

Reimplemented from Mat< M, N, ELT, CS, RS >.

References Rotation::row().

InverseRotation& operator~ (  )  [inline]

Reimplemented from Mat< M, N, ELT, CS, RS >.

References Rotation::updInvert().

const InverseRotation& operator~ (  )  const [inline]

Reimplemented from Mat< M, N, ELT, CS, RS >.

References Rotation::invert().

const RowType& row ( int  i  )  const [inline]

Reimplemented from Mat< M, N, ELT, CS, RS >.

References Rotation::asMat33().

Referenced by Rotation::operator[]().

Rotation& setRotationColFromUnitVecTrustMe ( int  coli,
const UnitVec3 uveci 
) [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

Referenced by Rotation::setRotationFromUnitVecsTrustMe().

Rotation& setRotationFromAngleAboutAxis ( Real  angle,
const CoordinateAxis axis 
) [inline]

Set this Rotation object to a right-handed rotation of an angle (in radians) about a coordinate axis.

References CoordinateAxis::isXAxis(), CoordinateAxis::isYAxis(), Rotation::setRotationFromAngleAboutX(), Rotation::setRotationFromAngleAboutY(), and Rotation::setRotationFromAngleAboutZ().

Referenced by Rotation::Rotation().

Rotation& setRotationFromAngleAboutNonUnitVector ( Real  angle,
const Vec3 nonUnitVector 
) [inline]

Set this Rotation object to a right-handed rotation of an angle (in radians) about an arbitrary vector.

References Rotation::setRotationFromAngleAboutUnitVector().

Referenced by Rotation::Rotation().

Rotation& setRotationFromAngleAboutUnitVector ( Real  angle,
const UnitVec3 unitVector 
)

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

Referenced by Rotation::Rotation(), and Rotation::setRotationFromAngleAboutNonUnitVector().

Rotation& setRotationFromAngleAboutX ( Real  cosAngle,
Real  sinAngle 
) [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

Rotation& setRotationFromAngleAboutX ( Real  angle  )  [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

References Rotation::setRotationFromAngleAboutX().

Referenced by Rotation::Rotation(), Rotation::setRotationFromAngleAboutAxis(), and Rotation::setRotationFromAngleAboutX().

Rotation& setRotationFromAngleAboutY ( Real  cosAngle,
Real  sinAngle 
) [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

Rotation& setRotationFromAngleAboutY ( Real  angle  )  [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

References Rotation::setRotationFromAngleAboutY().

Referenced by Rotation::Rotation(), Rotation::setRotationFromAngleAboutAxis(), and Rotation::setRotationFromAngleAboutY().

Rotation& setRotationFromAngleAboutZ ( Real  cosAngle,
Real  sinAngle 
) [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

Rotation& setRotationFromAngleAboutZ ( Real  angle  )  [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

References Rotation::setRotationFromAngleAboutZ().

Referenced by Rotation::Rotation(), Rotation::setRotationFromAngleAboutAxis(), and Rotation::setRotationFromAngleAboutZ().

Rotation& setRotationFromApproximateMat33 ( const Mat33 m  ) 

Set this Rotation object to an (hopefully nearby) orthogonal rotation matrix from a generic Mat33.

Referenced by Rotation::Rotation().

Rotation& setRotationFromMat33TrustMe ( const Mat33 m  )  [inline]

Set the Rotation matrix directly - but you had better know what you are doing!

Rotation& setRotationFromOneAxis ( const UnitVec3 uvec,
const CoordinateAxis  axis 
)

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

Referenced by Rotation::Rotation().

Rotation& setRotationFromQuaternion ( const Quaternion q  ) 

Method for relating a rotation matrix to a quaternion.

Referenced by Rotation::Rotation().

Rotation& setRotationFromThreeAnglesThreeAxes ( BodyOrSpaceType  bodyOrSpace,
Real  angle1,
const CoordinateAxis axis1,
Real  angle2,
const CoordinateAxis axis2,
Real  angle3,
const CoordinateAxis axis3 
)

Set this Rotation object to a three-angle Body-fixed or Space-fixed rotation sequences (angles are in radians).

Referenced by Rotation::Rotation(), and Rotation::setRotationToBodyFixedXYZ().

Rotation& setRotationFromTwoAnglesTwoAxes ( BodyOrSpaceType  bodyOrSpace,
Real  angle1,
const CoordinateAxis axis1,
Real  angle2,
const CoordinateAxis axis2 
)

Set this Rotation object to a two-angle, two-axes, Body-fixed or Space-fixed rotation sequences (angles are in radians).

Referenced by Rotation::Rotation(), and Rotation::setRotationToBodyFixedXY().

Rotation& setRotationFromTwoAxes ( const UnitVec3 uveci,
const CoordinateAxis axisi,
const Vec3 vecjApprox,
const CoordinateAxis axisjApprox 
)

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

Referenced by Rotation::Rotation().

Rotation& setRotationFromUnitVecsTrustMe ( const UnitVec3 colA,
const UnitVec3 colB,
const UnitVec3 colC 
) [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

References Rotation::setRotationColFromUnitVecTrustMe().

void setRotationToBodyFixedXY ( const Vec2 v  )  [inline]

Set this Rotation to represent a rotation characterized by subsequent rotations of: +v[0] about the body frame's X axis, followed by a rotation of +v[1] about the body frame's NEW Y axis.

See Kane, Spacecraft Dynamics, pg. 423, body-three: 1-2-3.

References SimTK::BodyRotationSequence, and Rotation::setRotationFromTwoAnglesTwoAxes().

void setRotationToBodyFixedXYZ ( const Vec3 v  )  [inline]

Set this Rotation to represent a rotation characterized by subsequent rotations of: +v[0] about the body frame's X axis, followed by a rotation of +v[1] about the body frame's NEW Y axis, followed by a rotation of +v[2] about the body frame's NEW Z axis.

See Kane, Spacecraft Dynamics, pg. 423, body-three: 1-2-3.

References SimTK::BodyRotationSequence, and Rotation::setRotationFromThreeAnglesThreeAxes().

Rotation& setRotationToIdentityMatrix (  )  [inline]

References Rotation::operator=().

Referenced by Transform::setToZero().

Rotation& setRotationToNaN (  )  [inline]

References Mat< 3, 3 >::setToNaN().

Referenced by Transform::setToNaN().

Mat33 toMat33 (  )  const [inline]

Constructor for right-handed rotation of an angle (in radians) about a coordinate axis.

References Rotation::asMat33().

const InverseRotation& transpose (  )  const [inline]

Reimplemented from Mat< M, N, ELT, CS, RS >.

References Rotation::invert().

InverseRotation& updInvert (  )  [inline]

Convert from Rotation to InverseRotation (no cost).

Referenced by Rotation::operator~(), and Rotation::updTranspose().

InverseRotation& updTranspose (  )  [inline]

Reimplemented from Mat< M, N, ELT, CS, RS >.

References Rotation::updInvert().

const ColType& x (  )  const [inline]

References Rotation::col().

Referenced by Transform::x().

const ColType& y (  )  const [inline]

References Rotation::col().

Referenced by Transform::y().

const ColType& z (  )  const [inline]

References Rotation::col().

Referenced by Transform::z().


The documentation for this class was generated from the following file:

Generated on Wed Dec 30 11:05:28 2009 for SimTKcore by  doxygen 1.6.1