Skip to content

Commit 729cadc

Browse files
committed
Add modelingtoolkitize
1 parent aabee57 commit 729cadc

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
@@ -84,3 +84,27 @@ function vars!(vars, O)
8484

8585
return vars
8686
end
87+
88+
"""
89+
$(SIGNATURES)
90+
91+
Generate `ODESystem`, dependent variables, and parameters from an `ODEProblem`.
92+
"""
93+
function modelingtoolkitize(prob::DiffEqBase.ODEProblem)
94+
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)];
95+
D, = @derivatives D'~t
96+
97+
rhs = [D(var) for var in vars]
98+
99+
if DiffEqBase.isinplace(prob)
100+
lhs = similar(vars, Any)
101+
prob.f(lhs, vars, params, t)
102+
else
103+
lhs = prob.f(vars, params, t)
104+
end
105+
106+
eqs = vcat([rhs[i] ~ lhs[i] for i in eachindex(prob.u0)]...)
107+
de = ODESystem(eqs)
108+
109+
de, vars, params
110+
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
@@ -197,3 +198,26 @@ ns = NonlinearSystem(eqs, [x,y,z])
197198
nlsys_func = generate_function(ns, [x,y,z], [σ,ρ,β])
198199
jac = calculate_jacobian(ns)
199200
jac = generate_jacobian(ns)
201+
202+
function lotka(u,p,t)
203+
x = u[1]
204+
y = u[2]
205+
[p[1]*x - p[2]*x*y,
206+
-p[3]*y + p[4]*x*y]
207+
end
208+
209+
prob = ODEProblem(ODEFunction{false}(lotka),[1.0,1.0],(0.0,1.0),[1.5,1.0,3.0,1.0])
210+
de, vars, params = modelingtoolkitize(prob)
211+
ODEFunction(de, vars, params)(similar(prob.u0), prob.u0, prob.p, 0.1)
212+
213+
function lotka(du,u,p,t)
214+
x = u[1]
215+
y = u[2]
216+
du[1] = p[1]*x - p[2]*x*y
217+
du[2] = -p[3]*y + p[4]*x*y
218+
end
219+
220+
prob = ODEProblem(lotka,[1.0,1.0],(0.0,1.0),[1.5,1.0,3.0,1.0])
221+
222+
de, vars, params = modelingtoolkitize(prob)
223+
ODEFunction(de, vars, params)(similar(prob.u0), prob.u0, prob.p, 0.1)

0 commit comments

Comments
 (0)