Skip to content
Draft
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 2 additions & 0 deletions Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -17,6 +17,7 @@ NonlinearSolveBase = "be0214bd-f91f-a760-ac4e-3421ce2b2da0"
PrecompileTools = "aea7be01-6a6a-4083-8856-8a6e6704d82a"
Reexport = "189a3867-3050-52da-a836-e630ba90ab69"
SciMLBase = "0bca4576-84f4-4d90-8ffe-ffa030f20462"
SciMLLogging = "a6db7da4-7206-11f0-1eab-35f2a5dbe1d1"
SparseArrays = "2f01184e-e22b-5df5-ae63-d93ebab69eaf"
Sundials_jll = "fb77eaff-e24c-56d4-86b1-d163f2edb164"
SymbolicIndexingInterface = "2efcf032-c050-4f8e-a9bb-153293bab1f5"
Expand Down Expand Up @@ -45,6 +46,7 @@ PrecompileTools = "1"
Reexport = "1.0"
SafeTestsets = "0.1"
SciMLBase = "2.119.0"
SciMLLogging = "1.9.1"
SparseArrays = "1"
SparseConnectivityTracer = "1"
Sundials_jll = "7.4.1"
Expand Down
3 changes: 2 additions & 1 deletion src/Sundials.jl
Original file line number Diff line number Diff line change
Expand Up @@ -11,7 +11,8 @@ using DiffEqBase: DiffEqBase, NonlinearFunction, ODEFunction, add_saveat!,
reeval_internals_due_to_modification!, reinit!, savevalues!,
set_proposed_dt!, solve, solve!, step!, terminate!, u_modified!,
update_coefficients!, warn_compat, DefaultInit, BrownFullBasicInit,
ShampineCollocationInit
ShampineCollocationInit, DEVerbosity
using SciMLLogging: SciMLLogging, Standard, @SciMLMessage
using SciMLBase: AbstractSciMLOperator, DAEProblem, ODEProblem, ReturnCode,
SciMLBase, SplitODEProblem, VectorContinuousCallback
import Accessors: @reset
Expand Down
18 changes: 9 additions & 9 deletions src/common_interface/integrator_types.jl
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,7 @@ mutable struct DEOptions{
callback::CType
abstol::abstolType
reltol::reltolType
verbose::Bool
verbose::DEVerbosity
advance_to_tstop::Bool
stop_at_next_tstop::Bool
progress::Bool
Expand Down Expand Up @@ -82,7 +82,7 @@ function (integrator::CVODEIntegrator)(
) where {T}
out = similar(integrator.u)
out_nvec = NVector(vec(out), integrator.ctx_handle.ctx)
integrator.flag = @checkflag CVodeGetDky(integrator.mem, t, Cint(T), out_nvec)
integrator.flag = @checkflag CVodeGetDky(integrator.mem, t, Cint(T), out_nvec) false integrator.opts.verbose
copyto!(out, out_nvec.v)
return idxs === nothing ? out : out[idxs]
end
Expand All @@ -94,7 +94,7 @@ function (integrator::CVODEIntegrator)(
idxs = nothing
) where {T}
out_nvec = NVector(vec(out), integrator.ctx_handle.ctx)
integrator.flag = @checkflag CVodeGetDky(integrator.mem, t, Cint(T), out_nvec)
integrator.flag = @checkflag CVodeGetDky(integrator.mem, t, Cint(T), out_nvec) false integrator.opts.verbose
copyto!(out, out_nvec.v)
return idxs === nothing ? out : @view out[idxs]
end
Expand Down Expand Up @@ -164,7 +164,7 @@ function (integrator::ARKODEIntegrator{
}
out = similar(integrator.u)
out_nvec = NVector(vec(out), integrator.ctx_handle.ctx)
integrator.flag = @checkflag ARKStepGetDky(integrator.mem, t, Cint(T), out_nvec)
integrator.flag = @checkflag ARKStepGetDky(integrator.mem, t, Cint(T), out_nvec) false integrator.opts.verbose
copyto!(out, out_nvec.v)
return idxs === nothing ? out : out[idxs]
end
Expand All @@ -182,7 +182,7 @@ function (integrator::ARKODEIntegrator{
}
out = similar(integrator.u)
out_nvec = NVector(vec(out), integrator.ctx_handle.ctx)
integrator.flag = @checkflag ERKStepGetDky(integrator.mem, t, Cint(T), out_nvec)
integrator.flag = @checkflag ERKStepGetDky(integrator.mem, t, Cint(T), out_nvec) false integrator.opts.verbose
copyto!(out, out_nvec.v)
return idxs === nothing ? out : out[idxs]
end
Expand All @@ -200,7 +200,7 @@ function (integrator::ARKODEIntegrator{
LStype, Atype, MLStype, Mtype, CallbackCacheType, IA, T,
}
out_nvec = NVector(vec(out), integrator.ctx_handle.ctx)
integrator.flag = @checkflag ARKStepGetDky(integrator.mem, t, Cint(T), out_nvec)
integrator.flag = @checkflag ARKStepGetDky(integrator.mem, t, Cint(T), out_nvec) false integrator.opts.verbose
copyto!(out, out_nvec.v)
return idxs === nothing ? out : @view out[idxs]
end
Expand All @@ -218,7 +218,7 @@ function (integrator::ARKODEIntegrator{
LStype, Atype, MLStype, Mtype, CallbackCacheType, IA, T,
}
out_nvec = NVector(vec(out), integrator.ctx_handle.ctx)
integrator.flag = @checkflag ERKStepGetDky(integrator.mem, t, Cint(T), out_nvec)
integrator.flag = @checkflag ERKStepGetDky(integrator.mem, t, Cint(T), out_nvec) false integrator.opts.verbose
copyto!(out, out_nvec.v)
return idxs === nothing ? out : @view out[idxs]
end
Expand Down Expand Up @@ -278,7 +278,7 @@ function (integrator::IDAIntegrator)(
) where {T}
out = similar(integrator.u)
out_nvec = NVector(vec(out), integrator.ctx_handle.ctx)
integrator.flag = @checkflag IDAGetDky(integrator.mem, t, Cint(T), out_nvec)
integrator.flag = @checkflag IDAGetDky(integrator.mem, t, Cint(T), out_nvec) false integrator.opts.verbose
copyto!(out, out_nvec.v)
return idxs === nothing ? out : out[idxs]
end
Expand All @@ -290,7 +290,7 @@ function (integrator::IDAIntegrator)(
idxs = nothing
) where {T}
out_nvec = NVector(vec(out), integrator.ctx_handle.ctx)
integrator.flag = @checkflag IDAGetDky(integrator.mem, t, Cint(T), out_nvec)
integrator.flag = @checkflag IDAGetDky(integrator.mem, t, Cint(T), out_nvec) false integrator.opts.verbose
copyto!(out, out_nvec.v)
return idxs === nothing ? out : @view out[idxs]
end
Expand Down
48 changes: 33 additions & 15 deletions src/common_interface/solve.jl
Original file line number Diff line number Diff line change
Expand Up @@ -115,7 +115,7 @@ function DiffEqBase.__init(
timeseries = [],
ts = [],
ks = [];
verbose = true,
verbose = Standard(),
callback = nothing,
abstol = 1 / 10^6,
reltol = 1 / 10^3,
Expand Down Expand Up @@ -153,10 +153,16 @@ function DiffEqBase.__init(
uType, tupType, isinplace, Method, LinearSolver,
}
tType = eltype(tupType)

if verbose
warned = !isempty(kwargs) && DiffEqBase.check_keywords(alg, kwargs, warnlist)
warned && DiffEqBase.warn_compat()
verbose = _process_verbose_param(verbose)

if !isempty(kwargs)
warned = DiffEqBase.check_keywords(alg, kwargs, warnlist)
if warned
@SciMLMessage(
"Unrecognized keyword arguments passed to Sundials solver.",
verbose, :warn_compat
)
end
end

if prob.f.mass_matrix != LinearAlgebra.I
Expand Down Expand Up @@ -527,7 +533,7 @@ function DiffEqBase.__init(
timeseries = [],
ts = [],
ks = [];
verbose = true,
verbose = Standard(),
callback = nothing,
abstol = 1 / 10^6,
reltol = 1 / 10^3,
Expand Down Expand Up @@ -567,10 +573,16 @@ function DiffEqBase.__init(
MassLinearSolver,
}
tType = eltype(tupType)

if verbose
warned = !isempty(kwargs) && DiffEqBase.check_keywords(alg, kwargs, warnlist)
warned && DiffEqBase.warn_compat()
verbose = _process_verbose_param(verbose)

if !isempty(kwargs)
warned = DiffEqBase.check_keywords(alg, kwargs, warnlist)
if warned
@SciMLMessage(
"Unrecognized keyword arguments passed to Sundials solver.",
verbose, :warn_compat
)
end
end

if reltol isa AbstractArray
Expand Down Expand Up @@ -1143,7 +1155,7 @@ function DiffEqBase.__init(
timeseries = [],
ts = [],
ks = [];
verbose = true,
verbose = Standard(),
dt = nothing,
dtmax = 0.0,
save_on = true,
Expand Down Expand Up @@ -1175,10 +1187,16 @@ function DiffEqBase.__init(
uType, duType, tupType, isinplace, LinearSolver,
}
tType = eltype(tupType)

if verbose
warned = !isempty(kwargs) && DiffEqBase.check_keywords(alg, kwargs, warnida)
warned && DiffEqBase.warn_compat()
verbose = _process_verbose_param(verbose)

if !isempty(kwargs)
warned = DiffEqBase.check_keywords(alg, kwargs, warnida)
if warned
@SciMLMessage(
"Unrecognized keyword arguments passed to Sundials solver.",
verbose, :warn_compat
)
end
end

if reltol isa AbstractArray
Expand Down
5 changes: 5 additions & 0 deletions src/common_interface/verbosity.jl
Original file line number Diff line number Diff line change
@@ -1,3 +1,8 @@
function null_error_handler(error_code::Cint, mod::Char, func::Char, eh_data::Ptr{Cvoid})
return nothing
end

# Normalize the verbose parameter to DEVerbosity
_process_verbose_param(v::SciMLLogging.AbstractVerbosityPreset) = DEVerbosity(v)
_process_verbose_param(v::Bool) = v ? DEVerbosity() : DEVerbosity(SciMLLogging.None())
_process_verbose_param(v::DEVerbosity) = v
8 changes: 5 additions & 3 deletions src/simple.jl
Original file line number Diff line number Diff line change
Expand Up @@ -8,16 +8,18 @@
Insert a check that the given function call returns 0,
throw an error otherwise. Only apply directly to function calls.
"""
macro checkflag(ex, throw_error = false)
macro checkflag(ex, throw_error = false, verbose = nothing)
@assert Base.Meta.isexpr(ex, :call)
fname = ex.args[1]
fname_str = string(fname)
msg_expr = :($fname_str * " failed with error code = " * string(flag))
return quote
flag = $(esc(ex))
if flag < 0
if $(esc(throw_error))
@error($(string(fname, " failed with error code = ")), flag)
else
@warn($(string(fname, " failed with error code = ")), flag)
elseif $(esc(verbose)) !== nothing
@SciMLMessage(() -> $msg_expr, $(esc(verbose)), :checkflag)
end
end
flag
Expand Down
49 changes: 49 additions & 0 deletions test/common_interface/verbosity.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,49 @@
using Sundials, SciMLLogging, Test
import ODEProblemLibrary: prob_ode_linear

prob = prob_ode_linear

# Default (Standard preset) — should work without warnings
sol = solve(prob, CVODE_BDF())
@test sol.retcode == ReturnCode.Success

# Explicit Standard preset
sol = solve(prob, CVODE_BDF(), verbose = SciMLLogging.Standard())
@test sol.retcode == ReturnCode.Success

# Explicit None preset — no output
sol = solve(prob, CVODE_BDF(), verbose = SciMLLogging.None())
@test sol.retcode == ReturnCode.Success

# Bool true — maps to DEVerbosity(Standard)
sol = solve(prob, CVODE_BDF(), verbose = true)
@test sol.retcode == ReturnCode.Success

# Bool false — maps to DEVerbosity(None)
sol = solve(prob, CVODE_BDF(), verbose = false)
@test sol.retcode == ReturnCode.Success

# DEVerbosity with individual toggle
sol = solve(prob, CVODE_BDF(), verbose = DEVerbosity(checkflag = SciMLLogging.WarnLevel()))
@test sol.retcode == ReturnCode.Success

# DEVerbosity with sundials_specific group
sol = solve(prob, CVODE_BDF(), verbose = DEVerbosity(sundials_specific = SciMLLogging.WarnLevel()))
@test sol.retcode == ReturnCode.Success

# ARKODE
sol = solve(prob, ARKODE(), verbose = false)
@test sol.retcode == ReturnCode.Success

# IDA
function robertson(out, du, u, _, _)
out[1] = -0.04u[1] + 1.0e4 * u[2] * u[3] - du[1]
out[2] = +0.04u[1] - 3.0e7 * u[2]^2 - 1.0e4 * u[2] * u[3] - du[2]
return out[3] = u[1] + u[2] + u[3] - 1.0
end
dae_prob = DAEProblem(
robertson, [-0.04, 0.04, 0.0], [1.0, 0.0, 0.0], (0.0, 1.0e5),
differential_vars = [true, true, false]
)
sol = solve(dae_prob, IDA(), verbose = false)
@test sol.retcode == ReturnCode.Success
3 changes: 3 additions & 0 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -90,6 +90,9 @@ end
@safetestset "Initialization" begin
include("common_interface/initialization.jl")
end
@safetestset "Verbosity" begin
include("common_interface/verbosity.jl")
end
end

@testset "Interpolation" begin
Expand Down
Loading