Drake
Drake C++ Documentation
Analysis

Detailed Description

Functions

template<typename T >
std::unique_ptr< LinearSystem< T > > DiscreteTimeApproximation (const LinearSystem< T > &linear_system, double time_period)
 Converts a continuous-time linear system to a discrete-time linear system using the zero-order hold method. More...
 
template<typename T >
std::unique_ptr< AffineSystem< T > > DiscreteTimeApproximation (const AffineSystem< T > &affine_system, double time_period)
 Converts a continuous-time affine system to a discrete-time affine system using the zero-order hold method. More...
 
template<typename T >
std::unique_ptr< System< T > > DiscreteTimeApproximation (std::shared_ptr< const System< T >> system, double time_period, const SimulatorConfig &integrator_config=SimulatorConfig())
 Converts a general continuous-time system \( \dot{x} = f(t,x(t),u(t)) \) to a discrete-time system with zero-order hold on the input. More...
 
template<typename T >
std::unique_ptr< System< T > > DiscreteTimeApproximation (const System< T > &system, double time_period, const SimulatorConfig &integrator_config=SimulatorConfig())
 Same as above, without claiming ownership of system. More...
 
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 RandomSimulatorFactory &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 RandomSimulatorFactory &make_simulator, const ScalarSystemFunction &output, double final_time, int num_samples, RandomGenerator *generator=nullptr, Parallelism parallelism=false)
 Generates samples of a scalar random variable output by running many random simulations drawn from independent samples of the distributions governing the stochastic simulation. More...
 
symbolic::Expression RegionOfAttraction (const System< double > &system, const Context< double > &context, const RegionOfAttractionOptions &options=RegionOfAttractionOptions())
 Estimates the region of attraction of the time-invariant system at the fixed point defined by context. More...
 

Function Documentation

◆ DiscreteTimeApproximation() [1/4]

std::unique_ptr<LinearSystem<T> > drake::systems::DiscreteTimeApproximation ( const LinearSystem< T > &  linear_system,
double  time_period 
)

Converts a continuous-time linear system to a discrete-time linear system using the zero-order hold method.

Parameters
linear_systemThe continuous-time LinearSystem.
time_periodThe discrete time period.
Returns
A discrete-time LinearSystem.
Exceptions
ifthe linear_system is not continuous or time_period <= 0.
Template Parameters
TThe scalar type, which must be one of the default scalars.

◆ DiscreteTimeApproximation() [2/4]

std::unique_ptr<AffineSystem<T> > drake::systems::DiscreteTimeApproximation ( const AffineSystem< T > &  affine_system,
double  time_period 
)

Converts a continuous-time affine system to a discrete-time affine system using the zero-order hold method.

Parameters
affine_systemThe continuous-time AffineSystem.
time_periodThe discrete time period.
Returns
A discrete-time AffineSystem.
Exceptions
ifthe affine_system is not continuous or time_period <= 0.
Template Parameters
TThe scalar type, which must be one of the default scalars.

◆ DiscreteTimeApproximation() [3/4]

std::unique_ptr<System<T> > drake::systems::DiscreteTimeApproximation ( std::shared_ptr< const System< T >>  system,
double  time_period,
const SimulatorConfig integrator_config = SimulatorConfig() 
)

Converts a general continuous-time system \( \dot{x} = f(t,x(t),u(t)) \) to a discrete-time system with zero-order hold on the input.

The approximate discrete-time dynamics is given by \( x[n+1] = f_d(n,x[n],u[n]) = x[n] + \int_{t[n]}^{t[n+1]} f(t,x(t),u[n]) \, dt \), where the integration is performed using numerical integration via an IntegratorBase.

Parameters
systemThe continuous-time System.
time_periodThe discrete time period.
integrator_configUse this parameter to configure the integrator (e.g. choose non-default integration_scheme, max_step_size, use_error_control, or accuracy).
Returns
A discrete-time System.
Exceptions
ifthe system is not continuous or time_period <= 0.
std::exceptionif the integration scheme does not support the scalar type T.
Template Parameters
TThe scalar type, which must be one of the default scalars.

◆ DiscreteTimeApproximation() [4/4]

std::unique_ptr<System<T> > drake::systems::DiscreteTimeApproximation ( const System< T > &  system,
double  time_period,
const SimulatorConfig integrator_config = SimulatorConfig() 
)

Same as above, without claiming ownership of system.

Warning
The system reference must remain valid for the lifetime of the returned system.
Template Parameters
TThe scalar type, which must be one of the default scalars.

◆ MonteCarloSimulation()

std::vector<RandomSimulationResult> drake::systems::analysis::MonteCarloSimulation ( const RandomSimulatorFactory make_simulator,
const ScalarSystemFunction output,
double  final_time,
int  num_samples,
RandomGenerator generator = nullptr,
Parallelism  parallelism = false 
)

Generates 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.
parallelismSpecify number of parallel executions to use while performing num_samples simulations. The default value (false) specifies that simulations should be executed in serial. To use the concurrency available on your hardware, specify either Parallellism::Max() or its terse abbreviation true.
Returns
a list of RandomSimulationResult's.

Thread safety when parallel execution is specified:

  • make_simulator and generator are only accessed from the main thread.
  • Each simulator created by make_simulator and its context are only accessed from within a single worker thread; however, any resource shared between these simulators must be safe for concurrent use.
  • output is called from within worker threads performing simulation with the simulator and context belonging to each worker thread. It must be safe to make concurrent calls to output (i.e. any mutable state inside the function must be safe for concurrent use).

◆ RandomSimulation()

double drake::systems::analysis::RandomSimulation ( const RandomSimulatorFactory 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.AdvanceTo(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.

◆ RegionOfAttraction()

symbolic::Expression drake::systems::analysis::RegionOfAttraction ( const System< double > &  system,
const Context< double > &  context,
const RegionOfAttractionOptions options = RegionOfAttractionOptions() 
)

Estimates the region of attraction of the time-invariant system at the fixed point defined by context.

This implementation only searches for the largest level set of the lyapunov_candidate function from options (or a candidate obtained from solving the Lyapunov equation on the linearization).

Parameters
systema time-invariant continuous-time System that supports scalar-type conversion to symbolic::Expression. The dynamics of the system must be polynomial.
contexta Context that defines the parameters of the system and the fixed-point about which we are analyzing the regional stability.
optionsprovides a variety of configuration options.
See also
RegionOfAttractionOptions.
Returns
a symbolic::Expression representing a Lyapunov function using the symbolic Variables named x0, x1..., where the order matches the continuous state vector in the context, or the vector state_variables passed in through the options structure (if it is non-empty). The level set {x | V(x)<=1} containing the fixed-point in context represents the region of attraction.
Precondition
For the given system and context, any required input ports on system must be "defined", i.e., connected to other systems in a larger diagram or holding fixed values; see System::FixInputPortsFrom for possible caveats. Analyzing a closed-loop system would typically be accomplished by having both the plant and the controller in a diagram (which then has no input ports), and passing the diagram into this method as system.

Note: There are more numerical recipes for region of attraction analysis that could extend the current implementation. Do report an issue if you discover a system for which this code does not perform well.

◆ SampleBasedLyapunovAnalysis()

Eigen::VectorXd drake::systems::analysis::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),