Skip to content
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

fix: handle nothing passed as u0 to ODEProblem #3060

Open
wants to merge 3 commits into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
5 changes: 3 additions & 2 deletions src/systems/diffeqs/abstractodesystem.jl
Original file line number Diff line number Diff line change
Expand Up @@ -817,11 +817,12 @@ function process_DEProblem(constructor, sys::AbstractODESystem, u0map, parammap;
ModelingToolkit.get_tearing_state(sys) !== nothing) ||
!isempty(initialization_equations(sys))) && t !== nothing
if eltype(u0map) <: Number
u0map = unknowns(sys) .=> u0map
u0map = unknowns(sys) .=> vec(u0map)
end
if isempty(u0map)
if u0map === nothing || isempty(u0map)
u0map = Dict()
end

initializeprob = ModelingToolkit.InitializationProblem(
sys, t, u0map, parammap; guesses, warn_initialize_determined,
initialization_eqs, eval_expression, eval_module, fully_determined, check_units)
Expand Down
7 changes: 7 additions & 0 deletions src/systems/discrete_system/discrete_system.jl
Original file line number Diff line number Diff line change
Expand Up @@ -247,6 +247,13 @@ function process_DiscreteProblem(constructor, sys::DiscreteSystem, u0map, paramm
dvs = unknowns(sys)
ps = parameters(sys)

if eltype(u0map) <: Number
u0map = unknowns(sys) .=> vec(u0map)
end
if u0map === nothing || isempty(u0map)
u0map = Dict()
end

trueu0map = Dict()
for (k, v) in u0map
k = unwrap(k)
Expand Down
8 changes: 8 additions & 0 deletions test/discrete_system.jl
Original file line number Diff line number Diff line change
Expand Up @@ -271,3 +271,11 @@ end
k = ShiftIndex(t)
@named sys = DiscreteSystem([x ~ x^2 + y^2, y ~ x(k - 1) + y(k - 1)], t)
@test_throws ["algebraic equations", "not yet supported"] structural_simplify(sys)

@testset "Passing `nothing` to `u0`" begin
@variables x(t) = 1
k = ShiftIndex()
@mtkbuild sys = DiscreteSystem([x(k) ~ x(k - 1) + 1], t)
prob = @test_nowarn DiscreteProblem(sys, nothing, (0.0, 1.0))
@test_nowarn solve(prob, FunctionMap())
end
7 changes: 7 additions & 0 deletions test/nonlinearsystem.jl
Original file line number Diff line number Diff line change
Expand Up @@ -318,3 +318,10 @@ sys = structural_simplify(ns; conservative = true)
sol = solve(prob, NewtonRaphson())
@test sol[x] ≈ sol[y] ≈ sol[z] ≈ -3
end

@testset "Passing `nothing` to `u0`" begin
@variables x = 1
@mtkbuild sys = NonlinearSystem([0 ~ x^2 - x^3 + 3])
prob = @test_nowarn NonlinearProblem(sys, nothing)
@test_nowarn solve(prob)
end
7 changes: 7 additions & 0 deletions test/odesystem.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1387,3 +1387,10 @@ end

@test obsfn(ones(2), 2ones(2), 3ones(4), 4.0) == 6ones(2)
end

@testset "Passing `nothing` to `u0`" begin
@variables x(t) = 1
@mtkbuild sys = ODESystem(D(x) ~ t, t)
prob = @test_nowarn ODEProblem(sys, nothing, (0.0, 1.0))
@test_nowarn solve(prob)
end
7 changes: 7 additions & 0 deletions test/optimizationsystem.jl
Original file line number Diff line number Diff line change
Expand Up @@ -340,3 +340,10 @@ end
prob.f.cons_h(H3, [1.0, 1.0], [1.0, 100.0])
@test prob.f.cons_h([1.0, 1.0], [1.0, 100.0]) == H3
end

@testset "Passing `nothing` to `u0`" begin
@variables x = 1.0
@mtkbuild sys = OptimizationSystem((x - 3)^2, [x], [])
prob = @test_nowarn OptimizationProblem(sys, nothing)
@test_nowarn solve(prob, NelderMead())
end
8 changes: 8 additions & 0 deletions test/sdesystem.jl
Original file line number Diff line number Diff line change
Expand Up @@ -776,3 +776,11 @@ end
prob = SDEProblem(de, u0map, (0.0, 100.0), parammap)
@test solve(prob, SOSRI()).retcode == ReturnCode.Success
end

@testset "Passing `nothing` to `u0`" begin
@variables x(t) = 1
@brownian b
@mtkbuild sys = System([D(x) ~ x + b], t)
prob = @test_nowarn SDEProblem(sys, nothing, (0.0, 1.0))
@test_nowarn solve(prob, ImplicitEM())
end
Loading