Skip to content

Commit 2f5c718

Browse files
Merge pull request #2861 from hersle/init_defaults_flipping
Left-to-right expand observed equations into defaults during initialization
2 parents 6d8e86e + 8862b97 commit 2f5c718

File tree

3 files changed

+36
-10
lines changed

3 files changed

+36
-10
lines changed

src/systems/diffeqs/abstractodesystem.jl

Lines changed: 10 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -725,10 +725,16 @@ function get_u0(
725725
defs = mergedefaults(defs, parammap, ps)
726726
end
727727

728-
obs = filter!(x -> !(x[1] isa Number),
729-
map(x -> isparameter(x.rhs) ? x.lhs => x.rhs : x.rhs => x.lhs, observed(sys)))
730-
observedmap = isempty(obs) ? Dict() : todict(obs)
731-
defs = mergedefaults(defs, observedmap, u0map, dvs)
728+
# Convert observed equations "lhs ~ rhs" into defaults.
729+
# Use the order "lhs => rhs" by default, but flip it to "rhs => lhs"
730+
# if "lhs" is known by other means (parameter, another default, ...)
731+
# TODO: Is there a better way to determine which equations to flip?
732+
obs = map(x -> x.lhs => x.rhs, observed(sys))
733+
obs = map(x -> x[1] in keys(defs) ? reverse(x) : x, obs)
734+
obs = filter!(x -> !(x[1] isa Number), obs) # exclude e.g. "0 => x^2 + y^2 - 25"
735+
obsmap = isempty(obs) ? Dict() : todict(obs)
736+
737+
defs = mergedefaults(defs, obsmap, u0map, dvs)
732738
if symbolic_u0
733739
u0 = varmap_to_vars(
734740
u0map, dvs; defaults = defs, tofloat = false, use_union = false, toterm)

test/guess_propagation.jl

Lines changed: 6 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -90,21 +90,21 @@ prob = ODEProblem(sys, [], (0.0, 1.0), [x0 => 1.0])
9090
@variables y(t) = x0
9191
@mtkbuild sys = ODESystem([x ~ x0, D(y) ~ x], t)
9292
prob = ODEProblem(sys, [], (0.0, 1.0), [x0 => 1.0])
93-
prob[x] == 1.0
94-
prob[y] == 1.0
93+
@test prob[x] == 1.0
94+
@test prob[y] == 1.0
9595

9696
@parameters x0
9797
@variables x(t)
9898
@variables y(t) = x0
9999
@mtkbuild sys = ODESystem([x ~ y, D(y) ~ x], t)
100100
prob = ODEProblem(sys, [], (0.0, 1.0), [x0 => 1.0])
101-
prob[x] == 1.0
102-
prob[y] == 1.0
101+
@test prob[x] == 1.0
102+
@test prob[y] == 1.0
103103

104104
@parameters x0
105105
@variables x(t) = x0
106106
@variables y(t) = x
107107
@mtkbuild sys = ODESystem([x ~ y, D(y) ~ x], t)
108108
prob = ODEProblem(sys, [], (0.0, 1.0), [x0 => 1.0])
109-
prob[x] == 1.0
110-
prob[y] == 1.0
109+
@test prob[x] == 1.0
110+
@test prob[y] == 1.0

test/odesystem.jl

Lines changed: 20 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1194,3 +1194,23 @@ end
11941194
@test_nowarn obsfn(buffer, [1.0], ps..., 3.0)
11951195
@test buffer [2.0, 3.0, 4.0]
11961196
end
1197+
1198+
# https://github.com/SciML/ModelingToolkit.jl/issues/2859
1199+
@testset "Initialization with defaults from observed equations (edge case)" begin
1200+
@variables x(t) y(t) z(t)
1201+
eqs = [D(x) ~ 0, y ~ x, D(z) ~ 0]
1202+
defaults = [x => 1, z => y]
1203+
@named sys = ODESystem(eqs, t; defaults)
1204+
ssys = structural_simplify(sys)
1205+
prob = ODEProblem(ssys, [], (0.0, 1.0), [])
1206+
@test prob[x] == prob[y] == prob[z] == 1.0
1207+
1208+
@parameters y0
1209+
@variables x(t) y(t) z(t)
1210+
eqs = [D(x) ~ 0, y ~ y0 / x, D(z) ~ y]
1211+
defaults = [y0 => 1, x => 1, z => y]
1212+
@named sys = ODESystem(eqs, t; defaults)
1213+
ssys = structural_simplify(sys)
1214+
prob = ODEProblem(ssys, [], (0.0, 1.0), [])
1215+
@test prob[x] == prob[y] == prob[z] == 1.0
1216+
end

0 commit comments

Comments
 (0)