Simbody  3.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
SimTK::System Class Reference

This is the base class that serves as the parent of all SimTK System objects; most commonly Simbody's MultibodySystem class. More...

#include <System.h>

+ Inheritance diagram for SimTK::System:

Classes

class  Guts
 This is the declaration for the System::Guts class, the abstract object to which a System handle points. More...
 

Public Member Functions

System-wide properties

These methods set parameters that provide System-level defaults that can be used by various algorithms.

This includes nominal time and length scales, and hints that can be useful when visualizing this System (for example, does it make sense to show a ground plane?). There are defaults for all of these parameters so you don't need to set them unless the defaults are not satisfactory.

SystemsetUpDirection (const CoordinateDirection &up)
 This is a hint to other software as to which way this System's designer considers to be "up". More...
 
SystemsetUseUniformBackground (bool useUniformBackground)
 This is a hint to visualization software that this System is best viewed against a uniform background (e.g. all white) rather than against a ground plane. More...
 
SystemsetDefaultTimeScale (Real tc)
 (Advanced) This is a hint used for some default behaviors, such as determining an initial step size for an integrator, or the default unit error for a constraint error derivative from the original constraint. More...
 
SystemsetDefaultLengthScale (Real lc)
 (Advanced) This is a hint that can be used to get a sense of what a "unit length" looks like for this System in the units it uses. More...
 
void setHasTimeAdvancedEvents (bool)
 This determines whether this System wants to be notified whenever time advances irreversibly. More...
 
CoordinateDirection getUpDirection () const
 Get the current setting of the "up" direction hint. More...
 
bool getUseUniformBackground () const
 Get the current setting of the "use uniform background" visualization hint. More...
 
Real getDefaultTimeScale () const
 Get the currently-set value for the default time scale, 0.1 time units if it has never been set. More...
 
Real getDefaultLengthScale () const
 Get the currently-set value for the default length scale, 1 length unit if it has never been set. More...
 
bool hasTimeAdvancedEvents () const
 Return the current value of the flag indicating whether this System wants an event generated whenever time advances irreversibly. More...
 
Event handlers and reporters

These methods are used to attach event handlers and reporters to this System.

These are actually added to the DefaultSystemSubsystem that is contained in every System.

void addEventHandler (ScheduledEventHandler *handler)
 Add a ScheduledEventHandler to this System, which takes over ownership of the event handler object. More...
 
void addEventHandler (TriggeredEventHandler *handler)
 Add a TriggeredEventHandler to this System, which takes over ownership of the event handler object. More...
 
void addEventReporter (ScheduledEventReporter *handler) const
 Add a ScheduledEventReporter to this System, which takes over ownership of the event reporter object. More...
 
void addEventReporter (TriggeredEventReporter *handler) const
 Add a TriggeredEventReporter to this System, which takes over ownership of the event reporter object. More...
 
Realization

These methods provide the ability to compute the consequences that follow from the values of state variables, leaving the results in the appropriate cache.

Usually this refers to values in a given State object, but we also use the same concept for values assigned to a System during its (extended) construction, which we call Topological variables. The realizeTopology() method is called at the end of construction, allocating resources, performing computations, and creating a default State, with all results stored in a cache that is maintained by the System itself. The other realize() methods are given a State object to work with and store their results in that State object without making any changes to the System.

const StaterealizeTopology () const
 The following call must be made after any topological change has been made to this System, before the System can be used to perform any computations. More...
 
const StategetDefaultState () const
 This is available after realizeTopology(), and will throw an exception if realizeTopology() has not been called since the most recent topological change to this System. More...
 
StateupdDefaultState ()
 Don't use this; make a copy of the default state instead and modify your copy. More...
 
void realizeModel (State &state) const
 Realize the model to be used for subsequent computations with the given state. More...
 
void realize (const State &state, Stage stage=Stage::HighestRuntime) const
 Realize the given state to the indicated stage. More...
 
The Constrained System

These methods deal with prescribed motion and constraints on the allowable values of q and u.

They are primarly for use by numerical integrators. Prescribed motion is always resolved exactly because it involves explicit computation of q=q(t) and u=u(t,q). Constraints are given by implicit equations like g(q)=0 and are solved to a specified tolerance. Also, prescribed motion has precedence over constraints in the sense that no changes will be made to prescribed variables in order to satisfy constraints; only the free variables are modified by constraint solvers.

void project (State &state, Real accuracy=-1) const
 Apply prescribed motion and attempt to satisfy all position and velocity constraints by projecting the generalized coordinates q and generalized speeds u onto the position and velocity constraint manifolds. More...
 
void projectQ (State &state, Real accuracy=-1) const
 Set prescribed positions q=q(t) and attempt to satisfy position constraints by modifying the remaining free q's. More...
 
void projectU (State &state, Real accuracy=-1) const
 Set prescribed velocities u=u(t,q) and attempt to satisfy velocity constraints by modifying the remaining free u's. More...
 
void projectQ (State &state, Vector &qErrEst, const ProjectOptions &options, ProjectResults &results) const
 (Advanced) Project free q's so that position constraints are satisfied and remove the corresponding error from the supplied error estimate. More...
 
void projectU (State &state, Vector &uErrEst, const ProjectOptions &options, ProjectResults &results) const
 (Advanced) Project free u's so that velocity constraints are satisfied and remove the corresponding error from the supplied error estimate. More...
 
void prescribe (State &state) const
 Set values for prescribed positions q and velocities u. More...
 
bool prescribeQ (State &state) const
 This solver sets prescribed q's to their known values q=q(t) as a function of time (and possibly earlier-stage discrete state variables). More...
 
bool prescribeU (State &state) const
 This solver sets prescribed u's to their known values u=u(t,q) as a function of time and position variables q (and possibly earlier-stage discrete state variables). More...
 
The Discrete System

These methods deal with the discrete (event-driven) aspects of this System.

These methods are primarily for the use of Integrator objects that detect events and TimeStepper objects that deal with them by calling appropriate handlers. In general, events interrupt the continuous evolution of a system, dividing a simulation into continuous segments interrupted by discontinuous changes of state. There are also discrete variables that can be updated without disrupting the continuous flow. These are called "auto update variables" and are updated at the beginning of every integration step, right after the initial state derivatives have been recorded. Thus the step proceeds with those derivatives even if the update could cause them to change. These updates are initiated by Integrator objects at the start of every continuous step, while trajectory-affecting discrete variable updates are initiated only by TimeStepper objects via calls to event handlers.

See Also
State::allocateAutoUpdateDiscreteVariable()
void handleEvents (State &state, Event::Cause cause, const Array_< EventId > &eventIds, const HandleEventsOptions &options, HandleEventsResults &results) const
 This solver handles a set of events which a TimeStepper has denoted as having occurred at the given time and state. More...
 
void reportEvents (const State &state, Event::Cause cause, const Array_< EventId > &eventIds) const
 This method is similar to handleEvents(), but does not allow the State to be modified. More...
 
void calcEventTriggerInfo (const State &state, Array_< EventTriggerInfo > &triggerInfo) const
 This routine provides the Integrator with information it needs about the individual event trigger functions, such as which sign transitions are relevant and how tightly we need to localize. More...
 
void calcTimeOfNextScheduledEvent (const State &state, Real &tNextEvent, Array_< EventId > &eventIds, bool includeCurrentTime) const
 This routine should be called to determine if and when there is an event scheduled to occur at a particular time. More...
 
void calcTimeOfNextScheduledReport (const State &state, Real &tNextEvent, Array_< EventId > &eventIds, bool includeCurrentTime) const
 This routine is similar to calcTimeOfNextScheduledEvent(), but is used for "reporting events" which do not modify the state. More...
 
The Fast System

(NOT IMPLEMENTED YET) These methods are for state variables whose dynamics are considered very fast compared to the rest of the system, meaning that their values can be determined by solving for algebraic equilibrium conditions rather than integrating differential equations.

void relax (State &state, Stage stage, Real accuracy=-1) const
 This optional method should modify fast variables at the given stage until they satisfy some relaxation criteria. More...
 
Kinematic differential equations

These methods are primarily for use by numerical integration methods.

Generalized coordinates q are not independent of generalized speeds u; they are related by the kinematic differential equation qdot=N(q)*u, where N is an nqXnu block diagonal, invertible matrix in the sense that u=pinv(N)*qdot, where pinv(N) is the pseudoinverse of N. N has full column rank, so pinv(N)*N = I, but N*pinv(N) != I.

Just as N provides the relation between velocities expressed in u-space and the equivalent in q-space, its transpose ~N relates forces in q-space to their equivalent in u-space: fu=~N*fq, and fq=~pinv(N)*fu. (Note that ~pinv(X)==pinv(~X).) This satisfies power=~fq*qdot==~fu*u as it must.

We provide fast O(n) operators for multiplication by N, pinv(N), and their transposes. N is often mostly an identity matrix so very little computation is required.

void multiplyByN (const State &state, const Vector &u, Vector &dq) const
 Calculate dq=N*u in O(n) time (very fast). More...
 
void multiplyByNTranspose (const State &state, const Vector &fq, Vector &fu) const
 Calculate fu=~N*fq in O(n) time (very fast). More...
 
void multiplyByNPInv (const State &state, const Vector &dq, Vector &u) const
 Calculate u=pinv(N)*dq in O(n) time (very fast). More...
 
void multiplyByNPInvTranspose (const State &state, const Vector &fu, Vector &fq) const
 Calculate fq=~pinv(N)*fu in O(n) time (very fast). More...
 
Statistics

The System keeps mutable statistics internally, initialized to zero at construction.

These must not affect results in any way.

void resetAllCountersToZero ()
 Zero out the statistics for this System. More...
 
int getNumRealizationsOfThisStage (Stage) const
 Whenever the system was realized from Stage-1 to the indicated Stage, this counter is bumped. More...
 
int getNumRealizeCalls () const
 Return the total number of calls to realizeTopology(), realizeModel(), or realize(), regardless of whether these routines actually did anything when called. More...
 
int getNumPrescribeQCalls () const
 Return the total number of calls to the System's prescribeQ() method. More...
 
int getNumPrescribeUCalls () const
 Return the total number of calls to the System's prescribeU() method. More...
 
int getNumProjectQCalls () const
 Return the total number of calls to projectQ(), regardless of whether the call did anything. More...
 
int getNumFailedProjectQCalls () const
 Return the total number of calls to projectQ() that failed. More...
 
int getNumQProjections () const
 How many of the successful projectQ() calls actually did a constraint projection, rather than returning quickly? More...
 
int getNumQErrorEstimateProjections () const
 How many of the projectQ() calls that did a constraint projection also projected an error estimate? More...
 
int getNumProjectUCalls () const
 Return the total number of calls to projectU(), regardless of whether the call did anything. More...
 
int getNumFailedProjectUCalls () const
 Return the total number of calls to projectU() that failed. More...
 
int getNumUProjections () const
 How many of the successful projectU() calls actually did a constraint projection, rather than returning quickly? More...
 
int getNumUErrorEstimateProjections () const
 How many of the projectU() calls that did a constraint projection also projected an error estimate? More...
 
int getNumHandlerCallsThatChangedStage (Stage) const
 handleEvents() reports the lowest Stage it modified and we bump the counter for that Stage. More...
 
int getNumHandleEventCalls () const
 This is the total number of calls to handleEvents() regardless of the outcome. More...
 
int getNumReportEventCalls () const
 This is the total number of calls to reportEvents() regardless of the outcome. More...
 
Construction and bookkeeping

These methods are mostly for use by concrete Systems and will not typically be of interest to users.

 System ()
 Default constructor creates an empty handle. More...
 
 System (const System &)
 Copy constructor (untested). More...
 
Systemoperator= (const System &)
 Copy assignment (untested). More...
 
 ~System ()
 Destructor here will invoke the virtual destructor in the System::Guts class and thus clean up all unneeded detritus owned by this System. More...
 
const StringgetName () const
 Return the name assigned to this System (arbitrary). More...
 
const StringgetVersion () const
 Return the version string assigned to this System (arbitrary). More...
 
SubsystemIndex adoptSubsystem (Subsystem &child)
 Take over ownership of the supplied subsystem and install it into the next free subsystem slot. More...
 
int getNumSubsystems () const
 How may Subsystems are in here? More...
 
const SubsystemgetSubsystem (SubsystemIndex) const
 Obtain read-only access to a particular subsystem by its index. More...
 
SubsystemupdSubsystem (SubsystemIndex)
 Obtain writable access to a particular subsystem by its index. More...
 
const DefaultSystemSubsystemgetDefaultSubsystem () const
 Get read-only access to the default subsystem which is present in every system. More...
 
DefaultSystemSubsystemupdDefaultSubsystem ()
 Get writable access to the default subsystem which is present in every system. More...
 
 operator const Subsystem & () const
 Implicitly convert this System into a const Subsystem reference; this actually returns a reference to the DefaultSystemSubsystem contained in this System. More...
 
 operator Subsystem & ()
 Implicitly convert this System into a writable Subsystem reference; this actually returns a reference to the DefaultSystemSubsystem contained in this System. More...
 
bool systemTopologyHasBeenRealized () const
 (Advanced) You can check whether realizeTopology() has been called since the last topological change to this Syatem. More...
 
StageVersion getSystemTopologyCacheVersion () const
 (Advanced) Return the current version number of this system's Topology cache information. More...
 
void setSystemTopologyCacheVersion (StageVersion topoVersion) const
 (Really advanced) Set the current version number of this system's Topology cache information. More...
 
void invalidateSystemTopologyCache () const
 (Advanced) Mark the Topology stage of this system and all its subsystems "not realized." This is normally handled automatically by whenever you make a Topology-stage change to any subsystem. More...
 
void calcDecorativeGeometryAndAppend (const State &, Stage, Array_< DecorativeGeometry > &) const
 (Advanced) Generate all decorative geometry computable at a specific stage; this is useful for visualizers. More...
 
bool isSameSystem (const System &otherSystem) const
 There can be multiple handles referring to the same System::Guts object; they are considered to be the same System. More...
 
const GutsgetSystemGuts () const
 Obtain a const reference to the System::Guts object to which this handle refers. More...
 
GutsupdSystemGuts ()
 Obtain a writable reference to the System::Guts object to which this handle refers. More...
 
void adoptSystemGuts (System::Guts *g)
 Put new unowned Guts into this *empty* handle and take over ownership. More...
 
 System (System::Guts *g)
 Constructor for internal use only. More...
 
bool hasGuts () const
 Return true if this System handle is not empty. More...
 
bool isOwnerHandle () const
 Internal use only. More...
 
bool isEmptyHandle () const
 Internal use only. More...
 

Friends

class Guts
 

Detailed Description

This is the base class that serves as the parent of all SimTK System objects; most commonly Simbody's MultibodySystem class.

System is also the object type on which the Simbody numerical integrators operate, and many of the methods provided here are primarily for use by integrators.

A System serves as a mediator for a group of interacting Subsystem objects. All will share a single system State, and typically subsystems will need access to content in the state which is managed by other subsystems. A System provides a unique SubsystemIndex (a small positive integer) for each of its subsystems, and the subsystems are constructed knowing their indices. The indices are used subsequently by the subsystems to find their own entries in the system state, and by each subsystem to refer to others within the same system. Index 0 is reserved for use by the System itself, e.g. for system-global state variables, which are maintained in subsystem 0.

Concrete Systems understand the kinds of subsystems they contain. For example, a MultibodySystem might contain a mechanical subsystem, some force subsystems, and a geometry subsystem. At each computation stage, a subsystem is realized in a single operation. That operation can refer to computations from already-realized subsystems, but cannot initiate computation in other subsystems. The System must know the proper order in which to realize the subsystems at each stage, and that ordering is likely to vary with stage. For example, at Position stage the mechanical positions must be realized before the configuration-dependent force elements. However, at Acceleration stage, the force elements must be realized before the mechanical accelerations can be calculated.

There are three distinct users of this class:

  • System Users: people who are making use of a concrete System (which will inherit methods from this class).
  • Numerical integrators: TimeStepper and Integrator objects requires a System object whose State can be advanced through time.
  • System Developers: people who are writing concrete System classes.

Note that System Users may include people who are writing studies, reporters, modelers and so on as well as end users who are accessing the System directly.

Methods intended for System Users, numerical integrators, and a few bookkeeping methods are in the main System class, which is a SimTK Handle class. A Handle class consists only of a single pointer, which in this case points to a System::Guts object. The System::Guts class is abstract, and virtual methods to be implemented by System Developers in the concrete System are defined there, along with other utilities of use to the concrete System Developer but not to the end user or numerical integrator. The System::Guts class is declared in a separate header file, and only people who are writing their own System classes need look there.

Constructor & Destructor Documentation

SimTK::System::System ( )
inline

Default constructor creates an empty handle.

SimTK::System::System ( const System )

Copy constructor (untested).

SimTK::System::~System ( )

Destructor here will invoke the virtual destructor in the System::Guts class and thus clean up all unneeded detritus owned by this System.

SimTK::System::System ( System::Guts g)
inlineexplicit

Constructor for internal use only.

Member Function Documentation

System& SimTK::System::setUpDirection ( const CoordinateDirection up)

This is a hint to other software as to which way this System's designer considers to be "up".

For visualization, this is the best direction to use as the default up direction for the camera, and the opposite direction is likely to be a good direction in which to apply gravitational forces. The default up direction is +YAxis, which is the same as the OpenGL convention for the camera up direction. You can set this to any of the coordinate axes in the positive or negative direction. For example, use setUpDirection(ZAxis) for the "virtual world" convention where ground is the x-y plane and setUpDirection(-ZAxis) might be appropriate for an aviation convention in which Z is directed downward. A visualizer that is showing a ground plane should make the ground plane normal be this up direction.

See Also
setUseUniformBackground()
System& SimTK::System::setUseUniformBackground ( bool  useUniformBackground)

This is a hint to visualization software that this System is best viewed against a uniform background (e.g. all white) rather than against a ground plane.

A molecular system will typically set this flag so that the visualizer will not attempt to place the molecule on the ground. The default is to consider this system best viewed with a ground plane displayed, perpendicular to the "up" direction and located at a height of zero.

See Also
setUpDirection()
System& SimTK::System::setDefaultTimeScale ( Real  tc)

(Advanced) This is a hint used for some default behaviors, such as determining an initial step size for an integrator, or the default unit error for a constraint error derivative from the original constraint.

Most users can ignore this and just take the default.

This should be set to roughly the time scale at which you expect to see interesting things happen, that is the scale at which you might choose to view reporting output. An orbital simulation using seconds as time units might set this to 10 or 100s, for example, while a biomechanical simulation could use 0.1s. This will affect the time scale on which velocity constraints are stabilized, with longer time scales being more demanding since there is more time to drift. By default this is 0.1 time units, so 100ms for systems measuring time in seconds and 100fs for systems measuring time in ps.

System& SimTK::System::setDefaultLengthScale ( Real  lc)

(Advanced) This is a hint that can be used to get a sense of what a "unit length" looks like for this System in the units it uses.

Most users can ignore this and just take the default.

For example, a model of a small toy car expressed in MKS units might set this to 0.01 (1 cm). The default for this is 1 length unit, meaning 1 meter in MKS and 1 nm in MD units.

void SimTK::System::setHasTimeAdvancedEvents ( bool  )

This determines whether this System wants to be notified whenever time advances irreversibly.

If set true, time advancement is treated as an event, and the handleEvents() method is invoked with its cause argument set to indicate a time-advanced event occurred. By default, time advancement proceeds silently.

CoordinateDirection SimTK::System::getUpDirection ( ) const

Get the current setting of the "up" direction hint.

bool SimTK::System::getUseUniformBackground ( ) const

Get the current setting of the "use uniform background" visualization hint.

Real SimTK::System::getDefaultTimeScale ( ) const

Get the currently-set value for the default time scale, 0.1 time units if it has never been set.

Real SimTK::System::getDefaultLengthScale ( ) const

Get the currently-set value for the default length scale, 1 length unit if it has never been set.

bool SimTK::System::hasTimeAdvancedEvents ( ) const

Return the current value of the flag indicating whether this System wants an event generated whenever time advances irreversibly.

void SimTK::System::addEventHandler ( ScheduledEventHandler handler)
inline

Add a ScheduledEventHandler to this System, which takes over ownership of the event handler object.

void SimTK::System::addEventHandler ( TriggeredEventHandler handler)
inline

Add a TriggeredEventHandler to this System, which takes over ownership of the event handler object.

void SimTK::System::addEventReporter ( ScheduledEventReporter handler) const
inline

Add a ScheduledEventReporter to this System, which takes over ownership of the event reporter object.

void SimTK::System::addEventReporter ( TriggeredEventReporter handler) const
inline

Add a TriggeredEventReporter to this System, which takes over ownership of the event reporter object.

const State& SimTK::System::realizeTopology ( ) const

The following call must be made after any topological change has been made to this System, before the System can be used to perform any computations.

Perhaps surprisingly, the method is const. That's because the topology cannot be changed by this method. Various mutable "cache" entries will get calculated, including the default State, a reference to which is returned. The returned State has a Topology stage version number that matches the System, and will have already been realized through the Model Stage, using the defaults for the Model-stage variables, meaning that all later stage variables have been allocated and set to their default values as well. You can access this same default State again using getDefaultState(). If the current topology has already been realized, this call does nothing but return a reference to the already-built default State.

See Also
getDefaultState(), realizeModel(), realize()
const State& SimTK::System::getDefaultState ( ) const

This is available after realizeTopology(), and will throw an exception if realizeTopology() has not been called since the most recent topological change to this System.

This method returns the same reference returned by realizeTopology(). The State to which a reference is returned was created by the most recent realizeTopology() call. It has a Topology version number that matches the one currently in this System, and has already been realized through the Model Stage, using default values for all the Model-stage variables. All later-stage variables have been allocated and set to their default values. You can use this state directly to obtain information about the System in its default state or you can use this state to initialize other States to which you have write access. Those States are then suitable for further computation with this System.

See Also
realizeTopology()
State& SimTK::System::updDefaultState ( )

Don't use this; make a copy of the default state instead and modify your copy.

void SimTK::System::realizeModel ( State state) const

Realize the model to be used for subsequent computations with the given state.

This call is required if Model-stage variables are changed from their default values. The System topology must already have been realized (that is, realizeTopology() must have been called since the last topological change made to the System). Also, the supplied state must already have been initialized to work with this System either by copying the default state or some other State of this System. If it has already been realized to Stage::Model or higher, nothing happens here. Otherwise, all the state variables at Stage::Instance or higher are allocated or reallocated (if necessary), and reinitialized to their default values.

Note
Any state information at Stage::Instance or higher in the passed-in state is destroyed here. The number, types and memory locations of those state variables will change, so any existing references or pointers to them are invalid after this call. Note that this routine modifies its state argument, but makes no changes at all to the System itself (not even to the System's mutable cache) and is hence const.
See Also
realizeTopology(), realize()
void SimTK::System::realize ( const State state,
Stage  stage = Stage::HighestRuntime 
) const

Realize the given state to the indicated stage.

The passed-in State must have been initialized to work with this System, and it must already have been realized through Stage::Model, since this realize() method doesn't have write access to the State. If the state has already been realized to the requested stage or higher, nothing happens here. Otherwise, the state is realized one stage at a time until it reaches the requested stage.

See Also
realizeTopology(), realizeModel()
void SimTK::System::project ( State state,
Real  accuracy = -1 
) const

Apply prescribed motion and attempt to satisfy all position and velocity constraints by projecting the generalized coordinates q and generalized speeds u onto the position and velocity constraint manifolds.

This method expects to be given a state that is close to satisfying the constraints already. It will apply prescribed motion and then attempt a least-squares projection to satisfy constraints, working first on q's and then on u's. If you are trying to assemble a system from an arbitrary configuration, would like control over which states get modified, or are attempting to track markers, use an Assembler solver instead.

Parameters
[in,out]stateA State whose q's and u's are to be modified. Must already have been realized to at least Stage::Model. On return will be realized to Stage::Velocity.
[in]accuracyIf provided, sets the required accuracy to which the constraints must be satisfied; an exception will be thrown if that accuracy cannot be achieved. If accuracy is left off or is nonpositive, the accuracy used is the default as provided by the ProjectOptions class.

This method is provided for convenience and involves several calls to realize(), prescribeQ(), prescribeU(), projectQ(), and projectU(); call those directly if you want more control.

On return the given state's q and u variables will have been modified and the state will be realized through Velocity stage.

Theory

This solver projects the given state back on to the position or velocity constraint manifold, by the shortest path possible in the scaled norm defined in that state, and modifying only free (non-prescribed) variables. Constraint errors are scaled by their unit error weightings, then satisfied to a given accuracy using an RMS norm, or optionally using the stricter infinity norm. You can specify accuracy here explicitly; otherwise it uses the default as provided by the ProjectOptions class.

Position constraints are satisfied as follows:

    solve |Tp*perr(t;q+dq)|_n <= accuracy for dq (n==rms or inf)
    such that |Wq*dq|_2 is minimized

Here Tp=diag(1./unit_p) scales each position constraint error to a fraction of its unit error. Wq=N*Wu*pinv(N) weights dq to include both the "unit change" weightings on u and the artifactual configuration-dependent weightings on q generated by choice of orientation coordinates such as quaternions or rotation angles. (N as in qdot=N*u; pinv() is pseudoinverse.) We do not allow relative weighting on dq based on the current values of q; Simbody always solves q to absolute accuracy since arbitrary translations and rotations by 2pi should not affect physically-significant results.

Velocity constraints are satisfied as follows:

    solve |Tpv*pverr(t,q;u+du)|_n <= accuracy for du (n==rms or inf)
    such that |Eu*du|_2 is minimized

Here Tpv=diag(ts./unit_p 1./unit_v) where ts is the system's time scale used to scale the holonomic constraint's unit errors, and unit_v is the unit error for the holonomic constraints. The error weighting matrix Eu combines relative and absolute accuracy requirements as follows:

    Eu_i={ min(Wu_i, 1/u_i), relative accuracy OK for u_i
         {       Wu_i,       otherwise
See Also
Assembler, projectQ(), projectU(), prescribeQ(), prescribeU()
void SimTK::System::projectQ ( State state,
Real  accuracy = -1 
) const

Set prescribed positions q=q(t) and attempt to satisfy position constraints by modifying the remaining free q's.

The given state will be realized as necessary and on return will have been realized through Position stage. See project() for more information.

Parameters
[in,out]stateA State whose q's are to be modified. Must already have been realized to at least Stage::Model. On return will be realized to Stage::Position.
[in]accuracyIf provided, sets the required accuracy to which the constraints must be satisfied; an exception will be thrown if that accuracy cannot be achieved. If accuracy is left off or is nonpositive, the accuracy used is the default as provided by the ProjectOptions class.
See Also
project(), projectU(), ProjectOptions
void SimTK::System::projectU ( State state,
Real  accuracy = -1 
) const

Set prescribed velocities u=u(t,q) and attempt to satisfy velocity constraints by modifying the remaining free u's.

The given state will be realized as necessary and on return will have been realized through Velocity stage. Note that no q's are modified by this method; you should already have ensured that they satisfy position constraints. See project() for more information.

Parameters
[in,out]stateA State whose u's are to be modified. Must already have been realized to at least Stage::Model. On return will be realized to Stage::Velocity.
[in]accuracyIf provided, sets the required accuracy to which the constraints must be satisfied; an exception will be thrown if that accuracy cannot be achieved. If accuracy is left off or is nonpositive, the accuracy used is the default as provided by the ProjectOptions class.
See Also
project(), projectQ(), ProjectOptions
void SimTK::System::projectQ ( State state,
Vector qErrEst,
const ProjectOptions options,
ProjectResults results 
) const

(Advanced) Project free q's so that position constraints are satisfied and remove the corresponding error from the supplied error estimate.

This is primarily intended for use by numerical integration algorithms. You must already have set prescribed q's; this method will not modify them but may depend on their current values. State must be realized to Position stage on entry and will still be realized through Position stage on return.

If the norm of the position constraint errors is already less than or equal to the accuracy requested in options on entry, nothing will happen unless you have selected the "force projection" option. You can find out what actually happened by looking in the returned results.

Optionally, this method can also project out the constraint-normal component of the passed-in error estimate vector qErrEst. This is part of the integration of the continuous DAE system and thus should never require an integrator restart. Just pass a 0-length (default constructed) Vector in place of the error estimate if you don't want error projection.

Some options are available, see ProjectOptions for complete list. For example,

  • use infinity norm
  • local projection only
  • force projection
  • don't throw an exception
  • force full Newton

This method and projectU() are not for satisfying acceleration constraints, which does not involve modifications to state variables. Acceleration constraints are satisfied automatically when you realize a state to Acceleration stage using realize(state); the resulting udots will satisfy the acceleration constraints (if possible), even if the position and velocity constraints are not satisfied.

See the Theory section in the documentation for project() for more information.

See Also
ProjectOptions, ProjectResults, projectU(), project()
void SimTK::System::projectU ( State state,
Vector uErrEst,
const ProjectOptions options,
ProjectResults results 
) const

(Advanced) Project free u's so that velocity constraints are satisfied and remove the corresponding error from the supplied error estimate.

This is primarily intended for use by numerical integration algorithms. You must already have dealt with q's and already set prescribed u's; this method will not modify them but may depend on their current values. State must be realized to Velocity stage on entry and will still be realized through Velocity stage on return.

If the norm of the velocity constraint errors is already less than or equal to the accuracy requested in options on entry, nothing will happen unless you have selected the "force projection" option. You can find out what actually happened by looking in the returned results.

Optionally, this method can also project out the constraint-normal component of the passed-in error estimate vector uErrEst. This is part of the integration of the continuous DAE system and thus should never require an integrator restart. Just pass a 0-length (default constructed) Vector in place of the error estimate if you don't want error projection.

See the Theory section in the documentation for project(), and the comments for projectQ(), for more information.

See Also
ProjectOptions, ProjectResults, projectQ(), project()
void SimTK::System::prescribe ( State state) const
inline

Set values for prescribed positions q and velocities u.

Prescribed positions are functions of time q(t) and prescribed velocities are functions of time and position u(t,q). Both can also depend on earlier-stage discrete variables such as modeling and instance parameters.

Parameters
[in,out]stateThe State to be modified. Time and the values of non-prescribed q's are obtained from state and prescribed q's and u's are modified on return. The state will be realized as needed and on return will have been realized through Position stage. The prescribed velocities will have been set but not yet realized.

This is a convenience method that performs several realize() calls and uses prescribeQ() and prescribeU(). If you want more control, use those methods directly instead of this one.

See Also
prescribeQ(), prescribeU(), project(), realize()
bool SimTK::System::prescribeQ ( State state) const

This solver sets prescribed q's to their known values q=q(t) as a function of time (and possibly earlier-stage discrete state variables).

We expect the supplied State already to have been realized to Stage::Time. Note that the derivatives qdot of prescribed q's (which are of necessity also prescribed but are not themselves state variables) are set in the subsequent realize() call. For example, prescribeU sets the prescribed u's, then the next realize(Velocity) call will use them to calculate the prescribed qdots=N*u. realize(Dynamics) calculates known forces and the prescribed udoti=udoti(t,q,u). realize(Acceleration) calculates the remaining udots, lambdas, taus, and all the zdots.

Note that this method is not used to set prescribed udots, because those are not state variables. Instead, prescribed udots (which depend on time, positions, and velocities) are set as part of realize(Dynamics).

Returns
true if any change was made to state.
See Also
prescribe(), prescribeU()
bool SimTK::System::prescribeU ( State state) const

This solver sets prescribed u's to their known values u=u(t,q) as a function of time and position variables q (and possibly earlier-stage discrete state variables).

We expect the supplied State already to have been realized to Stage::Position. Note that the derivatives of prescribed variables (which are of necessity also prescribed but are not themselves state variables) are set in the subsequent realize() calls. For example, prescribeU() sets the prescribed u's, then the next realize(Velocity) call will use them to calculate the prescribed qdots=N*u. realize(Dynamics) calculates known forces and the prescribed udoti=udoti(t,q,u). realize(Acceleration) calculates the remaining udots, lambdas, taus, and all the zdots.

Returns
true if any change was made to state.
See Also
prescribe(), prescribeQ()
void SimTK::System::handleEvents ( State state,
Event::Cause  cause,
const Array_< EventId > &  eventIds,
const HandleEventsOptions options,
HandleEventsResults results 
) const

This solver handles a set of events which a TimeStepper has denoted as having occurred at the given time and state.

The event handler may make discontinuous changes in the State, in general both to discrete and continuous variables, but not to time or topological information. If changes are made to continuous variables, the handler is required to make sure the returned state satisfies the constraints to the accuracy level specified in options.

On return, the handleEvents() method will set the output variable results to indicate what happened. If any invoked handler determines that the occurrence of some event requires that the simulation be terminated, that information is returned in results and a well-behaved TimeStepper will stop when it sees that.

Simbody will automatically set a field in results that says how much of the state was changed by the handler so that the calling TimeStepper will be able to determine how much reinitialization is required.

See Also
HandleEventsOptions, HandleEventsResults
void SimTK::System::reportEvents ( const State state,
Event::Cause  cause,
const Array_< EventId > &  eventIds 
) const

This method is similar to handleEvents(), but does not allow the State to be modified.

It is used for scheduled events that were marked as being reports.

void SimTK::System::calcEventTriggerInfo ( const State state,
Array_< EventTriggerInfo > &  triggerInfo 
) const

This routine provides the Integrator with information it needs about the individual event trigger functions, such as which sign transitions are relevant and how tightly we need to localize.

This is considered Instance stage information so cannot change during a continuous integration interval (so an Integrator can process it upon restart(Instance)), however it can be updated whenever a discrete update is made to the State. A default implementation is provided which returns default EventTriggerInfo for each event trigger in state. The state must already be realized to Stage::Instance.

void SimTK::System::calcTimeOfNextScheduledEvent ( const State state,
Real &  tNextEvent,
Array_< EventId > &  eventIds,
bool  includeCurrentTime 
) const

This routine should be called to determine if and when there is an event scheduled to occur at a particular time.

This is a lot cheaper than making the Integrator hunt these down like ordinary state-dependent events. The returned time can be passed to the Integrator's stepping function as the advance time limit.

void SimTK::System::calcTimeOfNextScheduledReport ( const State state,
Real &  tNextEvent,
Array_< EventId > &  eventIds,
bool  includeCurrentTime 
) const

This routine is similar to calcTimeOfNextScheduledEvent(), but is used for "reporting events" which do not modify the state.

Events returned by this method should be handled by invoking reportEvents() instead of handleEvents().

void SimTK::System::relax ( State state,
Stage  stage,
Real  accuracy = -1 
) const

This optional method should modify fast variables at the given stage until they satisfy some relaxation criteria.

The criteria may involve anything in the State *except* fast variables at higher stages. Anything that can be calculated from the state is also fair game provided that those calculations do not depend on higher-stage fast variables. "Relaxation" criteria may require that fast variables satisfy some implicit or explicit algebraic conditions (constraints), or reach some minimization or maximization condition. A common criterion is that fast q's are moved to minimize potential energy; that can be achieved by driving the fast mobilities' lambdas (calculated joint torques) to zero since they are the potential energy gradient. That may require repeated realization to Acceleration stage.

Note that when q's are fast, the corresponding u's and udots are prescribed (to zero), they are not fast. And when u's are fast, their udots are also zero (their q's are regular integrated variables). Fast z's have zero zdots. Any other variables (that is, x-partition variables) can also be fast but don't have derivatives.

TODO: should take options and return results.

void SimTK::System::multiplyByN ( const State state,
const Vector u,
Vector dq 
) const

Calculate dq=N*u in O(n) time (very fast).

void SimTK::System::multiplyByNTranspose ( const State state,
const Vector fq,
Vector fu 
) const

Calculate fu=~N*fq in O(n) time (very fast).

void SimTK::System::multiplyByNPInv ( const State state,
const Vector dq,
Vector u 
) const

Calculate u=pinv(N)*dq in O(n) time (very fast).

void SimTK::System::multiplyByNPInvTranspose ( const State state,
const Vector fu,
Vector fq 
) const

Calculate fq=~pinv(N)*fu in O(n) time (very fast).

void SimTK::System::resetAllCountersToZero ( )

Zero out the statistics for this System.

Although the statistics are mutable, we only allow them to be reset by a caller who has write access since setting the stats to zero is associated with construction.

int SimTK::System::getNumRealizationsOfThisStage ( Stage  ) const

Whenever the system was realized from Stage-1 to the indicated Stage, this counter is bumped.

Note that a single call to realize() can cause several counters to get bumped.

int SimTK::System::getNumRealizeCalls ( ) const

Return the total number of calls to realizeTopology(), realizeModel(), or realize(), regardless of whether these routines actually did anything when called.

int SimTK::System::getNumPrescribeQCalls ( ) const

Return the total number of calls to the System's prescribeQ() method.

int SimTK::System::getNumPrescribeUCalls ( ) const

Return the total number of calls to the System's prescribeU() method.

int SimTK::System::getNumProjectQCalls ( ) const

Return the total number of calls to projectQ(), regardless of whether the call did anything.

int SimTK::System::getNumFailedProjectQCalls ( ) const

Return the total number of calls to projectQ() that failed.

int SimTK::System::getNumQProjections ( ) const

How many of the successful projectQ() calls actually did a constraint projection, rather than returning quickly?

int SimTK::System::getNumQErrorEstimateProjections ( ) const

How many of the projectQ() calls that did a constraint projection also projected an error estimate?

int SimTK::System::getNumProjectUCalls ( ) const

Return the total number of calls to projectU(), regardless of whether the call did anything.

int SimTK::System::getNumFailedProjectUCalls ( ) const

Return the total number of calls to projectU() that failed.

int SimTK::System::getNumUProjections ( ) const

How many of the successful projectU() calls actually did a constraint projection, rather than returning quickly?

int SimTK::System::getNumUErrorEstimateProjections ( ) const

How many of the projectU() calls that did a constraint projection also projected an error estimate?

int SimTK::System::getNumHandlerCallsThatChangedStage ( Stage  ) const

handleEvents() reports the lowest Stage it modified and we bump the counter for that Stage.

We also count reportEvents() calls here as having "changed" Stage::Report.

int SimTK::System::getNumHandleEventCalls ( ) const

This is the total number of calls to handleEvents() regardless of the outcome.

int SimTK::System::getNumReportEventCalls ( ) const

This is the total number of calls to reportEvents() regardless of the outcome.

System& SimTK::System::operator= ( const System )

Copy assignment (untested).

const String& SimTK::System::getName ( ) const

Return the name assigned to this System (arbitrary).

const String& SimTK::System::getVersion ( ) const

Return the version string assigned to this System (arbitrary).

SubsystemIndex SimTK::System::adoptSubsystem ( Subsystem child)

Take over ownership of the supplied subsystem and install it into the next free subsystem slot.

The new slot index is returned.

int SimTK::System::getNumSubsystems ( ) const

How may Subsystems are in here?

const Subsystem& SimTK::System::getSubsystem ( SubsystemIndex  ) const

Obtain read-only access to a particular subsystem by its index.

Subsystem& SimTK::System::updSubsystem ( SubsystemIndex  )

Obtain writable access to a particular subsystem by its index.

const DefaultSystemSubsystem& SimTK::System::getDefaultSubsystem ( ) const

Get read-only access to the default subsystem which is present in every system.

DefaultSystemSubsystem& SimTK::System::updDefaultSubsystem ( )

Get writable access to the default subsystem which is present in every system.

SimTK::System::operator const Subsystem & ( ) const
inline

Implicitly convert this System into a const Subsystem reference; this actually returns a reference to the DefaultSystemSubsystem contained in this System.

SimTK::System::operator Subsystem & ( )
inline

Implicitly convert this System into a writable Subsystem reference; this actually returns a reference to the DefaultSystemSubsystem contained in this System.

bool SimTK::System::systemTopologyHasBeenRealized ( ) const

(Advanced) You can check whether realizeTopology() has been called since the last topological change to this Syatem.

If you don't check and just plunge ahead you are likely to encounter an exception since very few things will work without topology having been realized.

StageVersion SimTK::System::getSystemTopologyCacheVersion ( ) const

(Advanced) Return the current version number of this system's Topology cache information.

This is a counter that is incremented each time the Topology is invalidated. Any State to be used with this System must have the same Topology version number as the System does. The version number is returned regardless of whether topology has been realized; you can check that with systemTopologyHasBeenRealized().

See Also
State::getSystemTopologyStageVersion()
void SimTK::System::setSystemTopologyCacheVersion ( StageVersion  topoVersion) const

(Really advanced) Set the current version number of this system's Topology cache information.

Don't use this method unless you really know what you're doing! This has no effect on realization status; if topology has not yet been realized this is the version number it will have after the next realizeTopology() call.

void SimTK::System::invalidateSystemTopologyCache ( ) const

(Advanced) Mark the Topology stage of this system and all its subsystems "not realized." This is normally handled automatically by whenever you make a Topology-stage change to any subsystem.

Occasionally you may want to force recomputation of the Topology-stage cache, for example during testing. After this call the method systemTopologyHasBeenRealized() will return false and you will not be able to call getDefaultState(). A subsequent call to realizeTopology() will invoke all the subsystems' realizeTopology() methods. The Topology stage version number will have changed, so all previously-created State objects will be invalid.

void SimTK::System::calcDecorativeGeometryAndAppend ( const State ,
Stage  ,
Array_< DecorativeGeometry > &   
) const

(Advanced) Generate all decorative geometry computable at a specific stage; this is useful for visualizers.

This will throw an exception if the state hasn't already been realized to the given stage. Note that the list is not inclusive – you have to request geometry from each stage to get all of it. This routine asks each subsystem in succession to generate its decorative geometry and append it to the end of the array. If the stage is Stage::Topology, realizeTopology() must already have been called but the State is ignored.

bool SimTK::System::isSameSystem ( const System otherSystem) const

There can be multiple handles referring to the same System::Guts object; they are considered to be the same System.

const Guts& SimTK::System::getSystemGuts ( ) const
inline

Obtain a const reference to the System::Guts object to which this handle refers.

You should then dynamic_cast the returned reference to a reference to your concrete Guts class.

Guts& SimTK::System::updSystemGuts ( )
inline

Obtain a writable reference to the System::Guts object to which this handle refers.

You should then dynamic_cast the returned reference to a reference to your concrete Guts class.

void SimTK::System::adoptSystemGuts ( System::Guts g)

Put new unowned Guts into this *empty* handle and take over ownership.

If this handle is already in use, or if Guts is already owned this routine will throw an exception.

bool SimTK::System::hasGuts ( ) const
inline

Return true if this System handle is not empty.

bool SimTK::System::isOwnerHandle ( ) const

Internal use only.

bool SimTK::System::isEmptyHandle ( ) const

Internal use only.

Friends And Related Function Documentation

friend class Guts
friend

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