diff --git a/lib/OrdinaryDiffEqCore/src/solve.jl b/lib/OrdinaryDiffEqCore/src/solve.jl index ef77e4c6b4..61e3f9185b 100644 --- a/lib/OrdinaryDiffEqCore/src/solve.jl +++ b/lib/OrdinaryDiffEqCore/src/solve.jl @@ -242,6 +242,12 @@ function DiffEqBase.__init( resType = typeof(res_prototype) end + if tstops isa AbstractArray || tstops isa Tuple + _tstops = nothing + else + _tstops = tstops + tstops = () + end tstops_internal = initialize_tstops(tType, tstops, d_discontinuities, tspan) saveat_internal = initialize_saveat(tType, saveat, tspan) d_discontinuities_internal = initialize_d_discontinuities(tType, d_discontinuities, @@ -564,6 +570,13 @@ function DiffEqBase.__init( end end + if _tstops !== nothing + tstops = _tstops(parameter_values(integrator), prob.tspan) + for tstop in tstops + add_tstop!(integrator, tstop) + end + end + handle_dt!(integrator) integrator end diff --git a/test/interface/ode_tstops_tests.jl b/test/interface/ode_tstops_tests.jl index ce85f0e859..a911ecfe8a 100644 --- a/test/interface/ode_tstops_tests.jl +++ b/test/interface/ode_tstops_tests.jl @@ -76,3 +76,15 @@ end prob = ODEProblem(ff, [0.0], (0.0f0, 1.0f0)) sol = solve(prob, Tsit5(), tstops = [tval], callback = cb) end + +@testset "Late binding tstops" begin + function rhs(u, p, t) + u * p + t + end + prob = ODEProblem(rhs, 1.0, (0.0, 1.0), 0.1; tstops = (p, tspan) -> tspan[1]:p:tspan[2]) + sol = solve(prob, Tsit5()) + @test 0.0:0.1:1.0 ⊆ sol.t + prob2 = remake(prob; p = 0.07) + sol2 = solve(prob2, Tsit5()) + @test 0.0:0.07:1.0 ⊆ sol2.t +end