Skip to content

Commit 015ce60

Browse files
feat: preserve the u0map passed to ODEProblem and use it in remake_initializeprob
1 parent 8836abe commit 015ce60

File tree

2 files changed

+32
-0
lines changed

2 files changed

+32
-0
lines changed

src/systems/nonlinear/initializesystem.jl

Lines changed: 16 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -157,15 +157,22 @@ function generate_initializesystem(sys::ODESystem;
157157
for k in keys(defs)
158158
defs[k] = substitute(defs[k], paramsubs)
159159
end
160+
meta = InitializationSystemMetadata(Dict{Any, Any}(u0map), Dict{Any, Any}(pmap))
160161
return NonlinearSystem(eqs_ics,
161162
vars,
162163
pars;
163164
defaults = defs,
164165
checks = check_units,
165166
name,
167+
metadata = meta,
166168
kwargs...)
167169
end
168170

171+
struct InitializationSystemMetadata
172+
u0map::Dict{Any, Any}
173+
pmap::Dict{Any, Any}
174+
end
175+
169176
function is_parameter_solvable(p, pmap, defs, guesses)
170177
_val1 = pmap isa AbstractDict ? get(pmap, p, nothing) : nothing
171178
_val2 = get(defs, p, nothing)
@@ -253,6 +260,15 @@ function SciMLBase.remake_initializeprob(sys::ODESystem, odefn, u0, t0, p)
253260
!isempty(setobserved) || !isempty(setparobserved)) &&
254261
ModelingToolkit.get_tearing_state(sys) !== nothing) ||
255262
!isempty(initialization_equations(sys)))
263+
if SciMLBase.has_initializeprob(odefn)
264+
oldsys = odefn.initializeprob.f.sys
265+
meta = get_metadata(oldsys)
266+
if meta isa InitializationSystemMetadata
267+
u0 = merge(meta.u0map, u0)
268+
p = merge(meta.pmap, p)
269+
end
270+
end
271+
256272
initprob = InitializationProblem(sys, t0, u0, p)
257273
initprobmap = getu(initprob, unknowns(sys))
258274
punknowns = [p for p in all_variable_symbols(initprob) if is_parameter(sys, p)]

test/initializationsystem.jl

Lines changed: 16 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -780,3 +780,19 @@ end
780780
@test eltype(prob2.f.initializeprob.p.tunable) <: ForwardDiff.Dual
781781
@test prob2.f.initializeprob.u0 prob.f.initializeprob.u0
782782
end
783+
784+
@testset "`remake` preserves old u0map and pmap" begin
785+
@variables x(t) y(t)
786+
@parameters p
787+
@mtkbuild sys = ODESystem(
788+
[D(x) ~ x + p * y, y^2 + 4y * p^2 ~ x], t; guesses = [y => 1.0, p => 1.0])
789+
prob = ODEProblem(sys, [x => 1.0], (0.0, 1.0), [p => 1.0])
790+
@test is_variable(prob.f.initializeprob, y)
791+
prob2 = @test_nowarn remake(prob; p = [p => 3.0]) # ensure no over/under-determined warning
792+
@test is_variable(prob.f.initializeprob, y)
793+
794+
prob = ODEProblem(sys, [y => 1.0, x => 2.0], (0.0, 1.0), [p => missing])
795+
@test is_variable(prob.f.initializeprob, p)
796+
prob2 = @test_nowarn remake(prob; u0 = [y => 0.5])
797+
@test is_variable(prob.f.initializeprob, p)
798+
end

0 commit comments

Comments
 (0)