Drake
Analysis

Functions

Eigen::VectorXd SampleBasedLyapunovAnalysis (const System< double > &system, const Context< double > &context, const std::function< VectorX< AutoDiffXd >(const VectorX< AutoDiffXd > &state)> &basis_functions, const Eigen::Ref< const Eigen::MatrixXd > &state_samples, const Eigen::Ref< const Eigen::VectorXd > &V_zero_state)
 Sets up a linear program to search for the coefficients of a Lyapunov function that satisfies the Lyapunov conditions at a set of sample points. More...
 
double RandomSimulation (const SimulatorFactory &make_simulator, const ScalarSystemFunction &output, double final_time, RandomGenerator *generator)
 Run a deterministic simulation of a (stochastic) System using the generator to instantiate all "random" quantities. More...
 
std::vector< RandomSimulationResultMonteCarloSimulation (const SimulatorFactory &make_simulator, const ScalarSystemFunction &output, double final_time, int num_samples, RandomGenerator *generator=nullptr)
 Generate samples of a scalar random variable output by running many random simulations drawn from independent samples of the distributions governing the stochastic simulation. More...
 

Detailed Description

Function Documentation

◆ MonteCarloSimulation()

std::vector< RandomSimulationResult > MonteCarloSimulation ( const SimulatorFactory make_simulator,
const ScalarSystemFunction output,
double  final_time,
int  num_samples,
RandomGenerator generator = nullptr 
)

Generate samples of a scalar random variable output by running many random simulations drawn from independent samples of the distributions governing the stochastic simulation.

In pseudo-code, this algorithm implements:

for i=1:num_samples
const generator_snapshot = deepcopy(generator)
output = RandomSimulation(..., generator)
data(i) = std::pair(generator_snapshot, output)
return data
See also
RandomSimulation() for details about make_simulator, output, and final_time.
Parameters
num_samplesNumber of independent samples to draw from the distribution (and equivalently, the number of simulations to run).
generatorRandom number generator to be used to generate the random samples. If null, then a new RandomGenerator will be allocated and used internally (and repeated calls to this method will return identical results). To produce statistically "independent" samples on a future call to MonteCarloSimulation, you should make repeated uses of the same RandomGenerator object.
Returns
a list of RandomSimulationResult's.

◆ RandomSimulation()

double RandomSimulation ( const SimulatorFactory make_simulator,
const ScalarSystemFunction output,
double  final_time,
RandomGenerator generator 
)

Run a deterministic simulation of a (stochastic) System using the generator to instantiate all "random" quantities.

In pseudo-code, this algorithm implements:

simulator = make_simulator(generator)
simulator.get_system().SetRandomContext(generator)
simulator.StepTo(final_time)
return output(simulator.get_context())
Parameters
make_simulatorCallers to this method define a stochastic simulation by providing the make_simulator factory method to return a Simulator using the supplied RandomGenerator as the only source of randomness. This interface was designed to support cases where the System/Diagram is random (not only the Context), e.g. in the case where are variable number of objects are added to a multibody simulation.
outputThe scalar random variable output, denoted output, is defined as a function of the Simulator's System's Context, evaluated at the final_time. Monte-Carlo investigations that studying the details of an entire trajectory can still use this interface, e.g. by including a "runtime monitor" System that latches the worst-case deviation of a specification into it's Context to be queried at the final time.
final_timeThe time that each instance of the Simulator is stepped to. In many cases, this will be equivalent to the duration of the simulation, but it need not be because SetRandomContext() could initialize the time to a non-zero value, or an event could trigger premature termination of the simulation (see #4447).
generatorRandom number generator to be used to generate the random samples.
Returns
the output evaluated from the Context at final_time.

◆ SampleBasedLyapunovAnalysis()

Eigen::VectorXd SampleBasedLyapunovAnalysis ( const System< double > &  system,
const Context< double > &  context,
const std::function< VectorX< AutoDiffXd >(const VectorX< AutoDiffXd > &state)> &  basis_functions,
const Eigen::Ref< const Eigen::MatrixXd > &  state_samples,
const Eigen::Ref< const Eigen::VectorXd > &  V_zero_state 
)

Sets up a linear program to search for the coefficients of a Lyapunov function that satisfies the Lyapunov conditions at a set of sample points.

∀xᵢ, V(xᵢ) ≥ 0, ∀xᵢ, V̇(xᵢ) = ∂V/∂x f(xᵢ) ≤ 0. In order to provide boundary conditions to the problem, and improve numerical conditioning, we additionally impose the constraint V(x₀) = 0, and add an objective that pushes V̇(xᵢ) towards -1 (time-to-go): min ∑ |V̇(xᵢ) + 1|.

For background, and a description of this algorithm, see http://underactuated.csail.mit.edu/underactuated.html?chapter=lyapunov . It currently requires that the system to be optimized has only continuous state and it is assumed to be time invariant.

Parameters
systemto be verified. We currently require that the system has only continuous state, and it is assumed to be time invariant. Unlike many analysis algorithms, the system does not need to support conversion to other ScalarTypes (double is sufficient).
contextis used only to specify any parameters of the system, and to fix any input ports. The system/context must have all inputs assigned.
basis_functionsmust define an AutoDiffXd function that takes the state vector as an input argument and returns the vector of values of the basis functions at that state. The Lyapunov function will then have the form V(x) = ∑ pᵢ φᵢ(x), where p is the vector to be solved for and φ(x) is the vector of basis function evaluations returned by this function.
state_samplesis a list of sample states (one per column) at which to apply the optimization constraints and the objective.
V_zero_stateis a particular state, x₀, where we impose the condition: V(x₀) = 0.
Returns
params the VectorXd of parameters, p, that satisfies the Lyapunov conditions described above. The resulting Lyapunov function is V(x) = ∑ pᵢ φᵢ(x),