Split ODE Problems
SciMLBase.SplitODEProblem
— TypeDefines a split ordinary differential equation (ODE) problem. Documentation Page: https://docs.sciml.ai/DiffEqDocs/stable/types/splitodetypes/
Mathematical Specification of a Split ODE Problem
To define a SplitODEProblem
, you simply need to give two functions $f_1$ and $f_2$ along with an initial condition $u_0$ which define an ODE:
\[\frac{du}{dt} = f_1(u,p,t) + f_2(u,p,t)\]
f
should be specified as f(u,p,t)
(or in-place as f(du,u,p,t)
), and u₀
should be an AbstractArray (or number) whose geometry matches the desired geometry of u
. Note that we are not limited to numbers or vectors for u₀
; one is allowed to provide u₀
as arbitrary matrices / higher dimension tensors as well.
Many splits are at least partially linear. That is the equation:
\[\frac{du}{dt} = Au + f_2(u,p,t)\]
For how to define a linear function A
, see the documentation for the AbstractSciMLOperator.
Constructors
SplitODEProblem(f::SplitFunction,u0,tspan,p=NullParameters();kwargs...)
SplitODEProblem{isinplace}(f1,f2,u0,tspan,p=NullParameters();kwargs...)
The isinplace
parameter can be omitted and will be determined using the signature of f2
. Note that both f1
and f2
should support the in-place style if isinplace
is true
or they should both support the out-of-place style if isinplace
is false
. You cannot mix up the two styles.
Parameters are optional, and if not given, then a NullParameters()
singleton will be used which will throw nice errors if you try to index non-existent parameters. Any extra keyword arguments are passed on to the solvers. For example, if you set a callback
in the problem, then that callback
will be added in every solve call.
Under the hood, a SplitODEProblem
is just a regular ODEProblem
whose f
is a SplitFunction
. Therefore, you can solve a SplitODEProblem
using the same solvers for ODEProblem
. For solvers dedicated to split problems, see Split ODE Solvers.
For specifying Jacobians and mass matrices, see the DiffEqFunctions page.
Fields
f1
,f2
: The functions in the ODE.u0
: The initial condition.tspan
: The timespan for the problem.p
: The parameters for the problem. Defaults toNullParameters
kwargs
: The keyword arguments passed onto the solves.
SciMLBase.SplitFunction
— TypeDocStringExtensions.TypeDefinition()
A representation of a split ODE function f
, defined by:
\[M \frac{du}{dt} = f_1(u,p,t) + f_2(u,p,t)\]
and all of its related functions, such as the Jacobian of f
, its gradient with respect to time, and more. For all cases, u0
is the initial condition, p
are the parameters, and t
is the independent variable.
Generally, for ODE integrators the f_1
portion should be considered the "stiff portion of the model" with larger timescale separation, while the f_2
portion should be considered the "non-stiff portion". This interpretation is directly used in integrators like IMEX (implicit-explicit integrators) and exponential integrators.
Constructor
SplitFunction{iip,specialize}(f1,f2;
mass_matrix = __has_mass_matrix(f1) ? f1.mass_matrix : I,
analytic = __has_analytic(f1) ? f1.analytic : nothing,
tgrad= __has_tgrad(f1) ? f1.tgrad : nothing,
jac = __has_jac(f1) ? f1.jac : nothing,
jvp = __has_jvp(f1) ? f1.jvp : nothing,
vjp = __has_vjp(f1) ? f1.vjp : nothing,
jac_prototype = __has_jac_prototype(f1) ? f1.jac_prototype : nothing,
W_prototype = __has_W_prototype(f1) ? f1.W_prototype : nothing,
sparsity = __has_sparsity(f1) ? f1.sparsity : jac_prototype,
paramjac = __has_paramjac(f1) ? f1.paramjac : nothing,
colorvec = __has_colorvec(f1) ? f1.colorvec : nothing,
sys = __has_sys(f1) ? f1.sys : nothing)
Note that only the functions f_i
themselves are required. These functions should be given as f_i!(du,u,p,t)
or du = f_i(u,p,t)
. See the section on iip
for more details on in-place vs out-of-place handling.
All of the remaining functions are optional for improving or accelerating the usage of the SplitFunction
. These include:
mass_matrix
: the mass matrixM
represented in the ODE function. Can be used to determine that the equation is actually a differential-algebraic equation (DAE) ifM
is singular. Note that in this case special solvers are required, see the DAE solver page for more details: https://docs.sciml.ai/DiffEqDocs/stable/solvers/dae_solve/. Must be an AbstractArray or an AbstractSciMLOperator.analytic(u0,p,t)
: used to pass an analytical solution function for the analytical solution of the ODE. Generally only used for testing and development of the solvers.tgrad(dT,u,p,t)
or dT=tgrad(u,p,t): returns $\frac{\partial f_1(u,p,t)}{\partial t}$jac(J,u,p,t)
orJ=jac(u,p,t)
: returns $\frac{df_1}{du}$jvp(Jv,v,u,p,t)
orJv=jvp(v,u,p,t)
: returns the directional derivative $\frac{df_1}{du} v$vjp(Jv,v,u,p,t)
orJv=vjp(v,u,p,t)
: returns the adjoint derivative $\frac{df_1}{du}^\ast v$jac_prototype
: a prototype matrix matching the type that matches the Jacobian. For example, if the Jacobian is tridiagonal, then an appropriately sizedTridiagonal
matrix can be used as the prototype and integrators will specialize on this structure where possible. Non-structured sparsity patterns should use aSparseMatrixCSC
with a correct sparsity pattern for the Jacobian. The default isnothing
, which means a dense Jacobian.W_prototype
: a prototype matrix matching the type that matches the W matrix. For example, if the Jacobian is tridiagonal, and the mass_matrix is diagonal, then an appropriately sizedTridiagonal
matrix can be used as the prototype and integrators will specialize on this structure where possible. Non-structured sparsity patterns should use aSparseMatrixCSC
with a correct sparsity pattern for the W matrix. The default isnothing
, which means a W of appropriate type for the jacobian and linear solverparamjac(pJ,u,p,t)
: returns the parameter Jacobian $\frac{df_1}{dp}$.colorvec
: a color vector according to the SparseDiffTools.jl definition for the sparsity pattern of thejac_prototype
. This specializes the Jacobian construction when using finite differences and automatic differentiation to be computed in an accelerated manner based on the sparsity pattern. Defaults tonothing
, which means a color vector will be internally computed on demand when required. The cost of this operation is highly dependent on the sparsity pattern.
Note on the Derivative Definition
The derivatives, such as the Jacobian, are only defined on the f1
portion of the split ODE. This is used to treat the f1
implicit while keeping the f2
portion explicit.
iip: In-Place vs Out-Of-Place
For more details on this argument, see the ODEFunction documentation.
specialize: Controlling Compilation and Specialization
For more details on this argument, see the ODEFunction documentation.
Fields
The fields of the SplitFunction type directly match the names of the inputs.
Symbolically Generating the Functions
See the modelingtoolkitize
function from ModelingToolkit.jl for automatically symbolically generating the Jacobian and more from the numerically-defined functions. See ModelingToolkit.SplitODEProblem
for information on generating the SplitFunction from this symbolic engine.
Solution Type
SplitODEProblem
solutions return an ODESolution
. For more information, see the ODE problem definition page for the ODESolution
docstring.