Basis

Many of the methods require the definition of a Basis on observables or functional forms. A Basis is generated via:

Basis(h, u, parameters = [], iv = nothing)

where h is either a vector of ModelingToolkit Operations for the valid functional forms or a general function with the typical DiffEq signature h(u,p,t), which can be used with an Operation or vector of Operation. u are the ModelingToolkit Variables used to describe the Basis, and parameters are the optional ModelingToolkit Variables used to describe the parameters in the basis elements. iv represents the independent variable of the system - the time.

DataDrivenDiffEq.BasisType
Basis(f, u; p, iv, eval_expression)

A basis over the variables u with parameters p and independent variable iv. f can either be a Julia function which is able to use ModelingToolkit variables or a vector of Operation. It can be called with the typical DiffEq signature, meaning out of place with f(u,p,t) or in place with f(du, u, p, t).

Example

using ModelingToolkit
using DataDrivenDiffEq

@parameters w[1:2] t
@variables u[1:2]

Ψ = Basis([u; sin.(w.*u)], u, parameters = p, iv = t)

Note

The keyword argument eval_expression controls the function creation behavior. eval_expression=true means that eval is used, so normal world-age behavior applies (i.e. the functions cannot be called from the function that generates them). If eval_expression=false, then construction via GeneralizedGenerated.jl is utilized to allow for same world-age evaluation. However, this can cause Julia to segfault on sufficiently large basis functions. By default eval_expression=false.

source

Example

We start by crearting some Variables and Parameters using ModelingToolkit.

using LinearAlgebra
using DataDrivenDiffEq
using Plots
using ModelingToolkit

@variables u[1:3]
@parameters w[1:2]
(ModelingToolkit.Operation[w₁, w₂],)

To define a basis, simply write down the equations you want to be included as a Vector{Operation}. Possible used parameters have to be given to the constructor.

h = [u[1]; u[2]; cos(w[1]*u[2]+w[2]*u[3])]
b = Basis(h, u, parameters = w)
3 dimensional basis in ["u₁", "u₂", "u₃"]

What can a Basis do? Can it do stuff? Let's find out!

Basis are callable with the signature of functions to be used in DifferentialEquations. So, the function value at a single point looks like:

x = b([1;2;3])
3-element Array{ModelingToolkit.Expression,1}:
 ModelingToolkit.Constant(1)
 ModelingToolkit.Constant(2)
              cos(2w₁ + 3w₂)

Or, in place

dx = similar(x)
b(dx, [1;2;3])

Notice that since we did not use any numerical values for the parameters, the basis uses the symbolic values in the result.

To use numerical values, simply pass this on in the function call. Here, we evaluate over a trajectory with two parameters and 40 timestamps.

X = randn(3, 40)
Y = b(X, [2;4], 0:39)

Suppose we want to add another equation, say sin(u[1]). A Basis behaves like an array, so we can simply

push!(b, sin(u[1]))
size(b)
(4,)

To ensure that a basis is well-behaved, functions already present are not included again.

push!(b, sin(u[1]))
size(b)
(4,)

We can also define functions of time and add them

t = independent_variable(b)
push!(b, cos(t*π))
println(b)
5 dimensional basis in ["u₁", "u₂", "u₃"]
f_1 = u₁
f_2 = u₂
f_3 = cos(u₂ * w₁ + u₃ * w₂)
f_4 = sin(u₁)
f_5 = cos(πt)

Additionally, we can iterate over a Basis using [eq for eq in basis] or index specific equations, like basis[2].

We can also chain Basis via just using it in the constructor

@variables x[1:2]
y = [sin(x[1]); cos(x[1]); x[2]]
t = independent_variable(b)
b2 = Basis(b(y, parameters(b), t), x, parameters = w, iv = t)
println(b2)
5 dimensional basis in ["x₁", "x₂"]
f_1 = sin(x₁)
f_2 = cos(x₁)
f_3 = cos(cos(x₁) * w₁ + w₂ * x₂)
f_4 = sin(sin(x₁))
f_5 = cos(π * getindex(t, 1))

You can also use merge to create the union of two Basis:

b3 = merge(b, b2)
println(b3)
10 dimensional basis in ["u₁", "u₂", "u₃", "x₁", "x₂"]
f_1 = u₁
f_2 = u₂
f_3 = cos(u₂ * w₁ + u₃ * w₂)
f_4 = sin(u₁)
f_5 = cos(πt)
f_6 = sin(x₁)
f_7 = cos(x₁)
f_8 = cos(cos(x₁) * w₁ + w₂ * x₂)
f_9 = sin(sin(x₁))
f_10 = cos(π * getindex(t, 1))

which combines all the used variables and parameters ( and assumes the same independent_variable ):

variables(b)
3-element Array{ModelingToolkit.Operation,1}:
 u₁
 u₂
 u₃
parameters(b)
2-element Array{ModelingToolkit.Operation,1}:
 w₁
 w₂

If you have a function already defined as pure code, you can use this also to create a Basis. Only the signature has to be consistent, so use f(u,p,t).

f(u, p, t) = [u[1]; u[2]; cos(p[1]*u[2]+p[2]*u[3])]
b_f = Basis(f, u, parameters = w)
println(b_f)
3 dimensional basis in ["u₁", "u₂", "u₃"]
du₁ = u₁
du₂ = u₂
du₃ = cos(u₂ * w₁ + u₃ * w₂)

This works for every function defined over Operations. So to create a Basis from a Flux model, simply extend the activations used:

using Flux
NNlib.σ(x::Operation) = 1 / (1+exp(-x))

c = Chain(Dense(3,2,σ), Dense(2, 1, σ))
ps, re = Flux.destructure(c)

@parameters p[1:length(ps)]

g(u, p, t) = re(p)(u)
b = Basis(g, u, parameters = p)

Functions

DataDrivenDiffEq.jacobianFunction
jacobian(basis)

Returns a function representing the jacobian matrix / gradient of the `Basis` with respect to the
dependent variables as a function with the common signature `f(u,p,t)` for out of place and `f(du, u, p, t)` for in place computation.
source
Base.push!Function
push!(basis, op)

Push the operation(s) in `op` into the basis and update all internal fields accordingly.
`op` can either be a single `Operation` or an Array of `Operation`s.
source
Base.deleteat!Function
deleteat!(basis, inds)

Delete the entries specified by `inds` and update the `Basis` accordingly.
source
Base.mergeFunction
merge(basis_a, basis_b)

Return a new `Basis`, which is defined via the union of both bases.
source
Base.merge!Function
merge!(basis_a, basis_b)

Updates the `Basis` to include the union of both bases.
source