Simbody
|
The physical meaning of an inertia is the distribution of a rigid body's mass about a particular point. More...
#include <MassProperties.h>
Public Member Functions | |
Inertia_ () | |
Default is a NaN-ed out mess to avoid accidents, even in Release mode. | |
Inertia_ (const RealP &moment) | |
Create a principal inertia matrix with identical diagonal elements, like a sphere where moment=2/5 m r^2, or a cube where moment=1/6 m s^2, with m the total mass, r the sphere's radius and s the length of a side of the cube. | |
Inertia_ (const Vec3P &p, const RealP &mass) | |
Create an Inertia matrix for a point mass at a given location, measured from the origin OF of the implicit frame F, and expressed in F. | |
Inertia_ (const Vec3P &moments, const Vec3P &products=Vec3P(0)) | |
Create an inertia matrix from a vector of the moments of inertia (the inertia matrix diagonal) and optionally a vector of the products of inertia (the off-diagonals). | |
Inertia_ (const RealP &xx, const RealP &yy, const RealP &zz) | |
Create a principal inertia matrix (only non-zero on diagonal). | |
Inertia_ (const RealP &xx, const RealP &yy, const RealP &zz, const RealP &xy, const RealP &xz, const RealP &yz) | |
This is a general inertia matrix. | |
Inertia_ (const SymMat33P &I) | |
Construct an Inertia from a symmetric 3x3 matrix. | |
Inertia_ (const Mat33P &m) | |
Construct an Inertia matrix from a 3x3 symmetric matrix. | |
Inertia_ & | setInertia (const RealP &xx, const RealP &yy, const RealP &zz) |
Set an inertia matrix to have only principal moments (that is, it will be diagonal). | |
Inertia_ & | setInertia (const Vec3P &moments, const Vec3P &products=Vec3P(0)) |
Set principal moments and optionally off-diagonal terms. | |
Inertia_ & | setInertia (const RealP &xx, const RealP &yy, const RealP &zz, const RealP &xy, const RealP &xz, const RealP &yz) |
Set this Inertia to a general matrix. | |
Inertia_ & | operator+= (const Inertia_ &I) |
Add in another inertia matrix. | |
Inertia_ & | operator-= (const Inertia_ &I) |
Subtract off another inertia matrix. | |
Inertia_ & | operator*= (const P &s) |
Multiply this inertia matrix by a scalar. Cost is 6 flops. | |
Inertia_ & | operator/= (const P &s) |
Divide this inertia matrix by a scalar. | |
Inertia_ | shiftToMassCenter (const Vec3P &CF, const RealP &mass) const |
Assume that the current inertia is about the F frame's origin OF, and expressed in F. | |
Inertia_ & | shiftToMassCenterInPlace (const Vec3P &CF, const RealP &mass) |
Assume that the current inertia is about the F frame's origin OF, and expressed in F. | |
Inertia_ | shiftFromMassCenter (const Vec3P &p, const RealP &mass) const |
Assuming that the current inertia I is a central inertia (that is, it is inertia about the body center of mass CF), shift it to some other point p measured from the center of mass. | |
Inertia_ & | shiftFromMassCenterInPlace (const Vec3P &p, const RealP &mass) |
Assuming that the current inertia I is a central inertia (that is, it is inertia about the body center of mass CF), shift it to some other point p measured from the center of mass. | |
Inertia_ | reexpress (const Rotation_< P > &R_FB) const |
Return a new inertia matrix like this one but re-expressed in another frame (leaving the origin point unchanged). | |
Inertia_ | reexpress (const InverseRotation_< P > &R_FB) const |
Rexpress using an inverse rotation to avoid having to convert it. | |
Inertia_ & | reexpressInPlace (const Rotation_< P > &R_FB) |
Re-express this inertia matrix in another frame, changing the object in place; see reexpress() if you want to leave this object unmolested and get a new one instead. | |
Inertia_ & | reexpressInPlace (const InverseRotation_< P > &R_FB) |
Rexpress in place using an inverse rotation to avoid having to convert it. | |
RealP | trace () const |
operator const SymMat33P & () const | |
This is an implicit conversion to a const SymMat33. | |
const SymMat33P & | asSymMat33 () const |
Obtain a reference to the underlying symmetric matrix type. | |
Mat33P | toMat33 () const |
Expand the internal packed representation into a full 3x3 symmetric matrix with all elements set. | |
const Vec3P & | getMoments () const |
Obtain the inertia moments (diagonal of the Inertia matrix) as a Vec3. | |
const Vec3P & | getProducts () const |
Obtain the inertia products (off-diagonals of the Inertia matrix) as a Vec3 with elements ordered xx, xy, yz. | |
bool | isNaN () const |
bool | isInf () const |
bool | isFinite () const |
bool | isNumericallyEqual (const Inertia_< P > &other) const |
Compare this inertia matrix with another one and return true if they are close to within a default numerical tolerance. | |
bool | isNumericallyEqual (const Inertia_< P > &other, double tol) const |
Compare this inertia matrix with another one and return true if they are close to within a specified numerical tolerance. | |
Static Public Member Functions | |
static bool | isValidInertiaMatrix (const SymMat33P &m) |
Test some conditions that must hold for a valid Inertia matrix. | |
static Inertia_ | pointMassAtOrigin () |
Create an Inertia matrix for a point located at the origin -- that is, an all-zero matrix. | |
static Inertia_ | pointMassAt (const Vec3P &p, const RealP &m) |
Create an Inertia matrix for a point of a given mass, located at a given location measured from the origin of the implicit F frame. | |
Unit inertia matrix factories | |
These return UnitInertia matrices (inertias of unit-mass objects) converted to Inertias. Multiply the result by the actual mass to get the Inertia of an actual object of this shape. See the UnitInertia class for more information. | |
static Inertia_ | sphere (const RealP &r) |
Create a UnitInertia matrix for a unit mass sphere of radius r centered at the origin. | |
static Inertia_ | cylinderAlongZ (const RealP &r, const RealP &hz) |
Unit-mass cylinder aligned along z axis; use radius and half-length. | |
static Inertia_ | cylinderAlongY (const RealP &r, const RealP &hy) |
Unit-mass cylinder aligned along y axis; use radius and half-length. | |
static Inertia_ | cylinderAlongX (const RealP &r, const RealP &hx) |
Unit-mass cylinder aligned along x axis; use radius and half-length. | |
static Inertia_ | brick (const RealP &hx, const RealP &hy, const RealP &hz) |
Unit-mass brick given by half-lengths in each direction. | |
static Inertia_ | brick (const Vec3P &halfLengths) |
Alternate interface to brick() that takes a Vec3 for the half lengths. | |
static Inertia_ | ellipsoid (const RealP &hx, const RealP &hy, const RealP &hz) |
Unit-mass ellipsoid given by half-lengths in each direction. | |
static Inertia_ | ellipsoid (const Vec3P &halfLengths) |
Alternate interface to ellipsoid() that takes a Vec3 for the half lengths. | |
Protected Member Functions | |
const UnitInertia_< P > & | getAsUnitInertia () const |
UnitInertia_< P > & | updAsUnitInertia () |
void | errChk (const char *methodName) const |
Protected Attributes | |
SymMat33P | I_OF_F |
Related Functions | |
(Note that these are not member functions.) | |
template<class P > | |
Inertia_< P > | operator+ (const Inertia_< P > &l, const Inertia_< P > &r) |
Add two compatible inertia matrices, meaning they must be taken about the same point and expressed in the same frame. | |
template<class P > | |
Inertia_< P > | operator- (const Inertia_< P > &l, const Inertia_< P > &r) |
Subtract from one inertia matrix another one which is compatible, meaning that both must be taken about the same point and expressed in the same frame. | |
template<class P > | |
Inertia_< P > | operator* (const Inertia_< P > &i, const P &r) |
Multiply an inertia matrix by a scalar. | |
template<class P > | |
Inertia_< P > | operator* (const P &r, const Inertia_< P > &i) |
Multiply an inertia matrix by a scalar. | |
template<class P > | |
Inertia_< P > | operator* (const Inertia_< P > &i, int r) |
Multiply an inertia matrix by a scalar given as an int. | |
template<class P > | |
Inertia_< P > | operator* (int r, const Inertia_< P > &i) |
Multiply an inertia matrix by a scalar given as an int. | |
template<class P > | |
Inertia_< P > | operator/ (const Inertia_< P > &i, const P &r) |
Divide an inertia matrix by a scalar. | |
template<class P > | |
Inertia_< P > | operator/ (const Inertia_< P > &i, int r) |
Divide an inertia matrix by a scalar provided as an int. | |
template<class P > | |
Vec< 3, P > | operator* (const Inertia_< P > &I, const Vec< 3, P > &w) |
Multiply an inertia matrix I on the right by a vector w giving the vector result I*w. | |
template<class P > | |
bool | operator== (const Inertia_< P > &i1, const Inertia_< P > &i2) |
Compare two inertia matrices for exact (bitwise) equality. | |
template<class P > | |
std::ostream & | operator<< (std::ostream &o, const Inertia_< P > &inertia) |
Output a human-readable representation of an inertia matrix to the indicated stream. |
The physical meaning of an inertia is the distribution of a rigid body's mass about a particular point.
If that point is the center of mass of the body, then the measured inertia is called the "central inertia" of that body. To write down the inertia, we need to calculate the six scalars of the inertia tensor, which is a symmetric 3x3 matrix. These scalars must be expressed in an arbitrary but specified coordinate system. So an Inertia is meaningful only in conjunction with a particular set of axes, fixed to the body, whose origin is the point about which the inertia is being measured, and in whose coordinate system this measurement is being expressed. Note that changing the reference point results in a new physical quantity, but changing the reference axes only affects the measure numbers of that quantity. For any reference point, there is a unique set of reference axes in which the inertia tensor is diagonal; those are called the "principal axes" of the body at that point, and the resulting diagonal elements are the "principal moments of inertia". When we speak of an inertia being "in" a frame, we mean the physical quantity measured about the frame's origin and then expressed in the frame's axes.
This low-level Inertia class does not attempt to keep track of which frame it is in. It provides construction and operations involving inertia that can proceed using only an implicit frame F. Clients of this class are responsible for keeping track of that frame. In particular, in order to shift the inertia's "measured-about" point one must know whether either the starting or final inertia is central, because we must always shift inertias by passing through the central inertia. So this class provides operations for doing the shifting, but expects to be told by the client where to find the center of mass.
Re-expressing an Inertia in a different coordinate system does not entail a change of physical meaning in the way that shifting it to a different point does. Note that because inertia is a tensor, there is a "left frame" and "right frame". For our purposes, these will always be the same so we'll only indicate the frame once, as in 'I_pt_frame'. This should be understood to mean 'frame_I_pt_frame' and re-expressing an Inertia requires both a left and right multiply by the rotation matrix. So I_OB_B is the inertia about body B's origin point OB, expressed in B, while I_OB_G is the same physical quantity but expressed in Ground (the latter is a component of the Spatial Inertia which we usually want in the Ground frame). Frame conversion is done logically like this:
I_OB_G = R_GB * I_OB_B * R_BG (R_BG=~R_GB)
but we can save computation time by performing this as a single operation.
The central inertia would be I_CB_B for body B.
A Inertia is a symmetric matrix and is positive definite for nonsingular bodies (that is, a body composed of at least three noncollinear point masses).
Some attempt is made to check the validity of an Inertia matrix, at least when running in Debug mode. Some conditions it must satisfy are:
Typedefs exist for the most common invocations of Inertia_<P>:
SimTK::Inertia_< P >::Inertia_ | ( | ) | [inline] |
Default is a NaN-ed out mess to avoid accidents, even in Release mode.
Other than this value, an Inertia matrix should always be valid.
SimTK::Inertia_< P >::Inertia_ | ( | const RealP & | moment | ) | [inline, explicit] |
Create a principal inertia matrix with identical diagonal elements, like a sphere where moment=2/5 m r^2, or a cube where moment=1/6 m s^2, with m the total mass, r the sphere's radius and s the length of a side of the cube.
Note that many rigid bodies of different shapes and masses can have the same inertia matrix.
SimTK::Inertia_< P >::Inertia_ | ( | const Vec3P & | p, |
const RealP & | mass | ||
) | [inline] |
Create an Inertia matrix for a point mass at a given location, measured from the origin OF of the implicit frame F, and expressed in F.
Cost is 14 flops.
SimTK::Inertia_< P >::Inertia_ | ( | const Vec3P & | moments, |
const Vec3P & | products = Vec3P(0) |
||
) | [inline, explicit] |
Create an inertia matrix from a vector of the moments of inertia (the inertia matrix diagonal) and optionally a vector of the products of inertia (the off-diagonals).
Moments are in the order xx,yy,zz; products are xy,xz,yz.
SimTK::Inertia_< P >::Inertia_ | ( | const RealP & | xx, |
const RealP & | yy, | ||
const RealP & | zz | ||
) | [inline] |
Create a principal inertia matrix (only non-zero on diagonal).
SimTK::Inertia_< P >::Inertia_ | ( | const RealP & | xx, |
const RealP & | yy, | ||
const RealP & | zz, | ||
const RealP & | xy, | ||
const RealP & | xz, | ||
const RealP & | yz | ||
) | [inline] |
This is a general inertia matrix.
Note the order of these arguments: moments of inertia first, then products of inertia.
SimTK::Inertia_< P >::Inertia_ | ( | const SymMat33P & | I | ) | [inline, explicit] |
Construct an Inertia from a symmetric 3x3 matrix.
The diagonals must be nonnegative and satisfy the triangle inequality.
SimTK::Inertia_< P >::Inertia_ | ( | const Mat33P & | m | ) | [inline, explicit] |
Construct an Inertia matrix from a 3x3 symmetric matrix.
In Debug mode we'll test that the supplied matrix is numerically close to symmetric, and that it satisfies other requirements of an Inertia matrix.
Inertia_& SimTK::Inertia_< P >::setInertia | ( | const RealP & | xx, |
const RealP & | yy, | ||
const RealP & | zz | ||
) | [inline] |
Set an inertia matrix to have only principal moments (that is, it will be diagonal).
Returns a reference to "this" like an assignment operator.
Inertia_& SimTK::Inertia_< P >::setInertia | ( | const Vec3P & | moments, |
const Vec3P & | products = Vec3P(0) |
||
) | [inline] |
Set principal moments and optionally off-diagonal terms.
Returns a reference to "this" like an assignment operator.
Inertia_& SimTK::Inertia_< P >::setInertia | ( | const RealP & | xx, |
const RealP & | yy, | ||
const RealP & | zz, | ||
const RealP & | xy, | ||
const RealP & | xz, | ||
const RealP & | yz | ||
) | [inline] |
Set this Inertia to a general matrix.
Note the order of these arguments: moments of inertia first, then products of inertia. Behaves like an assignment statement. Will throw an error message in Debug mode if the supplied elements do not constitute a valid Inertia matrix.
Inertia_& SimTK::Inertia_< P >::operator+= | ( | const Inertia_< P > & | I | ) | [inline] |
Add in another inertia matrix.
Frames and reference point must be the same but we can't check. (6 flops)
Inertia_& SimTK::Inertia_< P >::operator-= | ( | const Inertia_< P > & | I | ) | [inline] |
Subtract off another inertia matrix.
Frames and reference point must be the same but we can't check. (6 flops)
Inertia_& SimTK::Inertia_< P >::operator*= | ( | const P & | s | ) | [inline] |
Multiply this inertia matrix by a scalar. Cost is 6 flops.
Inertia_& SimTK::Inertia_< P >::operator/= | ( | const P & | s | ) | [inline] |
Divide this inertia matrix by a scalar.
Cost is about 20 flops (a divide and 6 multiplies).
Inertia_ SimTK::Inertia_< P >::shiftToMassCenter | ( | const Vec3P & | CF, |
const RealP & | mass | ||
) | const [inline] |
Assume that the current inertia is about the F frame's origin OF, and expressed in F.
Given the vector from OF to the body center of mass CF, and the mass m of the body, we can shift the inertia to the center of mass. This produces a new Inertia I' whose (implicit) frame F' is aligned with F but has origin CF (an inertia like that is called a "central inertia". I' = I - Icom where Icom is the inertia of a fictitious point mass of mass m (that is, the same as the body mass) located at CF (measured in F) about OF. Cost is 20 flops.
Inertia_& SimTK::Inertia_< P >::shiftToMassCenterInPlace | ( | const Vec3P & | CF, |
const RealP & | mass | ||
) | [inline] |
Assume that the current inertia is about the F frame's origin OF, and expressed in F.
Given the vector from OF to the body center of mass CF, and the mass m of the body, we can shift the inertia to the center of mass. This produces a new Inertia I' whose (implicit) frame F' is aligned with F but has origin CF (an inertia like that is called a "central inertia". I' = I - Icom where Icom is the inertia of a fictitious point mass of mass m (that is, the same as the body mass) located at CF (measured in F) about OF. Cost is 20 flops.
Inertia_ SimTK::Inertia_< P >::shiftFromMassCenter | ( | const Vec3P & | p, |
const RealP & | mass | ||
) | const [inline] |
Assuming that the current inertia I is a central inertia (that is, it is inertia about the body center of mass CF), shift it to some other point p measured from the center of mass.
This produces a new inertia I' about the point p given by I' = I + Ip where Ip is the inertia of a fictitious point mass of mass mtot (the total body mass) located at p, about CF. Cost is 20 flops.
Inertia_& SimTK::Inertia_< P >::shiftFromMassCenterInPlace | ( | const Vec3P & | p, |
const RealP & | mass | ||
) | [inline] |
Assuming that the current inertia I is a central inertia (that is, it is inertia about the body center of mass CF), shift it to some other point p measured from the center of mass.
This produces a new inertia I' about the point p given by I' = I + Ip where Ip is the inertia of a fictitious point mass of mass mtot (the total body mass) located at p, about CF. Cost is 20 flops.
Inertia_ SimTK::Inertia_< P >::reexpress | ( | const Rotation_< P > & | R_FB | ) | const [inline] |
Return a new inertia matrix like this one but re-expressed in another frame (leaving the origin point unchanged).
Call this inertia matrix I_OF_F, that is, it is taken about the origin of some frame F, and expressed in F. We want to return I_OF_B, the same inertia matrix, still taken about the origin of F, but expressed in the B frame, given by I_OF_B=R_BF*I_OF_F*R_FB where R_FB is the rotation matrix giving the orientation of frame B in F. This is handled here by a special method of the Rotation class which rotates a symmetric tensor at a cost of 57 flops.
Reimplemented in SimTK::UnitInertia_< P >.
Inertia_ SimTK::Inertia_< P >::reexpress | ( | const InverseRotation_< P > & | R_FB | ) | const [inline] |
Rexpress using an inverse rotation to avoid having to convert it.
Reimplemented in SimTK::UnitInertia_< P >.
Inertia_& SimTK::Inertia_< P >::reexpressInPlace | ( | const Rotation_< P > & | R_FB | ) | [inline] |
Re-express this inertia matrix in another frame, changing the object in place; see reexpress() if you want to leave this object unmolested and get a new one instead.
Cost is 57 flops.
Reimplemented in SimTK::UnitInertia_< P >.
Inertia_& SimTK::Inertia_< P >::reexpressInPlace | ( | const InverseRotation_< P > & | R_FB | ) | [inline] |
Rexpress in place using an inverse rotation to avoid having to convert it.
Reimplemented in SimTK::UnitInertia_< P >.
RealP SimTK::Inertia_< P >::trace | ( | ) | const [inline] |
SimTK::Inertia_< P >::operator const SymMat33P & | ( | ) | const [inline] |
This is an implicit conversion to a const SymMat33.
Reimplemented in SimTK::UnitInertia_< P >.
const SymMat33P& SimTK::Inertia_< P >::asSymMat33 | ( | ) | const [inline] |
Obtain a reference to the underlying symmetric matrix type.
Mat33P SimTK::Inertia_< P >::toMat33 | ( | ) | const [inline] |
Expand the internal packed representation into a full 3x3 symmetric matrix with all elements set.
const Vec3P& SimTK::Inertia_< P >::getMoments | ( | ) | const [inline] |
Obtain the inertia moments (diagonal of the Inertia matrix) as a Vec3.
const Vec3P& SimTK::Inertia_< P >::getProducts | ( | ) | const [inline] |
Obtain the inertia products (off-diagonals of the Inertia matrix) as a Vec3 with elements ordered xx, xy, yz.
bool SimTK::Inertia_< P >::isNaN | ( | ) | const [inline] |
bool SimTK::Inertia_< P >::isInf | ( | ) | const [inline] |
bool SimTK::Inertia_< P >::isFinite | ( | ) | const [inline] |
bool SimTK::Inertia_< P >::isNumericallyEqual | ( | const Inertia_< P > & | other | ) | const [inline] |
Compare this inertia matrix with another one and return true if they are close to within a default numerical tolerance.
Cost is about 30 flops.
bool SimTK::Inertia_< P >::isNumericallyEqual | ( | const Inertia_< P > & | other, |
double | tol | ||
) | const [inline] |
Compare this inertia matrix with another one and return true if they are close to within a specified numerical tolerance.
Cost is about 30 flops.
static bool SimTK::Inertia_< P >::isValidInertiaMatrix | ( | const SymMat33P & | m | ) | [inline, static] |
Test some conditions that must hold for a valid Inertia matrix.
Cost is about 12 flops. TODO: this may not be comprehensive.
static Inertia_ SimTK::Inertia_< P >::pointMassAtOrigin | ( | ) | [inline, static] |
Create an Inertia matrix for a point located at the origin -- that is, an all-zero matrix.
Reimplemented in SimTK::UnitInertia_< P >.
static Inertia_ SimTK::Inertia_< P >::pointMassAt | ( | const Vec3P & | p, |
const RealP & | m | ||
) | [inline, static] |
Create an Inertia matrix for a point of a given mass, located at a given location measured from the origin of the implicit F frame.
This is equivalent to m*crossMatSq(p) but is implemented elementwise here for speed, giving a cost of 14 flops.
Inertia_< P > SimTK::Inertia_< P >::sphere | ( | const RealP & | r | ) | [inline, static] |
Create a UnitInertia matrix for a unit mass sphere of radius r centered at the origin.
Reimplemented in SimTK::UnitInertia_< P >.
Inertia_< P > SimTK::Inertia_< P >::cylinderAlongZ | ( | const RealP & | r, |
const RealP & | hz | ||
) | [inline, static] |
Unit-mass cylinder aligned along z axis; use radius and half-length.
If r==0 this is a thin rod; hz=0 it is a thin disk.
Reimplemented in SimTK::UnitInertia_< P >.
Inertia_< P > SimTK::Inertia_< P >::cylinderAlongY | ( | const RealP & | r, |
const RealP & | hy | ||
) | [inline, static] |
Unit-mass cylinder aligned along y axis; use radius and half-length.
If r==0 this is a thin rod; hy=0 it is a thin disk.
Reimplemented in SimTK::UnitInertia_< P >.
Inertia_< P > SimTK::Inertia_< P >::cylinderAlongX | ( | const RealP & | r, |
const RealP & | hx | ||
) | [inline, static] |
Unit-mass cylinder aligned along x axis; use radius and half-length.
If r==0 this is a thin rod; hx=0 it is a thin disk.
Reimplemented in SimTK::UnitInertia_< P >.
Inertia_< P > SimTK::Inertia_< P >::brick | ( | const RealP & | hx, |
const RealP & | hy, | ||
const RealP & | hz | ||
) | [inline, static] |
Unit-mass brick given by half-lengths in each direction.
One dimension zero gives inertia of a thin rectangular sheet; two zero gives inertia of a thin rod in the remaining direction.
Reimplemented in SimTK::UnitInertia_< P >.
Inertia_< P > SimTK::Inertia_< P >::brick | ( | const Vec3P & | halfLengths | ) | [inline, static] |
Alternate interface to brick() that takes a Vec3 for the half lengths.
Reimplemented in SimTK::UnitInertia_< P >.
Inertia_< P > SimTK::Inertia_< P >::ellipsoid | ( | const RealP & | hx, |
const RealP & | hy, | ||
const RealP & | hz | ||
) | [inline, static] |
Unit-mass ellipsoid given by half-lengths in each direction.
Reimplemented in SimTK::UnitInertia_< P >.
Inertia_< P > SimTK::Inertia_< P >::ellipsoid | ( | const Vec3P & | halfLengths | ) | [inline, static] |
Alternate interface to ellipsoid() that takes a Vec3 for the half lengths.
Reimplemented in SimTK::UnitInertia_< P >.
const UnitInertia_<P>& SimTK::Inertia_< P >::getAsUnitInertia | ( | ) | const [inline, protected] |
UnitInertia_<P>& SimTK::Inertia_< P >::updAsUnitInertia | ( | ) | [inline, protected] |
void SimTK::Inertia_< P >::errChk | ( | const char * | methodName | ) | const [inline, protected] |
Inertia_< P > operator+ | ( | const Inertia_< P > & | l, |
const Inertia_< P > & | r | ||
) | [related] |
Add two compatible inertia matrices, meaning they must be taken about the same point and expressed in the same frame.
There is no way to verify compatibility; make sure you know what you're doing. Cost is 6 flops.
Inertia_< P > operator- | ( | const Inertia_< P > & | l, |
const Inertia_< P > & | r | ||
) | [related] |
Subtract from one inertia matrix another one which is compatible, meaning that both must be taken about the same point and expressed in the same frame.
There is no way to verify compatibility; make sure you know what you're doing. Cost is 6 flops.
Multiply an inertia matrix by a scalar.
Cost is 6 flops.
Multiply an inertia matrix by a scalar.
Cost is 6 flops.
Multiply an inertia matrix by a scalar given as an int.
Cost is 6 flops.
Multiply an inertia matrix by a scalar given as an int.
Cost is 6 flops.
Divide an inertia matrix by a scalar.
Cost is about 20 flops (one divide and six multiplies).
Divide an inertia matrix by a scalar provided as an int.
Cost is about 20 flops (one divide and six multiplies).
Vec< 3, P > operator* | ( | const Inertia_< P > & | I, |
const Vec< 3, P > & | w | ||
) | [related] |
Multiply an inertia matrix I on the right by a vector w giving the vector result I*w.
Compare two inertia matrices for exact (bitwise) equality.
This is too strict for most purposes; use Inertia::isNumericallyEqual() instead to test for approximate equality. Cost here is 6 flops.
std::ostream & operator<< | ( | std::ostream & | o, |
const Inertia_< P > & | inertia | ||
) | [related] |
Output a human-readable representation of an inertia matrix to the indicated stream.
SymMat33P SimTK::Inertia_< P >::I_OF_F [protected] |