Skip to content

Commit 93fe6bd

Browse files
Merge pull request #146 from YingboMa/myb/modelingtoolkitize
Add modelingtoolkitize
2 parents e2f0d08 + 729cadc commit 93fe6bd

File tree

3 files changed

+49
-0
lines changed

3 files changed

+49
-0
lines changed

src/ModelingToolkit.jl

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -4,6 +4,7 @@ export Operation, Expression
44
export calculate_jacobian, generate_jacobian, generate_function
55
export independent_variables, dependent_variables, parameters
66
export @register
7+
export modelingtoolkitize
78

89

910
using DiffEqBase

src/utils.jl

Lines changed: 24 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -90,3 +90,27 @@ function vars!(vars, O)
9090

9191
return vars
9292
end
93+
94+
"""
95+
$(SIGNATURES)
96+
97+
Generate `ODESystem`, dependent variables, and parameters from an `ODEProblem`.
98+
"""
99+
function modelingtoolkitize(prob::DiffEqBase.ODEProblem)
100+
t, = @parameters t; vars = [Variable(Symbol(:x, i))(t) for i in eachindex(prob.u0)]; params = [Variable(Symbol(, i); known = true)() for i in eachindex(prob.p)];
101+
D, = @derivatives D'~t
102+
103+
rhs = [D(var) for var in vars]
104+
105+
if DiffEqBase.isinplace(prob)
106+
lhs = similar(vars, Any)
107+
prob.f(lhs, vars, params, t)
108+
else
109+
lhs = prob.f(vars, params, t)
110+
end
111+
112+
eqs = vcat([rhs[i] ~ lhs[i] for i in eachindex(prob.u0)]...)
113+
de = ODESystem(eqs)
114+
115+
de, vars, params
116+
end

test/system_construction.jl

Lines changed: 24 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,5 @@
11
using ModelingToolkit, StaticArrays, LinearAlgebra
2+
using DiffEqBase
23
using Test
34

45
# Define some variables
@@ -198,3 +199,26 @@ ns = NonlinearSystem(eqs, [x,y,z])
198199
nlsys_func = generate_function(ns, [x,y,z], [σ,ρ,β])
199200
jac = calculate_jacobian(ns)
200201
jac = generate_jacobian(ns)
202+
203+
function lotka(u,p,t)
204+
x = u[1]
205+
y = u[2]
206+
[p[1]*x - p[2]*x*y,
207+
-p[3]*y + p[4]*x*y]
208+
end
209+
210+
prob = ODEProblem(ODEFunction{false}(lotka),[1.0,1.0],(0.0,1.0),[1.5,1.0,3.0,1.0])
211+
de, vars, params = modelingtoolkitize(prob)
212+
ODEFunction(de, vars, params)(similar(prob.u0), prob.u0, prob.p, 0.1)
213+
214+
function lotka(du,u,p,t)
215+
x = u[1]
216+
y = u[2]
217+
du[1] = p[1]*x - p[2]*x*y
218+
du[2] = -p[3]*y + p[4]*x*y
219+
end
220+
221+
prob = ODEProblem(lotka,[1.0,1.0],(0.0,1.0),[1.5,1.0,3.0,1.0])
222+
223+
de, vars, params = modelingtoolkitize(prob)
224+
ODEFunction(de, vars, params)(similar(prob.u0), prob.u0, prob.p, 0.1)

0 commit comments

Comments
 (0)