API  4.5
For C++ developers
OpenSim::Manager Class Reference

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...
 
StoragegetStateStorage () const
 
TimeSeriesTable getStatesTable () const
 
void halt ()
 
void clearHalt ()
 
bool checkHalt ()
 

Configure the Integrator

Note
Call these functions before calling Manager::initialize().
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...
 

Detailed Description

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.

Examples:
exampleCustomImplicitAuxiliaryDynamics.cpp.

Member Enumeration Documentation

◆ IntegratorMethod

Supported integrator methods.

For MATLAB, int's must be used rather than enum's (see example in setIntegratorMethod(IntegratorMethod)).

Enumerator
ExplicitEuler 

0 : For details, see SimTK::ExplicitEulerIntegrator.

RungeKutta2 

1 : For details, see SimTK::RungeKutta2Integrator.

RungeKutta3 

2 : For details, see SimTK::RungeKutta3Integrator.

RungeKuttaFeldberg 

3 : For details, see SimTK::RungeKuttaFeldbergIntegrator.

RungeKuttaMerson 

4 : For details, see SimTK::RungeKuttaMersonIntegrator.

SemiExplicitEuler2 

5 : For details, see SimTK::SemiExplicitEuler2Integrator.

Verlet 

6 : For details, see SimTK::VerletIntegrator.

Constructor & Destructor Documentation

◆ Manager() [1/4]

OpenSim::Manager::Manager ( Model model)

Constructor that takes a model only and internally uses a SimTK::RungeKuttaMersonIntegrator with default settings (accuracy, constraint tolerance, etc.).

◆ Manager() [2/4]

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.

◆ Manager() [3/4]

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.

◆ Manager() [4/4]

OpenSim::Manager::Manager ( const Manager )
delete

Member Function Documentation

◆ checkHalt()

bool OpenSim::Manager::checkHalt ( )

◆ clearHalt()

void OpenSim::Manager::clearHalt ( )

◆ getDTArray()

const Array<double>& OpenSim::Manager::getDTArray ( )

◆ getDTArrayDT()

double OpenSim::Manager::getDTArrayDT ( int  aStep)

◆ getFixedStepSize()

double OpenSim::Manager::getFixedStepSize ( int  tArrayStep) const

◆ getIntegrator()

SimTK::Integrator& OpenSim::Manager::getIntegrator ( ) const

◆ getNextTimeArrayTime()

double OpenSim::Manager::getNextTimeArrayTime ( double  aTime)

◆ getSessionName()

const std::string& OpenSim::Manager::getSessionName ( ) const

◆ getState()

const SimTK::State& OpenSim::Manager::getState ( ) const

Get the current State from the Integrator associated with this Manager.

◆ getStatesTable()

TimeSeriesTable OpenSim::Manager::getStatesTable ( ) const

◆ getStateStorage()

Storage& OpenSim::Manager::getStateStorage ( ) const

◆ getTimeArray()

const Array<double>& OpenSim::Manager::getTimeArray ( )

◆ getTimeArrayStep()

int OpenSim::Manager::getTimeArrayStep ( double  aTime)

◆ getTimeArrayTime()

double OpenSim::Manager::getTimeArrayTime ( int  aStep)

◆ getUseConstantDT()

bool OpenSim::Manager::getUseConstantDT ( ) const

◆ getUseSpecifiedDT()

bool OpenSim::Manager::getUseSpecifiedDT ( ) const

◆ halt()

void OpenSim::Manager::halt ( )

◆ hasStateStorage()

bool OpenSim::Manager::hasStateStorage ( ) const

◆ initialize()

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.

◆ integrate()

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

SimTK::State state = model.initSystem();
Manager manager(model);
state.setTime(1.0);
manager.initialize(state);
state = manager.integrate(2.0);

Example: Integrating from time = 1s to time = 2s using the convenience constructor

SimTK::State state = model.initSystem();
state.setTime(1.0);
Manager manager(model, state);
state = manager.integrate(2.0);

Example: Integrate from time = 0s to time = 10s, in 2s increments

dTime = 2.0;
finalTime = 10.0;
int n = int(round(finalTime/dTime));
state.setTime(0.0);
manager.initialize(state);
for (int i = 1; i <= n; ++i) {
state = manager.integrate(i*dTime);
}

Example: Integrate from time = 0s to time = 10s, updating the state (e.g., the model's first coordinate value) every 2s

dTime = 2.0;
finalTime = 10.0;
int n = int(round(finalTime/dTime));
state.setTime(0.0);
for (int i = 0; i < n; ++i) {
model.getCoordinateSet().get(0).setValue(state, 0.1*i);
Manager manager(model);
state.setTime(i*dTime);
manager.initialize(state);
state = manager.integrate((i+1)*dTime);
}

◆ operator=()

void OpenSim::Manager::operator= ( const Manager )
delete

◆ printDTArray()

void OpenSim::Manager::printDTArray ( const char *  aFileName = NULL)

◆ printTimeArray()

void OpenSim::Manager::printTimeArray ( const char *  aFileName = NULL)

◆ resetTimeAndDTArrays()

void OpenSim::Manager::resetTimeAndDTArrays ( double  aTime)

◆ setDTArray()

void OpenSim::Manager::setDTArray ( const SimTK::Vector_< double > &  aDT,
double  aTI = 0.0 
)

◆ setIntegratorAccuracy()

void OpenSim::Manager::setIntegratorAccuracy ( double  accuracy)

Sets the accuracy of the integrator.

For more details, see SimTK::Integrator::setAccuracy(SimTK::Real).

◆ setIntegratorInternalStepLimit()

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).

◆ setIntegratorMaximumStepSize()

void OpenSim::Manager::setIntegratorMaximumStepSize ( double  hmax)

Sets the maximum step size of the integrator.

For more details, see SimTK::Integrator::setMaximumStepSize(SimTK::Real).

◆ setIntegratorMethod()

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

auto manager = Manager(model);
manager.setIntegratorMethod(Manager::IntegratorMethod::SemiExplicitEuler2);

Python example

import opensim
manager = opensim.Manager(model)
manager.setIntegratorMethod(opensim.Manager.IntegratorMethod_SemiExplicitEuler2)

MATLAB example

import org.opensim.modeling.*
manager = Manager(model);
manager.setIntegratorMethod(5);

◆ setIntegratorMinimumStepSize()

void OpenSim::Manager::setIntegratorMinimumStepSize ( double  hmin)

Sets the minimum step size of the integrator.

For more details, see SimTK::Integrator::setMinimumStepSize(SimTK::Real).

◆ setModel()

void OpenSim::Manager::setModel ( Model model)

◆ setPerformAnalyses()

void OpenSim::Manager::setPerformAnalyses ( bool  performAnalyses)
inline

◆ setSessionName()

void OpenSim::Manager::setSessionName ( const std::string &  name)

◆ setStateStorage()

void OpenSim::Manager::setStateStorage ( Storage aStorage)

Set the Storage object to be used for storing states.

The Manager takes ownership of the passed-in Storage.

◆ setUseConstantDT()

void OpenSim::Manager::setUseConstantDT ( bool  aTrueFalse)

◆ setUseSpecifiedDT()

void OpenSim::Manager::setUseSpecifiedDT ( bool  aTrueFalse)

◆ setWriteToStorage()

void OpenSim::Manager::setWriteToStorage ( bool  writeToStorage)
inline

◆ toString()

const std::string& OpenSim::Manager::toString ( ) const

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