-
-
Notifications
You must be signed in to change notification settings - Fork 208
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
feat: initial implementation of new DiscreteSystem
#2507
Changes from all commits
05952e1
fd308cb
62dee1b
5d7693e
3ac5867
63e9383
f73d8b1
6bfd3a9
27e1145
3f1c06b
88bcd3d
5a19c6a
e9cc50e
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,51 @@ | ||
# (Experimental) Modeling Discrete Systems | ||
|
||
In this example, we will use the new [`DiscreteSystem`](@ref) API | ||
to create an SIR model. | ||
|
||
```@example discrete | ||
using ModelingToolkit | ||
using ModelingToolkit: t_nounits as t | ||
using OrdinaryDiffEq: solve, FunctionMap | ||
|
||
@inline function rate_to_proportion(r, t) | ||
1 - exp(-r * t) | ||
end | ||
@parameters c δt β γ | ||
@constants h = 1 | ||
@variables S(t) I(t) R(t) | ||
k = ShiftIndex(t) | ||
infection = rate_to_proportion( | ||
β * c * I(k - 1) / (S(k - 1) * h + I(k - 1) + R(k - 1)), δt * h) * S(k - 1) | ||
recovery = rate_to_proportion(γ * h, δt) * I(k - 1) | ||
|
||
# Equations | ||
eqs = [S(k) ~ S(k - 1) - infection * h, | ||
I(k) ~ I(k - 1) + infection - recovery, | ||
R(k) ~ R(k - 1) + recovery] | ||
@mtkbuild sys = DiscreteSystem(eqs, t) | ||
|
||
u0 = [S(k - 1) => 990.0, I(k - 1) => 10.0, R(k - 1) => 0.0] | ||
p = [β => 0.05, c => 10.0, γ => 0.25, δt => 0.1] | ||
tspan = (0.0, 100.0) | ||
prob = DiscreteProblem(sys, u0, tspan, p) | ||
sol = solve(prob, FunctionMap()) | ||
``` | ||
|
||
All shifts must be non-positive, i.e., discrete-time variables may only be indexed at index | ||
`k, k-1, k-2, ...`. If default values are provided, they are treated as the value of the | ||
variable at the previous timestep. For example, consider the following system to generate | ||
the Fibonacci series: | ||
|
||
```@example discrete | ||
@variables x(t) = 1.0 | ||
@mtkbuild sys = DiscreteSystem([x ~ x(k - 1) + x(k - 2)], t) | ||
``` | ||
|
||
The "default value" here should be interpreted as the value of `x` at all past timesteps. | ||
For example, here `x(k-1)` and `x(k-2)` will be `1.0`, and the inital value of `x(k)` will | ||
thus be `2.0`. During problem construction, the _past_ value of a variable should be | ||
provided. For example, providing `[x => 1.0]` while constructing this problem will error. | ||
Provide `[x(k-1) => 1.0]` instead. Note that values provided during problem construction | ||
_do not_ apply to the entire history. Hence, if `[x(k-1) => 2.0]` is provided, the value of | ||
`x(k-2)` will still be `1.0`. |
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -382,8 +382,8 @@ function tearing_reassemble(state::TearingState, var_eq_matching; | |
dx = fullvars[dv] | ||
# add `x_t` | ||
order, lv = var_order(dv) | ||
x_t = lower_varname(fullvars[lv], iv, order) | ||
push!(fullvars, x_t) | ||
x_t = lower_varname_withshift(fullvars[lv], iv, order) | ||
push!(fullvars, simplify_shifts(x_t)) | ||
v_t = length(fullvars) | ||
v_t_idx = add_vertex!(var_to_diff) | ||
add_vertex!(graph, DST) | ||
|
@@ -437,11 +437,12 @@ function tearing_reassemble(state::TearingState, var_eq_matching; | |
# We cannot solve the differential variable like D(x) | ||
if isdervar(iv) | ||
order, lv = var_order(iv) | ||
dx = D(lower_varname(fullvars[lv], idep, order - 1)) | ||
eq = dx ~ ModelingToolkit.fixpoint_sub( | ||
dx = D(simplify_shifts(lower_varname_withshift( | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Why is the There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Ah, that's the structural info. |
||
fullvars[lv], idep, order - 1))) | ||
eq = dx ~ simplify_shifts(ModelingToolkit.fixpoint_sub( | ||
Symbolics.solve_for(neweqs[ieq], | ||
fullvars[iv]), | ||
total_sub) | ||
total_sub; operator = ModelingToolkit.Shift)) | ||
for e in 𝑑neighbors(graph, iv) | ||
e == ieq && continue | ||
for v in 𝑠neighbors(graph, e) | ||
|
@@ -450,7 +451,7 @@ function tearing_reassemble(state::TearingState, var_eq_matching; | |
rem_edge!(graph, e, iv) | ||
end | ||
push!(diff_eqs, eq) | ||
total_sub[eq.lhs] = eq.rhs | ||
total_sub[simplify_shifts(eq.lhs)] = eq.rhs | ||
push!(diffeq_idxs, ieq) | ||
push!(diff_vars, diff_to_var[iv]) | ||
continue | ||
|
@@ -469,7 +470,7 @@ function tearing_reassemble(state::TearingState, var_eq_matching; | |
neweq = var ~ ModelingToolkit.fixpoint_sub( | ||
simplify ? | ||
Symbolics.simplify(rhs) : rhs, | ||
total_sub) | ||
total_sub; operator = ModelingToolkit.Shift) | ||
push!(subeqs, neweq) | ||
push!(solved_equations, ieq) | ||
push!(solved_variables, iv) | ||
|
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -412,3 +412,40 @@ function numerical_nlsolve(f, u0, p) | |
# TODO: robust initial guess, better debugging info, and residual check | ||
sol.u | ||
end | ||
|
||
### | ||
### Misc | ||
### | ||
|
||
function lower_varname_withshift(var, iv, order) | ||
order == 0 && return var | ||
if ModelingToolkit.isoperator(var, ModelingToolkit.Shift) | ||
op = operation(var) | ||
return Shift(op.t, order)(var) | ||
end | ||
return lower_varname(var, iv, order) | ||
end | ||
|
||
function isdoubleshift(var) | ||
return ModelingToolkit.isoperator(var, ModelingToolkit.Shift) && | ||
ModelingToolkit.isoperator(arguments(var)[1], ModelingToolkit.Shift) | ||
end | ||
|
||
function simplify_shifts(var) | ||
ModelingToolkit.hasshift(var) || return var | ||
var isa Equation && return simplify_shifts(var.lhs) ~ simplify_shifts(var.rhs) | ||
if isdoubleshift(var) | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. What about triple shifts? Should we recurse on that, too? There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Thanks for pointing that out There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Fixed |
||
op1 = operation(var) | ||
vv1 = arguments(var)[1] | ||
op2 = operation(vv1) | ||
vv2 = arguments(vv1)[1] | ||
s1 = op1.steps | ||
s2 = op2.steps | ||
t1 = op1.t | ||
t2 = op2.t | ||
return simplify_shifts(ModelingToolkit.Shift(t1 === nothing ? t2 : t1, s1 + s2)(vv2)) | ||
else | ||
return similarterm(var, operation(var), simplify_shifts.(arguments(var)), | ||
Symbolics.symtype(var); metadata = unwrap(var).metadata) | ||
end | ||
end |
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Is it still experimental?
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
kind of, unless observed variables are working as intended now?
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Yeah, since shifting the entire system by one step is basically a workaround and eventually we'll need proper compiler support for this. I've discussed with @baggepinnen and that will end up changing things such as the state realization and observed equations in the simplified system.