forked from SciML/OrdinaryDiffEq.jl
-
Notifications
You must be signed in to change notification settings - Fork 0
Commit
This commit does not belong to any branch on this repository, and may belong to a fork outside of the repository.
Apply callbacks with a type-stable generated function.
Currently, as referenced in SciML/DifferentialEquations.jl#971, the old implementation of `handle_callbacks!` directly calls. `apply_callback!` on `continuous_callbacks[idx]`, which is inherently type-unstable because `apply_callback!` is specialized on the callback type. This commit adds a generated function `apply_ith_callback!` which generates type-stable code to do the same thing, where for each callback tuple type, the generated function unrolls the tuple by checking the callback index against static indicies. As a nice bonus, this generated function seems to often be converted into a switch statement at the LLVM level: ``` switch i64 %4, label %L46 [ i64 9, label %L3 i64 8, label %L8 i64 7, label %L13 i64 6, label %L18 i64 5, label %L23 i64 4, label %L28 i64 3, label %L33 i64 2, label %L38 i64 1, label %L43 ] ``` For testing, I added an allocation test which sets up a simple ODE problem, steps the integrator manually to before the first callback, then manipulates integrator state past the first callback point. This way, we can directly call `handle_callbacks!` and write a test on the allocation count. I confirm that (at least testing against commit SciML/DiffEqBase.jl@1799fc3, the current master branch tip in DiffEqBase.jl), the new method does not allocate, whereas the old one allocates. This may not be the case until a new release is cut of DiffEqBase.jl, because the old version of `find_first_continuous_callback` might allocate.
- Loading branch information
Showing
3 changed files
with
80 additions
and
6 deletions.
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,44 @@ | ||
using OrdinaryDiffEq, Test | ||
|
||
# Setup a simple ODE problem with several callbacks (to test LLVM code gen) | ||
# We will manually trigger the first callback and check its allocations. | ||
function f!(du, u, p, t) | ||
du .= .-u | ||
end | ||
|
||
cond_1(u, t, integrator) = u[1] - 0.5 | ||
cond_2(u, t, integrator) = u[2] + 0.5 | ||
cond_3(u, t, integrator) = u[2] + 0.6 | ||
cond_4(u, t, integrator) = u[2] + 0.7 | ||
cond_5(u, t, integrator) = u[2] + 0.8 | ||
cond_6(u, t, integrator) = u[2] + 0.9 | ||
cond_7(u, t, integrator) = u[2] + 0.1 | ||
cond_8(u, t, integrator) = u[2] + 0.11 | ||
cond_9(u, t, integrator) = u[2] + 0.12 | ||
|
||
function cb_affect!(integrator) | ||
integrator.p[1] += 1 | ||
end | ||
|
||
cbs = CallbackSet(ContinuousCallback(cond_1, cb_affect!), | ||
ContinuousCallback(cond_2, cb_affect!), | ||
ContinuousCallback(cond_3, cb_affect!), | ||
ContinuousCallback(cond_4, cb_affect!), | ||
ContinuousCallback(cond_5, cb_affect!), | ||
ContinuousCallback(cond_6, cb_affect!), | ||
ContinuousCallback(cond_7, cb_affect!), | ||
ContinuousCallback(cond_8, cb_affect!), | ||
ContinuousCallback(cond_9, cb_affect!)) | ||
|
||
integrator = init(ODEProblem(f!, [0.8, 1.0], (0.0, 100.0), [0, 0]), Tsit5(), callback = cbs, | ||
save_on = false); | ||
# Force a callback event to occur so we can call handle_callbacks! directly. | ||
# Step to a point where u[1] > 0.5, so we can force it below 0.5 and call handle callbacks | ||
step!(integrator, 0.1, true) | ||
|
||
function handle_allocs(integrator) | ||
integrator.u[1] = 0.4 | ||
@allocations OrdinaryDiffEq.handle_callbacks!(integrator) | ||
end | ||
handle_allocs(integrator); | ||
@test handle_allocs(integrator) == 0 |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters