Representing a system as a component
MTKHelpers.embed_system
— Functionembed_system(m;name)
Embeds system m
as the single component of a larger system. This helps to match the naming of the states, parameters, observables to the namespace of the system.
using ModelingToolkit, OrdinaryDiffEq
using MTKHelpers
using ModelingToolkit: t_nounits as t, D_nounits as D
function samplesystem(;name)
sts = @variables x(t) RHS(t) # RHS is observed
ps = @parameters τ # parameters
ODESystem([ RHS ~ (1 - x)/τ, D(x) ~ RHS ], t, sts, ps; name)
end
@named m = samplesystem()
@named sys = embed_system(m)
# note that keys are `m.x`,`m.τ` or `m.RHS`.
# Hence only m needs to be defined rather than all the states, parameters,
# and observables of the system.
prob = ODEProblem(sys, [m.x => 0.0], (0.0,10.0), [m.τ => 3.0])
sol = solve(prob, Tsit5());
dx1 = sol[m.RHS][1:3]
dx2 = sol[getproperty(m,:RHS)][1:3] # access by symbols
# using Plots
# plot(sol, vars=[m.x,m.RHS])
# plot(sol, vars=getproperty.(Ref(m),[:x, :RHS])) # access by symbols
length(dx2) == 3
# output
true
Overriding equations of an existing system
For debugging bigger systems, it is helpful to set some equations to zero or modify/simplify the system in other ways.
Function override_system
takes a set of equations and matches the right-hand site to the equations of the original system and replaces those equations.
# setting up a simple example composite system and problem
using ModelingToolkit, OrdinaryDiffEq, ComponentArrays
using MTKHelpers
using ModelingToolkit: t_nounits as t, D_nounits as D
function samplesystem(;name,τ = 3.0, p1=1.1, p2=1.2)
sts = @variables x(t) RHS(t) # RHS is observed
ps = @parameters τ=τ p1=p1 p2=p2 # parameters
ODESystem([ RHS ~ p1/p2 * (1 - x)/τ, D(x) ~ RHS ], t, sts, ps; name)
end
# simplify the system by setting RHS ~ RHS_0 * x
function samplesystem_const(RHS0; name)
# demonstrating override_system by setting the RHS to constant first order rate
m = samplesystem(;name)
@unpack RHS,x = m
@parameters t
ps = @parameters RHS_0=RHS0
D = Differential(t)
eqs = [
RHS ~ RHS_0 * x,
]
sys_ext = override_system(eqs, m; name, ps)
end
@named mc = samplesystem_const(-0.1)
@named sys = embed_system(mc)
prob = ODEProblem(sys, [mc.x => 1.0], (0.0,10.0))
sol = solve(prob, Tsit5())
isapprox(sol(8)[1], exp(-0.1*8), atol = 1e-5)
true
MTKHelpers.override_system
— Functionoverride_system(eqs, basesys::AbstractSystem;
name::Symbol=Symbol(string(nameof(basesys))*"_ext"),
ps=Term[],
obs=Equation[],
evs=ModelingToolkit.SymbolicContinuousCallback[],
defs=Dict()
)
Modify basesys
by replacing some equations matched by their left-hand-side. The keyword argument correspond to ODESystem.
Utilities
MTKHelpers.symbol_op
— Functionsymbol_op(t)
Extract the inner symbol_op from a Term, Num, or BasicSymbolic object.
MTKHelpers.strip_namespace
— Functionstrip_namespace(s)
Omit the part before the first dot.