Skip to content

Commit 89a643d

Browse files
committed
Fix handling of small timesteps
Sundials solver is allowed to warn and continue if it returns a timestep with no change in t. If the solver doesn't recover, it will fail with eg SciMLBase.ReturnCode.MaxIters Reverts PR #416, fixes #420 For CVode and ARKode, it looks like this is intended behaviour: the Sundial solver emits a warning message, with an API call CVodeSetMaxHnilWarns ARKStepSetMaxHnilWarns to set the maximum number of warning messages printed, which is exposed to Julia as max_hnil_warns see Sundials driver code https://github.com/LLNL/sundials/blob/e8a3e67e3883bc316c48bc534ee08319a5e8c620/src/cvode/cvode.c#L1339-L1347 For IDA, there is no API like this so it's not so clear what should happen, however the behaviour prior to #416 (restored by this PR) for the f_noconverge test added in test/common_interface/ida.jl is to return SciMLBase.ReturnCode.MaxIters which seems reasonable? (NB: @oscardssmith I can't reproduce the issue implied by the title of #416 so perhaps what is missing here is the original failure case that demonstrates the issue? The call to solve in the f_noconverge test test/common_interface/ida.jl, returns SciMLBase.ReturnCode.MaxIters, it doesn't return with no error? Also SciML/SciMLBase.jl#458 shows IDA printing [IDAS ERROR] IDASolve At t = 0 and h = 1.49012e-18, the error test failed repeatedly or with |h| = hmin. which suggests that IDA did flag an error in that case ? )
1 parent ca8ac7c commit 89a643d

File tree

2 files changed

+2
-4
lines changed

2 files changed

+2
-4
lines changed

src/common_interface/solve.jl

Lines changed: 1 addition & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -1406,9 +1406,7 @@ function DiffEqBase.solve!(integrator::AbstractSundialsIntegrator; early_free =
14061406
integrator.userfun.p = integrator.p
14071407
solver_step(integrator, tstop)
14081408
integrator.t = first(integrator.tout)
1409-
if integrator.t == integrator.tprev
1410-
integrator.flag = -3
1411-
end
1409+
# NB: CVode, ARKode may warn and then recover if integrator.t == integrator.tprev so don't flag this as an error
14121410
integrator.flag < 0 && break
14131411
handle_callbacks!(integrator) # this also updates the interpolation
14141412
integrator.flag < 0 && break

test/common_interface/ida.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -117,4 +117,4 @@ isapprox(only(sol.u[end]), exp(1), rtol = 1e-3)
117117
f_noconverge(out, du, u, p, t) = out .= [du[1] + u[1] / (t - 1)]
118118
prob = DAEProblem(f_noconverge, [1.0], [1.0], (0, 2); differential_vars = [true])
119119
sol = solve(prob, IDA())
120-
@test !(sol.retcode in (ReturnCode.Success, ReturnCode.MaxIters))
120+
@test !(sol.retcode in (ReturnCode.Success, ))

0 commit comments

Comments
 (0)