Skip to content

Fix _varmap_to_vars, re-enable dde.jl test #2545

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

Merged
merged 1 commit into from
Mar 28, 2024
Merged
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
3 changes: 2 additions & 1 deletion Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -115,6 +115,7 @@ AmplNLWriter = "7c4d4715-977e-5154-bfe0-e096adeac482"
BenchmarkTools = "6e4b80f9-dd63-53aa-95a3-0cdb28fa8baf"
ControlSystemsBase = "aaaaaaaa-a6ca-5380-bf3e-84a91bcd477e"
DeepDiffs = "ab62b9b5-e342-54a8-a765-a90f495de1a6"
DelayDiffEq = "bcd4f6db-9728-5f36-b5f7-82caef46ccdb"
ForwardDiff = "f6369f11-7733-5829-9624-2563aa707210"
Ipopt = "b6b21f68-93f8-5de0-b562-5493be1d77c9"
Ipopt_jll = "9cc047cb-c261-5740-88fc-0cf96f7bdcc7"
Expand All @@ -136,4 +137,4 @@ Sundials = "c3572dad-4567-51f8-b174-8c6c989267f4"
Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40"

[targets]
test = ["AmplNLWriter", "BenchmarkTools", "ControlSystemsBase", "NonlinearSolve", "ForwardDiff", "Ipopt", "Ipopt_jll", "ModelingToolkitStandardLibrary", "Optimization", "OptimizationOptimJL", "OptimizationMOI", "Random", "ReferenceTests", "SafeTestsets", "StableRNGs", "Statistics", "SteadyStateDiffEq", "Test", "StochasticDiffEq", "Sundials", "StochasticDelayDiffEq", "Pkg"]
test = ["AmplNLWriter", "BenchmarkTools", "ControlSystemsBase", "DelayDiffEq", "NonlinearSolve", "ForwardDiff", "Ipopt", "Ipopt_jll", "ModelingToolkitStandardLibrary", "Optimization", "OptimizationOptimJL", "OptimizationMOI", "Random", "ReferenceTests", "SafeTestsets", "StableRNGs", "Statistics", "SteadyStateDiffEq", "Test", "StochasticDiffEq", "Sundials", "StochasticDelayDiffEq", "Pkg"]
3 changes: 2 additions & 1 deletion src/systems/diffeqs/abstractodesystem.jl
Original file line number Diff line number Diff line change
Expand Up @@ -957,7 +957,6 @@ function process_DEProblem(constructor, sys::AbstractODESystem, u0map, parammap;
ddvs = nothing
end
check_eqs_u0(eqs, dvs, u0; kwargs...)

f = constructor(sys, dvs, ps, u0; ddvs = ddvs, tgrad = tgrad, jac = jac,
checkbounds = checkbounds, p = p,
linenumbers = linenumbers, parallel = parallel, simplify = simplify,
Expand Down Expand Up @@ -1243,6 +1242,8 @@ function DiffEqBase.SDDEProblem{iip}(sys::AbstractODESystem, u0map = [],
h_oop, h_iip = generate_history(sys, u0)
h(out, p, t) = h_iip(out, p, t)
h(p, t) = h_oop(p, t)
h(p::MTKParameters, t) = h_oop(p..., t)
h(out, p::MTKParameters, t) = h_iip(out, p..., t)
u0 = h(p, tspan[1])
cbs = process_events(sys; callback, kwargs...)
if has_discrete_subsystems(sys) && (dss = get_discrete_subsystems(sys)) !== nothing
Expand Down
2 changes: 1 addition & 1 deletion src/systems/parameter_buffer.jl
Original file line number Diff line number Diff line change
Expand Up @@ -90,7 +90,7 @@ function MTKParameters(
for (sym, val) in p
sym = unwrap(sym)
ctype = concrete_symtype(sym)
val = symconvert(ctype, fixpoint_sub(val, p))
val = symconvert(ctype, unwrap(fixpoint_sub(val, p)))
done = set_value(sym, val)
if !done && Symbolics.isarraysymbolic(sym)
done = all(set_value.(collect(sym), val))
Expand Down
2 changes: 1 addition & 1 deletion src/variables.jl
Original file line number Diff line number Diff line change
Expand Up @@ -197,7 +197,7 @@ function _varmap_to_vars(varmap::Dict, varlist; defaults = Dict(), check = false
for var in varlist
var = unwrap(var)
val = unwrap(fixpoint_sub(var, varmap; operator = Symbolics.Operator))
if symbolic_type(val) === NotSymbolic()
if !isequal(val, var)
values[var] = val
end
end
Expand Down
17 changes: 7 additions & 10 deletions test/dde.jl
Original file line number Diff line number Diff line change
@@ -1,4 +1,6 @@
using ModelingToolkit, DelayDiffEq, Test
using ModelingToolkit: t_nounits as t, D_nounits as D

p0 = 0.2;
q0 = 0.3;
v0 = 1;
Expand Down Expand Up @@ -29,14 +31,13 @@ prob2 = DDEProblem(bc_model, u0, h2, tspan, constant_lags = lags)
sol2 = solve(prob2, alg, reltol = 1e-7, abstol = 1e-10)

@parameters p0=0.2 p1=0.2 q0=0.3 q1=0.3 v0=1 v1=1 d0=5 d1=1 d2=1 beta0=1 beta1=1
@variables t x₀(t) x₁(t) x₂(..)
@variables x₀(t) x₁(t) x₂(..)
tau = 1
D = Differential(t)
eqs = [D(x₀) ~ (v0 / (1 + beta0 * (x₂(t - tau)^2))) * (p0 - q0) * x₀ - d0 * x₀
D(x₁) ~ (v0 / (1 + beta0 * (x₂(t - tau)^2))) * (1 - p0 + q0) * x₀ +
(v1 / (1 + beta1 * (x₂(t - tau)^2))) * (p1 - q1) * x₁ - d1 * x₁
D(x₂(t)) ~ (v1 / (1 + beta1 * (x₂(t - tau)^2))) * (1 - p1 + q1) * x₁ - d2 * x₂(t)]
@named sys = System(eqs)
@mtkbuild sys = System(eqs, t)
prob = DDEProblem(sys,
[x₀ => 1.0, x₁ => 1.0, x₂(t) => 1.0],
tspan,
Expand Down Expand Up @@ -70,21 +71,17 @@ prob = SDDEProblem(hayes_modelf, hayes_modelg, [1.0], h, tspan, pmul;
constant_lags = (pmul[1],));
sol = solve(prob, RKMil())

@variables t x(..)
@variables x(..)
@parameters a=-4.0 b=-2.0 c=10.0 α=-1.3 β=-1.2 γ=1.1
D = Differential(t)
@brownian η
τ = 1.0
eqs = [D(x(t)) ~ a * x(t) + b * x(t - τ) + c + (α * x(t) + γ) * η]
@named sys = System(eqs)
sys = structural_simplify(sys)
@mtkbuild sys = System(eqs, t)
@test equations(sys) == [D(x(t)) ~ a * x(t) + b * x(t - τ) + c]
@test isequal(ModelingToolkit.get_noiseeqs(sys), [α * x(t) + γ;;])
prob_mtk = SDDEProblem(sys, [x(t) => 1.0 + t], tspan; constant_lags = (τ,));
@test_nowarn sol_mtk = solve(prob_mtk, RKMil())

@variables t
D = Differential(t)
@parameters x(..) a

function oscillator(; name, k = 1.0, τ = 0.01)
Expand All @@ -93,7 +90,7 @@ function oscillator(; name, k = 1.0, τ = 0.01)
eqs = [D(x(t)) ~ y,
D(y) ~ -k * x(t - τ) + jcn,
delx ~ x(t - τ)]
return System(eqs; name = name)
return System(eqs, t; name = name)
end

systems = @named begin
Expand Down
1 change: 1 addition & 0 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -34,6 +34,7 @@ end
@safetestset "Mass Matrix Test" include("mass_matrix.jl")
@safetestset "SteadyStateSystem Test" include("steadystatesystems.jl")
@safetestset "SDESystem Test" include("sdesystem.jl")
@safetestset "DDESystem Test" include("dde.jl")
@safetestset "NonlinearSystem Test" include("nonlinearsystem.jl")
@safetestset "InitializationSystem Test" include("initializationsystem.jl")
@safetestset "PDE Construction Test" include("pde.jl")
Expand Down