DiffEqCallbacks.gauss_points
— Constantgauss_points::Vector{Vector{Float64}}
Precomputed Gaussian nodes up to degree 2*10-1 = 19. Computed using FastGaussQuadrature.jl with the command [gausslegendre(i)[1] for i in 1:10]
DiffEqCallbacks.gauss_weights
— Constantgauss_weights::Vector{Vector{Float64}}
Precomputed Gaussian node weights up to degree 2*10-1 = 19. Computed using FastGaussQuadrature.jl with the command [gausslegendre(i)[2] for i in 1:10]
DiffEqCallbacks.CachePool
— TypeCachePool
Simple memory-reusing cache that allows us to grow a cache and keep re-using those pieces of memory (in our case, typically u
vectors) until the solve is finished. Note that this datastructure is not thread-safe!
DiffEqCallbacks.IndependentlyLinearizedSolution
— TypeIndependentlyLinearizedSolution
Efficient datastructure that holds a set of independently linearized solutions (obtained via the LinearizingSavingCallback
) with related, but slightly different time vectors. Stores a single time vector with a packed BitMatrix
denoting which u
vectors are sampled at which timepoints. Provides an efficient iterate()
method that can be used to reconstruct coherent views of the state variables at all timepoints, as well as an efficient sample!()
method that can sample at arbitrary timesteps.
DiffEqCallbacks.IndependentlyLinearizedSolutionChunks
— TypeIndependentlyLinearizedSolutionChunks
When constructing an IndependentlyLinearizedSolution
via the IndependentlyLinearizingCallback
, we use this indermediate structure to reduce allocations and collect the unknown number of timesteps that the solve will generate.
DiffEqCallbacks.IntegrandValues
— TypeIntegrandValues{integrandType}
A struct used to save values of the integrand values in integrand::Vector{integrandType}
.
DiffEqCallbacks.IntegrandValues
— MethodIntegrandValues(integrandType::DataType)
Return IntegrandValues{integrandType}
with empty storage vectors.
DiffEqCallbacks.IntegrandValuesSum
— TypeIntegrandValuesSum{integrandType}
A struct used to save values of the integrand values in integrand::Vector{integrandType}
.
DiffEqCallbacks.IntegrandValuesSum
— MethodIntegrandValuesSum(integrandType::DataType)
Return IntegrandValuesSum{integrandType}
with empty storage vectors.
DiffEqCallbacks.ManifoldProjection
— TypeManifoldProjection(g; nlsolve = nothing, save = true, nlls = Val(true),
isinplace = Val(true), autonomous = nothing, resid_prototype = nothing, kwargs...)
In many cases, you may want to declare a manifold on which a solution lives. Mathematically, a manifold M
is defined by a function g
as the set of points where g(u) = 0
. An embedded manifold can be a lower dimensional object which constrains the solution. For example, g(u) = E(u) - C
where E
is the energy of the system in state u
, meaning that the energy must be constant (energy preservation). Thus by defining the manifold the solution should live on, you can retain desired properties of the solution.
ManifoldProjection
projects the solution of the differential equation to the chosen manifold g
, conserving a property while conserving the order. It is a consequence of convergence proofs both in the deterministic and stochastic cases that post-step projection to manifolds keep the same convergence rate, thus any algorithm can be easily extended to conserve properties. If the solution is supposed to live on a specific manifold or conserve such property, this guarantees the conservation law without modifying the convergence properties.
Arguments
g
: The residual function for the manifold.- This is an inplace function of form
g(resid, u, p)
org(resid, u, p, t)
which writes to the residualresid
the difference from the manifold components. Here, it is assumed thatresid
is of the same shape asu
. - If
isinplace = Val(false)
, theng
should be a function of the formg(u, p)
org(u, p, t)
which returns the residual.
- This is an inplace function of form
Keyword Arguments
nlsolve
: A nonlinear solver as defined in the NonlinearSolve.jl format. Defaults to a PolyAlgorithm.save
: Whether to do the standard saving (applied after the callback)nlls
: If the problem is a nonlinear least squares problem.nlls = Val(false)
generates aNonlinearProblem
which is typically faster thanNonlinearLeastSquaresProblem
, but is only applicable if the residual size is same as the state size.autonomous
: Whetherg
is an autonomous function of the formg(resid, u, p)
org(u, p)
. Specify it asVal(::Bool)
to ensure this function call is type stable.residual_prototype
: This needs to be specified ifnlls = Val(true)
andinplace = Val(true)
are specified together, else it is taken to be same asu
.kwargs
: All other keyword arguments are passed to NonlinearSolve.jl.
Saveat Warning
Note that the ManifoldProjection
callback modifies the endpoints of the integration intervals and thus breaks assumptions of internal interpolations. Because of this, the values for given by saveat will not be order-matching. However, the interpolation error can be proportional to the change by the projection, so if the projection is making small changes then one is still safe. However, if there are large changes from each projection, you should consider only saving at stopping/projection times. To do this, set tstops
to the same values as saveat
. There is a performance hit by doing so because now the integrator is forced to stop at every saving point, but this is guerenteed to match the order of the integrator even with the ManifoldProjection.
References
Ernst Hairer, Christian Lubich, Gerhard Wanner. Geometric Numerical Integration: Structure-Preserving Algorithms for Ordinary Differential Equations. Berlin ; New York :Springer, 2002.
DiffEqCallbacks.SavedValues
— TypeSavedValues{tType<:Real, savevalType}
A struct used to save values of the time in t::Vector{tType}
and additional values in saveval::Vector{savevalType}
.
DiffEqCallbacks.SavedValues
— MethodSavedValues(tType::DataType, savevalType::DataType)
Return SavedValues{tType, savevalType}
with empty storage vectors.
DiffEqCallbacks.AdaptiveProbIntsUncertainty
— FunctionAdaptiveProbIntsUncertainty(order, save = true)
The ProbInts method for uncertainty quantification involves the transformation of an ODE into an associated SDE where the noise is related to the timesteps and the order of the algorithm.
AdaptiveProbIntsUncertainty
is a more automated form of ProbIntsUncertainty
which uses the error estimate from within adaptive time stepping methods to estimate σ
at every step.
Arguments
order
is the order of the ODE solver algorithm.save
is for choosing whether this callback should control the saving behavior. Generally this is true unless one is stacking callbacks in aCallbackSet
.
References
Conrad P., Girolami M., Särkkä S., Stuart A., Zygalakis. K, Probability Measures for Numerical Solutions of Differential Equations, arXiv:1506.04592
DiffEqCallbacks.AutoAbstol
— FunctionAutoAbstol(save = true; init_curmax = 1e-6)
Provides a way to automatically adapt the absolute tolerance to the problem. This helps the solvers automatically “learn” what appropriate limits are. This callback set starts the absolute tolerance at init_curmax
(default 1e-6
), and at each iteration it is set to the maximum value that the state has thus far reached times the relative tolerance.
Keyword Arguments
save
determines whether this callback has saving enabledinit_curmax
is the initialabstol
.
If this callback is used in isolation, save=true
is required for normal saving behavior. Otherwise, save=false
should be set to ensure extra saves do not occur.
DiffEqCallbacks.FunctionCallingCallback
— MethodFunctionCallingCallback(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 calls 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 the initial condition.tdir
should besign(tspan[end]-tspan[1])
. It defaults to1
and should be adapted iftspan[1] > tspan[end]
.
DiffEqCallbacks.GeneralDomain
— FunctionGeneralDomain(
g, u = nothing; save = true, abstol = nothing, scalefactor = nothing,
autonomous = maximum(SciMLBase.numargs(g)) == 3, nlsolve_kwargs = (;
abstol = 10 * eps()), kwargs...)
A GeneralDomain
callback in DiffEqCallbacks.jl generalizes the concept of a PositiveDomain
callback to arbitrary domains. Domains are specified by in-place functions g(resid, u, p)
or g(resid, u, p, t)
that calculate residuals of a state vector u
at time t
relative to that domain, with p
the parameters of the corresponding integrator. As for PositiveDomain
, steps are accepted if residuals of the extrapolated values at the next time step are below a certain tolerance. Moreover, this callback is automatically coupled with a ManifoldProjection
that keeps all calculated state vectors close to the desired domain, but in contrast to a PositiveDomain
callback the nonlinear solver in a ManifoldProjection
cannot guarantee that all state vectors of the solution are actually inside the domain. Thus, a PositiveDomain
callback should generally be preferred.
Arguments
g
: the implicit definition of the domain as a functiong(resid, u, p)
org(resid, u, p, t)
which is zero when the value is in the domain.u
: A prototype of the state vector of the integrator. A copy of it is saved and extrapolated values are written to it. If it is not specified, every application of the callback allocates a new copy of the state vector.
Keyword Arguments
save
: Whether to do the standard saving (applied after the callback).abstol
: Tolerance up to, which residuals are accepted. Element-wise tolerances are allowed. If it is not specified, every application of the callback uses the current absolute tolerances of the integrator.scalefactor
: Factor by which an unaccepted time step is reduced. If it is not specified, time steps are halved.autonomous
: Whetherg
is an autonomous function of the formg(resid, u, p)
. If it is not specified, it is determined automatically.kwargs
: All other keyword arguments are passed toManifoldProjection
.nlsolve_kwargs
: All keyword arguments are passed to the nonlinear solver inManifoldProjection
. The default is(; abstol = 10 * eps())
.
References
Shampine, Lawrence F., Skip Thompson, Jacek Kierzenka and G. D. Byrne. Non-negative solutions of ODEs. Applied Mathematics and Computation 170 (2005): 556-569.
DiffEqCallbacks.IntegratingCallback
— MethodIntegratingCallback(integrand_func,
integrand_values::IntegrandValues, integrand_prototype)
Let one define a function integrand_func(u, t, integrator)::typeof(integrand_prototype)
which returns Integral(integrand_func(u(t),t)dt over the problem tspan.
Arguments
integrand_func(out, u, t, integrator)
for in-place problems andout = integrand_func(u, t, integrator)
for out-of-place problems. Returns the quantity in the integral for computing dG/dp. Note that for out-of-place problems, this should allocate the output (not as a view tou
).integrand_values::IntegrandValues
is the types thatintegrand_func
will return, i.e.integrand_func(t, u, integrator)::integrandType
. It's specified viaIntegrandValues(integrandType)
, i.e. give the type thatintegrand_func
will output (or higher compatible type).integrand_prototype
is a prototype of the output from the integrand.
The outputted values are saved into integrand_values
. The values are found via integrand_values.integrand
.
This method is currently limited to ODE solvers of order 10 or lower. Open an issue if other solvers are required.
If integrand_func
is in-place, you must use cache
to store the output of integrand_func
.
DiffEqCallbacks.IntegratingSumCallback
— MethodIntegratingCallback(integrand_func,
integrand_values::IntegrandValues,
cache = nothing)
Lets one define a function integrand_func(u, t, integrator)
which returns Integral(integrand_func(u(t),t)dt over the problem tspan.
Arguments
integrand_func(out, u, t, integrator)
for in-place problems andout = integrand_func(u, t, integrator)
for out-of-place problems. Returns the quantity in the integral for computing dG/dp. Note that for out-of-place problems, this should allocate the output (not as a view tou
).integrand_values::IntegrandValues
is the types thatintegrand_func
will return, i.e.integrand_func(t, u, integrator)::integrandType
. It's specified viaIntegrandValues(integrandType)
, i.e. give the type thatintegrand_func
will output (or higher compatible type).integrand_prototype
is a prototype of the output from the integrand.
The outputted values are saved into integrand_values
. The values are found via integrand_values.integrand
.
This method is currently limited to ODE solvers of order 10 or lower. Open an issue if other solvers are required.
DiffEqCallbacks.IterativeCallback
— FunctionIterativeCallback(time_choice, user_affect!, tType = Float64;
initial_affect = false, kwargs...)
A callback to be used to iteratively apply some affect. For example, if given the first effect at t₁
, you can define t₂
to apply the next effect.
Arguments
time_choice(integrator)
determines the time of the next callback. Ifnothing
is returned for the time choice, then the iterator ends.user_affect!
is the effect applied to the integrator at the stopping points.
Keyword Arguments
initial_affect
is whether to apply the affect att=0
which defaults tofalse
DiffEqCallbacks.LinearizingSavingCallback
— MethodLinearizingSavingCallback(ils::IndependentlyLinearizedSolution)
LinearizingSavingCallback(ilss::Vector{IndependentlyLinearizedSolution})
Provides a saving callback that inserts interpolation points into your signal such that a naive linear interpolation of the resultant saved values will be within abstol
/reltol
of the higher-order interpolation of your solution. This essentially makes a time/space tradeoff, where more points in time are saved, costing more memory, but interpolation is incredibly cheap and downstream algorithm complexity is reduced by not needing to bother with multiple interpolation types.
The algorithm internally checks 3 equidistant points between each time point to determine goodness of fit versus the linearly interpolated function; this should be sufficient for interpolations up to the 4th order, higher orders may need more points to ensure good fit. This has not been implemented yet.
This callback generator takes in an IndependentlyLinearizedSolution
object to store output into. An IndependentlyLinearizedSolution
object itself controls how many derivatives (if any) to linearize along with the primal states themselves.
Example usage:
ils = IndependentlyLinearizedSolution(prob)
solve(prob, solver; callback=LinearizingSavingCallback(ils))
Keyword Arguments
interpolate_mask::BitVector
: a set ofu
indices for which the integrator interpolant can be queried. Any false indices will be linearly-interpolated based on thesol.t
points instead (no subdivision). This is useful for when a solution needs to ignore certain indices due to badly-behaved interpolation.
DiffEqCallbacks.PeriodicCallback
— MethodPeriodicCallback(f, Δt::Number; initial_affect = false,
final_affect = false,
kwargs...)
PeriodicCallback
can be used when a function should be called periodically in terms of integration time (as opposed to wall time), i.e. at t = tspan[1]
, t = tspan[1] + Δt
, t = tspan[1] + 2Δt
, and so on. This callback can, for example, be used to model a discrete-time controller for a continuous-time system, running at a fixed rate.
Arguments
f
theaffect!(integrator)
function to be called periodicallyΔt
is the period
Keyword Arguments
initial_affect
is whether to apply the affect att=0
, which defaults tofalse
final_affect
is whether to apply the affect at the final time, which defaults tofalse
kwargs
are keyword arguments accepted by theDiscreteCallback
constructor.
DiffEqCallbacks.PositiveDomain
— FunctionPositiveDomain(u = nothing; save = true, abstol = nothing, scalefactor = nothing)
Especially in biology and other natural sciences, a desired property of dynamical systems is the positive invariance of the positive cone, i.e. non-negativity of variables at time $t_0$ ensures their non-negativity at times $t \geq t_0$ for which the solution is defined. However, even if a system satisfies this property mathematically it can be difficult for ODE solvers to ensure it numerically, as these MATLAB examples show.
To deal with this problem, one can specify isoutofdomain=(u,p,t) -> any(x -> x < 0, u)
as additional solver option, which will reject any step that leads to non-negative values and reduce the next time step. However, since this approach only rejects steps and hence calculations might be repeated multiple times until a step is accepted, it can be computationally expensive.
Another approach is taken by a PositiveDomain
callback in DiffEqCallbacks.jl, which is inspired by Shampine's et al. paper about non-negative ODE solutions. It reduces the next step by a certain scale factor until the extrapolated value at the next time point is non-negative with a certain tolerance. Extrapolations are cheap to compute but might be inaccurate, so if a time step is changed it is additionally reduced by a safety factor of 0.9. Since extrapolated values are only non-negative up to a certain tolerance and in addition actual calculations might lead to negative values, also any negative values at the current time point are set to 0. Hence, by this callback non-negative values at any time point are ensured in a computationally cheap way, but the quality of the solution depends on how accurately extrapolations approximate next time steps.
Please note, that the system should be defined also outside the positive domain, since even with these approaches, negative variables might occur during the calculations. Moreover, one should follow Shampine's et al. advice and set the derivative $x'_i$ of a negative component $x_i$ to $\max \{0, f_i(x, t)\}$, where $t$ denotes the current time point with state vector $x$ and $f_i$ is the $i$-th component of function $f$ in an ODE system $x' = f(x, t)$.
Arguments
u
: A prototype of the state vector of the integrator. A copy of it is saved and extrapolated values are written to it. If it is not specified, every application of the callback allocates a new copy of the state vector.
Keyword Arguments
save
: Whether to do the standard saving (applied after the callback).abstol
: Tolerance up to, which negative extrapolated values are accepted. Element-wise tolerances are allowed. If it is not specified, every application of the callback uses the current absolute tolerances of the integrator.scalefactor
: Factor by which an unaccepted time step is reduced. If it is not specified, time steps are halved.
References
Shampine, Lawrence F., Skip Thompson, Jacek Kierzenka and G. D. Byrne. Non-negative solutions of ODEs. Applied Mathematics and Computation 170 (2005): 556-569.
DiffEqCallbacks.PresetTimeCallback
— MethodPresetTimeCallback(tstops, user_affect!;
initialize = DiffEqBase.INITIALIZE_DEFAULT,
filter_tstops = true,
kwargs...)
A callback that adds callback affect!
calls at preset times. No playing around with tstops
or anything is required: this callback adds the triggers for you to make it automatic.
Arguments
tstops
: the times for theaffect!
to trigger at.user_affect!
: anaffect!(integrator)
function to use at the time points.
Keyword Arguments
filter_tstops
: Whether to filter out tstops beyond the end of the integration timespan. Defaults to true. If false, then tstops can extend the interval of integration.
DiffEqCallbacks.ProbIntsUncertainty
— FunctionProbIntsUncertainty(σ, order, save = true)
The ProbInts method for uncertainty quantification involves the transformation of an ODE into an associated SDE where the noise is related to the timesteps and the order of the algorithm.
Arguments
σ
is the noise scaling factor. It is recommended thatσ
is representative of the size of the errors in a single step of the equation. If such a value is unknown, it can be estimated automatically in adaptive time-stepping algorithms via AdaptiveProbIntsUncertaintyorder
is the order of the ODE solver algorithm.save
is for choosing whether this callback should control the saving behavior. Generally this is true unless one is stacking callbacks in aCallbackSet
.
References
Conrad P., Girolami M., Särkkä S., Stuart A., Zygalakis. K, Probability Measures for Numerical Solutions of Differential Equations, arXiv:1506.04592
DiffEqCallbacks.SavingCallback
— MethodSavingCallback(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.StepsizeLimiter
— MethodStepsizeLimiter(dtFE;safety_factor=9//10,max_step=false,cached_dtcache=0.0)
In many cases, there is a known maximal stepsize for which the computation is stable and produces correct results. For example, in hyperbolic PDEs one normally needs to ensure that the stepsize stays below some $\Delta t_{FE}$ determined by the CFL condition. For nonlinear hyperbolic PDEs this limit can be a function dtFE(u,p,t)
which changes throughout the computation. The stepsize limiter lets you pass a function which will adaptively limit the stepsizes to match these constraints.
Arguments
dtFE
is the maximal timestep and is calculated using the previoust
andu
.
Keyword Arguments
safety_factor
is the factor below the true maximum that will be stepped to which defaults to9//10
.max_step=true
makes every step equal tosafety_factor*dtFE(u,p,t)
when the solver is set toadaptive=false
.cached_dtcache
should be set to match the type for time when not using Float64 values.
DiffEqCallbacks.TerminateSteadyState
— MethodTerminateSteadyState(abstol = 1e-8, reltol = 1e-6, test = allDerivPass; min_t = nothing,
wrap_test::Val = Val(true))
TerminateSteadyState
can be used to solve the problem for the steady-state by running the solver until the derivatives of the problem converge to 0 or tspan[2]
is reached. This is an alternative approach to root finding (see the Steady State Solvers section).
Arguments
abstol
andreltol
are the absolute and relative tolerance, respectively. These tolerances may be specified as scalars or as arrays of the same length as the states of the problem.test
represents the function that evaluates the condition for termination. The default condition is that all derivatives should become smaller thanabstol
and the states timesreltol
. The user can pass any other function to implement a different termination condition. Such function should take four arguments:integrator
,abstol
,reltol
, andmin_t
.wrap_test
can be set toVal(false)
, in which casetest
must have the definitiontest(u, t, integrator)
. Otherwise,test
must have the definitiontest(integrator, abstol, reltol, min_t)
.
Keyword Arguments
min_t
specifies an optional minimumt
before the steady state calculations are allowed to terminate.
DiffEqCallbacks.affect!
— Methodaffect!(integrator, f::AbstractDomainAffect)
Apply domain callback f
to integrator
.
DiffEqCallbacks.allocate_vjp
— Methodallocate_vjp(λ, x)
allocate_vjp(x)
similar(λ, size(x))
for generic x
. This is used to handle non-array parameters!
DiffEqCallbacks.allocate_zeros
— Methodallocate_zeros(x)
zero.(x)
for generic x
. This is used to handle non-array parameters!
DiffEqCallbacks.isaccepted
— Methodisaccepted(u, abstol, f::AbstractDomainAffect, args...)
Return whether u
is an acceptable state vector at the next time point given absolute tolerance abstol
, callback f
, and other optional arguments.
DiffEqCallbacks.modify_u!
— Methodmodify_u!(integrator, f::AbstractDomainAffect)
Modify current state vector u
of integrator
if required, and return whether it actually was modified.
DiffEqCallbacks.recursive_add!
— Methodrecursive_add!(y, x)
y .+= x
for generic x
and y
. This is used to handle non-array parameters!
DiffEqCallbacks.recursive_adjoint
— Methodrecursive_adjoint(y)
adjoint(y)
for generic y
. This is used to handle non-array parameters!
DiffEqCallbacks.recursive_copy
— Methodrecursive_copy(y)
copy(y)
for generic y
. This is used to handle non-array parameters!
DiffEqCallbacks.recursive_copyto!
— Methodrecursive_copyto!(y, x)
y[:] .= vec(x)
for generic x
and y
. This is used to handle non-array parameters!
DiffEqCallbacks.recursive_neg!
— Methodneg!(x)
x .*= -1
for generic x
. This is used to handle non-array parameters!
DiffEqCallbacks.recursive_sub!
— Methodrecursive_sub!(y, x)
y .-= x
for generic x
and y
. This is used to handle non-array parameters!
DiffEqCallbacks.recursive_zero!
— Methodzero!(x)
x .= 0
for generic x
. This is used to handle non-array parameters!
DiffEqCallbacks.sample!
— Methodsample!(out::Matrix{S}, ils::IndependentlyLinearizedSolution, ts::Vector{T}, deriv_idx::Int = 0)
Batch-sample ils
at the given timepoints for the given derivative level, storing into out
.
DiffEqCallbacks.seek_forward
— Methodseek_forward(ils::IndependentlyLinearizedSolution, cursor::ILSStateCursor, t_target)
Seek the given cursor
forward until it contains t_target
. Does not seek backward, use seek()
for the more general formulation, this form is optimized for the inner loop of iterate()
.
DiffEqCallbacks.setup
— Methodsetup(f::AbstractDomainAffect, integrator)
Setup callback f
and return an arbitrary tuple whose elements are used as additional arguments in checking whether time step is accepted.
DiffEqCallbacks.store!
— Methodstore!(ilsc::IndependentlyLinearizedSolutionChunks, t, us, u_mask)
Store a new us
matrix (one row per derivative level) into our ilsc
, but only the values identified by the given u_mask
. The us
matrix should be of the size (num_us(ilsc), num_derivatives(ilsc))
.