Skip to content

Commit ebd3131

Browse files
Merge pull request #3593 from BenChung/bc/fix_ia_namespacing
Fix ImperativeAffect namespacing for arrays in non-root components
2 parents 81d8c72 + dc1e1c3 commit ebd3131

File tree

2 files changed

+40
-1
lines changed

2 files changed

+40
-1
lines changed

src/systems/imperative_affect.jl

+13-1
Original file line numberDiff line numberDiff line change
@@ -101,10 +101,22 @@ end
101101

102102
namespace_affects(af::ImperativeAffect, s) = namespace_affect(af, s)
103103
function namespace_affect(affect::ImperativeAffect, s)
104+
rmn = []
105+
for modded in modified(affect)
106+
if symbolic_type(modded) == NotSymbolic() && modded isa AbstractArray
107+
res = []
108+
for m in modded
109+
push!(res, renamespace(s, m))
110+
end
111+
push!(rmn, res)
112+
else
113+
push!(rmn, renamespace(s, modded))
114+
end
115+
end
104116
ImperativeAffect(func(affect),
105117
namespace_expr.(observed(affect), (s,)),
106118
observed_syms(affect),
107-
renamespace.((s,), modified(affect)),
119+
rmn,
108120
modified_syms(affect),
109121
context(affect),
110122
affect.skip_checks)

test/symbolic_events.jl

+27
Original file line numberDiff line numberDiff line change
@@ -1434,3 +1434,30 @@ end
14341434
sol2 = solve(ODEProblem(sys2, [], (0.0, 1.0)), Tsit5())
14351435
@test 100.0 sol2[sys2.wd2.θ]
14361436
end
1437+
1438+
@testset "Array parameter updates of parent components in ImperativeEffect" begin
1439+
function child(vals; name, max_time = 0.1)
1440+
vars = @variables begin
1441+
x(t) = 0.0
1442+
end
1443+
eqs = reduce(vcat, Symbolics.scalarize.([
1444+
D(x) ~ 1.0
1445+
]))
1446+
reset = ModelingToolkit.ImperativeAffect(
1447+
modified = (; vals = Symbolics.scalarize(ParentScope.(vals)), x)) do m, o, _, i
1448+
@set! m.vals = m.vals .+ 1
1449+
@set! m.x = 0.0
1450+
return m
1451+
end
1452+
return ODESystem(eqs, t, vars, []; name = name,
1453+
continuous_events = [[x ~ max_time] => reset])
1454+
end
1455+
shared_pars = @parameters begin
1456+
vals(t)[1:2] = 0.0
1457+
end
1458+
1459+
@named sys = ODESystem(Equation[], t, [], Symbolics.scalarize(vals);
1460+
systems = [child(vals; name = :child)])
1461+
sys = structural_simplify(sys)
1462+
sol = solve(ODEProblem(sys, [], (0.0, 1.0)), Tsit5())
1463+
end

0 commit comments

Comments
 (0)