Skip to content
This repository was archived by the owner on Mar 24, 2026. It is now read-only.

Commit e1c9d14

Browse files
Merge pull request #1232 from dcourteville/fix_1231
Fix first occuring event search for backward integration
2 parents 2ff29cc + 6b55c42 commit e1c9d14

File tree

2 files changed

+32
-1
lines changed

2 files changed

+32
-1
lines changed

src/callbacks.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -140,7 +140,7 @@ end
140140
integrator,
141141
callbacks[$i],
142142
$i)
143-
if event_occurred2 && (tmin2 < tmin || !event_occurred)
143+
if event_occurred2 && (!event_occurred || integrator.tdir * tmin2 < integrator.tdir * tmin)
144144
tmin = tmin2
145145
upcrossing = upcrossing2
146146
event_occurred = true

test/downstream/community_callback_tests.jl

Lines changed: 31 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -230,3 +230,34 @@ sol = solve(prob, DFBDF())
230230
# test that the callback flipping p caused u[2] to get flipped.
231231
first_t = findfirst(isequal(0.5), sol.t)
232232
@test sol.u[first_t][2] == -sol.u[first_t + 1][2]
233+
234+
# https://github.com/SciML/DiffEqBase.jl/issues/1231
235+
@testset "Successive callbacks in same integration step" begin
236+
cb = ContinuousCallback(
237+
(u, t, integrator) -> t - 0e-8,
238+
(integrator) -> push!(record, 0)
239+
)
240+
241+
vcb = VectorContinuousCallback(
242+
(out, u, t, integrator) -> out .= (t - 1e-8, t - 2e-8),
243+
(integrator, event_index) -> push!(record, event_index),
244+
2
245+
)
246+
247+
f(u, p, t) = 1.0
248+
u0 = 0.0
249+
250+
# Forward propagation with successive events
251+
record = []
252+
tspan = (-1.0, 1.0)
253+
prob = ODEProblem(f, u0, tspan)
254+
sol = solve(prob, Tsit5(), dt = 2.0, callback = CallbackSet(cb, vcb))
255+
@test record == [0, 1, 2]
256+
257+
# Backward propagation with successive events
258+
record = []
259+
tspan = (1.0, -1.0)
260+
prob = ODEProblem(f, u0, tspan)
261+
sol = solve(prob, Tsit5(), dt = 2.0, callback = CallbackSet(cb, vcb))
262+
@test record == [2, 1, 0]
263+
end

0 commit comments

Comments
 (0)