Solving a Neural Lyapunov Problem

NeuralLyapunov.jl represents neural Lyapunov problems as systems of partial differential equations, using the ModelingToolkit.PDESystem type. Such a PDESystem can then be solved using a physics-informed neural network through NeuralPDE.jl.

Candidate Lyapunov functions will be trained within a box domain subset of the state space.

NeuralLyapunov.NeuralLyapunovPDESystemFunction
NeuralLyapunovPDESystem(dynamics::ODESystem, bounds, spec; <keyword_arguments>)
NeuralLyapunovPDESystem(dynamics::Function, lb, ub, spec; <keyword_arguments>)

Construct a ModelingToolkit.PDESystem representing the specified neural Lyapunov problem.

Positional Arguments

  • dynamics: the dynamical system being analyzed, represented as an ODESystem or the function f such that ẋ = f(x[, u], p, t); either way, the ODE should not depend on time and only t = 0.0 will be used. (For an example of when f would have a u argument, see add_policy_search.)
  • bounds: an array of domains, defining the training domain by bounding the states (and derivatives, when applicable) of dynamics; only used when dynamics isa ODESystem, otherwise use lb and ub.
  • lb and ub: the training domain will be $[lb_1, ub_1]×[lb_2, ub_2]×...$; not used when dynamics isa ODESystem, then use bounds.
  • spec: a NeuralLyapunovSpecification defining the Lyapunov function structure, as well as the minimization and decrease conditions.

Keyword Arguments

  • fixed_point: the equilibrium being analyzed; defaults to the origin.
  • p: the values of the parameters of the dynamical system being analyzed; defaults to SciMLBase.NullParameters(); not used when dynamics isa ODESystem, then use the default parameter values of dynamics.
  • state_syms: an array of the Symbol representing each state; not used when dynamics isa ODESystem, then the symbols from dynamics are used; if dynamics isa ODEFunction, symbols stored there are used, unless overridden here; if not provided here and cannot be inferred, [:state1, :state2, ...] will be used.
  • parameter_syms: an array of the Symbol representing each parameter; not used when dynamics isa ODESystem, then the symbols from dynamics are used; if dynamics isa ODEFunction, symbols stored there are used, unless overridden here; if not provided here and cannot be inferred, [:param1, :param2, ...] will be used.
  • policy_search::Bool: whether or not to include a loss term enforcing fixed_point to actually be a fixed point; defaults to false; only used when dynamics isa Function && !(dynamics isa ODEFunction); when dynamics isa ODEFunction, policy_search should not be supplied (as it must be false); when dynamics isa ODESystem, value inferred by the presence of unbound inputs.
  • name: the name of the constructed PDESystem
source
Warning

When using NeuralLyapunovPDESystem, the Lyapuonv function structure, minimization and decrease conditions, and dynamics will all be symbolically traced to generate the resulting PDESystem equations. In some cases tracing results in more efficient code, but in others it can result in inefficiencies or even errors.

If the generated PDESystem is then used with NeuralPDE.jl, that library's parser will convert the equations into Julia functions representing the loss, which presents another opportunity for unexpected errors.

Extracting the numerical Lyapunov function

We provide the following convenience function for generating the Lyapunov function after the parameters have been found. If the PDESystem was solved using NeuralPDE.jl and Optimization.jl, then the argument phi is a field of the output of NeuralPDE.discretize and the argument θ is res.u.depvar where res is the result of the optimization.

NeuralLyapunov.get_numerical_lyapunov_functionFunction
get_numerical_lyapunov_function(phi, θ, structure, dynamics, fixed_point;
                                <keyword_arguments>)

Combine Lyapunov function structure, dynamics, and neural network weights to generate Julia functions representing the Lyapunov function and its time derivative: $V(x), V̇(x)$.

These functions can operate on a state vector or columnwise on a matrix of state vectors.

Positional Arguments

  • phi: the neural network, represented as phi(x, θ) if the neural network has a single output, or a Vector of the same with one entry per neural network output.
  • θ: the parameters of the neural network; If the neural network has multiple outputs, θ[:φ1] should be the parameters of the first neural network output, θ[:φ2] the parameters of the second (if there are multiple), and so on. If the neural network has a single output, θ should be the parameters of the network.
  • structure: a NeuralLyapunovStructure representing the structure of the neural Lyapunov function.
  • dynamics: the system dynamics, as a function to be used in conjunction with structure.f_call.
  • fixed_point: the equilibrium point being analyzed by the Lyapunov function.

Keyword Arguments

  • p: parameters to be passed into dynamics; defaults to SciMLBase.NullParameters().
  • use_V̇_structure: when true, $V̇(x)$ is calculated using structure.V̇; when false, $V̇(x)$ is calculated using deriv as $\frac{∂}{∂t} V(x + t f(x))$ at $t = 0$; defaults to false, as it is more efficient in many cases.
  • deriv: a function for calculating derivatives; defaults to (and expects same arguments as) ForwardDiff.derivative; only used when use_V̇_structure is false.
  • jac: a function for calculating Jacobians; defaults to (and expects same arguments as) ForwardDiff.jacobian; only used when use_V̇_structure is true.
  • J_net: the Jacobian of the neural network, specified as a function J_net(phi, θ, state); if isnothing(J_net) (as is the default), J_net will be calculated using jac; only used when use_V̇_structure is true.
source