SDE Problems

SciMLBase.SDEProblemType

Defines an stochastic differential equation (SDE) problem. Documentation Page: https://docs.sciml.ai/DiffEqDocs/stable/types/sde_types/

Mathematical Specification of a SDE Problem

To define an SDE Problem, you simply need to give the forcing function f, the noise function g, and the initial condition u₀ which define an SDE:

\[du = f(u,p,t)dt + Σgᵢ(u,p,t)dWⁱ\]

f and g should be specified as f(u,p,t) and g(u,p,t) respectively, and u₀ should be an AbstractArray 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. A vector of gs can also be defined to determine an SDE of higher Ito dimension.

Problem Type

Wraps the data which defines an SDE problem

\[u = f(u,p,t)dt + Σgᵢ(u,p,t)dWⁱ\]

with initial condition u0.

Constructors

  • SDEProblem(f::SDEFunction,u0,tspan,p=NullParameters();noise=WHITE_NOISE,noise_rate_prototype=nothing)
  • SDEProblem{isinplace,specialize}(f,g,u0,tspan,p=NullParameters();noise=WHITE_NOISE,noise_rate_prototype=nothing) : Defines the SDE with the specified functions. The default noise is WHITE_NOISE. isinplace optionally sets whether the function is inplace or not. This is determined automatically, but not inferred. specialize optionally controls the specialization level. See the specialization levels section of the SciMLBase documentation for more details. The default is `AutoSpecialize.

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.

For specifying Jacobians and mass matrices, see the DiffEqFunctions page.

Fields

  • f: The drift function in the SDE.
  • g: The noise function in the SDE.
  • u0: The initial condition.
  • tspan: The timespan for the problem.
  • p: The optional parameters for the problem. Defaults to NullParameters.
  • noise: The noise process applied to the noise upon generation. Defaults to Gaussian white noise. For information on defining different noise processes, see the noise process documentation page.
  • noise_rate_prototype: A prototype type instance for the noise rates, that is the output g. It can be any type which overloads A_mul_B! with itself being the middle argument. Commonly, this is a matrix or sparse matrix. If this is not given, it defaults to nothing, which means the problem should be interpreted as having diagonal noise.
  • kwargs: The keyword arguments passed onto the solves.

Example Problems

Examples problems can be found in DiffEqProblemLibrary.jl.

To use a sample problem, such as prob_sde_linear, you can do something like:

#] add SDEProblemLibrary
using SDEProblemLibrary
prob = SDEProblemLibrary.prob_sde_linear
sol = solve(prob)
source
SciMLBase.SDEFunctionType

DocStringExtensions.TypeDefinition()

A representation of an SDE function f, defined by:

\[M du = f(u,p,t)dt + g(u,p,t) dW\]

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.

Constructor

SDEFunction{iip,specialize}(f,g;
                           mass_matrix = __has_mass_matrix(f) ? f.mass_matrix : I,
                           analytic = __has_analytic(f) ? f.analytic : nothing,
                           tgrad= __has_tgrad(f) ? f.tgrad : nothing,
                           jac = __has_jac(f) ? f.jac : nothing,
                           jvp = __has_jvp(f) ? f.jvp : nothing,
                           vjp = __has_vjp(f) ? f.vjp : nothing,
                           ggprime = nothing,
                           jac_prototype = __has_jac_prototype(f) ? f.jac_prototype : nothing,
                           sparsity = __has_sparsity(f) ? f.sparsity : jac_prototype,
                           paramjac = __has_paramjac(f) ? f.paramjac : nothing,
                           colorvec = __has_colorvec(f) ? f.colorvec : nothing,
                           sys = __has_sys(f) ? f.sys : nothing)

Note that both the function f and g are required. This function should be given as f!(du,u,p,t) or du = f(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 f. These include:

  • mass_matrix: the mass matrix M represented in the ODE function. Can be used to determine that the equation is actually a differential-algebraic equation (DAE) if M 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(u,p,t)}{\partial t}$
  • jac(J,u,p,t) or J=jac(u,p,t): returns $\frac{df}{du}$
  • jvp(Jv,v,u,p,t) or Jv=jvp(v,u,p,t): returns the directional derivative $\frac{df}{du} v$
  • vjp(Jv,v,u,p,t) or Jv=vjp(v,u,p,t): returns the adjoint derivative $\frac{df}{du}^\ast v$
  • ggprime(J,u,p,t) or J = ggprime(u,p,t): returns the Milstein derivative $\frac{dg(u,p,t)}{du} g(u,p,t)$
  • jac_prototype: a prototype matrix matching the type that matches the Jacobian. For example, if the Jacobian is tridiagonal, then an appropriately sized Tridiagonal matrix can be used as the prototype and integrators will specialize on this structure where possible. Non-structured sparsity patterns should use a SparseMatrixCSC with a correct sparsity pattern for the Jacobian. The default is nothing, which means a dense Jacobian.
  • paramjac(pJ,u,p,t): returns the parameter Jacobian $\frac{df}{dp}$.
  • colorvec: a color vector according to the SparseDiffTools.jl definition for the sparsity pattern of the jac_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 to nothing, 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.

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 ODEFunction type directly match the names of the inputs.

source

Solution Type

SDEProblem solutions return an RODESolution. For more information, see the RODE problem definition page for the RODESolution docstring.

Alias Specifier

SciMLBase.SDEAliasSpecifierType

Holds information on what variables to alias when solving an SDEProblem. Conforms to the AbstractAliasSpecifier interface. SDEAliasSpecifier(;alias_p = nothing, alias_f = nothing, alias_u0 = nothing, alias_tstops = nothing, alias = nothing)

When a keyword argument is nothing, the default behaviour of the solver is used.

Keywords

  • alias_p::Union{Bool, Nothing}
  • alias_f::Union{Bool, Nothing}
  • alias_u0::Union{Bool, Nothing}: alias the u0 array. Defaults to false .
  • alias_tstops::Union{Bool, Nothing}: alias the tstops array
  • alias_jumps::Union{Bool, Nothing}: alias jump process if wrapped in a JumpProcess
  • alias::Union{Bool, Nothing}: sets all fields of the SDEAliasSpecifier to alias
source

Example Problems

Examples problems can be found in DiffEqProblemLibrary.jl.

To use a sample problem, such as prob_sde_linear, you can do something like:

#] add DiffEqProblemLibrary
using DiffEqProblemLibrary.SDEProblemLibrary
# load problems
SDEProblemLibrary.importsdeproblems()
prob = SDEProblemLibrary.prob_sde_linear
sol = solve(prob)
SDEProblemLibrary.prob_sde_linearConstant

\[du_t = αudt + βudW_t\]

where $α=1.01$, $β=0.87$, and initial condition $u_0=1/2$, with solution

\[u(u_0,p,t,W_t)=u_0\exp((α-\frac{β^2}{2})t+βW_t)\]

source
SDEProblemLibrary.prob_sde_2DlinearConstant

8 linear SDEs (as a 4x2 matrix):

\[du_t = αudt + βudW_t\]

where $α=1.01$, $β=0.87$, and initial condition $u_0=\frac{1}{2}$ with solution

\[u(u_0,p,t,W_t)=u_0\exp((α-\frac{β^2}{2})t+βW_t)\]

source
SDEProblemLibrary.prob_sde_waveConstant

\[du_t = -\frac{1}{100}\sin(u)\cos^3(u)dt + \frac{1}{10}\cos^{2}(u_t) dW_t\]

and initial condition $u_0=1$ with solution

\[u(u_0,p,t,W_t)=\arctan(\frac{W_t}{10} + \tan(u_0))\]

source
SDEProblemLibrary.prob_sde_lorenzConstant

Lorenz Attractor with additive noise

\[dx = σ(y-x)dt + αdW_t\]

\[dy = (x(ρ-z) - y)dt + αdW_t\]

\[dz = (xy - βz)dt + αdW_t\]

with $σ=10$, $ρ=28$, $β=8/3$, $α=3.0$ and initial condition $u_0=[1;1;1]$.

source
SDEProblemLibrary.prob_sde_cubicConstant

\[du_t = \frac{1}{4}u(1-u^2)dt + \frac{1}{2}(1-u^2)dW_t\]

and initial condition $u_0=\frac{1}{2}$, with solution

\[u(u0,p,t,W_t)=\frac{(1+u_0)\exp(W_t)+u)0-1}{(1+u_0)\exp(W_t)+1-u_0}\]

source
SDEProblemLibrary.prob_sde_additiveConstant

Additive noise problem

\[u_t = (\frac{β}{\sqrt{1+t}}-\frac{1}{2(1+t)}u_t)dt + \frac{αβ}{\sqrt{1+t}}dW_t\]

and initial condition $u_0=1$ with $α=0.1$ and $β=0.05$, with solution

\[u(u_0,p,t,W_t)=\frac{u_0}{\sqrt{1+t}} + \frac{β(t+αW_t)}{\sqrt{1+t}}\]

source
SDEProblemLibrary.oval2ModelExampleFunction

oval2ModelExample(;largeFluctuations=false,useBigs=false,noiseLevel=1)

A function which generates the Oval2 Epithelial-Mesenchymal Transition model from:

Rackauckas, C., & Nie, Q. (2017). Adaptive methods for stochastic differential equations via natural embeddings and rejection sampling with memory. Discrete and continuous dynamical systems. Series B, 22(7), 2731.

19 SDEs which are only stiff during transitions between biological states.

source
SDEProblemLibrary.prob_sde_stiffquadstratConstant

The composite Euler method for stiff stochastic differential equations

Kevin Burrage, Tianhai Tian

And

S-ROCK: CHEBYSHEV METHODS FOR STIFF STOCHASTIC DIFFERENTIAL EQUATIONS

ASSYR ABDULLE AND STEPHANE CIRILLI

Stiffness of Euler is determined by α+β²<1 Higher α or β is stiff, with α being deterministic stiffness and β being noise stiffness (and grows by square).

source
SDEProblemLibrary.prob_sde_stiffquaditoConstant

The composite Euler method for stiff stochastic differential equations

Kevin Burrage, Tianhai Tian

And

S-ROCK: CHEBYSHEV METHODS FOR STIFF STOCHASTIC DIFFERENTIAL EQUATIONS

ASSYR ABDULLE AND STEPHANE CIRILLI

Stiffness of Euler is determined by α+β²<1 Higher α or β is stiff, with α being deterministic stiffness and β being noise stiffness (and grows by square).

source
SDEProblemLibrary.generate_stiff_stoch_heatFunction

Stochastic Heat Equation with scalar multiplicative noise

S-ROCK: CHEBYSHEV METHODS FOR STIFF STOCHASTIC DIFFERENTIAL EQUATIONS

ASSYR ABDULLE AND STEPHANE CIRILLI

Raising D or k increases stiffness

source