Output and Saving Controls
These callbacks extend the output and saving controls available during time stepping.
DiffEqCallbacks.SavingCallback
— FunctionSavingCallback(save_func, saved_values::SavedValues;
saveat = Vector{eltype(saved_values.t)}(),
save_everystep = isempty(saveat),
save_start = true,
tdir = 1)
The saving callback lets you define a function save_func(u, t, integrator)
which returns quantities of interest that shall be saved.
Arguments
save_func(u, t, integrator)
returns the quantities which shall be saved. Note that this should allocate the output (not as a view tou
).saved_values::SavedValues
is the types thatsave_func
will return, i.e.save_func(t, u, integrator)::savevalType
. It's specified viaSavedValues(typeof(t),savevalType)
, i.e. give the type for time and the type thatsave_func
will output (or higher compatible type).
Keyword Arguments
saveat
mimicssaveat
insolve
fromsolve
.save_everystep
mimicssave_everystep
fromsolve
.save_start
mimicssave_start
fromsolve
.save_end
mimicssave_end
fromsolve
.tdir
should besign(tspan[end]-tspan[1])
. It defaults to1
and should be adapted iftspan[1] > tspan[end]
.
The outputted values are saved into saved_values
. Time points are found via saved_values.t
and the values are saved_values.saveval
.
DiffEqCallbacks.FunctionCallingCallback
— FunctionFunctionCallingCallback(func;
funcat = Vector{Float64}(),
func_everystep = isempty(funcat),
func_start = true,
tdir = 1)
The function calling callback lets you define a function func(u,t,integrator)
which gets called at the time points of interest. The constructor is:
func(u, t, integrator)
is the function to be called.funcat
values or interval that the function is sure to be evaluated at.func_everystep
whether to call the function after each integrator step.func_start
whether the function is called at the initial condition.tdir
should besign(tspan[end]-tspan[1])
. It defaults to1
and should be adapted iftspan[1] > tspan[end]
.
Saving Example
In this example, we will solve a matrix equation and at each step save a tuple of values which contains the current trace and the norm of the matrix. We build the SavedValues
cache to use Float64
for time and Tuple{Float64,Float64}
for the saved values, and then call the solver with the callback.
using DiffEqCallbacks, OrdinaryDiffEq, LinearAlgebra
prob = ODEProblem((du, u, p, t) -> du .= u, rand(4, 4), (0.0, 1.0))
saved_values = SavedValues(Float64, Tuple{Float64, Float64})
cb = SavingCallback((u, t, integrator) -> (tr(u), norm(u)), saved_values)
sol = solve(prob, Tsit5(), callback = cb)
print(saved_values.saveval)
[(2.6084247209807607, 2.3098485766342063), (2.8831073547743635, 2.5530893669820784), (3.695481302812439, 3.2724740563224386), (5.1676986164505685, 4.576172429923441), (7.090433208291806, 6.27881913637696)]
Note that the values are retrieved from the cache as .saveval
, and the time points are found as .t
. If we want to control the saved times, we use saveat
in the callback. The save controls like saveat
act analogously to how they act in the solve
function.
saved_values = SavedValues(Float64, Tuple{Float64, Float64})
cb = SavingCallback((u, t, integrator) -> (tr(u), norm(u)), saved_values,
saveat = 0.0:0.1:1.0)
sol = solve(prob, Tsit5(), callback = cb)
print(saved_values.saveval)
print(saved_values.t)
[(2.6084247209807607, 2.3098485766342063), (2.8827551436783425, 2.5527774721084118), (3.1859372623976996, 2.8212555231527747), (3.521004871686383, 3.1179692577552767), (3.8913131281219084, 3.4458897808830757), (4.300565098403878, 3.808296283718671), (4.752857879366456, 4.208816884495497), (5.252722343294279, 4.651463824323444), (5.8051564370635935, 5.1406629547079605), (6.415687495162425, 5.681308917843256), (7.090433208291806, 6.27881913637696)][0.0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0]