1 #ifndef SimTK_SimTKCOMMON_ROTATION_H_
2 #define SimTK_SimTKCOMMON_ROTATION_H_
46 template <
class P>
class Rotation_;
140 Rotation_&
setRotationFromAngleAboutX( RealP cosAngle, RealP sinAngle ) {
Mat33P& R = *
this; R[0][0] = 1; R[0][1] = R[0][2] = R[1][0] = R[2][0] = 0; R[1][1] = R[2][2] = cosAngle; R[1][2] = -(R[2][1] = sinAngle);
return *
this; }
141 Rotation_&
setRotationFromAngleAboutY( RealP cosAngle, RealP sinAngle ) {
Mat33P& R = *
this; R[1][1] = 1; R[0][1] = R[1][0] = R[1][2] = R[2][1] = 0; R[0][0] = R[2][2] = cosAngle; R[2][0] = -(R[0][2] = sinAngle);
return *
this; }
142 Rotation_&
setRotationFromAngleAboutZ( RealP cosAngle, RealP sinAngle ) {
Mat33P& R = *
this; R[2][2] = 1; R[0][2] = R[1][2] = R[2][0] = R[2][1] = 0; R[0][0] = R[1][1] = cosAngle; R[0][1] = -(R[1][0] = sinAngle);
return *
this; }
159 Rotation_(
BodyOrSpaceType bodyOrSpace, RealP angle1,
const CoordinateAxis& axis1, RealP angle2,
const CoordinateAxis& axis2, RealP angle3,
const CoordinateAxis& axis3 ) {
setRotationFromThreeAnglesThreeAxes(bodyOrSpace,angle1,axis1,angle2,axis2,angle3,axis3); }
249 for(
int i=0; i<=2; i++ )
for(
int j=0; j<=2; j++ )
250 { RealP absDiff = std::fabs(A[i][j] - B[i][j]);
251 if( absDiff > maxDiff ) maxDiff = absDiff; }
308 {
Mat33P& R = *
this; R=m;
return *
this; }
310 {
Mat33P& R = *
this; R(colj)=uvecj.
asVec3();
return *
this; }
328 const RealP s0s1=s[0]*s[1], s2c0=s[2]*c[0], c0c2=c[0]*c[2], nc1= -c[1];
330 R =
Mat33P( c[1]*c[2] , s[2]*nc1 , s[1] ,
331 s2c0 + s0s1*c[2] , c0c2 - s0s1*s[2] , s[0]*nc1 ,
332 s[0]*s[2] - s[1]*c0c2 , s[0]*c[2] + s[1]*s2c0 , c[0]*c[1] );
342 const RealP s1 = std::sin(q[1]), c1 = std::cos(q[1]);
343 const RealP s2 = std::sin(q[2]), c2 = std::cos(q[2]);
344 const RealP ooc1 = RealP(1)/c1;
345 const RealP s2oc1 = s2*ooc1, c2oc1 = c2*ooc1;
347 const Mat33P E( 0, s2oc1 , c2oc1 ,
349 1, s1*s2oc1 , s1*c2oc1 );
356 const RealP s1 = std::sin(q[1]), c1 = std::cos(q[1]);
357 const RealP s2 = std::sin(q[2]), c2 = std::cos(q[2]);
359 const Mat33P Einv( -s1 , 0 , 1 ,
372 const RealP s1 = std::sin(q[1]), c1 = std::cos(q[1]);
373 const RealP s2 = std::sin(q[2]), c2 = std::cos(q[2]);
374 const RealP ooc1 = 1/c1;
375 const RealP s2oc1 = s2*ooc1, c2oc1 = c2*ooc1, s1oc1 = s1*ooc1;
377 const Mat33P E( 0 , s2oc1 , c2oc1 ,
379 1 , s1*s2oc1 , s1*c2oc1 );
380 const Vec3P qdot = E * w_PB_B;
382 const RealP t = qdot[1]*s1oc1;
383 const RealP a = t*s2oc1 + qdot[2]*c2oc1;
384 const RealP b = t*c2oc1 - qdot[2]*s2oc1;
386 const Mat33P Edot( 0 , a , b ,
387 0 , -qdot[2]*s2 , -qdot[2]*c2 ,
388 0 , s1*a + qdot[1]*s2 , s1*b + qdot[1]*c2 );
390 return E*wdot_PB_B + Edot*w_PB_B;
408 (
Vec3P(0, std::cos(q[1]), std::cos(q[2])),
409 Vec3P(0, std::sin(q[1]), std::sin(q[2])));
418 const RealP s1 = sq[1], c1 = cq[1];
419 const RealP s2 = sq[2], c2 = cq[2];
420 const RealP ooc1 = 1/c1;
421 const RealP s2oc1 = s2*ooc1, c2oc1 = c2*ooc1;
423 return Mat33P( c2oc1 , -s2oc1 , 0,
425 -s1*c2oc1 , s1*s2oc1, 1 );
444 (
Vec3P(std::cos(q[0]), std::cos(q[1]), 0),
445 Vec3P(std::sin(q[0]), std::sin(q[1]), 0));
455 const RealP s0 = sq[0], c0 = cq[0];
456 const RealP s1 = sq[1], c1 = cq[1];
457 const RealP ooc1 = 1/c1;
458 const RealP s0oc1 = s0*ooc1, c0oc1 = c0*ooc1;
460 return Mat33P( 1 , s1*s0oc1 , -s1*c0oc1,
462 0 , -s0oc1 , c0oc1 );
476 const RealP s0 = sinxy[0], c0 = cosxy[0];
477 const RealP s1 = sinxy[1];
478 const RealP w0 = w_PB[0], w1 = w_PB[1], w2 = w_PB[2];
480 const RealP t = (s0*w1-c0*w2)*oocosy;
481 return Vec3P( w0 + t*s1, c0*w1 + s0*w2, -t );
492 const RealP s0 = sinxy[0], c0 = cosxy[0];
493 const RealP s1 = sinxy[1];
494 const RealP q0 = q[0], q1 = q[1], q2 = q[2];
496 const RealP t = (q0*s1-q2) * oocosy;
497 return Vec3P( q0, c0*q1 + t*s0, s0*q1 - t*c0 );
533 const RealP s1 = sinxy[1], c1 = cosxy[1];
534 const RealP q0 = qdot[0], q1 = qdot[1], q2 = qdot[2];
539 const RealP q1oc1 = q1*oocosy;
540 const Vec3P NDotw((q0*s1-q2)*q1oc1,
554 const RealP s0 = sinxy[0], c0 = cosxy[0];
555 const RealP s1 = sinxy[1], c1 = cosxy[1];
556 const RealP q0 = qdot[0], q1 = qdot[1], q2 = qdot[2];
557 const RealP c1q2 = c1*q2;
559 return Vec3P( q0 + s1*q2,
570 const RealP s0 = sinxy[0], c0 = cosxy[0];
571 const RealP s1 = sinxy[1], c1 = cosxy[1];
572 const RealP w0 = v_P[0], w1 = v_P[1], w2 = v_P[2];
576 s1*w0 - s0*c1*w1 + c0*c1*w2);
593 (
Vec3P(0, std::cos(q[1]), std::cos(q[2])),
594 Vec3P(0, std::sin(q[1]), std::sin(q[2])),
606 const RealP s1 = sq[1], c1 = cq[1];
607 const RealP s2 = sq[2], c2 = cq[2];
608 const RealP ooc1 = 1/c1;
609 const RealP s2oc1 = s2*ooc1, c2oc1 = c2*ooc1;
611 const RealP t = qdot[1]*s1*ooc1;
612 const RealP a = t*s2oc1 + qdot[2]*c2oc1;
613 const RealP b = t*c2oc1 - qdot[2]*s2oc1;
615 return Mat33P( b , -a , 0,
616 qdot[2]*c2 , -qdot[2]*s2 , 0,
617 -(s1*b + qdot[1]*c2) , s1*a + qdot[1]*s2 , 0 );
633 const RealP cy = std::cos(q[1]);
635 (
Vec2P(std::cos(q[0]), cy),
636 Vec2P(std::sin(q[0]), std::sin(q[1])),
646 const RealP s0 = sq[0], c0 = cq[0];
647 const RealP s1 = sq[1], c1 = cq[1];
648 const RealP s0oc1 = s0*ooc1, c0oc1 = c0*ooc1;
650 const RealP t = qdot[1]*s1*ooc1;
651 const RealP a = t*s0oc1 + qdot[0]*c0oc1;
652 const RealP b = t*c0oc1 - qdot[0]*s0oc1;
654 return Mat33P( 0, s1*a + qdot[1]*s0, -(s1*b + qdot[1]*c0),
655 0, -qdot[0]*s0 , qdot[0]*c0 ,
671 (
Vec3P(0, std::cos(q[1]), std::cos(q[2])),
672 Vec3P(0, std::sin(q[1]), std::sin(q[2])));
682 const RealP s1 = sq[1], c1 = cq[1];
683 const RealP s2 = sq[2], c2 = cq[2];
685 return Mat33P( c1*c2 , s2 , 0 ,
701 (
Vec3P(std::cos(q[0]), std::cos(q[1]), 0),
702 Vec3P(std::sin(q[0]), std::sin(q[1]), 0));
712 const RealP s0 = sq[0], c0 = cq[0];
713 const RealP s1 = sq[1], c1 = cq[1];
715 return Mat33P( 1 , 0 , s1 ,
731 (
Vec3P(0, std::cos(q[1]), std::cos(q[2])),
732 Vec3P(0, std::sin(q[1]), std::sin(q[2])),
754 (
Vec3P(0, std::cos(q[1]), std::cos(q[2])),
755 Vec3P(0, std::sin(q[1]), std::sin(q[2])),
780 (
Vec3P(0, std::cos(q[1]), std::cos(q[2])),
781 Vec3P(0, std::sin(q[1]), std::sin(q[2])),
796 const Vec3P qdot = N * w_PB_B;
799 return N*wdot_PB_B + NDot*w_PB_B;
809 const RealP ne1 = -e[1], ne2 = -e[2], ne3 = -e[3];
810 return Mat43P( ne1, ne2, ne3,
823 const Vec4P ed = qdot/2;
824 const RealP ned1 = -ed[1], ned2 = -ed[2], ned3 = -ed[3];
825 return Mat43P( ned1, ned2, ned3,
841 const RealP ne1 = -e[1], ne2 = -e[2], ne3 = -e[3];
842 return Mat34P(ne1, e[0], ne3, e[2],
843 ne2, e[3], e[0], ne1,
844 ne3, ne2, e[1], e[0]);
873 const Vec4P qdot = N*w_PB_P;
876 return N*wdot_PB_P + NDot*w_PB_P;
889 const Vec4P Nb = N*b_PB;
900 Rotation_(
const RealP& xx,
const RealP& xy,
const RealP& xz,
901 const RealP& yx,
const RealP& yy,
const RealP& yz,
902 const RealP& zx,
const RealP& zy,
const RealP& zz )
903 : Mat33P( xx,xy,xz, yx,yy,yz, zx,zy,zz ) {}
907 SimTK_SimTKCOMMON_EXPORT Rotation_& setTwoAngleTwoAxesBodyFixedForwardCyclicalRotation( RealP cosAngle1, RealP sinAngle1,
const CoordinateAxis& axis1, RealP cosAngle2, RealP sinAngle2,
const CoordinateAxis& axis2 );
908 SimTK_SimTKCOMMON_EXPORT Rotation_& setThreeAngleTwoAxesBodyFixedForwardCyclicalRotation( RealP cosAngle1, RealP sinAngle1,
const CoordinateAxis& axis1, RealP cosAngle2, RealP sinAngle2,
const CoordinateAxis& axis2, RealP cosAngle3, RealP sinAngle3 );
909 SimTK_SimTKCOMMON_EXPORT Rotation_& setThreeAngleThreeAxesBodyFixedForwardCyclicalRotation( RealP cosAngle1, RealP sinAngle1,
const CoordinateAxis& axis1, RealP cosAngle2, RealP sinAngle2,
const CoordinateAxis& axis2, RealP cosAngle3, RealP sinAngle3,
const CoordinateAxis& axis3 );
913 SimTK_SimTKCOMMON_EXPORT Vec2P convertTwoAxesBodyFixedRotationToTwoAngles(
const CoordinateAxis& axis1,
const CoordinateAxis& axis2 )
const;
914 SimTK_SimTKCOMMON_EXPORT Vec3P convertTwoAxesBodyFixedRotationToThreeAngles(
const CoordinateAxis& axis1,
const CoordinateAxis& axis2 )
const;
915 SimTK_SimTKCOMMON_EXPORT Vec3P convertThreeAxesBodyFixedRotationToThreeAngles(
const CoordinateAxis& axis1,
const CoordinateAxis& axis2,
const CoordinateAxis& axis3 )
const;
923 static Mat33P calcQBlockForBodyXYZInBodyFrame(
const Vec3P& a)
926 static Mat33P calcQInvBlockForBodyXYZInBodyFrame(
const Vec3P& a)
929 static Mat<4,3,P> calcUnnormalizedQBlockForQuaternion(
const Vec4P& q)
932 static Mat<3,4,P> calcUnnormalizedQInvBlockForQuaternion(
const Vec4P& q)
935 static Vec3P convertAngVelToBodyFixed123Dot(
const Vec3P& q,
const Vec3P& w_PB_B)
938 static Vec3P convertBodyFixed123DotToAngVel(
const Vec3P& q,
const Vec3P& qdot)
941 static Vec3P convertAngVelDotToBodyFixed123DotDot
942 (
const Vec3P& q,
const Vec3P& w_PB_B,
const Vec3P& wdot_PB_B)
964 static Rotation_ aboutAxis(
const RealP& angleInRad,
const UnitVec3P& axis ) {
return Rotation_(angleInRad,axis); }
965 static Rotation_ aboutAxis(
const RealP& angleInRad,
const Vec3P& axis ) {
return Rotation_(angleInRad,axis); }
978 static Rotation_ aboutYThenNewX(
const RealP& yInRad,
const RealP& xInRad) {
return aboutXThenOldY(xInRad, yInRad); }
979 static Rotation_ aboutXThenNewZ(
const RealP& xInRad,
const RealP& zInRad) {
return aboutZThenOldX(zInRad, xInRad); }
980 static Rotation_ aboutZThenNewX(
const RealP& zInRad,
const RealP& xInRad) {
return aboutXThenOldZ(xInRad, zInRad); }
981 static Rotation_ aboutYThenNewZ(
const RealP& yInRad,
const RealP& zInRad) {
return aboutZThenOldY(zInRad, yInRad); }
982 static Rotation_ aboutZThenNewY(
const RealP& zInRad,
const RealP& yInRad) {
return aboutYThenOldZ(yInRad, zInRad); }
1029 class InverseRotation_ :
public Mat<3,3,P>::TransposeType {
1031 typedef Rotation_<P> RotationP;
1032 typedef Mat<3,3,P> Mat33P;
1033 typedef SymMat<3,P> SymMat33P;
1034 typedef Mat<2,2,P> Mat22P;
1035 typedef Mat<3,2,P> Mat32P;
1036 typedef Vec<2,P> Vec2P;
1037 typedef Vec<3,P> Vec3P;
1038 typedef Vec<4,P> Vec4P;
1039 typedef Quaternion_<P> QuaternionP;
1063 { BaseMat::operator=(R.asMat33());
return *
this; }
1113 operator<<(std::ostream&, const Rotation_<P>&);
1116 operator<<(std::ostream&, const InverseRotation_<P>&);
1122 template <
class P,
int S>
inline UnitVec<P,1>
1124 template <
class P,
int S>
inline UnitRow<P,1>
1126 template <
class P,
int S>
inline UnitVec<P,1>
1128 template <
class P,
int S>
inline UnitRow<P,1>
1133 template <
class P>
inline
1150 template <
class P>
inline Rotation_<P>
1152 template <
class P>
inline Rotation_<P>
1154 template <
class P>
inline Rotation_<P>
1161 template <
class P>
inline Rotation_<P>
1163 template <
class P>
inline Rotation_<P>
1165 template <
class P>
inline Rotation_<P>
1167 template <
class P>
inline Rotation_<P>
1176 #endif // SimTK_SimTKCOMMON_ROTATION_H_
Rotation_ & setRotationFromAngleAboutZ(RealP cosAngle, RealP sinAngle)
Set this Rotation_ object to a right-handed rotation by an angle (in radians) about a coordinate axis...
Definition: Rotation.h:142
Rotation_(const Rotation_ &R)
Definition: Rotation.h:124
Matrix_< E > operator/(const MatrixBase< E > &l, const typename CNT< E >::StdNumber &r)
Definition: BigMatrix.h:2693
InverseRotation_< Real > InverseRotation
Definition: Rotation.h:53
Rotation_ & setRotationFromQuaternion(const QuaternionP &q)
Method for creating a rotation matrix from a quaternion.
#define SimTK_SimTKCOMMON_EXPORT
Definition: SimTKcommon/include/SimTKcommon/internal/common.h:202
static Vec3P multiplyByBodyXYZ_NInvT_P(const Vec2P &cosxy, const Vec2P &sinxy, const Vec3P &v_P)
Fastest way to form the product q=~NInv_P*v_P=~(~v_P*NInv_P).
Definition: Rotation.h:566
----------------------------------------------------------------------------- This InverseRotation cl...
Definition: Rotation.h:47
InverseRotation_< double > dInverseRotation
Definition: Rotation.h:55
Defines the CoordinateAxis and CoordinateDirection classes.
bool areAllRotationElementsSameToEpsilon(const Rotation_ &R, RealP epsilon) const
Definition: Rotation.h:255
Rotation_(RealP angle, const Vec3P &nonUnitVector)
Constructor for right-handed rotation by an angle (in radians) about an arbitrary vector...
Definition: Rotation.h:148
static Mat33P calcNDotForBodyXYZInParentFrame(const Vec3P &q, const Vec3P &qdot)
Given Euler angles forming a body-fixed X-Y-Z (123) sequence q, and their time derivatives qdot...
Definition: Rotation.h:630
This class is a Vec3 plus an ironclad guarantee either that:
Definition: UnitVec.h:40
const RowType & operator[](int i) const
Definition: Rotation.h:302
QuaternionP convertRotationToQuaternion() const
Converts rotation matrix to a quaternion.
static Vec3P convertAngVelDotInBodyFrameToBodyXYZDotDot(const Vec3P &q, const Vec3P &w_PB_B, const Vec3P &wdot_PB_B)
TODO: sherm: is this right? Warning: everything is measured in the *PARENT* frame, but has to be expressed in the *BODY* frame.
Definition: Rotation.h:775
const RotationP & transpose() const
Transpose, and transpose operators (override BaseMat versions of transpose).
Definition: Rotation.h:1082
static Vec3P convertBodyXYZDotToAngVelInBodyFrame(const Vec3P &q, const Vec3P &qdot)
Inverse of the above routine.
Definition: Rotation.h:752
Rotation_(BodyOrSpaceType bodyOrSpace, RealP angle1, const CoordinateAxis &axis1, RealP angle2, const CoordinateAxis &axis2, RealP angle3, const CoordinateAxis &axis3)
Constructor for three-angle Body-fixed or Space-fixed rotation sequences (angles are in radians) ...
Definition: Rotation.h:159
Vec3P convertRotationToBodyFixedXYZ() const
A convenient special case of convertThreeAxesRotationToThreeAngles().
Definition: Rotation.h:230
const RowType & row(int i) const
Definition: Rotation.h:297
Rotation_ & operator=(const Rotation_ &R)
Definition: Rotation.h:125
Rotation_< double > dRotation
Definition: Rotation.h:51
static Vec3P convertAngVelDotToBodyFixed321DotDot(const Vec3P &q, const Vec3P &w_PB_B, const Vec3P &wdot_PB_B)
Definition: Rotation.h:370
RotationP & operator~()
Transpose, and transpose operators (override BaseMat versions of transpose).
Definition: Rotation.h:1085
Definition: CoordinateAxis.h:196
Definition: CoordinateAxis.h:199
bool isSameRotationToWithinAngle(const Rotation_ &R, RealP okPointingAngleErrorRads) const
Return true if "this" Rotation is nearly identical to "R" within a specified pointing angle error...
static Vec4P OLDconvertAngVelDotToQuaternionDotDot(const Vec4P &q, const Vec3P &w_PB_P, const Vec3P &wdot_PB_P)
Given a quaternion q representing R_PB, angular velocity of B in P, and the time derivative of the an...
Definition: Rotation.h:870
static Mat34P calcUnnormalizedNInvForQuaternion(const Vec4P &q)
Given a (possibly unnormalized) quaternion q, calculate the 3x4 matrix NInv (= N^-1) which maps quate...
Definition: Rotation.h:839
Rotation_ & setRotationFromAngleAboutY(RealP cosAngle, RealP sinAngle)
Set this Rotation_ object to a right-handed rotation by an angle (in radians) about a coordinate axis...
Definition: Rotation.h:141
Definition: CoordinateAxis.h:202
Definition: Rotation.h:42
static Vec3P convertAngAccInParentToBodyXYZDotDot(const Vec2P &cosxy, const Vec2P &sinxy, RealP oocosy, const Vec3P &qdot, const Vec3P &b_PB)
Calculate second time derivative qdotdot of body-fixed XYZ Euler angles q given sines and cosines of ...
Definition: Rotation.h:527
static Vec4P convertAngVelToQuaternionDot(const Vec4P &q, const Vec3P &w_PB_P)
Given a possibly unnormalized quaternion (0th element is the scalar) and the relative angular velocit...
Definition: Rotation.h:852
InverseRotation_< P > & updInvert()
Convert from Rotation_ to writable InverseRotation_ (no cost).
Definition: Rotation.h:268
const InverseRotation_< P > & operator~() const
Transpose, and transpose operators.
Definition: Rotation.h:274
Rotation_ & setRotationFromAngleAboutUnitVector(RealP angle, const UnitVec3P &unitVector)
Set this Rotation_ object to a right-handed rotation of an angle (in radians) about an arbitrary vect...
InverseRotation_< P > & updTranspose()
Transpose, and transpose operators.
Definition: Rotation.h:275
const ColType & y() const
Access individual rows and columns of this InverseRotation; no cost or copying since suitably-cast re...
Definition: Rotation.h:1097
Rotation_ & setRotationFromApproximateMat33(const Mat33P &m)
Set this Rotation_ object to an (hopefully nearby) orthogonal rotation matrix from a generic Mat33P...
SymMat33P reexpressSymMat33(const SymMat33P &S_BB) const
Perform an efficient transform of a symmetric matrix that must be re-expressed with a multiply from b...
Rotation_(RealP angle, const CoordinateAxis::YCoordinateAxis)
Constructor for right-handed rotation by an angle (in radians) about a coordinate axis...
Definition: Rotation.h:131
InverseRotation_< P > & operator~()
Transpose, and transpose operators.
Definition: Rotation.h:276
const CoordinateAxis::ZCoordinateAxis ZAxis
Constant representing the Z coordinate axis; will implicitly convert to the integer 2 when used in a ...
const ColType & z() const
Access individual rows and columns of this InverseRotation; no cost or copying since suitably-cast re...
Definition: Rotation.h:1098
Vec4P convertQuaternionToAngleAxis() const
Returns [ a vx vy vz ] with (a,v) in canonical form, i.e., -180 < a <= 180 and |v|=1.
This class, along with its sister class CoordinateDirection, provides convenient manipulation of the ...
Definition: CoordinateAxis.h:53
RealP convertOneAxisRotationToOneAngle(const CoordinateAxis &axis1) const
Converts rotation matrix to a single orientation angle.
Rotation_ & setRotationFromThreeAnglesThreeAxes(BodyOrSpaceType bodyOrSpace, RealP angle1, const CoordinateAxis &axis1, RealP angle2, const CoordinateAxis &axis2, RealP angle3, const CoordinateAxis &axis3)
Set this Rotation_ object to a three-angle Body-fixed or Space-fixed rotation sequences (angles are i...
Rotation_ & setRotationFromAngleAboutAxis(RealP angle, const CoordinateAxis &axis)
Set this Rotation_ object to a right-handed rotation by an angle (in radians) about a coordinate axis...
Definition: Rotation.h:136
const ColType & y() const
Definition: Rotation.h:300
static Mat43P calcUnnormalizedNDotForQuaternion(const Vec4P &qdot)
Given the time derivative qdot of a possibly unnormalized quaternion q, calculate the 4x3 matrix NDot...
Definition: Rotation.h:822
static Mat33P calcNInvForBodyXYZInParentFrame(const Vec3P &q)
Inverse of the above routine.
Definition: Rotation.h:697
const RowType & row(int i) const
Access individual rows and columns of this InverseRotation; no cost or copying since suitably-cast re...
Definition: Rotation.h:1094
Declares and defines the UnitVec and UnitRow classes.
const CoordinateAxis::YCoordinateAxis YAxis
Constant representing the Y coordinate axis; will implicitly convert to the integer 1 when used in a ...
UnitVec< P, BaseMat::RowSpacing > ColType
Note that the unit vectors representing the rows and columns of this matrix do not necessarily have u...
Definition: Rotation.h:1049
Rotation_(const UnitVec3P &uvec, const CoordinateAxis axis)
Calculate R_AB by knowing one of B's unit vector expressed in A.
Definition: Rotation.h:191
Rotation_ & setRotationToIdentityMatrix()
Definition: Rotation.h:120
InverseRotation_()
You should not ever construct one of these as they should only occur as expression intermediates resu...
Definition: Rotation.h:1057
Rotation_< Real > Rotation
Definition: Rotation.h:47
static Mat33P calcNForBodyXYZInParentFrame(const Vec3P &q)
Given Euler angles q forming a body-fixed X-Y-Z (123) sequence return the block N_P of the system N m...
Definition: Rotation.h:440
bool isYAxis() const
Return true if this is the Y axis.
Definition: CoordinateAxis.h:92
static Vec3P multiplyByBodyXYZ_NInv_P(const Vec2P &cosxy, const Vec2P &sinxy, const Vec3P &qdot)
Fastest way to form the product w_PB=NInv_P*qdot.
Definition: Rotation.h:550
This is a fixed length column vector designed for no-overhead inline computation. ...
Definition: Vec.h:131
const BaseVec & asVec3() const
Return a reference to the underlying Vec3 (no copying here).
Definition: UnitVec.h:118
void setRotationToBodyFixedXY(const Vec2P &v)
Set this Rotation_ to represent a rotation characterized by subsequent rotations of: +v[0] about the ...
Definition: Rotation.h:168
BaseMat toMat33() const
Conversion from InverseRotation_ to BaseMat.
Definition: Rotation.h:1107
const ColType & x() const
Access individual rows and columns of this InverseRotation; no cost or copying since suitably-cast re...
Definition: Rotation.h:1096
const BaseRow & asRow3() const
Return a const reference to the Row3 underlying this UnitRow.
Definition: UnitVec.h:258
RotationP & updTranspose()
Transpose, and transpose operators (override BaseMat versions of transpose).
Definition: Rotation.h:1084
InverseRotation_ & operator=(const InverseRotation_ &R)
This is an explicit implementation of the default copy assignment operator.
Definition: Rotation.h:1062
Matrix_< E > operator*(const MatrixBase< E > &l, const typename CNT< E >::StdNumber &r)
Definition: BigMatrix.h:2685
RealP getMaxAbsDifferenceInRotationElements(const Rotation_ &R) const
Definition: Rotation.h:247
Rotation_ & setRotationToNaN()
Definition: Rotation.h:121
Definition: Rotation.h:42
Vec4P convertRotationToAngleAxis() const
Converts rotation matrix to angle-axis form.
Definition: Rotation.h:225
Mat< 3, 3, P >::TransposeType BaseMat
This is the type of the underlying 3x3 matrix; note that it will have unusual row and column spacing ...
Definition: Rotation.h:1043
Rotation_ & setRotationColFromUnitVecTrustMe(int colj, const UnitVec3P &uvecj)
Set the Rotation_ matrix directly - but you had better know what you are doing!
Definition: Rotation.h:309
bool areAllRotationElementsSameToMachinePrecision(const Rotation_ &R) const
Definition: Rotation.h:257
Rotation_(BodyOrSpaceType bodyOrSpace, RealP angle1, const CoordinateAxis &axis1, RealP angle2, const CoordinateAxis &axis2)
Constructor for two-angle, two-axes, Body-fixed or Space-fixed rotation sequences (angles are in radi...
Definition: Rotation.h:157
static Vec4P convertAngVelDotToQuaternionDotDot(const Vec4P &q, const Vec3P &w_PB, const Vec3P &b_PB)
We want to differentiate qdot=N(q)*w to get qdotdot=N*b+NDot*w where b is angular acceleration wdot...
Definition: Rotation.h:886
const RotationP & invert() const
We can invert an InverseRotation just by recasting it to a Rotation at zero cost. ...
Definition: Rotation.h:1075
static Vec3P multiplyByBodyXYZ_N_P(const Vec2P &cosxy, const Vec2P &sinxy, RealP oocosy, const Vec3P &w_PB)
This is the fastest way to form the product qdot=N_P*w_PB for a body-fixed XYZ sequence where angular...
Definition: Rotation.h:471
Rotation_ & setRotationFromTwoAnglesTwoAxes(BodyOrSpaceType bodyOrSpace, RealP angle1, const CoordinateAxis &axis1, RealP angle2, const CoordinateAxis &axis2)
Set this Rotation_ object to a two-angle, two-axes, Body-fixed or Space-fixed rotation sequences (ang...
Rotation_(RealP angle, const UnitVec3P &unitVector)
Constructor for right-handed rotation by an angle (in radians) about an arbitrary vector...
Definition: Rotation.h:147
Rotation_ & setRotationFromAngleAboutNonUnitVector(RealP angle, const Vec3P &nonUnitVector)
Set this Rotation_ object to a right-handed rotation of an angle (in radians) about an arbitrary vect...
Definition: Rotation.h:152
static Mat33P calcNDotForBodyXYZInBodyFrame(const Vec3P &q, const Vec3P &qdot)
Given Euler angles forming a body-fixed X-Y-Z (123) sequence q, and their time derivatives qdot...
Definition: Rotation.h:589
This file is the user-includeable header to be included in user programs to provide fixed-length Vec ...
Rotation_ & setRotationFromAngleAboutX(RealP cosAngle, RealP sinAngle)
Set this Rotation_ object to a right-handed rotation by an angle (in radians) about a coordinate axis...
Definition: Rotation.h:140
const ColType & z() const
Definition: Rotation.h:301
void setToNaN()
Definition: Mat.h:856
Rotation_(const Mat33P &m, bool)
Construct a Rotation_ directly from a Mat33P (we trust that m is a valid Rotation_!) ...
Definition: Rotation.h:181
Rotation_(const UnitVec3P &uveci, const CoordinateAxis &axisi, const Vec3P &vecjApprox, const CoordinateAxis &axisjApprox)
Calculate R_AB by knowing one of B's unit vectors u1 (could be Bx, By, or Bz) expressed in A and a ve...
Definition: Rotation.h:202
Rotation_(RealP angle, const CoordinateAxis::XCoordinateAxis)
Constructor for right-handed rotation by an angle (in radians) about a coordinate axis...
Definition: Rotation.h:130
const ColType & col(int j) const
Definition: Rotation.h:298
static Vec3P convertQuaternionDotToAngVel(const Vec4P &q, const Vec4P &qdot)
Inverse of the above routine.
Definition: Rotation.h:859
const ColType & operator()(int j) const
Definition: Rotation.h:303
BodyOrSpaceType
Definition: Rotation.h:42
const Mat33P & asMat33() const
Conversion from Rotation to its base class Mat33.
Definition: Rotation.h:290
Vec3P convertThreeAxesRotationToThreeAngles(BodyOrSpaceType bodyOrSpace, const CoordinateAxis &axis1, const CoordinateAxis &axis2, const CoordinateAxis &axis3) const
Converts rotation matrix to three orientation angles.
Rotation_()
Definition: Rotation.h:119
Rotation_ & setRotationFromUnitVecsTrustMe(const UnitVec3P &colA, const UnitVec3P &colB, const UnitVec3P &colC)
Set the Rotation_ matrix directly - but you had better know what you are doing!
Definition: Rotation.h:311
Rotation_ & setRotationFromAngleAboutX(RealP angle)
Set this Rotation_ object to a right-handed rotation by an angle (in radians) about a coordinate axis...
Definition: Rotation.h:137
UnitRow< P, Mat33P::ColSpacing > RowType
Definition: Rotation.h:296
Rotation_(RealP angle, const CoordinateAxis &axis)
Constructor for right-handed rotation by an angle (in radians) about a coordinate axis...
Definition: Rotation.h:129
bool isSameRotationToWithinAngleOfMachinePrecision(const Rotation_ &R) const
Return true if "this" Rotation is nearly identical to "R" within a specified pointing angle error...
Definition: Rotation.h:244
Mat33P toMat33() const
Conversion from Rotation to its base class Mat33.
Definition: Rotation.h:291
bool isXAxis() const
Return true if this is the X axis.
Definition: CoordinateAxis.h:90
UnitVec< P, Mat33P::RowSpacing > ColType
The column and row unit vector types do not necessarily have unit spacing.
Definition: Rotation.h:295
Mat & operator=(const Mat &src)
Copy assignment copies only the elements that are present and does not touch any unused memory space ...
Definition: Mat.h:258
const ColType & x() const
Definition: Rotation.h:299
InverseRotation_< float > fInverseRotation
Definition: Rotation.h:54
InverseRotation_(const InverseRotation_ &R)
This is an explicit implementation of the default copy constructor.
Definition: Rotation.h:1060
static Mat43P calcUnnormalizedNForQuaternion(const Vec4P &q)
Given a possibly unnormalized quaternion q, calculate the 4x3 matrix N which maps angular velocity w ...
Definition: Rotation.h:807
Rotation_ & setRotationFromMat33TrustMe(const Mat33P &m)
Set the Rotation_ matrix directly - but you had better know what you are doing!
Definition: Rotation.h:307
const RotationP & operator~() const
Transpose, and transpose operators (override BaseMat versions of transpose).
Definition: Rotation.h:1083
Vec2P convertRotationToBodyFixedXY() const
A convenient special case of convertTwoAxesRotationToTwoAngles().
Definition: Rotation.h:228
static Vec3P multiplyByBodyXYZ_NT_P(const Vec2P &cosxy, const Vec2P &sinxy, RealP oocosy, const Vec3P &q)
This is the fastest way to form the product v_P=~N_P*q=~(~q*N_P); see the untransposed method multipl...
Definition: Rotation.h:487
static Mat33P calcNForBodyXYZInBodyFrame(const Vec3P &q)
Given Euler angles q forming a body-fixed X-Y-Z sequence return the block N_B of the system N matrix ...
Definition: Rotation.h:404
ScalarNormSq normSqr() const
Definition: Vec.h:553
static Vec3P convertBodyFixed321DotToAngVel(const Vec3P &q, const Vec3P &qd)
Inverse of the above routine.
Definition: Rotation.h:355
Rotation_ & operator*=(const Rotation_< P > &R)
In-place composition of Rotation matrices.
Definition: Rotation.h:1138
Rotation_ & setRotationFromTwoAxes(const UnitVec3P &uveci, const CoordinateAxis &axisi, const Vec3P &vecjApprox, const CoordinateAxis &axisjApprox)
Calculate R_AB by knowing one of B's unit vectors u1 (could be Bx, By, or Bz) expressed in A and a ve...
static Vec3P convertAngVelInParentToBodyXYZDot(const Vec2P &cosxy, const Vec2P &sinxy, RealP oocosy, const Vec3P &w_PB)
Calculate first time derivative qdot of body-fixed XYZ Euler angles q given sines and cosines of the ...
Definition: Rotation.h:507
Rotation_< float > fRotation
Definition: Rotation.h:50
const BaseMat & asMat33() const
Conversion from InverseRotation_ to BaseMat.
Definition: Rotation.h:1106
A Quaternion is a Vec4 with the following behavior:
Definition: Quaternion.h:41
Rotation_(RealP angle, const CoordinateAxis::ZCoordinateAxis)
Constructor for right-handed rotation by an angle (in radians) about a coordinate axis...
Definition: Rotation.h:132
Rotation_ & setRotationFromAngleAboutZ(RealP angle)
Set this Rotation_ object to a right-handed rotation by an angle (in radians) about a coordinate axis...
Definition: Rotation.h:139
static Mat33P calcNForBodyXYZInBodyFrame(const Vec3P &cq, const Vec3P &sq)
This faster version of calcNForBodyXYZInBodyFrame() assumes you have already calculated the cosine an...
Definition: Rotation.h:417
UnitRow< P, BaseMat::ColSpacing > RowType
This is the type of a row of this InverseRotation.
Definition: Rotation.h:1051
static Vec3P convertAngVelToBodyFixed321Dot(const Vec3P &q, const Vec3P &w_PB_B)
Given Euler angles forming a body-fixed 3-2-1 sequence, and the relative angular velocity vector of B...
Definition: Rotation.h:341
void setRotationToBodyFixedXYZ(const Vec3P &c, const Vec3P &s)
Given cosines and sines (in that order) of three angles, set this Rotation matrix to the body-fixed 1...
Definition: Rotation.h:326
Rotation_(const Mat33P &m)
Constructs an (hopefully nearby) orthogonal rotation matrix from a generic Mat33P.
Definition: Rotation.h:184
Vec2P convertTwoAxesRotationToTwoAngles(BodyOrSpaceType bodyOrSpace, const CoordinateAxis &axis1, const CoordinateAxis &axis2) const
Converts rotation matrix to two orientation angles.
This class represents a small matrix whose size is known at compile time, containing elements of any ...
Definition: Mat.h:51
This type is used for the transpose of UnitVec, and as the returned row type of a Rotation...
Definition: UnitVec.h:41
void setRotationToBodyFixedXYZ(const Vec3P &v)
Set this Rotation_ to represent a rotation characterized by subsequent rotations of: +v[0] about the ...
Definition: Rotation.h:173
const InverseRotation_< P > & invert() const
Convert from Rotation_ to InverseRotation_ (no cost). Overrides base class invert().
Definition: Rotation.h:266
Rotation_(const QuaternionP &q)
Constructor for creating a rotation matrix from a quaternion.
Definition: Rotation.h:176
Rotation_ & operator/=(const Rotation_< P > &R)
In-place composition of Rotation matrices.
Definition: Rotation.h:1140
SymMat33P reexpressSymMat33(const SymMat33P &S_BB) const
Assuming this InverseRotation_ is R_AB, and given a symmetric dyadic matrix S_BB expressed in B...
const InverseRotation_< P > & transpose() const
Transpose, and transpose operators.
Definition: Rotation.h:273
const RowType & operator[](int i) const
Access individual rows and columns of this InverseRotation; no cost or copying since suitably-cast re...
Definition: Rotation.h:1099
static Mat33P calcNForBodyXYZInParentFrame(const Vec3P &cq, const Vec3P &sq)
This faster version of calcNForBodyXYZInParentFrame() assumes you have already calculated the cosine ...
Definition: Rotation.h:454
Rotation_ & setRotationFromAngleAboutY(RealP angle)
Set this Rotation_ object to a right-handed rotation by an angle (in radians) about a coordinate axis...
Definition: Rotation.h:138
const ColType & operator()(int j) const
Access individual rows and columns of this InverseRotation; no cost or copying since suitably-cast re...
Definition: Rotation.h:1100
static Vec3P convertAngVelInBodyFrameToBodyXYZDot(const Vec3P &q, const Vec3P &w_PB_B)
Given Euler angles forming a body-fixed X-Y-Z (123) sequence, and the relative angular velocity vecto...
Definition: Rotation.h:729
Rotation_ & setRotationFromOneAxis(const UnitVec3P &uvec, const CoordinateAxis axis)
Calculate R_AB by knowing one of B's unit vector expressed in A.
RotationP & updInvert()
We can invert an InverseRotation just by recasting it to a Rotation at zero cost. ...
Definition: Rotation.h:1076
const ColType & col(int j) const
Access individual rows and columns of this InverseRotation; no cost or copying since suitably-cast re...
Definition: Rotation.h:1095
static Mat33P calcNInvForBodyXYZInBodyFrame(const Vec3P &q)
Inverse of routine calcNForBodyXYZInBodyFrame().
Definition: Rotation.h:667
const CoordinateAxis::XCoordinateAxis XAxis
Constant representing the X coordinate axis; will implicitly convert to the integer 0 when used in a ...