-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathsimulate.jl
More file actions
192 lines (172 loc) · 7.99 KB
/
simulate.jl
File metadata and controls
192 lines (172 loc) · 7.99 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
# ── Phase 3: ODE simulation with DifferentialEquations / MTK ──────────────────
import DifferentialEquations
import Logging
import ModelingToolkit
import Printf: @sprintf
"""Module-level default simulation settings. Modify via `configure_simulate!`."""
const _SIM_SETTINGS = SimulateSettings(solver = DifferentialEquations.Rodas5P())
"""
configure_simulate!(; solver, saveat_n) → SimulateSettings
Update the module-level simulation settings in-place and return them.
# Keyword arguments
- `solver` — any SciML ODE/DAE algorithm instance (e.g. `Rodas5P`, `FBDF()`).
- `saveat_n` — number of uniform time points for purely algebraic systems.
# Example
```julia
using OrdinaryDiffEqBDF
configure_simulate!(solver = FBDF())
```
"""
function configure_simulate!(;
solver :: Union{Any,Nothing} = nothing,
saveat_n :: Union{Int,Nothing} = nothing,
)
isnothing(solver) || (_SIM_SETTINGS.solver = solver)
isnothing(saveat_n) || (_SIM_SETTINGS.saveat_n = saveat_n)
return _SIM_SETTINGS
end
"""
simulate_settings() → SimulateSettings
Return the current module-level simulation settings.
"""
simulate_settings() = _SIM_SETTINGS
"""
run_simulate(ode_prob, model_dir, model; settings, cmp_signals, csv_max_size_mb) → (success, time, error, sol)
Solve `ode_prob` using the algorithm in `settings.solver`. On success, also writes the
solution as a CSV file `<Short>_sim.csv` in `model_dir`.
Writes a `<model>_sim.log` file in `model_dir`.
Returns `nothing` as the fourth element on failure.
When `cmp_signals` is non-empty, only observed variables whose names appear in
that list are written to the CSV, keeping file sizes small when only a subset
of signals will be compared.
CSV files larger than `csv_max_size_mb` MiB are replaced with a
`<Short>_sim.csv.toobig` marker so that the report can note the omission.
"""
function run_simulate(ode_prob,
model_dir::String,
model::String;
settings ::SimulateSettings = _SIM_SETTINGS,
cmp_signals ::Vector{String} = String[],
csv_max_size_mb::Int = CSV_MAX_SIZE_MB)::Tuple{Bool,Float64,String,Any}
sim_success = false
sim_time = 0.0
sim_error = ""
sol = nothing
solver_settings_string = ""
log_file = open(joinpath(model_dir, "$(model)_sim.log"), "w")
println(log_file, "Model: $model")
logger = Logging.SimpleLogger(log_file, Logging.Debug)
t0 = time()
solver = settings.solver
try
# Redirect all library log output (including Symbolics/MTK warnings)
# to the log file so they don't clutter stdout.
sol = Logging.with_logger(logger) do
# Special handling for systems without states
sys = ode_prob.f.sys
n_unknowns = length(ModelingToolkit.unknowns(sys))
eqs = ModelingToolkit.equations(sys)
n_states = count(ModelingToolkitBase.isdiffeq, eqs)
kwargs, active_solver = if n_states == 0
# No state variables (e.g. BusUsage, Modelica.Electrical.Analog.Examples.OpAmps.Subtracter):
# Use Rodas5Pr on a fixed grid so observed variables can be evaluated at each step.
dt_s = (ode_prob.tspan[end] - ode_prob.tspan[1]) / (settings.saveat_n - 1)
saveat_s = collect(range(ode_prob.tspan[1], ode_prob.tspan[end]; length = settings.saveat_n))
(saveat = saveat_s, dt = dt_s, adaptive = false, dense = false),
DifferentialEquations.Rodas5Pr()
else
saveat_s = Float64[]
(saveat = saveat_s, dense = true), solver
end
# Log solver settings
integ = DifferentialEquations.init(ode_prob, active_solver; kwargs...)
saveat_s = kwargs.saveat
solver_settings_string = if hasproperty(integ, :opts)
sv_str = isempty(saveat_s) ? "[]" : "$(length(saveat_s)) points in [$(first(saveat_s)), $(last(saveat_s))]"
"""
Solver $(parentmodule(typeof(solver))).$(nameof(typeof(solver)))
saveat: $sv_str
abstol: $(@sprintf("%.2e", integ.opts.abstol))
reltol: $(@sprintf("%.2e", integ.opts.reltol))
adaptive: $(integ.opts.adaptive)
dense: $(integ.opts.dense)
"""
else
sv_str = isempty(saveat_s) ? "[]" : "$(length(saveat_s)) points in [$(first(saveat_s)), $(last(saveat_s))]"
"Solver (NullODEIntegrator — no unknowns)
saveat: $sv_str
dense: true"
end
# Solve
DifferentialEquations.solve(ode_prob, active_solver; kwargs...)
end
sim_time = time() - t0
if sol.retcode == DifferentialEquations.ReturnCode.Success
sys = sol.prob.f.sys
n_vars = length(ModelingToolkit.unknowns(sys))
n_obs = length(ModelingToolkit.observed(sys))
if isempty(sol.t)
sim_error = "Simulation produced no time points"
elseif n_vars == 0 && n_obs == 0
sim_error = "Simulation produced no output variables (no states or observed)"
else
sim_success = true
end
else
sim_error = "Solver returned: $(sol.retcode)"
end
catch e
sim_time = time() - t0
sim_error = sprint(showerror, e, catch_backtrace())
end
println(log_file, solver_settings_string)
println(log_file, "Time: $(round(sim_time; digits=3)) s")
println(log_file, "Success: $sim_success")
isempty(sim_error) || println(log_file, "\n--- Error ---\n$sim_error")
close(log_file)
# Write simulation results CSV (time + state variables + observed variables)
if sim_success && sol !== nothing
short_name = split(model, ".")[end]
sim_csv = joinpath(model_dir, "$(short_name)_sim.csv")
try
sys = sol.prob.f.sys
vars = ModelingToolkit.unknowns(sys)
obs_eqs = ModelingToolkit.observed(sys)
# Only save observed variables that appear in cmp_signals.
# This avoids writing thousands of algebraic variables to disk when
# only a handful are actually verified during comparison.
norm_cmp = Set(_normalize_var(s) for s in cmp_signals)
obs_eqs_filtered = isempty(norm_cmp) ? obs_eqs :
filter(eq -> _normalize_var(string(eq.lhs)) in norm_cmp, obs_eqs)
obs_syms = [eq.lhs for eq in obs_eqs_filtered]
col_names = vcat(
[_clean_var_name(string(v)) for v in vars],
[_clean_var_name(string(s)) for s in obs_syms],
)
open(sim_csv, "w") do f
println(f, join(["time"; col_names], ","))
for (ti, t) in enumerate(sol.t)
row = [@sprintf("%.10g", t)]
for vi in eachindex(vars)
push!(row, @sprintf("%.10g", sol[vi, ti]))
end
for sym in obs_syms
val = try Float64(sol(t; idxs = sym)) catch; NaN end
push!(row, @sprintf("%.10g", val))
end
println(f, join(row, ","))
end
end
csv_bytes = filesize(sim_csv)
if csv_bytes > csv_max_size_mb * 1024^2
csv_mb = round(csv_bytes / 1024^2; digits=1)
@warn "Simulation CSV for $model is $(csv_mb) MB (> $(csv_max_size_mb) MB limit); skipping."
rm(sim_csv)
write(sim_csv * ".toobig", string(csv_bytes))
end
catch e
@warn "Failed to write simulation CSV for $model: $(sprint(showerror, e))"
end
end
return sim_success, sim_time, sim_error, sol
end