Algorithms

Available Algorithms

This section of the manual describes the algorithms that are available in this package.

Continuous-time reachability

Decomposition-based approach

BFFPSV18 <: ContinuousPost

Implementation of the reachability algorithm using decompositions from [1].

This algorithm applies to continuous affine ordinary differential equations of the form

    x'(t) = Ax(t) + u(t),\qquad x(0) ∈ X0, u(t) ∈ U ∀ t ∈ [0, T],    (1)

where T is the (finite) time horizon. The result of the algorithm is a sequence of sets Xk whose union overapproximates the exact flowpipe of (1).

Fields

  • options – an Options structure that holds the algorithm-specific options

The following options are available:

option :discretization of type String has default value 'forward'; model for bloating/continuous time analysis
option :algorithm of type String has default value 'explicit'; algorithm backend
option :δ of type Float64 with alias :sampling_time has default value '0.01'; time step
option :vars of type AbstractArray{Int64,1} has default value 'Int64[]'; variables of interest; default: all variables
option :partition of type AbstractArray{#s117,1} where #s117<:AbstractArray{Int64,1} has default value 'Array{Int64,1}[[]]'; block partition; a block is represented by a vector containing its indices
option :sih_method of type String has default value 'concrete'; method to compute the symmetric interval hull in discretization
option :exp_method of type String has default value 'base'; method to compute the matrix exponential
option :assume_sparse of type Bool has default value 'false'; use an analysis for sparse discretized matrices?
option :lazy_inputs_interval of type Union{Int64, Function} has default value 'getfield(Reachability.ReachSets, Symbol("##33#34"))()'; length of interval in which the inputs are handled as a lazy set (``-1`` for 'never'); may generally also be a predicate over indices; the default corresponds to ``-1``
option :block_options of type Any has default value 'nothing'; short hand to set ':block_options_init' and ':block_options_iter'
option :block_options_init of type Any has default value 'nothing'; option for the approximation of the initial states (during decomposition)
option :block_options_iter of type Any has default value 'nothing'; option for the approximation of the states ``X_k``, ``k>0``
option :assume_homogeneous of type Bool has default value 'false'; ignore dynamic inputs during the analysis?
option :eager_checking of type Bool has default value 'true'; terminate as soon as property violation was detected?

See the Notes below for examples and an explanation of the options.

Algorithm

We refer to [1] for technical details.

Notes

The compositional approach offers a tradeoff between accuracy and performance, since one can choose different partitions and different set representations.

These notes illustrate the options that can be passed to the solver. We begin by explaining the default options used by the algorithm. Then we show how to specify the computation of a subset of the flowpipe, which helps to improve performance, specially for large systems.

We continue to show some lower-level options such as how to specify a custom block decomposition (the set types used and the sizes of the blocks).

We conclude showing that there are two different modes of operation: flowpipe computation and safety property checking, and examples of use.

Running example and default options

By default, BFFPSV18 uses a uniform partition of dimension one over all variables.

For example, consider the two-dimensional example from [2]:

julia> using Reachability, MathematicalSystems

julia> A = [-1. -4   0  0   0;
             4  -1   0  0   0;
             0   0  -3  1   0;
             0   0  -1  3   0;
             0   0   0  0  -2];

julia> X₀ = BallInf(ones(5), 0.1);

julia> problem = InitialValueProblem(LinearContinuousSystem(A), X₀);

julia> sol = solve(problem, Options(:T=>5.0), op=BFFPSV18(Options(:δ=>0.04)));

Let's check that the dimension of the first set computed is 5:

julia> X₁ = set(first(first(sol.flowpipes).reachsets));

julia> dim(X₁)
5

julia> sol.options[:partition]
5-element Array{Array{Int64,1},1}:
 [1]
 [2]
 [3]
 [4]
 [5]

Here, [1] corresponds to a block of size one for variable 5, etc., until variable 5. Let's check the set type used is Interval:

julia> all(Xi isa Interval for Xi in array(X₁))
true

Specifying a subset of all variables

Suppose that we are only interested in variables 1 and 2. Then we can specify it using the :vars option:

julia> sol = solve(problem, Options(:T=>5.0),
                   op=BFFPSV18(Options(:δ=>0.04, :vars=>[1, 2])));

julia> X₁ = set(first(first(sol.flowpipes).reachsets));

julia> dim(X₁)
2

The sets used are Intervals:

julia> all(Xi isa Interval for Xi in array(X₁))
true

Specifying the set type

We can as well specify using other set types for the overapproximation of the flowpipe, such as hyperrectangles. It can be specified using the :overapproximation option:

julia> sol = solve(problem, Options(:T=>5.0),
                   op=BFFPSV18(Options(:δ=>0.04, :vars=>[1, 2],
                                       :block_options=>Hyperrectangle)));

julia> X₁ = set(first(first(sol.flowpipes).reachsets));

julia> all(Xi isa Hyperrectangle for Xi in array(X₁))
true

To increase the accuracy, the overapproximation option offers other possibilities. For example, one can use template directions instead of a fixed set type. Below we use octagonal directions:

julia> sol = solve(problem, Options(:T=>5.0),
                   op=BFFPSV18(Options(:δ=>0.04, :vars=>[1, 2],
                                       :block_options=>OctDirections)));

julia> X₁ = set(first(first(sol.flowpipes).reachsets));

julia> all(Xi isa HPolytope for Xi in array(X₁))
true

References

[1] Reach Set Approximation through Decomposition with Low-dimensional Sets and High-dimensional Matrices. S. Bogomolov, M. Forets, G. Frehse, A. Podelski, C. Schilling, F. Viry. HSCC '18 Proceedings of the 21st International Conference on Hybrid Systems: Computation and Control (part of CPS Week).

[2] Althoff, Matthias. Reachability analysis and its application to the safety assessment of autonomous cars. Diss. Technische Universität München, 2010.

source