API
4.5
For C++ developers
|
A class that manages the execution of a simulation. More...
Public Member Functions | |
Manager (Model &model) | |
Constructor that takes a model only and internally uses a SimTK::RungeKuttaMersonIntegrator with default settings (accuracy, constraint tolerance, etc.). More... | |
Manager (Model &model, const SimTK::State &state) | |
Convenience constructor for creating and initializing a Manager (by calling initialize()). More... | |
Manager () | |
(Deprecated) A Constructor that does not take a model or controllerSet. More... | |
Manager (const Manager &)=delete | |
void | operator= (const Manager &)=delete |
void | setSessionName (const std::string &name) |
void | setModel (Model &model) |
const std::string & | getSessionName () const |
const std::string & | toString () const |
void | setPerformAnalyses (bool performAnalyses) |
void | setWriteToStorage (bool writeToStorage) |
void | setUseSpecifiedDT (bool aTrueFalse) |
bool | getUseSpecifiedDT () const |
void | setUseConstantDT (bool aTrueFalse) |
bool | getUseConstantDT () const |
const Array< double > & | getDTArray () |
void | setDTArray (const SimTK::Vector_< double > &aDT, double aTI=0.0) |
double | getDTArrayDT (int aStep) |
void | printDTArray (const char *aFileName=NULL) |
const Array< double > & | getTimeArray () |
double | getTimeArrayTime (int aStep) |
int | getTimeArrayStep (double aTime) |
void | printTimeArray (const char *aFileName=NULL) |
void | resetTimeAndDTArrays (double aTime) |
double | getNextTimeArrayTime (double aTime) |
void | initialize (const SimTK::State &s) |
Initializes the Manager by creating and initializing the underlying SimTK::TimeStepper. More... | |
const SimTK::State & | integrate (double finalTime) |
Integrate the equations of motion for the specified model, given the current state (at which the integration will start) and a finalTime. More... | |
const SimTK::State & | getState () const |
Get the current State from the Integrator associated with this Manager. More... | |
double | getFixedStepSize (int tArrayStep) const |
bool | hasStateStorage () const |
void | setStateStorage (Storage &aStorage) |
Set the Storage object to be used for storing states. More... | |
Storage & | getStateStorage () const |
TimeSeriesTable | getStatesTable () const |
void | halt () |
void | clearHalt () |
bool | checkHalt () |
Configure the Integrator | |
| |
enum | IntegratorMethod { IntegratorMethod::ExplicitEuler = 0, IntegratorMethod::RungeKutta2 = 1, IntegratorMethod::RungeKutta3 = 2, IntegratorMethod::RungeKuttaFeldberg = 3, IntegratorMethod::RungeKuttaMerson = 4, IntegratorMethod::SemiExplicitEuler2 = 5, IntegratorMethod::Verlet = 6 } |
Supported integrator methods. More... | |
void | setIntegratorMethod (IntegratorMethod integMethod) |
Sets the integrator method used via IntegratorMethod enum. More... | |
SimTK::Integrator & | getIntegrator () const |
void | setIntegratorAccuracy (double accuracy) |
Sets the accuracy of the integrator. More... | |
void | setIntegratorMinimumStepSize (double hmin) |
Sets the minimum step size of the integrator. More... | |
void | setIntegratorMaximumStepSize (double hmax) |
Sets the maximum step size of the integrator. More... | |
void | setIntegratorInternalStepLimit (int nSteps) |
Sets the limit of steps the integrator can take per call of stepTo() . More... | |
A class that manages the execution of a simulation.
This class uses a SimTK::Integrator and SimTK::TimeStepper to perform the simulation. By default, a Runge-Kutta-Merson integrator is used, but can be changed by using setIntegratorMethod().
In order to prevent an inconsistency between the Integrator and TimeStepper, we only create a TimeStepper once, specifically when we call initialize(). To ensure this, the Manager will throw an exception if initialize() is called more than once. Note that editing the SimTK::State after calling initialize() will not affect the simulation.
Note that this interface means that you cannot "reinitialize" a Manager. If you make changes to the SimTK::State, a new Manager must be created before integrating again.
|
strong |
Supported integrator methods.
For MATLAB, int's must be used rather than enum's (see example in setIntegratorMethod(IntegratorMethod)).
OpenSim::Manager::Manager | ( | Model & | model | ) |
Constructor that takes a model only and internally uses a SimTK::RungeKuttaMersonIntegrator with default settings (accuracy, constraint tolerance, etc.).
OpenSim::Manager::Manager | ( | Model & | model, |
const SimTK::State & | state | ||
) |
Convenience constructor for creating and initializing a Manager (by calling initialize()).
Do not use this constructor if you intend to change integrator settings; changes to the integrator may not take effect after initializing.
OpenSim::Manager::Manager | ( | ) |
(Deprecated) A Constructor that does not take a model or controllerSet.
This constructor also does not set an integrator; you must call setIntegrator() on your own. You should use one of the other constructors.
|
delete |
bool OpenSim::Manager::checkHalt | ( | ) |
void OpenSim::Manager::clearHalt | ( | ) |
const Array<double>& OpenSim::Manager::getDTArray | ( | ) |
double OpenSim::Manager::getDTArrayDT | ( | int | aStep | ) |
double OpenSim::Manager::getFixedStepSize | ( | int | tArrayStep | ) | const |
SimTK::Integrator& OpenSim::Manager::getIntegrator | ( | ) | const |
double OpenSim::Manager::getNextTimeArrayTime | ( | double | aTime | ) |
const std::string& OpenSim::Manager::getSessionName | ( | ) | const |
const SimTK::State& OpenSim::Manager::getState | ( | ) | const |
Get the current State from the Integrator associated with this Manager.
TimeSeriesTable OpenSim::Manager::getStatesTable | ( | ) | const |
Storage& OpenSim::Manager::getStateStorage | ( | ) | const |
const Array<double>& OpenSim::Manager::getTimeArray | ( | ) |
int OpenSim::Manager::getTimeArrayStep | ( | double | aTime | ) |
double OpenSim::Manager::getTimeArrayTime | ( | int | aStep | ) |
bool OpenSim::Manager::getUseConstantDT | ( | ) | const |
bool OpenSim::Manager::getUseSpecifiedDT | ( | ) | const |
void OpenSim::Manager::halt | ( | ) |
bool OpenSim::Manager::hasStateStorage | ( | ) | const |
void OpenSim::Manager::initialize | ( | const SimTK::State & | s | ) |
Initializes the Manager by creating and initializing the underlying SimTK::TimeStepper.
This must be called before calling Manager::integrate() Subsequent changes to the State object passed in here will not affect the simulation. Calling this function multiple times with the same Manager will trigger an Exception.
Changes to the integrator (e.g., setIntegratorAccuracy()) after calling initialize() may not have any effect.
const SimTK::State& OpenSim::Manager::integrate | ( | double | finalTime | ) |
Integrate the equations of motion for the specified model, given the current state (at which the integration will start) and a finalTime.
You must call Manager::initialize() before calling this function.
If you must update states or controls in a loop, you must recreate the manager within the loop (such discontinuous changes are considered "events" and cannot be handled during integration of the otherwise continuous system). The proper way to handle this situation is to create a SimTK::EventHandler.
Example: Integrating from time = 1s to time = 2s
Example: Integrating from time = 1s to time = 2s using the convenience constructor
Example: Integrate from time = 0s to time = 10s, in 2s increments
Example: Integrate from time = 0s to time = 10s, updating the state (e.g., the model's first coordinate value) every 2s
|
delete |
void OpenSim::Manager::printDTArray | ( | const char * | aFileName = NULL | ) |
void OpenSim::Manager::printTimeArray | ( | const char * | aFileName = NULL | ) |
void OpenSim::Manager::resetTimeAndDTArrays | ( | double | aTime | ) |
void OpenSim::Manager::setDTArray | ( | const SimTK::Vector_< double > & | aDT, |
double | aTI = 0.0 |
||
) |
void OpenSim::Manager::setIntegratorAccuracy | ( | double | accuracy | ) |
Sets the accuracy of the integrator.
For more details, see SimTK::Integrator::setAccuracy(SimTK::Real)
.
void OpenSim::Manager::setIntegratorInternalStepLimit | ( | int | nSteps | ) |
Sets the limit of steps the integrator can take per call of stepTo()
.
Note that Manager::integrate() calls stepTo()
for each interval when a fixed step size is used. For more details, see SimTK::Integrator::setInternalStepLimit(int).
void OpenSim::Manager::setIntegratorMaximumStepSize | ( | double | hmax | ) |
Sets the maximum step size of the integrator.
For more details, see SimTK::Integrator::setMaximumStepSize(SimTK::Real)
.
void OpenSim::Manager::setIntegratorMethod | ( | IntegratorMethod | integMethod | ) |
Sets the integrator method used via IntegratorMethod enum.
The integrator will be set to its default options, even if the caller requests the same integrator method. Note that this function must be called before Manager::initialize()
.
C++ example
Python example
MATLAB example
void OpenSim::Manager::setIntegratorMinimumStepSize | ( | double | hmin | ) |
Sets the minimum step size of the integrator.
For more details, see SimTK::Integrator::setMinimumStepSize(SimTK::Real)
.
void OpenSim::Manager::setModel | ( | Model & | model | ) |
|
inline |
void OpenSim::Manager::setSessionName | ( | const std::string & | name | ) |
void OpenSim::Manager::setStateStorage | ( | Storage & | aStorage | ) |
void OpenSim::Manager::setUseConstantDT | ( | bool | aTrueFalse | ) |
void OpenSim::Manager::setUseSpecifiedDT | ( | bool | aTrueFalse | ) |
|
inline |
const std::string& OpenSim::Manager::toString | ( | ) | const |