1 #ifndef SimTK_SIMMATRIX_MASS_PROPERTIES_H_
2 #define SimTK_SIMMATRIX_MASS_PROPERTIES_H_
212 explicit Inertia_(
const RealP& moment) : I_OF_F(moment)
213 { errChk(
"Inertia::Inertia(moment)"); }
218 Inertia_(
const Vec3P& p,
const RealP& mass) : I_OF_F(pointMassAt(p,mass)) {}
225 { I_OF_F.updDiag() = moments;
226 I_OF_F.updLower() = products;
227 errChk(
"Inertia::Inertia(moments,products)"); }
230 Inertia_(
const RealP& xx,
const RealP& yy,
const RealP& zz)
234 errChk(
"Inertia::setInertia(xx,yy,zz)"); }
238 Inertia_(
const RealP& xx,
const RealP& yy,
const RealP& zz,
239 const RealP& xy,
const RealP& xz,
const RealP& yz)
243 errChk(
"Inertia::setInertia(xx,yy,zz,xy,xz,yz)"); }
248 { errChk(
"Inertia::Inertia(SymMat33)"); }
255 "Inertia(Mat33)",
"The supplied matrix was not symmetric.");
257 errChk(
"Inertia(Mat33)"); }
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)");
271 I_OF_F.updDiag() = moments;
272 I_OF_F.updLower() = products;
273 errChk(
"Inertia::setInertia(moments,products)");
282 Inertia_& setInertia(
const RealP& xx,
const RealP& yy,
const RealP& zz,
283 const RealP& xy,
const RealP& xz,
const RealP& yz) {
285 errChk(
"Inertia::setInertia(xx,yy,zz,xy,xz,yz)");
294 errChk(
"Inertia::operator+=()");
301 errChk(
"Inertia::operator-=()");
321 {
Inertia_ I(*
this); I -= pointMassAt(CF, mass);
322 I.
errChk(
"Inertia::shiftToMassCenter()");
336 { (*this) -= pointMassAt(CF, mass);
337 errChk(
"Inertia::shiftToMassCenterInPlace()");
348 {
Inertia_ I(*
this); I += pointMassAt(p, mass);
349 I.
errChk(
"Inertia::shiftFromMassCenter()");
361 { (*this) += pointMassAt(p, mass);
362 errChk(
"Inertia::shiftFromMassCenterInPlace()");
376 {
return Inertia_((~R_FB).reexpressSymMat33(I_OF_F)); }
381 {
return Inertia_((~R_FB).reexpressSymMat33(I_OF_F)); }
388 { I_OF_F = (~R_FB).reexpressSymMat33(I_OF_F);
return *
this; }
393 { I_OF_F = (~R_FB).reexpressSymMat33(I_OF_F);
return *
this; }
395 RealP
trace()
const {
return I_OF_F.trace();}
414 bool isNaN()
const {
return I_OF_F.isNaN();}
415 bool isInf()
const {
return I_OF_F.isInf();}
422 {
return I_OF_F.isNumericallyEqual(other.
I_OF_F); }
428 {
return I_OF_F.isNumericallyEqual(other.
I_OF_F, tol); }
433 if (m.
isNaN())
return false;
436 if (!(d >= 0))
return false;
441 if (!(d[0]+d[1]+Slop>=d[2] && d[0]+d[2]+Slop>=d[1] && d[1]+d[2]+Slop>=d[0]))
463 const Vec3P mp = m*p;
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] );
482 inline static Inertia_ sphere(
const RealP& r);
486 inline static Inertia_ cylinderAlongZ(
const RealP& r,
const RealP& hz);
490 inline static Inertia_ cylinderAlongY(
const RealP& r,
const RealP& hy);
494 inline static Inertia_ cylinderAlongX(
const RealP& r,
const RealP& hx);
499 inline static Inertia_ brick(
const RealP& hx,
const RealP& hy,
const RealP& hz);
502 inline static Inertia_ brick(
const Vec3P& halfLengths);
505 inline static Inertia_ ellipsoid(
const RealP& hx,
const RealP& hy,
const RealP& hz);
508 inline static Inertia_ ellipsoid(
const Vec3P& halfLengths);
526 void errChk(
const char* methodName)
const {
529 "Inertia matrix contains a NaN.");
531 const Vec3P& d = I_OF_F.getDiag();
532 const Vec3P& p = I_OF_F.getLower();
533 const RealP Ixx = d[0], Iyy = d[1], Izz = d[2];
534 const RealP Ixy = p[0], Ixz = p[1], Iyz = p[2];
537 "Diagonals of an Inertia matrix must be nonnegative; got %g,%g,%g.",
538 (
double)Ixx,(
double)Iyy,(
double)Izz);
548 && Iyy+Izz+Slop>=Ixx,
550 "Diagonals of an Inertia matrix must satisfy the triangle "
551 "inequality; got %g,%g,%g.",
552 (
double)Ixx,(
double)Iyy,(
double)Izz);
559 "The magnitude of a product of inertia was too large to be physical.");
637 template <
class P>
inline bool
644 template <
class P>
inline std::ostream&
645 operator<<(std::ostream& o, const Inertia_<P>& inertia)
646 {
return o << inertia.toMat33(); }
706 const RealP& xy,
const RealP& xz,
const RealP& yz)
727 UnitInertia_& setUnitInertia(
const RealP& xx,
const RealP& yy,
const RealP& zz)
728 { InertiaP::setInertia(xx,yy,zz);
return *
this; }
733 { InertiaP::setInertia(moments,products);
return *
this; }
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; }
759 G.
Inertia_<P>::operator-=(pointMassAt(CF));
775 { InertiaP::operator-=(pointMassAt(CF));
787 G.
Inertia_<P>::operator+=(pointMassAt(p));
799 { InertiaP::operator+=(pointMassAt(p));
813 {
return UnitInertia_((~R_FB).reexpressSymMat33(this->I_OF_F)); }
818 {
return UnitInertia_((~R_FB).reexpressSymMat33(this->I_OF_F)); }
825 { InertiaP::reexpressInPlace(R_FB);
return *
this; }
830 { InertiaP::reexpressInPlace(R_FB);
return *
this; }
834 operator const SymMat33P&()
const {
return this->I_OF_F;}
851 static bool isValidUnitInertiaMatrix(
const SymMat33P& m)
880 const RealP Ixx = (r*r)/4 + (hz*hz)/3;
887 const RealP Ixx = (r*r)/4 + (hy*hy)/3;
894 const RealP Iyy = (r*r)/4 + (hx*hx)/3;
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));
909 {
return brick(halfLengths[0],halfLengths[1],halfLengths[2]); }
913 const RealP hx2=hx*hx, hy2=hy*hy, hz2=hz*hz;
914 return UnitInertia_((hy2+hz2)/5, (hx2+hz2)/5, (hx2+hy2)/5);
919 {
return ellipsoid(halfLengths[0],halfLengths[1],halfLengths[2]); }
928 void operator+=(
int) {}
929 void operator-=(
int) {}
930 void operator*=(
int) {}
931 void operator/=(
int) {}
949 brick(
const RealP& hx,
const RealP& hy,
const RealP& hz)
955 ellipsoid(
const RealP& hx,
const RealP& hy,
const RealP& hz)
1010 : m(nanP()), p(nanP()) {}
1012 : m(mass), p(com), G(gyration) {}
1018 "Negative mass %g is illegal.", (
double)mass);
1019 m=mass;
return *
this; }
1021 { p=com;
return *
this;}
1023 { G=gyration;
return *
this; }
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;
1046 G.setFromUnitInertia(oomtot*(calcInertia()+src.
calcInertia()));
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;
1060 G.setFromUnitInertia(oomtot*(calcInertia()-src.
calcInertia()));
1076 {
return m*
SpatialVecP(G*v[0]+p%v[1], v[1]-p%v[0]); }
1115 G.shiftToCentroidInPlace(p);
1116 G.shiftFromCentroidInPlace(S);
1147 shiftInPlace(X_FB.
p());
1148 reexpressInPlace(X_FB.
R());
1155 shiftInPlace(X_FB.
p());
1156 reexpressInPlace(X_FB.
R());
1176 template <
class P>
inline SpatialInertia_<P>
1255 : M(mass), J(inertia), F(massMoment) {}
1260 : M(rbi.
getMass()), J(rbi.calcInertia()), F(
crossMat(rbi.calcMassMoment())) {}
1281 { M+=src.M; J+=src.J; F+=src.F;
return *
this; }
1285 { M-=src.M; J-=src.J; F-=src.F;
return *
this; }
1289 {
return SpatialVecP(J*v[0]+F*v[1], ~F*v[0]+M*v[1]); }
1295 for (
int j=0; j < N; ++j)
1296 res.
col(j) = (*this) * m.
col(j);
1327 template <
class P>
inline ArticulatedInertia_<P>
1374 { setMassProperties(m,com,inertia); }
1378 { setMassProperties(m,com,gyration); }
1400 { mass=m; comInB=com; unitInertia_OB_B=gyration;
return *
this; }
1430 return mass*unitInertia_OB_B -
InertiaP(comInB, mass);
1437 return calcCentralInertia() +
InertiaP(newOriginB-comInB, mass);
1444 return calcShiftedInertia(X_BC.
p()).reexpress(X_BC.
R());
1453 calcShiftedInertia(newOriginB));
1466 return MassProperties_(mass, ~X_BC*comInB, calcTransformedInertia(X_BC));
1478 return MassProperties_(mass, ~R_BC*comInB, unitInertia_OB_B.reexpress(R_BC));
1504 return comInB.normSqr() <= tol*tol;
1515 if (
isNaN())
return false;
1516 return SimTK::isInf(mass) || comInB.isInf() || unitInertia_OB_B.isInf();
1522 return SimTK::isFinite(mass) && comInB.isFinite() && unitInertia_OB_B.isFinite();
1530 M(0,0) = mass*unitInertia_OB_B.toMat33();
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;
1554 UnitInertiaP unitInertia_OB_B;
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()
1570 #endif // SimTK_SIMMATRIX_MASS_PROPERTIES_H_