Nonlinear Time Continuous System

Similarly, we can use the Extended Dynamic Mode Decomposition via a nonlinear Basis of observables. Here, we will look at a rather famous example with a finite dimensional solution.

using DataDrivenDiffEq
using LinearAlgebra
using OrdinaryDiffEq
using DataDrivenDMD
using Plots

function slow_manifold(du, u, p, t)
    du[1] = p[1] * u[1]
    du[2] = p[2] * (u[2] - u[1]^2)
end

u0 = [3.0; -2.0]
tspan = (0.0, 5.0)
p = [-0.8; -0.7]

problem = ODEProblem{true, SciMLBase.NoSpecialize}(slow_manifold, u0, tspan, p)
solution = solve(problem, Tsit5(), saveat = 0.1);
plot(solution)
Example block output

Since we are dealing with a continuous system in time, we define the associated DataDrivenProblem accordingly using the measured states X, their derivatives DX and the time t.

prob = DataDrivenProblem(solution)
plot(prob)
Example block output

Additionally, we need to define the Basis for our lifting, before we solve the problem in the lifted space.

@parameters t
@variables u(t)[1:2]
Ψ = Basis([u; u[1]^2], u, independent_variable = t)
res = solve(prob, Ψ, DMDPINV(), digits = 2)
"DataDrivenSolution{Float64}" with 2 equations and 3 parameters.
Returncode: Success
Residual sum of squares: 9.030744649079049e-19

We can also use different metrics on the DataDrivenSolution like the aic

aic(res)
-4703.692384728004

The aicc

aicc(res)
-4703.44748676882

The bic

bic(res)
-4695.817466288151

The loglikelihood

loglikelihood(res)
2354.846192364002

And the number of parameters

dof(res)
3

Let's have a closer look at the Basis

basis = get_basis(res)
Model ##Basis#314 with 2 equations
States : (u(t))[1] (u(t))[2]
Parameters : p₁ p₂ p₃
Independent variable: t
Equations
Differential(t)((u(t))[1]) = p₁*(u(t))[1]
Differential(t)((u(t))[2]) = p₂*(u(t))[2] + p₃*((u(t))[1]^2)

And the connected parameters

get_parameter_map(basis)
3-element Vector{Pair{SymbolicUtils.BasicSymbolic{Real}, Float64}}:
 p₁ => -0.7999999999
 p₂ => -0.6999999999
 p₃ => 0.7

And plot the results

plot(res)
Example block output

Copy-Pasteable Code

using DataDrivenDiffEq
using LinearAlgebra
using OrdinaryDiffEq
using DataDrivenDMD

function slow_manifold(du, u, p, t)
    du[1] = p[1] * u[1]
    du[2] = p[2] * (u[2] - u[1]^2)
end

u0 = [3.0; -2.0]
tspan = (0.0, 5.0)
p = [-0.8; -0.7]

problem = ODEProblem{true, SciMLBase.NoSpecialize}(slow_manifold, u0, tspan, p)
solution = solve(problem, Tsit5(), saveat = 0.1);

prob = DataDrivenProblem(solution)

@parameters t
@variables u(t)[1:2]
Ψ = Basis([u; u[1]^2], u, independent_variable = t)
res = solve(prob, Ψ, DMDPINV(), digits = 2)

# This file was generated using Literate.jl, https://github.com/fredrikekre/Literate.jl

This page was generated using Literate.jl.