CoupledSDEs

DynamicalSystemsBase.CoupledSDEsType
CoupledSDEs <: ContinuousTimeDynamicalSystem
CoupledSDEs(f, u0 [, p]; kwargs...)

A stochastic continuous time dynamical system defined by a set of coupled stochastic differential equations (SDEs) as follows:

\[\text{d}\mathbf{u} = \mathbf{f}(\mathbf{u}, p, t) \text{d}t + \mathbf{g}(\mathbf{u}, p, t) \text{d}\mathcal{N}_t\]

where $\mathbf{u}(t)$ is the state vector at time $t$, $\mathbf{f}$ describes the deterministic dynamics, and the noise term $\mathbf{g}(\mathbf{u}, p, t) \text{d}\mathcal{N}_t$ describes the stochastic forcing in terms of a noise function (or diffusion function) $\mathbf{g}$ and a noise process $\mathcal{N}_t$. The parameters of the functions $\mathcal{f}$ and $\mathcal{g}$ are contained in the vector $p$.

There are multiple ways to construct a CoupledSDEs depending on the type of stochastic forcing.

The only required positional arguments are the deterministic dynamic rule f(u, p, t), the initial state u0, and optinally the parameter container p (by default p = nothing). For construction instructions regarding f, u0 see the DynamicalSystems.jl tutorial .

By default, the noise term is standard Brownian motion, i.e. additive Gaussian white noise with identity covariance matrix. To construct different noise structures, see below.

Noise term

The noise term can be specified via several keyword arguments. Based on these keyword arguments, the noise function g is constructed behind the scenes unless explicitly given.

  • The noise strength (i.e. the magnitude of the stochastic forcing) can be scaled with noise_strength (defaults to 1.0). This factor is multiplied with the whole noise term.
  • For non-diagonal and correlated noise, a covariance matrix can be provided via covariance (defaults to identity matrix of size length(u0).)
  • For more complicated noise structures, including state- and time-dependent noise, the noise function g can be provided explicitly as a keyword argument (defaults to nothing). For construction instructions, continue reading.

The function g interfaces to the diffusion function specified in an SDEProblem of DynamicalSystems.jl. g must follow the same syntax as f, i.e., g(u, p, t) for out-of-place (oop) and g!(du, u, p, t) for in-place (iip).

Unless g is of vector form and describes diagonal noise, a prototype type instance for the output of g must be specified via the keyword argument noise_prototype. It can be of any type A that has the method LinearAlgebra.mul!(Y, A, B) -> Y defined. Commonly, this is a matrix or sparse matrix. If this is not given, it defaults to nothing, which means the g should be interpreted as being diagonal.

The noise process can be specified via noise_process. It defaults to a standard Wiener process (Gaussian white noise). For details on defining noise processes, see the docs of DiffEqNoiseProcess.jl . A complete list of the pre-defined processes can be found here. Note that DiffEqNoiseProcess.jl also has an interface for defining custom noise processes.

By combining g and noise_process, you can define different types of stochastic systems. Examples of different types of stochastic systems are listed on the StochasticDiffEq.jl tutorial page. A quick overview of common types of stochastic systems can also be found in the online docs for CoupledSDEs.

Keyword arguments

  • g: noise function (default nothing)
  • noise_strength: scaling factor for noise strength (default 1.0)
  • covariance: noise covariance matrix (default nothing)
  • noise_prototype: prototype instance for the output of g (default nothing)
  • noise_process: stochastic process as provided by DiffEqNoiseProcess.jl (default nothing, i.e. standard Wiener process)
  • t0: initial time (default 0.0)
  • diffeq: DiffEq solver settings (see below)
  • seed: random number seed (default UInt64(0))

DifferentialEquations.jl interfacing

The CoupledSDEs is evolved using the solvers of DifferentialEquations.jl. To specify a solver via the diffeq keyword argument, use the flag alg, which can be accessed after loading StochasticDiffEq.jl (using StochasticDiffEq). The default diffeq is:

(alg = SOSRA(), abstol = 1.0e-2, reltol = 1.0e-2)

diffeq keywords can also include a callback for event handling .

Dev note: CoupledSDEs is a light wrapper of SDEIntegrator from StochasticDiffEq.jl. The integrator is available as the field integ, and the SDEProblem is integ.sol.prob. The convenience syntax SDEProblem(ds::CoupledSDEs, tspan = (t0, Inf)) is available to extract the problem.

Converting between CoupledSDEs and CoupledODEs

You can convert a CoupledSDEs system to CoupledODEs to analyze its deterministic part using the function CoupledODEs(ds::CoupledSDEs; diffeq, t0). Similarly, use CoupledSDEs(ds::CoupledODEs, p; kwargs...) to convert a CoupledODEs into a CoupledSDEs.

source

Examples defining stochastic dynamics

Let's look at some examples of the different types of stochastic systems that can be defined.

For simplicity, we choose a slow exponential growth in 2 dimensions as the deterministic dynamics f:

using DynamicalSystemsBase, StochasticDiffEq, DiffEqNoiseProcess
using CairoMakie
f!(du, u, p, t) = du .= 1.01u # deterministic part

function plot_trajectory(Y, t)
    fig = Figure()
    ax = Axis(fig[1,1]; xlabel = "time", ylabel = "variable")
    for var in columns(Y)
        lines!(ax, t, var)
    end
    fig
end;
plot_trajectory (generic function with 1 method)

Additive noise

When $g(u, p, t)$ is independent of the state $u$, the noise is called additive; otherwise, it is multiplicative. We can define a simple additive noise system as follows:

sde = CoupledSDEs(f!, zeros(2));
2-dimensional CoupledSDEs
 deterministic: false
 discrete time: false
 in-place:      true
 dynamic rule:  f!
 SDE solver:    SOSRA
 SDE kwargs:    (abstol = 0.01, reltol = 0.01, dt = 0.1)
 Noise type:    (additive = true, autonomous = true, linear = true, invertible = true)
 parameters:    SciMLBase.NullParameters()
 time:          0.0
 state:         [0.0, 0.0]

which is equivalent to

t0 = 0.0; W0 = zeros(2);
W = WienerProcess(t0, W0, 0.0)
sde = CoupledSDEs(f!, zeros(2);
    noise_process=W, covariance=[1 0; 0 1], noise_strength=1.0
    );
2-dimensional CoupledSDEs
 deterministic: false
 discrete time: false
 in-place:      true
 dynamic rule:  f!
 SDE solver:    SOSRA
 SDE kwargs:    (abstol = 0.01, reltol = 0.01, dt = 0.1)
 Noise type:    (additive = true, autonomous = true, linear = true, invertible = true)
 parameters:    SciMLBase.NullParameters()
 time:          0.0
 state:         [0.0, 0.0]

We defined a Wiener process W, whose increments are vectors of normally distributed random numbers of length matching the output of g. The noise is applied element-wise, i.e., g.*dW. Since the noise processes are uncorrelated, meaning the covariance matrix is diagonal, this type of noise is referred to as diagonal.

We can sample a trajectory from this system using the trajectory function also used for the deterministic systems:

tr = trajectory(sde, 1.0)
plot_trajectory(tr...)
Example block output

Correlated noise

In the case of correlated noise, the random numbers in a vector increment dW are correlated. This can be achieved by specifying the covariance matrix $\Sigma$ via the covariance keyword:

ρ = 0.3
Σ = [1 ρ; ρ 1]
diffeq = (alg = LambaEM(), dt=0.1)
sde = CoupledSDEs(f!, zeros(2); covariance=Σ, diffeq=diffeq)
2-dimensional CoupledSDEs
 deterministic: false
 discrete time: false
 in-place:      true
 dynamic rule:  f!
 SDE solver:    LambaEM
 SDE kwargs:    (dt = 0.1,)
 Noise type:    (additive = true, autonomous = true, linear = true, invertible = true)
 parameters:    SciMLBase.NullParameters()
 time:          0.0
 state:         [0.0, 0.0]

Alternatively, we can parametrise the covariance matrix by defining the diffusion function $g$ ourselves:

g!(du, u, p, t) = (du .= [1 p[1]; p[1] 1]; return nothing)
sde = CoupledSDEs(f!, zeros(2), (ρ); g=g!, noise_prototype=zeros(2, 2))
2-dimensional CoupledSDEs
 deterministic: false
 discrete time: false
 in-place:      true
 dynamic rule:  f!
 SDE solver:    SOSRA
 SDE kwargs:    (abstol = 0.01, reltol = 0.01, dt = 0.1)
 Noise type:    (additive = true, autonomous = true, linear = true, invertible = true)
 parameters:    0.3
 time:          0.0
 state:         [0.0, 0.0]

Here, we had to provide noise_prototype to indicate that the diffusion function g will output a 2x2 matrix.

Scalar noise

If all state variables are forced by the same single random variable, we have scalar noise. To define scalar noise, one has to give an one-dimensional noise process to the noise_process keyword of the CoupledSDEs constructor.

t0 = 0.0; W0 = 0.0;
noise = WienerProcess(t0, W0, 0.0)
sde = CoupledSDEs(f!, rand(2)/10; noise_process=noise)

tr = trajectory(sde, 1.0)
plot_trajectory(tr...)
Example block output

We can see that noise applied to each variable is the same.

Multiplicative and time-dependent noise

In the SciML ecosystem, multiplicative noise is defined through the condition $g_i(t, u)=a_i u$. However, in the literature the name is more broadly used for any situation where the noise is non-additive and depends on the state $u$, possibly also in a non-linear way. When defining a CoupledSDEs, we can make the noise term time- and state-dependent by specifying an explicit time- or state-dependence in the noise function g, just like we would define f. For example, we can define a system with temporally decreasing multiplicative noise as follows:

function g!(du, u, p, t)
    du .= u ./ (1+t)
    return nothing
end
sde = CoupledSDEs(f!, rand(2)./10; g=g!)
2-dimensional CoupledSDEs
 deterministic: false
 discrete time: false
 in-place:      true
 dynamic rule:  f!
 SDE solver:    SOSRA
 SDE kwargs:    (abstol = 0.01, reltol = 0.01, dt = 0.1)
 Noise type:    (additive = false, autonomous = false, linear = true, invertible = false)
 parameters:    SciMLBase.NullParameters()
 time:          0.0
 state:         [0.07752822361739639, 0.0597275434811622]

Non-diagonal noise

Non-diagonal noise allows for the terms to be linearly mixed (correlated) via g being a matrix. Suppose we have two Wiener processes and two state variables such that the output of g is a 2x2 matrix. Therefore, we have

\[du_1 = f_1(u,p,t)dt + g_{11}(u,p,t)dW_1 + g_{12}(u,p,t)dW_2 \\ du_2 = f_2(u,p,t)dt + g_{21}(u,p,t)dW_1 + g_{22}(u,p,t)dW_2\]

To indicate the structure that g should have, we must use the noise_prototype keyword. Let us define a special type of non-diagonal noise called commutative noise. For this we can utilize the RKMilCommute algorithm which is designed to utilize the structure of commutative noise.

σ = 0.25 # noise strength
function g!(du, u, p, t)
  du[1,1] = σ*u[1]
  du[2,1] = σ*u[2]
  du[1,2] = σ*u[1]
  du[2,2] = σ*u[2]
    return nothing
end
diffeq = (alg = RKMilCommute(), reltol = 1e-3, abstol = 1e-3, dt=0.1)
sde = CoupledSDEs(f!, rand(2)./10; g=g!, noise_prototype = zeros(2, 2), diffeq = diffeq)
2-dimensional CoupledSDEs
 deterministic: false
 discrete time: false
 in-place:      true
 dynamic rule:  f!
 SDE solver:    RKMilCommute
 SDE kwargs:    (reltol = 0.001, abstol = 0.001, dt = 0.1)
 Noise type:    (additive = false, autonomous = true, linear = true, invertible = false)
 parameters:    SciMLBase.NullParameters()
 time:          0.0
 state:         [0.03806281791672672, 0.09681978040051963]
Warning

Non-diagonal problems need specific solvers. See the SciML recommendations.