Simbody  3.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MassProperties.h
Go to the documentation of this file.
1 #ifndef SimTK_SIMMATRIX_MASS_PROPERTIES_H_
2 #define SimTK_SIMMATRIX_MASS_PROPERTIES_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: Michael Sherman *
14  * Contributors: Paul Mitiguy *
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 
33 #include "SimTKcommon/Scalar.h"
36 
37 #include <iostream>
38 
39 namespace SimTK {
57 
67 
79 
80 // These are templatized by precision (float or double).
81 template <class P> class UnitInertia_;
82 template <class P> class Inertia_;
83 template <class P> class SpatialInertia_;
84 template <class P> class ArticulatedInertia_;
85 template <class P> class MassProperties_;
86 
87 // The "no trailing underscore" typedefs use whatever the
88 // compile-time precision is set to.
89 
96 
103 
110 
117 
124 
127 
128 // -----------------------------------------------------------------------------
129 // INERTIA MATRIX
130 // -----------------------------------------------------------------------------
192 template <class P>
194  typedef P RealP;
195  typedef Vec<3,P> Vec3P;
196  typedef SymMat<3,P> SymMat33P;
197  typedef Mat<3,3,P> Mat33P;
198  typedef Rotation_<P> RotationP;
199 public:
202 Inertia_() : I_OF_F(NTraits<P>::getNaN()) {}
203 
204 // Default copy constructor, copy assignment, destructor.
205 
212 explicit Inertia_(const RealP& moment) : I_OF_F(moment)
213 { errChk("Inertia::Inertia(moment)"); }
214 
218 Inertia_(const Vec3P& p, const RealP& mass) : I_OF_F(pointMassAt(p,mass)) {}
219 
224 explicit Inertia_(const Vec3P& moments, const Vec3P& products=Vec3P(0))
225 { I_OF_F.updDiag() = moments;
226  I_OF_F.updLower() = products;
227  errChk("Inertia::Inertia(moments,products)"); }
228 
230 Inertia_(const RealP& xx, const RealP& yy, const RealP& zz)
231 { I_OF_F = SymMat33P(xx,
232  0, yy,
233  0, 0, zz);
234  errChk("Inertia::setInertia(xx,yy,zz)"); }
235 
238 Inertia_(const RealP& xx, const RealP& yy, const RealP& zz,
239  const RealP& xy, const RealP& xz, const RealP& yz)
240 { I_OF_F = SymMat33P(xx,
241  xy, yy,
242  xz, yz, zz);
243  errChk("Inertia::setInertia(xx,yy,zz,xy,xz,yz)"); }
244 
247 explicit Inertia_(const SymMat33P& I) : I_OF_F(I)
248 { errChk("Inertia::Inertia(SymMat33)"); }
249 
253 explicit Inertia_(const Mat33P& m)
255  "Inertia(Mat33)", "The supplied matrix was not symmetric.");
256  I_OF_F = SymMat33P(m);
257  errChk("Inertia(Mat33)"); }
258 
259 
262 Inertia_& setInertia(const RealP& xx, const RealP& yy, const RealP& zz) {
263  I_OF_F = RealP(0); I_OF_F(0,0) = xx; I_OF_F(1,1) = yy; I_OF_F(2,2) = zz;
264  errChk("Inertia::setInertia(xx,yy,zz)");
265  return *this;
266 }
267 
270 Inertia_& setInertia(const Vec3P& moments, const Vec3P& products=Vec3P(0)) {
271  I_OF_F.updDiag() = moments;
272  I_OF_F.updLower() = products;
273  errChk("Inertia::setInertia(moments,products)");
274  return *this;
275 }
276 
282 Inertia_& setInertia(const RealP& xx, const RealP& yy, const RealP& zz,
283  const RealP& xy, const RealP& xz, const RealP& yz) {
284  setInertia(Vec3P(xx,yy,zz), Vec3P(xy,xz,yz));
285  errChk("Inertia::setInertia(xx,yy,zz,xy,xz,yz)");
286  return *this;
287 }
288 
289 
292 Inertia_& operator+=(const Inertia_& I)
293 { I_OF_F += I.I_OF_F;
294  errChk("Inertia::operator+=()");
295  return *this; }
296 
299 Inertia_& operator-=(const Inertia_& I)
300 { I_OF_F -= I.I_OF_F;
301  errChk("Inertia::operator-=()");
302  return *this; }
303 
305 Inertia_& operator*=(const P& s) {I_OF_F *= s; return *this;}
306 
309 Inertia_& operator/=(const P& s) {I_OF_F /= s; return *this;}
310 
320 Inertia_ shiftToMassCenter(const Vec3P& CF, const RealP& mass) const
321 { Inertia_ I(*this); I -= pointMassAt(CF, mass);
322  I.errChk("Inertia::shiftToMassCenter()");
323  return I; }
324 
335 Inertia_& shiftToMassCenterInPlace(const Vec3P& CF, const RealP& mass)
336 { (*this) -= pointMassAt(CF, mass);
337  errChk("Inertia::shiftToMassCenterInPlace()");
338  return *this; }
339 
347 Inertia_ shiftFromMassCenter(const Vec3P& p, const RealP& mass) const
348 { Inertia_ I(*this); I += pointMassAt(p, mass);
349  I.errChk("Inertia::shiftFromMassCenter()");
350  return I; }
351 
360 Inertia_& shiftFromMassCenterInPlace(const Vec3P& p, const RealP& mass)
361 { (*this) += pointMassAt(p, mass);
362  errChk("Inertia::shiftFromMassCenterInPlace()");
363  return *this; }
364 
375 Inertia_ reexpress(const Rotation_<P>& R_FB) const
376 { return Inertia_((~R_FB).reexpressSymMat33(I_OF_F)); }
377 
380 Inertia_ reexpress(const InverseRotation_<P>& R_FB) const
381 { return Inertia_((~R_FB).reexpressSymMat33(I_OF_F)); }
382 
387 Inertia_& reexpressInPlace(const Rotation_<P>& R_FB)
388 { I_OF_F = (~R_FB).reexpressSymMat33(I_OF_F); return *this; }
389 
392 Inertia_& reexpressInPlace(const InverseRotation_<P>& R_FB)
393 { I_OF_F = (~R_FB).reexpressSymMat33(I_OF_F); return *this; }
394 
395 RealP trace() const {return I_OF_F.trace();}
396 
398 operator const SymMat33P&() const {return I_OF_F;}
399 
401 const SymMat33P& asSymMat33() const {return I_OF_F;}
402 
405 Mat33P toMat33() const {return Mat33P(I_OF_F);}
406 
409 const Vec3P& getMoments() const {return I_OF_F.getDiag();}
412 const Vec3P& getProducts() const {return I_OF_F.getLower();}
413 
414 bool isNaN() const {return I_OF_F.isNaN();}
415 bool isInf() const {return I_OF_F.isInf();}
416 bool isFinite() const {return I_OF_F.isFinite();}
417 
421 bool isNumericallyEqual(const Inertia_<P>& other) const
422 { return I_OF_F.isNumericallyEqual(other.I_OF_F); }
423 
427 bool isNumericallyEqual(const Inertia_<P>& other, double tol) const
428 { return I_OF_F.isNumericallyEqual(other.I_OF_F, tol); }
429 
432 static bool isValidInertiaMatrix(const SymMat33P& m) {
433  if (m.isNaN()) return false;
434 
435  const Vec3P& d = m.diag();
436  if (!(d >= 0)) return false; // diagonals must be nonnegative
437 
438  const RealP Slop = std::max(d.sum(),RealP(1))
440 
441  if (!(d[0]+d[1]+Slop>=d[2] && d[0]+d[2]+Slop>=d[1] && d[1]+d[2]+Slop>=d[0]))
442  return false; // must satisfy triangle inequality
443 
444  // Thanks to Paul Mitiguy for this condition on products of inertia.
445  const Vec3P& p = m.getLower();
446  if (!( d[0]+Slop>=std::abs(2*p[2])
447  && d[1]+Slop>=std::abs(2*p[1])
448  && d[2]+Slop>=std::abs(2*p[0])))
449  return false; // max products are limited by moments
450 
451  return true;
452 }
453 
456 static Inertia_ pointMassAtOrigin() {return Inertia_(0);}
457 
462 static Inertia_ pointMassAt(const Vec3P& p, const RealP& m) {
463  const Vec3P mp = m*p; // 3 flops
464  const RealP mxx = mp[0]*p[0];
465  const RealP myy = mp[1]*p[1];
466  const RealP mzz = mp[2]*p[2];
467  const RealP nmx = -mp[0];
468  const RealP nmy = -mp[1];
469  return Inertia_( myy+mzz, mxx+mzz, mxx+myy,
470  nmx*p[1], nmx*p[2], nmy*p[2] );
471 }
472 
478 
479 
482 inline static Inertia_ sphere(const RealP& r);
483 
486 inline static Inertia_ cylinderAlongZ(const RealP& r, const RealP& hz);
487 
490 inline static Inertia_ cylinderAlongY(const RealP& r, const RealP& hy);
491 
494 inline static Inertia_ cylinderAlongX(const RealP& r, const RealP& hx);
495 
499 inline static Inertia_ brick(const RealP& hx, const RealP& hy, const RealP& hz);
500 
502 inline static Inertia_ brick(const Vec3P& halfLengths);
503 
505 inline static Inertia_ ellipsoid(const RealP& hx, const RealP& hy, const RealP& hz);
506 
508 inline static Inertia_ ellipsoid(const Vec3P& halfLengths);
509 
511 
512 protected:
513 // Reinterpret this Inertia matrix as a UnitInertia matrix, that is, as the
514 // inertia of something with unit mass. This is useful in implementing
515 // methods of the UnitInertia class in terms of Inertia methods. Be sure you
516 // know that this is a unit-mass inertia!
517 const UnitInertia_<P>& getAsUnitInertia() const
518 { return *static_cast<const UnitInertia_<P>*>(this); }
519 UnitInertia_<P>& updAsUnitInertia()
520 { return *static_cast<UnitInertia_<P>*>(this); }
521 
522 // If error checking is enabled (only in Debug mode), this
523 // method will run some tests on the current contents of this Inertia
524 // matrix and throw an error message if it is not valid. This should be
525 // the same set of tests as run by the isValidInertiaMatrix() method above.
526 void errChk(const char* methodName) const {
527 #ifndef NDEBUG
528  SimTK_ERRCHK(!isNaN(), methodName,
529  "Inertia matrix contains a NaN.");
530 
531  const Vec3P& d = I_OF_F.getDiag(); // moments
532  const Vec3P& p = I_OF_F.getLower(); // products
533  const RealP Ixx = d[0], Iyy = d[1], Izz = d[2];
534  const RealP Ixy = p[0], Ixz = p[1], Iyz = p[2];
535 
536  SimTK_ERRCHK3(d >= -SignificantReal, methodName,
537  "Diagonals of an Inertia matrix must be nonnegative; got %g,%g,%g.",
538  (double)Ixx,(double)Iyy,(double)Izz);
539 
540  // TODO: This is looser than it should be as a workaround for distorted
541  // rotation matrices that were produced by an 11,000 body chain that
542  // Sam Flores encountered.
543  const RealP Slop = std::max(d.sum(),RealP(1))
544  * std::sqrt(NTraits<P>::getEps());
545 
546  SimTK_ERRCHK3( Ixx+Iyy+Slop>=Izz
547  && Ixx+Izz+Slop>=Iyy
548  && Iyy+Izz+Slop>=Ixx,
549  methodName,
550  "Diagonals of an Inertia matrix must satisfy the triangle "
551  "inequality; got %g,%g,%g.",
552  (double)Ixx,(double)Iyy,(double)Izz);
553 
554  // Thanks to Paul Mitiguy for this condition on products of inertia.
555  SimTK_ERRCHK( Ixx+Slop>=std::abs(2*Iyz)
556  && Iyy+Slop>=std::abs(2*Ixz)
557  && Izz+Slop>=std::abs(2*Ixy),
558  methodName,
559  "The magnitude of a product of inertia was too large to be physical.");
560 #endif
561 }
562 
563 // Inertia expressed in frame F and about F's origin OF. Note that frame F
564 // is implicit here; all we actually have are the inertia scalars.
566 };
567 
572 template <class P> inline Inertia_<P>
573 operator+(const Inertia_<P>& l, const Inertia_<P>& r)
574 { return Inertia_<P>(l) += r; }
575 
581 template <class P> inline Inertia_<P>
582 operator-(const Inertia_<P>& l, const Inertia_<P>& r)
583 { return Inertia_<P>(l) -= r; }
584 
587 template <class P> inline Inertia_<P>
588 operator*(const Inertia_<P>& i, const P& r)
589 { return Inertia_<P>(i) *= r; }
590 
593 template <class P> inline Inertia_<P>
594 operator*(const P& r, const Inertia_<P>& i)
595 { return Inertia_<P>(i) *= r; }
596 
597 
601 template <class P> inline Inertia_<P>
602 operator*(const Inertia_<P>& i, int r)
603 { return Inertia_<P>(i) *= P(r); }
604 
608 template <class P> inline Inertia_<P>
609 operator*(int r, const Inertia_<P>& i)
610 { return Inertia_<P>(i) *= P(r); }
611 
615 template <class P> inline Inertia_<P>
616 operator/(const Inertia_<P>& i, const P& r)
617 { return Inertia_<P>(i) /= r; }
618 
622 template <class P> inline Inertia_<P>
623 operator/(const Inertia_<P>& i, int r)
624 { return Inertia_<P>(i) /= P(r); }
625 
629 template <class P> inline Vec<3,P>
630 operator*(const Inertia_<P>& I, const Vec<3,P>& w)
631 { return I.asSymMat33() * w; }
632 
637 template <class P> inline bool
638 operator==(const Inertia_<P>& i1, const Inertia_<P>& i2)
639 { return i1.asSymMat33() == i2.asSymMat33(); }
640 
644 template <class P> inline std::ostream&
645 operator<<(std::ostream& o, const Inertia_<P>& inertia)
646 { return o << inertia.toMat33(); }
647 
648 
649 // -----------------------------------------------------------------------------
650 // UNIT INERTIA MATRIX
651 // -----------------------------------------------------------------------------
672 template <class P>
674  typedef P RealP;
675  typedef Vec<3,P> Vec3P;
676  typedef SymMat<3,P> SymMat33P;
677  typedef Mat<3,3,P> Mat33P;
678  typedef Rotation_<P> RotationP;
679  typedef Inertia_<P> InertiaP;
680 public:
684 
685 // Default copy constructor, copy assignment, destructor.
686 
690 explicit UnitInertia_(const RealP& moment) : InertiaP(moment) {}
691 
696 explicit UnitInertia_(const Vec3P& moments, const Vec3P& products=Vec3P(0))
697 : InertiaP(moments,products) {}
698 
700 UnitInertia_(const RealP& xx, const RealP& yy, const RealP& zz)
701 : InertiaP(xx,yy,zz) {}
702 
705 UnitInertia_(const RealP& xx, const RealP& yy, const RealP& zz,
706  const RealP& xy, const RealP& xz, const RealP& yz)
707 : InertiaP(xx,yy,zz,xy,xz,yz) {}
708 
711 explicit UnitInertia_(const SymMat33P& m) : InertiaP(m) {}
712 
716 explicit UnitInertia_(const Mat33P& m) : InertiaP(m) {}
717 
722 explicit UnitInertia_(const Inertia_<P>& I) : InertiaP(I) {}
723 
727 UnitInertia_& setUnitInertia(const RealP& xx, const RealP& yy, const RealP& zz)
728 { InertiaP::setInertia(xx,yy,zz); return *this; }
729 
732 UnitInertia_& setUnitInertia(const Vec3P& moments, const Vec3P& products=Vec3P(0))
733 { InertiaP::setInertia(moments,products); return *this; }
734 
740 UnitInertia_& setUnitInertia(const RealP& xx, const RealP& yy, const RealP& zz,
741  const RealP& xy, const RealP& xz, const RealP& yz)
742 { InertiaP::setInertia(xx,yy,zz,xy,xz,yz); return *this; }
743 
744 
745 // No +=, -=, etc. operators because those don't result in a UnitInertia
746 // matrix. The parent class ones are suppressed below.
747 
757 UnitInertia_ shiftToCentroid(const Vec3P& CF) const
758 { UnitInertia_ G(*this);
759  G.Inertia_<P>::operator-=(pointMassAt(CF));
760  return G; }
761 
774 UnitInertia_& shiftToCentroidInPlace(const Vec3P& CF)
775 { InertiaP::operator-=(pointMassAt(CF));
776  return *this; }
777 
785 UnitInertia_ shiftFromCentroid(const Vec3P& p) const
786 { UnitInertia_ G(*this);
787  G.Inertia_<P>::operator+=(pointMassAt(p));
788  return G; }
789 
798 UnitInertia_& shiftFromCentroidInPlace(const Vec3P& p)
799 { InertiaP::operator+=(pointMassAt(p));
800  return *this; }
801 
812 UnitInertia_ reexpress(const Rotation_<P>& R_FB) const
813 { return UnitInertia_((~R_FB).reexpressSymMat33(this->I_OF_F)); }
814 
817 UnitInertia_ reexpress(const InverseRotation_<P>& R_FB) const
818 { return UnitInertia_((~R_FB).reexpressSymMat33(this->I_OF_F)); }
819 
824 UnitInertia_& reexpressInPlace(const Rotation_<P>& R_FB)
825 { InertiaP::reexpressInPlace(R_FB); return *this; }
826 
829 UnitInertia_& reexpressInPlace(const InverseRotation_<P>& R_FB)
830 { InertiaP::reexpressInPlace(R_FB); return *this; }
831 
832 
834 operator const SymMat33P&() const {return this->I_OF_F;}
835 
839 const Inertia_<P>& asUnitInertia() const
840 { return *static_cast<const Inertia_<P>*>(this); }
841 
844 UnitInertia_& setFromUnitInertia(const Inertia_<P>& I)
846  return *this; }
847 
851 static bool isValidUnitInertiaMatrix(const SymMat33P& m)
853 
860 
861 
865 
870 static UnitInertia_ pointMassAt(const Vec3P& p)
871 { return UnitInertia_(crossMatSq(p)); }
872 
875 static UnitInertia_ sphere(const RealP& r) {return UnitInertia_(RealP(0.4)*r*r);}
876 
879 static UnitInertia_ cylinderAlongZ(const RealP& r, const RealP& hz) {
880  const RealP Ixx = (r*r)/4 + (hz*hz)/3;
881  return UnitInertia_(Ixx,Ixx,(r*r)/2);
882 }
883 
886 static UnitInertia_ cylinderAlongY(const RealP& r, const RealP& hy) {
887  const RealP Ixx = (r*r)/4 + (hy*hy)/3;
888  return UnitInertia_(Ixx,(r*r)/2,Ixx);
889 }
890 
893 static UnitInertia_ cylinderAlongX(const RealP& r, const RealP& hx) {
894  const RealP Iyy = (r*r)/4 + (hx*hx)/3;
895  return UnitInertia_((r*r)/2,Iyy,Iyy);
896 }
897 
901 static UnitInertia_ brick(const RealP& hx, const RealP& hy, const RealP& hz) {
902  const RealP oo3 = RealP(1)/RealP(3);
903  const RealP hx2=hx*hx, hy2=hy*hy, hz2=hz*hz;
904  return UnitInertia_(oo3*(hy2+hz2), oo3*(hx2+hz2), oo3*(hx2+hy2));
905 }
906 
908 static UnitInertia_ brick(const Vec3P& halfLengths)
909 { return brick(halfLengths[0],halfLengths[1],halfLengths[2]); }
910 
912 static UnitInertia_ ellipsoid(const RealP& hx, const RealP& hy, const RealP& hz) {
913  const RealP hx2=hx*hx, hy2=hy*hy, hz2=hz*hz;
914  return UnitInertia_((hy2+hz2)/5, (hx2+hz2)/5, (hx2+hy2)/5);
915 }
916 
918 static UnitInertia_ ellipsoid(const Vec3P& halfLengths)
919 { return ellipsoid(halfLengths[0],halfLengths[1],halfLengths[2]); }
920 
922 private:
923 // Suppress Inertia_ methods which are not allowed for UnitInertia_.
924 
925 // These kill all flavors of Inertia_::setInertia() and the
926 // Inertia_ assignment ops.
927 void setInertia() {}
928 void operator+=(int) {}
929 void operator-=(int) {}
930 void operator*=(int) {}
931 void operator/=(int) {}
932 };
933 
934 // Implement Inertia methods which are pass-throughs to UnitInertia methods.
935 
936 template <class P> inline Inertia_<P> Inertia_<P>::
937 sphere(const RealP& r)
938 { return UnitInertia_<P>::sphere(r); }
939 template <class P> inline Inertia_<P> Inertia_<P>::
940 cylinderAlongZ(const RealP& r, const RealP& hz)
941 { return UnitInertia_<P>::cylinderAlongZ(r,hz); }
942 template <class P> inline Inertia_<P> Inertia_<P>::
943 cylinderAlongY(const RealP& r, const RealP& hy)
944 { return UnitInertia_<P>::cylinderAlongY(r,hy); }
945 template <class P> inline Inertia_<P> Inertia_<P>::
946 cylinderAlongX(const RealP& r, const RealP& hx)
947 { return UnitInertia_<P>::cylinderAlongX(r,hx); }
948 template <class P> inline Inertia_<P> Inertia_<P>::
949 brick(const RealP& hx, const RealP& hy, const RealP& hz)
950 { return UnitInertia_<P>::brick(hx,hy,hz); }
951 template <class P> inline Inertia_<P> Inertia_<P>::
952 brick(const Vec3P& halfLengths)
953 { return UnitInertia_<P>::brick(halfLengths); }
954 template <class P> inline Inertia_<P> Inertia_<P>::
955 ellipsoid(const RealP& hx, const RealP& hy, const RealP& hz)
956 { return UnitInertia_<P>::ellipsoid(hx,hy,hz); }
957 template <class P> inline Inertia_<P> Inertia_<P>::
958 ellipsoid(const Vec3P& halfLengths)
959 { return UnitInertia_<P>::ellipsoid(halfLengths); }
960 
961 
962 // -----------------------------------------------------------------------------
963 // SPATIAL INERTIA MATRIX
964 // -----------------------------------------------------------------------------
996 template <class P>
998  typedef P RealP;
999  typedef Vec<3,P> Vec3P;
1000  typedef UnitInertia_<P> UnitInertiaP;
1001  typedef Mat<3,3,P> Mat33P;
1002  typedef Vec<2, Vec3P> SpatialVecP;
1003  typedef Mat<2,2,Mat33P> SpatialMatP;
1004  typedef Rotation_<P> RotationP;
1005  typedef Transform_<P> TransformP;
1006  typedef Inertia_<P> InertiaP;
1007 public:
1010 : m(nanP()), p(nanP()) {} // inertia is already NaN
1011 SpatialInertia_(RealP mass, const Vec3P& com, const UnitInertiaP& gyration)
1012 : m(mass), p(com), G(gyration) {}
1013 
1014 // default copy constructor, copy assignment, destructor
1015 
1016 SpatialInertia_& setMass(RealP mass)
1017 { SimTK_ERRCHK1(mass >= 0, "SpatialInertia::setMass()",
1018  "Negative mass %g is illegal.", (double)mass);
1019  m=mass; return *this; }
1020 SpatialInertia_& setMassCenter(const Vec3P& com)
1021 { p=com; return *this;}
1022 SpatialInertia_& setUnitInertia(const UnitInertiaP& gyration)
1023 { G=gyration; return *this; }
1024 
1025 RealP getMass() const {return m;}
1026 const Vec3P& getMassCenter() const {return p;}
1027 const UnitInertiaP& getUnitInertia() const {return G;}
1028 
1031 Vec3P calcMassMoment() const {return m*p;}
1032 
1035 InertiaP calcInertia() const {return m*G;}
1036 
1042  SimTK_ERRCHK(m+src.m != 0, "SpatialInertia::operator+=()",
1043  "The combined mass cannot be zero.");
1044  const RealP mtot = m+src.m, oomtot = 1/mtot; // ~11 flops
1045  p = oomtot*(calcMassMoment() + src.calcMassMoment()); // 10 flops
1046  G.setFromUnitInertia(oomtot*(calcInertia()+src.calcInertia())); // 19 flops
1047  m = mtot; // must do this last
1048  return *this;
1049 }
1050 
1056  SimTK_ERRCHK(m != src.m, "SpatialInertia::operator-=()",
1057  "The combined mass cannot be zero.");
1058  const RealP mtot = m-src.m, oomtot = 1/mtot; // ~11 flops
1059  p = oomtot*(calcMassMoment() - src.calcMassMoment()); // 10 flops
1060  G.setFromUnitInertia(oomtot*(calcInertia()-src.calcInertia())); // 19 flops
1061  m = mtot; // must do this last
1062  return *this;
1063 }
1064 
1067 SpatialInertia_& operator*=(const RealP& s) {m *= s; return *this;}
1068 
1071 SpatialInertia_& operator/=(const RealP& s) {m /= s; return *this;}
1072 
1076 { return m*SpatialVecP(G*v[0]+p%v[1], v[1]-p%v[0]); }
1077 
1083 SpatialInertia_ reexpress(const Rotation_<P>& R_FB) const
1084 { return SpatialInertia_(*this).reexpressInPlace(R_FB); }
1085 
1088 SpatialInertia_ reexpress(const InverseRotation_<P>& R_FB) const
1089 { return SpatialInertia_(*this).reexpressInPlace(R_FB); }
1090 
1095 SpatialInertia_& reexpressInPlace(const Rotation_<P>& R_FB)
1096 { p = (~R_FB)*p; G.reexpressInPlace(R_FB); return *this; }
1097 
1100 SpatialInertia_& reexpressInPlace(const InverseRotation_<P>& R_FB)
1101 { p = (~R_FB)*p; G.reexpressInPlace(R_FB); return *this; }
1102 
1107 SpatialInertia_ shift(const Vec3P& S) const
1108 { return SpatialInertia_(*this).shiftInPlace(S); }
1109 
1115  G.shiftToCentroidInPlace(p); // change to central inertia
1116  G.shiftFromCentroidInPlace(S); // now inertia is about S
1117  p -= S; // was p=com-OF, now want p'=com-(OF+S)=p-S
1118  return *this;
1119 }
1120 
1129 SpatialInertia_ transform(const Transform_<P>& X_FB) const
1130 { return SpatialInertia_(*this).transformInPlace(X_FB); }
1131 
1134 SpatialInertia_ transform(const InverseTransform_<P>& X_FB) const
1135 { return SpatialInertia_(*this).transformInPlace(X_FB); }
1136 
1147  shiftInPlace(X_FB.p()); // shift to the new origin OB.
1148  reexpressInPlace(X_FB.R()); // get everything in B
1149  return *this;
1150 }
1151 
1155  shiftInPlace(X_FB.p()); // shift to the new origin OB.
1156  reexpressInPlace(X_FB.R()); // get everything in B
1157  return *this;
1158 }
1159 
1160 const SpatialMatP toSpatialMat() const {
1161  Mat33P offDiag = crossMat(m*p);
1162  return SpatialMatP(m*G.toMat33(), offDiag,
1163  -offDiag, Mat33P(m));
1164 }
1165 
1166 private:
1167 RealP m; // mass of this rigid body F
1168 Vec3P p; // location of body's COM from OF, expressed in F
1169 UnitInertiaP G; // mass distribution; inertia is mass*gyration
1170 
1171 static P nanP() {return NTraits<P>::getNaN();}
1172 };
1173 
1176 template <class P> inline SpatialInertia_<P>
1178 { return SpatialInertia_<P>(l) += r; }
1179 
1183 template <class P> inline SpatialInertia_<P>
1185 { return SpatialInertia_<P>(l) -= r; }
1186 
1187 
1188 // -----------------------------------------------------------------------------
1189 // ARTICULATED BODY INERTIA MATRIX
1190 // -----------------------------------------------------------------------------
1237 template <class P>
1239  typedef P RealP;
1240  typedef Vec<3,P> Vec3P;
1241  typedef UnitInertia_<P> UnitInertiaP;
1242  typedef Mat<3,3,P> Mat33P;
1243  typedef SymMat<3,P> SymMat33P;
1244  typedef Vec<2, Vec3P> SpatialVecP;
1245  typedef Mat<2,2,Mat33P> SpatialMatP;
1246  typedef Rotation_<P> RotationP;
1247  typedef Transform_<P> TransformP;
1248  typedef Inertia_<P> InertiaP;
1249 public:
1254 ArticulatedInertia_(const SymMat33P& mass, const Mat33P& massMoment, const SymMat33P& inertia)
1255 : M(mass), J(inertia), F(massMoment) {}
1256 
1260 : M(rbi.getMass()), J(rbi.calcInertia()), F(crossMat(rbi.calcMassMoment())) {}
1261 
1263 ArticulatedInertia_& setMass (const SymMat33P& mass) {M=mass; return *this;}
1265 ArticulatedInertia_& setMassMoment(const Mat33P& massMoment) {F=massMoment; return *this;}
1267 ArticulatedInertia_& setInertia (const SymMat33P& inertia) {J=inertia; return *this;}
1268 
1270 const SymMat33P& getMass() const {return M;}
1272 const Mat33P& getMassMoment() const {return F;}
1274 const SymMat33P& getInertia() const {return J;}
1275 
1276 // default destructor, copy constructor, copy assignment
1277 
1281 { M+=src.M; J+=src.J; F+=src.F; return *this; }
1285 { M-=src.M; J-=src.J; F-=src.F; return *this; }
1286 
1289 { return SpatialVecP(J*v[0]+F*v[1], ~F*v[0]+M*v[1]); }
1290 
1292 template <int N>
1294  Mat<2,N,Vec3P> res;
1295  for (int j=0; j < N; ++j)
1296  res.col(j) = (*this) * m.col(j); // above this*SpatialVec operator
1297  return res;
1298 }
1299 
1307 SimTK_SimTKCOMMON_EXPORT ArticulatedInertia_ shift(const Vec3P& s) const;
1308 
1312 
1315 const SpatialMatP toSpatialMat() const {
1316  return SpatialMatP( Mat33P(J), F,
1317  ~F, Mat33P(M) );
1318 }
1319 private:
1320 SymMat33P M;
1321 SymMat33P J;
1322 Mat33P F;
1323 };
1324 
1327 template <class P> inline ArticulatedInertia_<P>
1329 { return ArticulatedInertia_<P>(l) += r; }
1330 
1334 template <class P> inline ArticulatedInertia_<P>
1336 { return ArticulatedInertia_<P>(l) -= r; }
1337 
1338 
1339 // -----------------------------------------------------------------------------
1340 // MASS PROPERTIES
1341 // -----------------------------------------------------------------------------
1354 template <class P>
1356  typedef P RealP;
1357  typedef Vec<3,P> Vec3P;
1358  typedef UnitInertia_<P> UnitInertiaP;
1359  typedef Mat<3,3,P> Mat33P;
1360  typedef Mat<6,6,P> Mat66P;
1361  typedef SymMat<3,P> SymMat33P;
1362  typedef Mat<2,2,Mat33P> SpatialMatP;
1363  typedef Rotation_<P> RotationP;
1364  typedef Transform_<P> TransformP;
1365  typedef Inertia_<P> InertiaP;
1366 public:
1369 MassProperties_() { setMassProperties(0,Vec3P(0),UnitInertiaP()); }
1373 MassProperties_(const RealP& m, const Vec3P& com, const InertiaP& inertia)
1374  { setMassProperties(m,com,inertia); }
1377 MassProperties_(const RealP& m, const Vec3P& com, const UnitInertiaP& gyration)
1378  { setMassProperties(m,com,gyration); }
1379 
1383 MassProperties_& setMassProperties(const RealP& m, const Vec3P& com, const InertiaP& inertia) {
1384  mass = m;
1385  comInB = com;
1386  if (m == 0) {
1387  SimTK_ASSERT(inertia.asSymMat33().getAsVec() == 0, "Mass is zero but inertia is nonzero");
1388  unitInertia_OB_B = UnitInertiaP(0);
1389  }
1390  else {
1391  unitInertia_OB_B = UnitInertiaP(inertia*(1/m));
1392  }
1393  return *this;
1394 }
1395 
1398 MassProperties_& setMassProperties
1399  (const RealP& m, const Vec3P& com, const UnitInertiaP& gyration)
1400 { mass=m; comInB=com; unitInertia_OB_B=gyration; return *this; }
1401 
1403 const RealP& getMass() const {return mass;}
1408 const Vec3P& getMassCenter() const {return comInB;}
1413 const UnitInertiaP& getUnitInertia() const {return unitInertia_OB_B;}
1419 const InertiaP calcInertia() const {return mass*unitInertia_OB_B;}
1424 const InertiaP getInertia() const {return calcInertia();}
1425 
1430  return mass*unitInertia_OB_B - InertiaP(comInB, mass);
1431 }
1436 InertiaP calcShiftedInertia(const Vec3P& newOriginB) const {
1437  return calcCentralInertia() + InertiaP(newOriginB-comInB, mass);
1438 }
1444  return calcShiftedInertia(X_BC.p()).reexpress(X_BC.R());
1445 }
1451 MassProperties_ calcShiftedMassProps(const Vec3P& newOriginB) const {
1452  return MassProperties_(mass, comInB-newOriginB,
1453  calcShiftedInertia(newOriginB));
1454 }
1455 
1466  return MassProperties_(mass, ~X_BC*comInB, calcTransformedInertia(X_BC));
1467 }
1468 
1478  return MassProperties_(mass, ~R_BC*comInB, unitInertia_OB_B.reexpress(R_BC));
1479 }
1480 
1484 bool isExactlyMassless() const { return mass==0; }
1490 bool isNearlyMassless(const RealP& tol=SignificantReal) const {
1491  return mass <= tol;
1492 }
1493 
1497 bool isExactlyCentral() const { return comInB==Vec3P(0); }
1503 bool isNearlyCentral(const RealP& tol=SignificantReal) const {
1504  return comInB.normSqr() <= tol*tol;
1505 }
1506 
1509 bool isNaN() const {return SimTK::isNaN(mass) || comInB.isNaN() || unitInertia_OB_B.isNaN();}
1514 bool isInf() const {
1515  if (isNaN()) return false;
1516  return SimTK::isInf(mass) || comInB.isInf() || unitInertia_OB_B.isInf();
1517 }
1521 bool isFinite() const {
1522  return SimTK::isFinite(mass) && comInB.isFinite() && unitInertia_OB_B.isFinite();
1523 }
1524 
1529  SpatialMatP M;
1530  M(0,0) = mass*unitInertia_OB_B.toMat33();
1531  M(0,1) = mass*crossMat(comInB);
1532  M(1,0) = ~M(0,1);
1533  M(1,1) = mass; // a diagonal matrix
1534  return M;
1535 }
1536 
1542 Mat66P toMat66() const {
1543  Mat66P M;
1544  M.template updSubMat<3,3>(0,0) = mass*unitInertia_OB_B.toMat33();
1545  M.template updSubMat<3,3>(0,3) = mass*crossMat(comInB);
1546  M.template updSubMat<3,3>(3,0) = ~M.template getSubMat<3,3>(0,3);
1547  M.template updSubMat<3,3>(3,3) = mass; // a diagonal matrix
1548  return M;
1549 }
1550 
1551 private:
1552 RealP mass;
1553 Vec3P comInB; // meas. from B origin, expr. in B
1554 UnitInertiaP unitInertia_OB_B; // about B origin, expr. in B
1555 };
1556 
1559 template <class P> static inline std::ostream&
1560 operator<<(std::ostream& o, const MassProperties_<P>& mp) {
1561  return o << "{ mass=" << mp.getMass()
1562  << "\n com=" << mp.getMassCenter()
1563  << "\n Ixx,yy,zz=" << mp.getUnitInertia().getMoments()
1564  << "\n Ixy,xz,yz=" << mp.getUnitInertia().getProducts()
1565  << "\n}\n";
1566 }
1567 
1568 } // namespace SimTK
1569 
1570 #endif // SimTK_SIMMATRIX_MASS_PROPERTIES_H_