Simbody  3.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Rotation.h
Go to the documentation of this file.
1 #ifndef SimTK_SimTKCOMMON_ROTATION_H_
2 #define SimTK_SimTKCOMMON_ROTATION_H_
3 
4 /* -------------------------------------------------------------------------- *
5  * Simbody(tm): SimTKcommon *
6  * -------------------------------------------------------------------------- *
7  * This is part of the SimTK biosimulation toolkit originating from *
8  * Simbios, the NIH National Center for Physics-Based Simulation of *
9  * Biological Structures at Stanford, funded under the NIH Roadmap for *
10  * Medical Research, grant U54 GM072970. See https://simtk.org/home/simbody. *
11  * *
12  * Portions copyright (c) 2005-12 Stanford University and the Authors. *
13  * Authors: Paul Mitiguy, Michael Sherman *
14  * Contributors: *
15  * *
16  * Licensed under the Apache License, Version 2.0 (the "License"); you may *
17  * not use this file except in compliance with the License. You may obtain a *
18  * copy of the License at http://www.apache.org/licenses/LICENSE-2.0. *
19  * *
20  * Unless required by applicable law or agreed to in writing, software *
21  * distributed under the License is distributed on an "AS IS" BASIS, *
22  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. *
23  * See the License for the specific language governing permissions and *
24  * limitations under the License. *
25  * -------------------------------------------------------------------------- */
26 
27 //------------------------------------------------------------------------------
28 
33 
34 //------------------------------------------------------------------------------
35 #include <iosfwd> // Forward declaration of iostream
36 //------------------------------------------------------------------------------
37 
38 //------------------------------------------------------------------------------
39 namespace SimTK {
40 
41 
43 
44 //------------------------------------------------------------------------------
45 // Forward declarations
46 template <class P> class Rotation_;
47 template <class P> class InverseRotation_;
48 
52 
56 
57 //------------------------------------------------------------------------------
102 //------------------------------------------------------------------------------
103 template <class P> // templatized by precision
104 class Rotation_ : public Mat<3,3,P> {
105  typedef P RealP;
106  typedef Mat<2,2,P> Mat22P;
107  typedef Mat<3,2,P> Mat32P;
108  typedef Mat<3,3,P> Mat33P;
109  typedef Mat<4,3,P> Mat43P;
110  typedef Mat<3,4,P> Mat34P;
111  typedef Vec<2,P> Vec2P;
112  typedef Vec<3,P> Vec3P;
113  typedef Vec<4,P> Vec4P;
114  typedef UnitVec<P,1> UnitVec3P; // stride is 1 here, length is always 3
115  typedef SymMat<3,P> SymMat33P;
116  typedef Quaternion_<P> QuaternionP;
117 public:
118  // Default constructor and constructor-like methods
119  Rotation_() : Mat33P(1) {}
121  Rotation_& setRotationToNaN() { Mat33P::setToNaN(); return *this; }
122 
123  // Default copy constructor and assignment operator
124  Rotation_( const Rotation_& R ) : Mat33P(R) {}
125  Rotation_& operator=( const Rotation_& R ) { Mat33P::operator=( R.asMat33() ); return *this; }
126 
128 
129  Rotation_( RealP angle, const CoordinateAxis& axis ) { setRotationFromAngleAboutAxis( angle, axis ); }
130  Rotation_( RealP angle, const CoordinateAxis::XCoordinateAxis ) { setRotationFromAngleAboutX( std::cos(angle), std::sin(angle) ); }
131  Rotation_( RealP angle, const CoordinateAxis::YCoordinateAxis ) { setRotationFromAngleAboutY( std::cos(angle), std::sin(angle) ); }
132  Rotation_( RealP angle, const CoordinateAxis::ZCoordinateAxis ) { setRotationFromAngleAboutZ( std::cos(angle), std::sin(angle) ); }
134 
137  Rotation_& setRotationFromAngleAboutX( RealP angle ) { return setRotationFromAngleAboutX( std::cos(angle), std::sin(angle) ); }
138  Rotation_& setRotationFromAngleAboutY( RealP angle ) { return setRotationFromAngleAboutY( std::cos(angle), std::sin(angle) ); }
139  Rotation_& setRotationFromAngleAboutZ( RealP angle ) { return setRotationFromAngleAboutZ( std::cos(angle), std::sin(angle) ); }
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; }
144 
146 
147  Rotation_( RealP angle, const UnitVec3P& unitVector ) { setRotationFromAngleAboutUnitVector(angle,unitVector); }
148  Rotation_( RealP angle, const Vec3P& nonUnitVector ) { setRotationFromAngleAboutNonUnitVector(angle,nonUnitVector); }
150 
152  Rotation_& setRotationFromAngleAboutNonUnitVector( RealP angle, const Vec3P& nonUnitVector ) { return setRotationFromAngleAboutUnitVector( angle, UnitVec3P(nonUnitVector) ); }
153  SimTK_SimTKCOMMON_EXPORT Rotation_& setRotationFromAngleAboutUnitVector( RealP angle, const UnitVec3P& unitVector );
155 
157  Rotation_( BodyOrSpaceType bodyOrSpace, RealP angle1, const CoordinateAxis& axis1, RealP angle2, const CoordinateAxis& axis2 ) { setRotationFromTwoAnglesTwoAxes( bodyOrSpace,angle1,axis1,angle2,axis2); }
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); }
161  SimTK_SimTKCOMMON_EXPORT Rotation_& setRotationFromTwoAnglesTwoAxes( BodyOrSpaceType bodyOrSpace, RealP angle1, const CoordinateAxis& axis1, RealP angle2, const CoordinateAxis& axis2 );
163  SimTK_SimTKCOMMON_EXPORT Rotation_& setRotationFromThreeAnglesThreeAxes( BodyOrSpaceType bodyOrSpace, RealP angle1, const CoordinateAxis& axis1, RealP angle2, const CoordinateAxis& axis2, RealP angle3, const CoordinateAxis& axis3 );
164 
174 
176  explicit Rotation_( const QuaternionP& q ) { setRotationFromQuaternion(q); }
179 
181  Rotation_( const Mat33P& m, bool ) : Mat33P(m) {}
182 
184  explicit Rotation_( const Mat33P& m ) { setRotationFromApproximateMat33(m); }
187 
190 
191  Rotation_( const UnitVec3P& uvec, const CoordinateAxis axis ) { setRotationFromOneAxis(uvec,axis); }
192  SimTK_SimTKCOMMON_EXPORT Rotation_& setRotationFromOneAxis( const UnitVec3P& uvec, const CoordinateAxis axis );
194 
201 
202  Rotation_( const UnitVec3P& uveci, const CoordinateAxis& axisi, const Vec3P& vecjApprox, const CoordinateAxis& axisjApprox ) { setRotationFromTwoAxes(uveci,axisi,vecjApprox,axisjApprox); }
203  SimTK_SimTKCOMMON_EXPORT Rotation_& setRotationFromTwoAxes( const UnitVec3P& uveci, const CoordinateAxis& axisi, const Vec3P& vecjApprox, const CoordinateAxis& axisjApprox );
205 
206  // Converts rotation matrix to one or two or three orientation angles.
207  // Note: The result is most meaningful if the Rotation_ matrix is one that can be produced by such a sequence.
208  // Use1: someRotation.convertOneAxisRotationToOneAngle( XAxis );
209  // Use2: someRotation.convertTwoAxesRotationToTwoAngles( SpaceRotationSequence, YAxis, ZAxis );
210  // Use3: someRotation.convertThreeAxesRotationToThreeAngles( SpaceRotationSequence, ZAxis, YAxis, XAxis );
211  // Use4: someRotation.convertRotationToAngleAxis(); Return: [angleInRadians, unitVectorX, unitVectorY, unitVectorZ].
212 
218  SimTK_SimTKCOMMON_EXPORT Vec2P convertTwoAxesRotationToTwoAngles( BodyOrSpaceType bodyOrSpace, const CoordinateAxis& axis1, const CoordinateAxis& axis2 ) const;
221  SimTK_SimTKCOMMON_EXPORT Vec3P convertThreeAxesRotationToThreeAngles( BodyOrSpaceType bodyOrSpace, const CoordinateAxis& axis1, const CoordinateAxis& axis2, const CoordinateAxis& axis3 ) const;
226 
231 
239  SimTK_SimTKCOMMON_EXPORT SymMat33P reexpressSymMat33(const SymMat33P& S_BB) const;
240 
242 
243  SimTK_SimTKCOMMON_EXPORT bool isSameRotationToWithinAngle( const Rotation_& R, RealP okPointingAngleErrorRads ) const;
248  const Mat33P& A = asMat33(); const Mat33P& B = R.asMat33(); RealP maxDiff = 0;
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; }
252  return maxDiff;
253  }
254 
255  bool areAllRotationElementsSameToEpsilon( const Rotation_& R, RealP epsilon ) const
256  { return getMaxAbsDifferenceInRotationElements(R) <= epsilon ; }
259 
261  inline Rotation_( const InverseRotation_<P>& );
263  inline Rotation_& operator=( const InverseRotation_<P>& );
264 
266  const InverseRotation_<P>& invert() const { return *reinterpret_cast<const InverseRotation_<P>*>(this); }
268  InverseRotation_<P>& updInvert() { return *reinterpret_cast<InverseRotation_<P>*>(this); }
269 
272 
273  const InverseRotation_<P>& transpose() const { return invert(); }
274  const InverseRotation_<P>& operator~() const { return invert(); }
278 
280 
281  inline Rotation_& operator*=( const Rotation_<P>& R );
282  inline Rotation_& operator/=( const Rotation_<P>& R );
283  inline Rotation_& operator*=( const InverseRotation_<P>& );
284  inline Rotation_& operator/=( const InverseRotation_<P>& );
286 
289 
290  const Mat33P& asMat33() const { return *static_cast<const Mat33P*>(this); }
291  Mat33P toMat33() const { return asMat33(); }
293 
297  const RowType& row( int i ) const { return reinterpret_cast<const RowType&>(asMat33()[i]); }
298  const ColType& col( int j ) const { return reinterpret_cast<const ColType&>(asMat33()(j)); }
299  const ColType& x() const { return col(0); }
300  const ColType& y() const { return col(1); }
301  const ColType& z() const { return col(2); }
302  const RowType& operator[]( int i ) const { return row(i); }
303  const ColType& operator()( int j ) const { return col(j); }
304 
306 
308  { Mat33P& R = *this; R=m; return *this; }
310  { Mat33P& R = *this; R(colj)=uvecj.asVec3(); return *this; }
311  Rotation_& setRotationFromUnitVecsTrustMe( const UnitVec3P& colA, const UnitVec3P& colB, const UnitVec3P& colC )
312  { Mat33P& R = *this; R(0)=colA.asVec3(); R(1)=colB.asVec3(); R(2)=colC.asVec3(); return *this; }
314 
315 //--------------------------- PAUL CONTINUE FROM HERE --------------------------
316 public:
317 //------------------------------------------------------------------------------
318 
319  // These are ad hoc routines that don't match the nice API Paul Mitiguy
320  // implemented above.
321 
322 
326  void setRotationToBodyFixedXYZ(const Vec3P& c, const Vec3P& s) {
327  Mat33P& R = *this;
328  const RealP s0s1=s[0]*s[1], s2c0=s[2]*c[0], c0c2=c[0]*c[2], nc1= -c[1];
329 
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] );
333  }
334 
335 
341  static Vec3P convertAngVelToBodyFixed321Dot(const Vec3P& q, const Vec3P& w_PB_B) {
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;
346 
347  const Mat33P E( 0, s2oc1 , c2oc1 ,
348  0, c2 , -s2 ,
349  1, s1*s2oc1 , s1*c2oc1 );
350  return E * w_PB_B;
351  }
352 
355  static Vec3P convertBodyFixed321DotToAngVel(const Vec3P& q, const Vec3P& qd) {
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]);
358 
359  const Mat33P Einv( -s1 , 0 , 1 ,
360  c1*s2 , c2 , 0 ,
361  c1*c2 , -s2 , 0 );
362  return Einv*qd;
363  }
364 
365  // TODO: sherm: is this right? Warning: everything is measured in the
366  // *PARENT* frame, but angular velocities and accelerations are
367  // expressed in the *BODY* frame.
368  // TODO: this is not an efficient way to do this computation.
370  (const Vec3P& q, const Vec3P& w_PB_B, const Vec3P& wdot_PB_B)
371  {
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;
376 
377  const Mat33P E( 0 , s2oc1 , c2oc1 ,
378  0 , c2 , -s2 ,
379  1 , s1*s2oc1 , s1*c2oc1 );
380  const Vec3P qdot = E * w_PB_B;
381 
382  const RealP t = qdot[1]*s1oc1;
383  const RealP a = t*s2oc1 + qdot[2]*c2oc1; // d/dt s2oc1
384  const RealP b = t*c2oc1 - qdot[2]*s2oc1; // d/dt c2oc1
385 
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 );
389 
390  return E*wdot_PB_B + Edot*w_PB_B;
391  }
392 
405  // Note: q[0] is not referenced so we won't waste time calculating
406  // its cosine and sine here.
408  (Vec3P(0, std::cos(q[1]), std::cos(q[2])),
409  Vec3P(0, std::sin(q[1]), std::sin(q[2])));
410  }
411 
417  static Mat33P calcNForBodyXYZInBodyFrame(const Vec3P& cq, const Vec3P& sq) {
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;
422 
423  return Mat33P( c2oc1 , -s2oc1 , 0,
424  s2 , c2 , 0,
425  -s1*c2oc1 , s1*s2oc1, 1 );
426  }
427 
441  // Note: q[2] is not referenced so we won't waste time calculating
442  // its cosine and sine here.
444  (Vec3P(std::cos(q[0]), std::cos(q[1]), 0),
445  Vec3P(std::sin(q[0]), std::sin(q[1]), 0));
446  }
447 
454  static Mat33P calcNForBodyXYZInParentFrame(const Vec3P& cq, const Vec3P& sq) {
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;
459 
460  return Mat33P( 1 , s1*s0oc1 , -s1*c0oc1,
461  0 , c0 , s0,
462  0 , -s0oc1 , c0oc1 );
463  }
464 
471  static Vec3P multiplyByBodyXYZ_N_P(const Vec2P& cosxy,
472  const Vec2P& sinxy,
473  RealP oocosy,
474  const Vec3P& w_PB)
475  {
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];
479 
480  const RealP t = (s0*w1-c0*w2)*oocosy;
481  return Vec3P( w0 + t*s1, c0*w1 + s0*w2, -t ); // qdot
482  }
483 
487  static Vec3P multiplyByBodyXYZ_NT_P(const Vec2P& cosxy,
488  const Vec2P& sinxy,
489  RealP oocosy,
490  const Vec3P& q)
491  {
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];
495 
496  const RealP t = (q0*s1-q2) * oocosy;
497  return Vec3P( q0, c0*q1 + t*s0, s0*q1 - t*c0 ); // v_P
498  }
499 
507  (const Vec2P& cosxy,
508  const Vec2P& sinxy,
509  RealP oocosy,
510  const Vec3P& w_PB)
511  {
512  return multiplyByBodyXYZ_N_P(cosxy,sinxy,oocosy,w_PB);
513  }
514 
527  (const Vec2P& cosxy,
528  const Vec2P& sinxy,
529  RealP oocosy,
530  const Vec3P& qdot,
531  const Vec3P& b_PB)
532  {
533  const RealP s1 = sinxy[1], c1 = cosxy[1];
534  const RealP q0 = qdot[0], q1 = qdot[1], q2 = qdot[2];
535 
536  // 10 flops
537  const Vec3P Nb = multiplyByBodyXYZ_N_P(cosxy,sinxy,oocosy,b_PB);
538 
539  const RealP q1oc1 = q1*oocosy;
540  const Vec3P NDotw((q0*s1-q2)*q1oc1, // NDot_P*w_PB
541  q0*q2*c1, // = NDot_P*(NInv_P*qdot)
542  (q2*s1-q0)*q1oc1 ); // (9 flops)
543 
544  return Nb + NDotw; // 3 flops
545  }
546 
547 
550  static Vec3P multiplyByBodyXYZ_NInv_P(const Vec2P& cosxy,
551  const Vec2P& sinxy,
552  const Vec3P& qdot)
553  {
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;
558 
559  return Vec3P( q0 + s1*q2, // w_PB
560  c0*q1 - s0*c1q2,
561  s0*q1 + c0*c1q2 );
562  }
563 
566  static Vec3P multiplyByBodyXYZ_NInvT_P(const Vec2P& cosxy,
567  const Vec2P& sinxy,
568  const Vec3P& v_P)
569  {
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];
573 
574  return Vec3P( w0, // qdot-like
575  c0*w1 + s0*w2,
576  s1*w0 - s0*c1*w1 + c0*c1*w2);
577  }
578 
588  static Mat33P calcNDotForBodyXYZInBodyFrame
589  (const Vec3P& q, const Vec3P& qdot) {
590  // Note: q[0] is not referenced so we won't waste time calculating
591  // its cosine and sine here.
593  (Vec3P(0, std::cos(q[1]), std::cos(q[2])),
594  Vec3P(0, std::sin(q[1]), std::sin(q[2])),
595  qdot);
596  }
597 
603  static Mat33P calcNDotForBodyXYZInBodyFrame
604  (const Vec3P& cq, const Vec3P& sq, const Vec3P& qdot)
605  {
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;
610 
611  const RealP t = qdot[1]*s1*ooc1;
612  const RealP a = t*s2oc1 + qdot[2]*c2oc1; // d/dt s2oc1
613  const RealP b = t*c2oc1 - qdot[2]*s2oc1; // d/dt c2oc1
614 
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 );
618  }
619 
629  static Mat33P calcNDotForBodyXYZInParentFrame
630  (const Vec3P& q, const Vec3P& qdot) {
631  // Note: q[2] is not referenced so we won't waste time calculating
632  // its cosine and sine here.
633  const RealP cy = std::cos(q[1]); // cos(y)
635  (Vec2P(std::cos(q[0]), cy),
636  Vec2P(std::sin(q[0]), std::sin(q[1])),
637  1/cy, qdot);
638  }
639 
644  static Mat33P calcNDotForBodyXYZInParentFrame
645  (const Vec2P& cq, const Vec2P& sq, RealP ooc1, const Vec3P& qdot) {
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;
649 
650  const RealP t = qdot[1]*s1*ooc1;
651  const RealP a = t*s0oc1 + qdot[0]*c0oc1; // d/dt s0oc1
652  const RealP b = t*c0oc1 - qdot[0]*s0oc1; // d/dt c0oc1
653 
654  return Mat33P( 0, s1*a + qdot[1]*s0, -(s1*b + qdot[1]*c0),
655  0, -qdot[0]*s0 , qdot[0]*c0 ,
656  0, -a , b );
657  }
658 
668  // Note: q[0] is not referenced so we won't waste time calculating
669  // its cosine and sine here.
671  (Vec3P(0, std::cos(q[1]), std::cos(q[2])),
672  Vec3P(0, std::sin(q[1]), std::sin(q[2])));
673  }
674 
680  static Mat33P calcNInvForBodyXYZInBodyFrame
681  (const Vec3P& cq, const Vec3P& sq) {
682  const RealP s1 = sq[1], c1 = cq[1];
683  const RealP s2 = sq[2], c2 = cq[2];
684 
685  return Mat33P( c1*c2 , s2 , 0 ,
686  -c1*s2 , c2 , 0 ,
687  s1 , 0 , 1 );
688  }
689 
698  // Note: q[0] is not referenced so we won't waste time calculating
699  // its cosine and sine here.
701  (Vec3P(std::cos(q[0]), std::cos(q[1]), 0),
702  Vec3P(std::sin(q[0]), std::sin(q[1]), 0));
703  }
704 
710  static Mat33P calcNInvForBodyXYZInParentFrame
711  (const Vec3P& cq, const Vec3P& sq) {
712  const RealP s0 = sq[0], c0 = cq[0];
713  const RealP s1 = sq[1], c1 = cq[1];
714 
715  return Mat33P( 1 , 0 , s1 ,
716  0 , c0 , -s0*c1 ,
717  0 , s0 , c0*c1 );
718  }
719 
729  (const Vec3P& q, const Vec3P& w_PB_B) {
731  (Vec3P(0, std::cos(q[1]), std::cos(q[2])),
732  Vec3P(0, std::sin(q[1]), std::sin(q[2])),
733  w_PB_B);
734  }
735 
741  //TODO: reimplement
743  (const Vec3P& cq, const Vec3P& sq, const Vec3P& w_PB_B)
744  { return calcNForBodyXYZInBodyFrame(cq,sq)*w_PB_B; }
745 
752  (const Vec3P& q, const Vec3P& qdot) {
754  (Vec3P(0, std::cos(q[1]), std::cos(q[2])),
755  Vec3P(0, std::sin(q[1]), std::sin(q[2])),
756  qdot);
757  }
758 
764  // TODO: reimplement
766  (const Vec3P& cq, const Vec3P& sq, const Vec3P& qdot)
767  { return calcNInvForBodyXYZInBodyFrame(cq,sq)*qdot; }
768 
775  (const Vec3P& q, const Vec3P& w_PB_B, const Vec3P& wdot_PB_B)
776  {
777  // Note: q[0] is not referenced so we won't waste time calculating
778  // its cosine and sine here.
780  (Vec3P(0, std::cos(q[1]), std::cos(q[2])),
781  Vec3P(0, std::sin(q[1]), std::sin(q[2])),
782  w_PB_B, wdot_PB_B);
783  }
784 
791  // TODO: reimplement
793  (const Vec3P& cq, const Vec3P& sq, const Vec3P& w_PB_B, const Vec3P& wdot_PB_B)
794  {
795  const Mat33P N = calcNForBodyXYZInBodyFrame(cq,sq);
796  const Vec3P qdot = N * w_PB_B; // 15 flops
797  const Mat33P NDot = calcNDotForBodyXYZInBodyFrame(cq,sq,qdot);
798 
799  return N*wdot_PB_B + NDot*w_PB_B; // 33 flops
800  }
801 
808  const Vec4P e = q/2;
809  const RealP ne1 = -e[1], ne2 = -e[2], ne3 = -e[3];
810  return Mat43P( ne1, ne2, ne3,
811  e[0], e[3], ne2,
812  ne3, e[0], e[1],
813  e[2], ne1, e[0]);
814  }
815 
823  const Vec4P ed = qdot/2;
824  const RealP ned1 = -ed[1], ned2 = -ed[2], ned3 = -ed[3];
825  return Mat43P( ned1, ned2, ned3,
826  ed[0], ed[3], ned2,
827  ned3, ed[0], ed[1],
828  ed[2], ned1, ed[0]);
829  }
830 
840  const Vec4P e = 2*q;
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]);
845  }
846 
847 
852  static Vec4P convertAngVelToQuaternionDot(const Vec4P& q, const Vec3P& w_PB_P) {
853  return calcUnnormalizedNForQuaternion(q)*w_PB_P;
854  }
855 
859  static Vec3P convertQuaternionDotToAngVel(const Vec4P& q, const Vec4P& qdot) {
860  return calcUnnormalizedNInvForQuaternion(q)*qdot;
861  }
862 
868  // TODO: reimplement
870  (const Vec4P& q, const Vec3P& w_PB_P, const Vec3P& wdot_PB_P)
871  {
873  const Vec4P qdot = N*w_PB_P; // 20 flops
874  const Mat43P NDot = calcUnnormalizedNDotForQuaternion(qdot);
875 
876  return N*wdot_PB_P + NDot*w_PB_P; // 44 flops
877  }
878 
886  (const Vec4P& q, const Vec3P& w_PB, const Vec3P& b_PB)
887  {
888  const Mat43P N = calcUnnormalizedNForQuaternion(q); // 7 flops
889  const Vec4P Nb = N*b_PB; // 20 flops
890  const Vec4P NDotw = RealP(-.25)*w_PB.normSqr()*q; // 10 flops
891  return Nb + NDotw; // 4 flops
892  }
893 
894 
895 private:
896  // This is only for the most trustworthy of callers, that is, methods of
897  // the Rotation_ class. There are a lot of ways for this NOT to be a
898  // legitimate rotation matrix -- be careful!!
899  // Note that these are supplied in rows.
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 ) {}
904 
905  // These next methods are highly-efficient power-user methods. Read the
906  // code to understand them.
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 );
910 
911  // These next methods are highly-efficient power-user methods to convert
912  // Rotation matrices to orientation angles. Read the code to understand them.
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;
916 
917 //------------------------------------------------------------------------------
918 // These are obsolete names from a previous release, listed here so that
919 // users will get a decipherable compilation error. (sherm 091101)
920 //------------------------------------------------------------------------------
921 private:
922  // REPLACED BY: calcNForBodyXYZInBodyFrame()
923  static Mat33P calcQBlockForBodyXYZInBodyFrame(const Vec3P& a)
924  { return calcNForBodyXYZInBodyFrame(a); }
925  // REPLACED BY: calcNInvForBodyXYZInBodyFrame()
926  static Mat33P calcQInvBlockForBodyXYZInBodyFrame(const Vec3P& a)
927  { return calcNInvForBodyXYZInBodyFrame(a); }
928  // REPLACED BY: calcUnnormalizedNForQuaternion()
929  static Mat<4,3,P> calcUnnormalizedQBlockForQuaternion(const Vec4P& q)
930  { return calcUnnormalizedNForQuaternion(q); }
931  // REPLACED BY: calcUnnormalizedNInvForQuaternion()
932  static Mat<3,4,P> calcUnnormalizedQInvBlockForQuaternion(const Vec4P& q)
933  { return calcUnnormalizedNInvForQuaternion(q); }
934  // REPLACED BY: convertAngVelInBodyFrameToBodyXYZDot
935  static Vec3P convertAngVelToBodyFixed123Dot(const Vec3P& q, const Vec3P& w_PB_B)
936  { return convertAngVelInBodyFrameToBodyXYZDot(q,w_PB_B); }
937  // REPLACED BY: convertBodyXYZDotToAngVelInBodyFrame
938  static Vec3P convertBodyFixed123DotToAngVel(const Vec3P& q, const Vec3P& qdot)
939  { return convertBodyXYZDotToAngVelInBodyFrame(q,qdot); }
940  // REPLACED BY: convertAngVelDotInBodyFrameToBodyXYZDotDot
941  static Vec3P convertAngVelDotToBodyFixed123DotDot
942  (const Vec3P& q, const Vec3P& w_PB_B, const Vec3P& wdot_PB_B)
943  { return convertAngVelDotInBodyFrameToBodyXYZDotDot(q,w_PB_B,wdot_PB_B); }
944 
945 //------------------------------------------------------------------------------
946 // The following code is obsolete - it is here temporarily for backward
947 // compatibility (Mitiguy 9/5/2007)
948 //------------------------------------------------------------------------------
949 private:
950  // These static methods are like constructors with friendlier names.
951  static Rotation_ zero() { return Rotation_(); }
952  static Rotation_ NaN() { Rotation_ r; r.setRotationToNaN(); return r; }
953 
955  Rotation_& setToZero() { return setRotationToIdentityMatrix(); }
956  Rotation_& setToIdentityMatrix() { return setRotationToIdentityMatrix(); }
957  Rotation_& setToNaN() { return setRotationToNaN(); }
958  static Rotation_ trustMe( const Mat33P& m ) { return Rotation_(m,true); }
959 
960  // One-angle rotations.
961  static Rotation_ aboutX( const RealP& angleInRad ) { return Rotation_( angleInRad, XAxis ); }
962  static Rotation_ aboutY( const RealP& angleInRad ) { return Rotation_( angleInRad, YAxis ); }
963  static Rotation_ aboutZ( const RealP& angleInRad ) { return Rotation_( angleInRad, ZAxis ); }
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); }
966  void setToRotationAboutZ( const RealP& q ) { setRotationFromAngleAboutZ( q ); }
967 
968  // Two-angle space-fixed rotations.
969  static Rotation_ aboutXThenOldY(const RealP& xInRad, const RealP& yInRad) { return Rotation_( SpaceRotationSequence, xInRad, XAxis, yInRad, YAxis ); }
970  static Rotation_ aboutYThenOldX(const RealP& yInRad, const RealP& xInRad) { return Rotation_( SpaceRotationSequence, yInRad, YAxis, xInRad, XAxis ); }
971  static Rotation_ aboutXThenOldZ(const RealP& xInRad, const RealP& zInRad) { return Rotation_( SpaceRotationSequence, xInRad, XAxis, zInRad, ZAxis ); }
972  static Rotation_ aboutZThenOldX(const RealP& zInRad, const RealP& xInRad) { return Rotation_( SpaceRotationSequence, zInRad, ZAxis, xInRad, XAxis ); }
973  static Rotation_ aboutYThenOldZ(const RealP& yInRad, const RealP& zInRad) { return Rotation_( SpaceRotationSequence, yInRad, YAxis, zInRad, ZAxis ); }
974  static Rotation_ aboutZThenOldY(const RealP& zInRad, const RealP& yInRad) { return Rotation_( SpaceRotationSequence, zInRad, ZAxis, yInRad, YAxis ); }
975 
976  // Two-angle body fixed rotations (reversed space-fixed ones).
977  static Rotation_ aboutXThenNewY(const RealP& xInRad, const RealP& yInRad) { return Rotation_( BodyRotationSequence, xInRad, XAxis, yInRad, YAxis ); }
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); }
983 
986  explicit Rotation_( const UnitVec3P& uvecZ ) { setRotationFromOneAxis(uvecZ,ZAxis); }
987 
989  // We will take x seriously after normalizing, but use the y only to create z = normalize(x X y),
990  // then y = z X x. Bad things happen if x and y are aligned but we may not catch it.
991  Rotation_( const Vec3P& x, const Vec3P& yish ) { setRotationFromTwoAxes( UnitVec3P(x), XAxis, yish, YAxis ); }
992 
994  void setToQuaternion( const QuaternionP& q ) { setRotationFromQuaternion(q); }
995 
1000  // Similarly for BodyFixed123.
1001  void setToBodyFixed321( const Vec3P& v) { setRotationFromThreeAnglesThreeAxes( BodyRotationSequence, v[0], ZAxis, v[1], YAxis, v[2], XAxis ); }
1002  void setToBodyFixed123( const Vec3P& v) { setRotationToBodyFixedXYZ(v); }
1003 
1006  Vec4P convertToAngleAxis() const { return convertRotationToAngleAxis(); }
1007 
1009  QuaternionP convertToQuaternion() const { return convertRotationToQuaternion(); }
1010 
1013  void setToSpaceFixed12( const Vec2P& q ) { setRotationFromTwoAnglesTwoAxes( SpaceRotationSequence, q[0], XAxis, q[1], YAxis ); }
1014 
1018  Vec3P convertToBodyFixed123() const { return convertRotationToBodyFixedXYZ(); }
1019  Vec2P convertToBodyFixed12() const { return convertRotationToBodyFixedXY(); }
1020  Vec2P convertToSpaceFixed12() const { return convertTwoAxesRotationToTwoAngles( SpaceRotationSequence, XAxis, YAxis ); }
1021 };
1022 
1023 
1028 template <class P>
1029 class InverseRotation_ : public Mat<3,3,P>::TransposeType {
1030  typedef P RealP;
1031  typedef Rotation_<P> RotationP;
1032  typedef Mat<3,3,P> Mat33P; // not the base type!
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;
1040 public:
1044 
1047 
1053 
1058 
1063  { BaseMat::operator=(R.asMat33()); return *this; }
1064 
1071  SimTK_SimTKCOMMON_EXPORT SymMat33P reexpressSymMat33(const SymMat33P& S_BB) const;
1072 
1074 
1075  const RotationP& invert() const {return *reinterpret_cast<const RotationP*>(this);}
1076  RotationP& updInvert() {return *reinterpret_cast<RotationP*>(this);}
1078 
1081 
1082  const RotationP& transpose() const { return invert(); }
1083  const RotationP& operator~() const { return invert(); }
1085  RotationP& operator~() { return updInvert(); }
1087 
1093 
1094  const RowType& row( int i ) const { return reinterpret_cast<const RowType&>(asMat33()[i]); }
1095  const ColType& col( int j ) const { return reinterpret_cast<const ColType&>(asMat33()(j)); }
1096  const ColType& x() const { return col(0); }
1097  const ColType& y() const { return col(1); }
1098  const ColType& z() const { return col(2); }
1099  const RowType& operator[]( int i ) const { return row(i); }
1100  const ColType& operator()( int j ) const { return col(j); }
1102 
1105 
1106  const BaseMat& asMat33() const { return *static_cast<const BaseMat*>(this); }
1107  BaseMat toMat33() const { return asMat33(); }
1109 };
1110 
1112 template <class P> SimTK_SimTKCOMMON_EXPORT std::ostream&
1113 operator<<(std::ostream&, const Rotation_<P>&);
1115 template <class P> SimTK_SimTKCOMMON_EXPORT std::ostream&
1116 operator<<(std::ostream&, const InverseRotation_<P>&);
1117 
1121 
1122 template <class P, int S> inline UnitVec<P,1>
1123 operator*(const Rotation_<P>& R, const UnitVec<P,S>& v) {return UnitVec<P,1>(R.asMat33()* v.asVec3(), true);}
1124 template <class P, int S> inline UnitRow<P,1>
1125 operator*(const UnitRow<P,S>& r, const Rotation_<P>& R) {return UnitRow<P,1>(r.asRow3() * R.asMat33(), true);}
1126 template <class P, int S> inline UnitVec<P,1>
1127 operator*(const InverseRotation_<P>& R, const UnitVec<P,S>& v) {return UnitVec<P,1>(R.asMat33()* v.asVec3(), true);}
1128 template <class P, int S> inline UnitRow<P,1>
1129 operator*(const UnitRow<P,S>& r, const InverseRotation_<P>& R) {return UnitRow<P,1>(r.asRow3() * R.asMat33(), true);}
1131 
1132 // Couldn't implement these Rotation_ methods until InverseRotation_ was defined.
1133 template <class P> inline
1134 Rotation_<P>::Rotation_(const InverseRotation_<P>& R) : Mat<3,3,P>( R.asMat33() ) {}
1135 template <class P> inline Rotation_<P>&
1136 Rotation_<P>::operator=(const InverseRotation_<P>& R) {static_cast<Mat<3,3,P>&>(*this) = R.asMat33(); return *this;}
1137 template <class P> inline Rotation_<P>&
1138 Rotation_<P>::operator*=(const Rotation_<P>& R) {static_cast<Mat<3,3,P>&>(*this) *= R.asMat33(); return *this;}
1139 template <class P> inline Rotation_<P>&
1140 Rotation_<P>::operator/=(const Rotation_<P>& R) {static_cast<Mat<3,3,P>&>(*this) *= (~R).asMat33(); return *this;}
1141 template <class P> inline Rotation_<P>&
1142 Rotation_<P>::operator*=(const InverseRotation_<P>& R) {static_cast<Mat<3,3,P>&>(*this) *= R.asMat33(); return *this;}
1143 template <class P> inline Rotation_<P>&
1144 Rotation_<P>::operator/=(const InverseRotation_<P>& R) {static_cast<Mat<3,3,P>&>(*this) *= (~R).asMat33(); return *this;}
1145 
1147 
1148 template <class P> inline Rotation_<P>
1149 operator*(const Rotation_<P>& R1, const Rotation_<P>& R2) {return Rotation_<P>(R1) *= R2;}
1150 template <class P> inline Rotation_<P>
1151 operator*(const Rotation_<P>& R1, const InverseRotation_<P>& R2) {return Rotation_<P>(R1) *= R2;}
1152 template <class P> inline Rotation_<P>
1153 operator*(const InverseRotation_<P>& R1, const Rotation_<P>& R2) {return Rotation_<P>(R1) *= R2;}
1154 template <class P> inline Rotation_<P>
1155 operator*(const InverseRotation_<P>& R1, const InverseRotation_<P>& R2) {return Rotation_<P>(R1) *= R2;}
1157 
1160 
1161 template <class P> inline Rotation_<P>
1162 operator/( const Rotation_<P>& R1, const Rotation_<P>& R2 ) {return Rotation_<P>(R1) /= R2;}
1163 template <class P> inline Rotation_<P>
1164 operator/( const Rotation_<P>& R1, const InverseRotation& R2 ) {return Rotation_<P>(R1) /= R2;}
1165 template <class P> inline Rotation_<P>
1166 operator/( const InverseRotation_<P>& R1, const Rotation_<P>& R2 ) {return Rotation_<P>(R1) /= R2;}
1167 template <class P> inline Rotation_<P>
1168 operator/( const InverseRotation_<P>& R1, const InverseRotation_<P>& R2 ) {return Rotation_<P>(R1) /= R2;}
1170 
1171 
1172 //------------------------------------------------------------------------------
1173 } // End of namespace SimTK
1174 
1175 //--------------------------------------------------------------------------
1176 #endif // SimTK_SimTKCOMMON_ROTATION_H_
1177 //--------------------------------------------------------------------------
1178 
1179 
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
Definition: negator.h:64
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 ...