OpenSim::SimbodyEngine Class Reference

A wrapper class to use the SimTK Simbody dynamics engine as the underlying engine for OpenSim. More...

#include <SimbodyEngine.h>

Inheritance diagram for OpenSim::SimbodyEngine:
OpenSim::Object

List of all members.

Public Member Functions

virtual ~SimbodyEngine ()
 Destructor.
 SimbodyEngine ()
 Default constructor.
 SimbodyEngine (const std::string &aFileName)
 SimbodyEngine (const SimbodyEngine &aEngine)
 Copy constructor.
virtual Objectcopy () const
 Copy this engine and return a pointer to the copy.
SimbodyEngineoperator= (const SimbodyEngine &aEngine)
 Assignment operator.
const ModelgetModel () const
ModelgetModel ()
void setModel (Model &aModel)
void initializeState (SimTK::State &s)
 Go through the coordinate set and make sure the locking flags and coordinate values are in the state.
virtual void setup (Model &aModel)
 Perform some set up functions that happen after the object has been deserialized or copied.
virtual void updateCoordinateSet (CoordinateSet &aCoordinateSet)
 Update all coordinates in the model with the ones in the passed-in coordinate set.
virtual void getUnlockedCoordinates (const SimTK::State &s, CoordinateSet &rUnlockedCoordinates) const
 Get the set of coordinates that are not locked.
virtual bool projectConfigurationToSatisfyConstraints (SimTK::State &s, const double cTol) const
 Project (change) the current configuration such that it satisfies the constraints acting on the system.
virtual bool scale (SimTK::State &s, const ScaleSet &aScaleSet, double aFinalMass=-1.0, bool aPreserveMassDist=false)
 Scale the dynamics engine.
virtual BodygetGroundBody () const
 Get the body that is being used as ground.
virtual WrapObjectgetWrapObject (const std::string &aName) const
virtual double getMass () const
 Get the total mass of the model.
virtual void getSystemInertia (const SimTK::State &s, double &rM, SimTK::Vec3 &rCOM, double rI[3][3]) const
virtual void getSystemInertia (const SimTK::State &s, double &rM, double *rCOM, double *rI) const
 getSystemInertia
virtual void getPosition (const SimTK::State &s, const OpenSim::Body &aBody, const SimTK::Vec3 &aPoint, SimTK::Vec3 &rPos) const
virtual void getVelocity (const SimTK::State &s, const OpenSim::Body &aBody, const SimTK::Vec3 &aPoint, SimTK::Vec3 &rVel) const
virtual void getAcceleration (const SimTK::State &s, const OpenSim::Body &aBody, const SimTK::Vec3 &aPoint, SimTK::Vec3 &rAcc) const
virtual void getDirectionCosines (const SimTK::State &s, const OpenSim::Body &aBody, double rDirCos[3][3]) const
 Get the body orientation with respect to the ground.
virtual void getDirectionCosines (const SimTK::State &s, const OpenSim::Body &aBody, double *rDirCos) const
 Get the body orientation with respect to the ground.
virtual void getAngularVelocity (const SimTK::State &s, const OpenSim::Body &aBody, SimTK::Vec3 &rAngVel) const
virtual void getAngularVelocityBodyLocal (const SimTK::State &s, const OpenSim::Body &aBody, SimTK::Vec3 &rAngVel) const
virtual void getAngularAcceleration (const SimTK::State &s, const OpenSim::Body &aBody, SimTK::Vec3 &rAngAcc) const
virtual void getAngularAccelerationBodyLocal (const SimTK::State &s, const OpenSim::Body &aBody, SimTK::Vec3 &rAngAcc) const
virtual Transform getTransform (const SimTK::State &s, const OpenSim::Body &aBody) const
 get a copy of the transform from the inertial frame to a body
virtual void computeReactions (const SimTK::State &s, SimTK::Vector_< SimTK::Vec3 > &rForces, SimTK::Vector_< SimTK::Vec3 > &rTorques) const
virtual void formCompleteStorages (const SimTK::State &s, const OpenSim::Storage &aQIn, OpenSim::Storage *&rQComplete, OpenSim::Storage *&rUComplete) const
 From a potentially partial specification of the generalized coordinates, form a complete storage of the generalized coordinates (q's) and generalized speeds (u's).
virtual void formEulerTransform (const SimTK::State &s, const OpenSim::Body &aBody, double *rE) const
virtual void formMassMatrix (double *rI)
virtual void formJacobianTranslation (const OpenSim::Body &aBody, const SimTK::Vec3 &aPoint, double *rJ, const OpenSim::Body *aRefBody=NULL) const
virtual void formJacobianOrientation (const OpenSim::Body &aBody, double *rJ0, const OpenSim::Body *aRefBody=NULL) const
virtual void formJacobianEuler (const OpenSim::Body &aBody, double *rJE, const OpenSim::Body *aRefBody=NULL) const
virtual void computeDerivatives (const SimTK::State &s, double *dqdt, double *dudt)
 Compute the derivatives of the generalized coordinates and speeds.
virtual void transform (const SimTK::State &s, const OpenSim::Body &aBodyFrom, const double aVec[3], const OpenSim::Body &aBodyTo, double rVec[3]) const
 Transform a vector from one body to another.
virtual void transform (const SimTK::State &s, const OpenSim::Body &aBodyFrom, const SimTK::Vec3 &aVec, const OpenSim::Body &aBodyTo, SimTK::Vec3 &rVec) const
virtual void transformPosition (const SimTK::State &s, const OpenSim::Body &aBodyFrom, const double aPos[3], const OpenSim::Body &aBodyTo, double rPos[3]) const
 Transform a point from one body to another.
virtual void transformPosition (const SimTK::State &s, const OpenSim::Body &aBodyFrom, const SimTK::Vec3 &aPos, const OpenSim::Body &aBodyTo, SimTK::Vec3 &rPos) const
virtual void transformPosition (const SimTK::State &s, const OpenSim::Body &aBodyFrom, const double aPos[3], double rPos[3]) const
 Transform a point from one body to the ground body.
virtual void transformPosition (const SimTK::State &s, const OpenSim::Body &aBodyFrom, const SimTK::Vec3 &aPos, SimTK::Vec3 &rPos) const
virtual double calcDistance (const SimTK::State &s, const OpenSim::Body &aBody1, const double aPoint1[3], const OpenSim::Body &aBody2, const double aPoint2[3]) const
 Calculate the distance between a point on one body and a point on another body.
virtual double calcDistance (const SimTK::State &s, const OpenSim::Body &aBody1, const SimTK::Vec3 &aPoint1, const OpenSim::Body &aBody2, const SimTK::Vec3 &aPoint2) const
void convertRadiansToDegrees (Storage &rStorage) const
 Convert the rotational generalized coordinates or speeds from units of radians to units of degrees for all the state-vectors in an Storage object.
void convertDegreesToRadians (Storage &rStorage) const
 Convert the rotational generalized coordinates or speeds from units of degrees to units of radians for all the state-vectors in an Storage object.
void convertDegreesToRadians (double *aQDeg, double *rQRad) const
 Convert an array of Q/U values from degrees to radians.
void convertRadiansToDegrees (double *aQRad, double *rQDeg) const
 Convert an array of Q/U values from radians to degrees.
virtual void convertAnglesToDirectionCosines (double aE1, double aE2, double aE3, double rDirCos[3][3]) const
 Convert angles to direction cosines.
virtual void convertAnglesToDirectionCosines (double aE1, double aE2, double aE3, double *rDirCos) const
 Convert angles to direction cosines.
virtual void convertDirectionCosinesToAngles (double aDirCos[3][3], double *rE1, double *rE2, double *rE3) const
 Convert direction cosines to angles.
virtual void convertDirectionCosinesToAngles (double *aDirCos, double *rE1, double *rE2, double *rE3) const
 Convert direction cosines to angles.
virtual void convertDirectionCosinesToQuaternions (double aDirCos[3][3], double *rQ1, double *rQ2, double *rQ3, double *rQ4) const
 Convert direction cosines to quaternions.
virtual void convertDirectionCosinesToQuaternions (double *aDirCos, double *rQ1, double *rQ2, double *rQ3, double *rQ4) const
 Convert direction cosines to quaternions.
virtual void convertQuaternionsToDirectionCosines (double aQ1, double aQ2, double aQ3, double aQ4, double rDirCos[3][3]) const
 Convert quaternions to direction cosines.
virtual void convertQuaternionsToDirectionCosines (double aQ1, double aQ2, double aQ3, double aQ4, double *rDirCos) const
 Convert quaternions to direction cosines.

Public Attributes

Model_model
 Pointer to the model that owns this dynamics engine.

Protected Attributes

Body_groundBody
 Body used for ground, the inertial frame.

Friends

class Body
class Coordinate
class Joint
class Constraint
class WeldConstraint
class CoordinateCouplerConstraint

Detailed Description

A wrapper class to use the SimTK Simbody dynamics engine as the underlying engine for OpenSim.

Authors:
Frank C. Anderson, Ajay Seth
Version:
1.0

Constructor & Destructor Documentation

SimbodyEngine::~SimbodyEngine (  )  [virtual]

Destructor.

SimbodyEngine::SimbodyEngine (  ) 

Default constructor.

This constructor constructs a dynamic model of a simple pendulum.

OpenSim::SimbodyEngine::SimbodyEngine ( const std::string &  aFileName  ) 
SimbodyEngine::SimbodyEngine ( const SimbodyEngine aEngine  ) 

Copy constructor.


Member Function Documentation

virtual double OpenSim::SimbodyEngine::calcDistance ( const SimTK::State &  s,
const OpenSim::Body aBody1,
const SimTK::Vec3 &  aPoint1,
const OpenSim::Body aBody2,
const SimTK::Vec3 &  aPoint2 
) const [virtual]
double SimbodyEngine::calcDistance ( const SimTK::State &  s,
const OpenSim::Body aBody1,
const double  aPoint1[3],
const OpenSim::Body aBody2,
const double  aPoint2[3] 
) const [virtual]

Calculate the distance between a point on one body and a point on another body.

Parameters:
aBody1 the body that the first point is expressed in
aPoint1 the XYZ coordinates of the first point
aBody2 the body that the second point is expressed in
aPoint2 the XYZ coordinates of the second point
Returns:
the distance between aPoint1 and aPoint2
void SimbodyEngine::computeDerivatives ( const SimTK::State &  s,
double *  dqdt,
double *  dudt 
) [virtual]

Compute the derivatives of the generalized coordinates and speeds.

Parameters:
dqdt Derivatives of generalized coordinates.
dudt Derivatives of generalized speeds.
virtual void OpenSim::SimbodyEngine::computeReactions ( const SimTK::State &  s,
SimTK::Vector_< SimTK::Vec3 > &  rForces,
SimTK::Vector_< SimTK::Vec3 > &  rTorques 
) const [virtual]
void SimbodyEngine::convertAnglesToDirectionCosines ( double  aE1,
double  aE2,
double  aE3,
double *  rDirCos 
) const [virtual]

Convert angles to direction cosines.

Parameters:
aE1 1st Euler angle.
aE2 2nd Euler angle.
aE3 3rd Euler angle.
rDirCos Vector of direction cosines.
void SimbodyEngine::convertAnglesToDirectionCosines ( double  aE1,
double  aE2,
double  aE3,
double  rDirCos[3][3] 
) const [virtual]

Convert angles to direction cosines.

Parameters:
aE1 1st Euler angle.
aE2 2nd Euler angle.
aE3 3rd Euler angle.
rDirCos Vector of direction cosines.
void SimbodyEngine::convertDegreesToRadians ( double *  aQDeg,
double *  rQRad 
) const

Convert an array of Q/U values from degrees to radians.

The sizes of the arrays are assumed to be equal to the number of Us.

Parameters:
aQDeg Array of values, in degrees
rQRad Array of values, in radians
void SimbodyEngine::convertDegreesToRadians ( Storage rStorage  )  const

Convert the rotational generalized coordinates or speeds from units of degrees to units of radians for all the state-vectors in an Storage object.

Coordinates/speeds are identified by column names.

Parameters:
rStorage Storage object.
void SimbodyEngine::convertDirectionCosinesToAngles ( double *  aDirCos,
double *  rE1,
double *  rE2,
double *  rE3 
) const [virtual]

Convert direction cosines to angles.

Parameters:
aDirCos Vector of direction cosines.
rE1 1st Euler angle.
rE2 2nd Euler angle.
rE3 3rd Euler angle.
void SimbodyEngine::convertDirectionCosinesToAngles ( double  aDirCos[3][3],
double *  rE1,
double *  rE2,
double *  rE3 
) const [virtual]

Convert direction cosines to angles.

Parameters:
aDirCos Vector of direction cosines.
rE1 1st Euler angle.
rE2 2nd Euler angle.
rE3 3rd Euler angle.
void SimbodyEngine::convertDirectionCosinesToQuaternions ( double *  aDirCos,
double *  rQ1,
double *  rQ2,
double *  rQ3,
double *  rQ4 
) const [virtual]

Convert direction cosines to quaternions.

Parameters:
aDirCos Vector of direction cosines.
rQ1 1st Quaternion.
rQ2 2nd Quaternion.
rQ3 3rd Quaternion.
rQ4 4th Quaternion.
void SimbodyEngine::convertDirectionCosinesToQuaternions ( double  aDirCos[3][3],
double *  rQ1,
double *  rQ2,
double *  rQ3,
double *  rQ4 
) const [virtual]

Convert direction cosines to quaternions.

Parameters:
aDirCos Vector of direction cosines.
rQ1 1st Quaternion.
rQ2 2nd Quaternion.
rQ3 3rd Quaternion.
rQ4 4th Quaternion.
void SimbodyEngine::convertQuaternionsToDirectionCosines ( double  aQ1,
double  aQ2,
double  aQ3,
double  aQ4,
double *  rDirCos 
) const [virtual]

Convert quaternions to direction cosines.

Parameters:
aQ1 1st Quaternion.
aQ2 2nd Quaternion.
aQ3 3rd Quaternion.
aQ4 4th Quaternion.
rDirCos Vector of direction cosines.
void SimbodyEngine::convertQuaternionsToDirectionCosines ( double  aQ1,
double  aQ2,
double  aQ3,
double  aQ4,
double  rDirCos[3][3] 
) const [virtual]

Convert quaternions to direction cosines.

Parameters:
aQ1 1st Quaternion.
aQ2 2nd Quaternion.
aQ3 3rd Quaternion.
aQ4 4th Quaternion.
rDirCos Vector of direction cosines.
void SimbodyEngine::convertRadiansToDegrees ( double *  aQRad,
double *  rQDeg 
) const

Convert an array of Q/U values from radians to degrees.

The sizes of the arrays are assumed to be equal to the number of Us.

Parameters:
aQRad Array of values, in radians
rQDeg Array of values, in degrees
void SimbodyEngine::convertRadiansToDegrees ( Storage rStorage  )  const

Convert the rotational generalized coordinates or speeds from units of radians to units of degrees for all the state-vectors in an Storage object.

Coordinates/speeds are identified by column names.

Parameters:
rStorage Storage object.
Object * SimbodyEngine::copy (  )  const [virtual]

Copy this engine and return a pointer to the copy.

The copy constructor for this class is used.

Returns:
Pointer to a copy of this SimbodyEngine.

Reimplemented from OpenSim::Object.

void SimbodyEngine::formCompleteStorages ( const SimTK::State &  s,
const OpenSim::Storage aQIn,
OpenSim::Storage *&  rQComplete,
OpenSim::Storage *&  rUComplete 
) const [virtual]

From a potentially partial specification of the generalized coordinates, form a complete storage of the generalized coordinates (q's) and generalized speeds (u's).

Parameters:
aQIn Storage containing the q's or a subset of the q's. Rotational q's should be in degrees.
rQComplete Storage containing all the q's. If q's were not in aQIn, the values are set to 0.0. When a q is constrained, its value is altered to be consistent with the constraint. The caller is responsible for deleting the memory associated with this storage.
rUComplete Storage containing all the u's. The generalized speeds are obtained by spline fitting the q's and differentiating the splines. When a u is constrained, its value is altered to be consisten with the constraint. The caller is responsible for deleting the memory associated with this storage.
void SimbodyEngine::formEulerTransform ( const SimTK::State &  s,
const OpenSim::Body aBody,
double *  rE 
) const [virtual]
virtual void OpenSim::SimbodyEngine::formJacobianEuler ( const OpenSim::Body aBody,
double *  rJE,
const OpenSim::Body aRefBody = NULL 
) const [inline, virtual]
virtual void OpenSim::SimbodyEngine::formJacobianOrientation ( const OpenSim::Body aBody,
double *  rJ0,
const OpenSim::Body aRefBody = NULL 
) const [inline, virtual]
virtual void OpenSim::SimbodyEngine::formJacobianTranslation ( const OpenSim::Body aBody,
const SimTK::Vec3 &  aPoint,
double *  rJ,
const OpenSim::Body aRefBody = NULL 
) const [inline, virtual]
virtual void OpenSim::SimbodyEngine::formMassMatrix ( double *  rI  )  [inline, virtual]
virtual void OpenSim::SimbodyEngine::getAcceleration ( const SimTK::State &  s,
const OpenSim::Body aBody,
const SimTK::Vec3 &  aPoint,
SimTK::Vec3 &  rAcc 
) const [virtual]
virtual void OpenSim::SimbodyEngine::getAngularAcceleration ( const SimTK::State &  s,
const OpenSim::Body aBody,
SimTK::Vec3 &  rAngAcc 
) const [virtual]
virtual void OpenSim::SimbodyEngine::getAngularAccelerationBodyLocal ( const SimTK::State &  s,
const OpenSim::Body aBody,
SimTK::Vec3 &  rAngAcc 
) const [virtual]
virtual void OpenSim::SimbodyEngine::getAngularVelocity ( const SimTK::State &  s,
const OpenSim::Body aBody,
SimTK::Vec3 &  rAngVel 
) const [virtual]
virtual void OpenSim::SimbodyEngine::getAngularVelocityBodyLocal ( const SimTK::State &  s,
const OpenSim::Body aBody,
SimTK::Vec3 &  rAngVel 
) const [virtual]
void SimbodyEngine::getDirectionCosines ( const SimTK::State &  s,
const OpenSim::Body aBody,
double *  rDirCos 
) const [virtual]

Get the body orientation with respect to the ground.

Parameters:
aBody Pointer to body.
rDirCos Orientation of the body with respect to the ground frame.
void SimbodyEngine::getDirectionCosines ( const SimTK::State &  s,
const OpenSim::Body aBody,
double  rDirCos[3][3] 
) const [virtual]

Get the body orientation with respect to the ground.

Parameters:
aBody Pointer to body.
rDirCos Orientation of the body with respect to the ground frame.
OpenSim::Body & SimbodyEngine::getGroundBody (  )  const [virtual]

Get the body that is being used as ground.

Returns:
Pointer to the ground body.
double SimbodyEngine::getMass (  )  const [virtual]

Get the total mass of the model.

Returns:
the mass of the model
Model& OpenSim::SimbodyEngine::getModel (  )  [inline]
const Model& OpenSim::SimbodyEngine::getModel (  )  const [inline]
virtual void OpenSim::SimbodyEngine::getPosition ( const SimTK::State &  s,
const OpenSim::Body aBody,
const SimTK::Vec3 &  aPoint,
SimTK::Vec3 &  rPos 
) const [virtual]
void SimbodyEngine::getSystemInertia ( const SimTK::State &  s,
double &  rM,
double *  rCOM,
double *  rI 
) const [virtual]

getSystemInertia

Parameters:
rM 
rCOM 
rI 
virtual void OpenSim::SimbodyEngine::getSystemInertia ( const SimTK::State &  s,
double &  rM,
SimTK::Vec3 &  rCOM,
double  rI[3][3] 
) const [virtual]
SimTK::Transform SimbodyEngine::getTransform ( const SimTK::State &  s,
const OpenSim::Body aBody 
) const [virtual]

get a copy of the transform from the inertial frame to a body

Parameters:
aBody 
Returns:
Transform from inertial frame to body
void SimbodyEngine::getUnlockedCoordinates ( const SimTK::State &  s,
CoordinateSet rUnlockedCoordinates 
) const [virtual]

Get the set of coordinates that are not locked.

Parameters:
rUnlockedCoordinates set of unlocked coordinates is returned here
virtual void OpenSim::SimbodyEngine::getVelocity ( const SimTK::State &  s,
const OpenSim::Body aBody,
const SimTK::Vec3 &  aPoint,
SimTK::Vec3 &  rVel 
) const [virtual]
virtual WrapObject* OpenSim::SimbodyEngine::getWrapObject ( const std::string &  aName  )  const [virtual]
void SimbodyEngine::initializeState ( SimTK::State &  s  ) 

Go through the coordinate set and make sure the locking flags and coordinate values are in the state.

Parameters:
s current SimTK state
SimbodyEngine & SimbodyEngine::operator= ( const SimbodyEngine aEngine  ) 

Assignment operator.

Returns:
Reference to this object.

Reimplemented from OpenSim::Object.

bool SimbodyEngine::projectConfigurationToSatisfyConstraints ( SimTK::State &  s,
const double  cTol 
) const [virtual]

Project (change) the current configuration such that it satisfies the constraints acting on the system.

This alters the configuration (state) down in Simbody and then returns the projected configuration. This method is intended for corrections to the state during integration and should be called after a successful integrations step. It can also be called to pose the model (i.e. in coord->setValue() after ) so that the configuration satisfies the constraints but WARNING this may produce upredictable results when the current configuration is far from satisfying the constraints.

returns true iff the engine configuration (state) was changed by the projection

Parameters:
cTol constraint tolerance. (input)
bool SimbodyEngine::scale ( SimTK::State &  s,
const ScaleSet aScaleSet,
double  aFinalMass = -1.0,
bool  aPreserveMassDist = false 
) [virtual]

Scale the dynamics engine.

Parameters:
aScaleSet the set of XYZ scale factors for the bodies
aFinalMass the mass that the scaled model should have
aPreserveMassDist whether or not the masses of the individual bodies should be scaled with the body scale factors.
Returns:
Whether or not scaling was successful.
void OpenSim::SimbodyEngine::setModel ( Model aModel  )  [inline]
void SimbodyEngine::setup ( Model aModel  )  [virtual]

Perform some set up functions that happen after the object has been deserialized or copied.

Parameters:
aModel model containing this SimbodyEngine.
virtual void OpenSim::SimbodyEngine::transform ( const SimTK::State &  s,
const OpenSim::Body aBodyFrom,
const SimTK::Vec3 &  aVec,
const OpenSim::Body aBodyTo,
SimTK::Vec3 &  rVec 
) const [virtual]
void SimbodyEngine::transform ( const SimTK::State &  s,
const OpenSim::Body aBodyFrom,
const double  aVec[3],
const OpenSim::Body aBodyTo,
double  rVec[3] 
) const [virtual]

Transform a vector from one body to another.

Parameters:
aBodyFrom the body in which the vector is currently expressed
aPos the vector to be transformed
aBodyTo the body the vector will be transformed into
rPos the vector in the aBodyTo frame is returned here
virtual void OpenSim::SimbodyEngine::transformPosition ( const SimTK::State &  s,
const OpenSim::Body aBodyFrom,
const SimTK::Vec3 &  aPos,
SimTK::Vec3 &  rPos 
) const [virtual]
void SimbodyEngine::transformPosition ( const SimTK::State &  s,
const OpenSim::Body aBodyFrom,
const double  aPos[3],
double  rPos[3] 
) const [virtual]

Transform a point from one body to the ground body.

Parameters:
aBodyFrom the body in which the point is currently expressed
aPos the XYZ coordinates of the point
rPos the XYZ coordinates of the point in the ground frame are returned here
virtual void OpenSim::SimbodyEngine::transformPosition ( const SimTK::State &  s,
const OpenSim::Body aBodyFrom,
const SimTK::Vec3 &  aPos,
const OpenSim::Body aBodyTo,
SimTK::Vec3 &  rPos 
) const [virtual]
void SimbodyEngine::transformPosition ( const SimTK::State &  s,
const OpenSim::Body aBodyFrom,
const double  aPos[3],
const OpenSim::Body aBodyTo,
double  rPos[3] 
) const [virtual]

Transform a point from one body to another.

Parameters:
aBodyFrom the body in which the point is currently expressed
aPos the XYZ coordinates of the point
aBodyTo the body the point will be transformed into
rPos the XYZ coordinates of the point in the aBodyTo frame are returned here
void SimbodyEngine::updateCoordinateSet ( CoordinateSet aCoordinateSet  )  [virtual]

Update all coordinates in the model with the ones in the passed-in coordinate set.

If the coordinate does not exist in the model, it is not added.

Parameters:
aCoordinateSet set of coordinates to be updated/added

Friends And Related Function Documentation

friend class Body [friend]
friend class Constraint [friend]
friend class Coordinate [friend]
friend class CoordinateCouplerConstraint [friend]
friend class Joint [friend]
friend class WeldConstraint [friend]

Member Data Documentation

Body used for ground, the inertial frame.

Pointer to the model that owns this dynamics engine.


The documentation for this class was generated from the following files:

Generated on Wed Dec 16 15:03:48 2009 for OpenSim by  doxygen 1.6.1