Skip to content

Commit

Permalink
feat: add late binding for tstops
Browse files Browse the repository at this point in the history
  • Loading branch information
AayushSabharwal committed Nov 15, 2024
1 parent 6ec790a commit c0c112f
Show file tree
Hide file tree
Showing 2 changed files with 25 additions and 0 deletions.
13 changes: 13 additions & 0 deletions lib/OrdinaryDiffEqCore/src/solve.jl
Original file line number Diff line number Diff line change
Expand Up @@ -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,
Expand Down Expand Up @@ -542,6 +548,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
Expand Down
12 changes: 12 additions & 0 deletions test/interface/ode_tstops_tests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -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

0 comments on commit c0c112f

Please sign in to comment.