Noise Process Interface

This page documents the interface functions for working with noise processes.

Step Management

DiffEqNoiseProcess.accept_step!Function
accept_step!(W::NoiseProcess, dt, u, p, setup_next = true)

Accept a noise step, updating the current time and noise values. This function is called by the SDE solvers after a successful integration step.

Arguments

  • W: The noise process
  • dt: Time step size (may differ from W.dt due to adaptivity)
  • u: Current solution value (for state-dependent noise)
  • p: Parameters (for state-dependent noise)
  • setup_next: Whether to prepare for the next step
source
DiffEqNoiseProcess.reject_step!Function
reject_step!(W::NoiseProcess, dtnew, u, p)

Handle rejection of a time step in adaptive algorithms.

When an adaptive SDE solver rejects a step, this function uses bridge interpolation to generate appropriate noise values for the smaller time step, maintaining distributional correctness.

Arguments

  • W: The noise process
  • dtnew: New (smaller) time step after rejection
  • u: Current solution value (for state-dependent noise)
  • p: Parameters (for state-dependent noise)

Details

The function stores unused portions of noise in the RSWM stacks for later use, ensuring the noise process remains distributionally exact despite adaptivity.

source
DiffEqNoiseProcess.calculate_step!Function
calculate_step!(W::NoiseProcess, dt, u, p)

Calculate noise increments for the given time step.

This function generates new random values for the noise process according to its distribution function.

Arguments

  • W: The noise process
  • dt: Time step size
  • u: Current solution value (for state-dependent noise)
  • p: Parameters (for state-dependent noise)

Effects

Updates W.dW (and W.dZ if auxiliary process exists) with new noise increments

source
DiffEqNoiseProcess.setup_next_step!Function
setup_next_step!(W::NoiseProcess, u, p)

Prepare the noise process for the next integration step.

This function manages the Rejection Sampling with Memory (RSWM) algorithm, handling the stacks of pre-computed noise values for adaptive time stepping.

Arguments

  • W: The noise process
  • u: Current solution value (for state-dependent noise)
  • p: Parameters (for state-dependent noise)

Details

The function implements different variants of the RSWM algorithm:

  • RSwM1: Basic rejection sampling
  • RSwM2: Improved version with better memory management
  • RSwM3: Most advanced version with two-stack system
source
DiffEqNoiseProcess.save_noise!Function
save_noise!(W::NoiseProcess)

Save the current noise value and time to the process history if it hasn't been saved already. This function is called automatically by the stepping algorithms to maintain the noise trajectory.

source

Noise Process Types

Abstract Types

Missing docstring.

Missing docstring for AbstractNoiseProcess. Check Documenter's build log for details.

Core Types

Missing docstring.

Missing docstring for NoiseProcess. Check Documenter's build log for details.

Missing docstring.

Missing docstring for SimpleNoiseProcess. Check Documenter's build log for details.

Configuration

Rejection Sampling with Memory (RSWM)

DiffEqNoiseProcess.RSWMType
RSWM(; discard_length = 1e-15, adaptivealg = :RSwM3)

Rejection Sampling with Memory (RSWM) algorithm configuration for noise processes.

RSWM ensures distributional exactness when adaptive time stepping is used with noise processes. It maintains memory of rejected values to avoid biasing the noise distribution.

Fields

  • discard_length: Threshold for discarding stored values to save memory. Smaller values use more memory but are more accurate.
  • adaptivealg: The adaptive algorithm variant to use (:RSwM3 is recommended)

Examples

# Conservative (high accuracy)
rswm_accurate = RSWM(discard_length = 1e-12)

# Aggressive (lower memory usage)  
rswm_fast = RSWM(discard_length = 1e-6)

# Use in noise process
W = WienerProcess(0.0, 0.0, 1.0; rswm = rswm_accurate)
source

Internal Functions

These functions are used internally by the noise process implementations:

Distribution Functions

DiffEqNoiseProcess.WHITE_NOISE_DISTFunction
WHITE_NOISE_DIST(dW, W, dt, u, p, t, rng)

Generate white noise distributed according to N(0, dt) for use in Wiener processes.

Arguments

  • dW: Noise increment container
  • W: Current noise value (unused for white noise)
  • dt: Time step
  • u: Current state (for state-dependent noise)
  • p: Parameters
  • t: Current time
  • rng: Random number generator

Returns

Random values distributed as N(0, dt), scaled by √dt

source
DiffEqNoiseProcess.WHITE_NOISE_BRIDGEFunction
WHITE_NOISE_BRIDGE(dW, W, W0, Wh, q, h, u, p, t, rng)

Generate white noise for Brownian bridge interpolation between two points.

Arguments

  • dW: Noise increment container
  • W: Current noise value
  • W0: Starting noise value
  • Wh: Target noise value at end of interval
  • q: Interpolation parameter (0 to 1)
  • h: Total time interval
  • u, p, t: State, parameters, and time (for compatibility)
  • rng: Random number generator

Returns

Interpolated noise value that maintains correct distribution

source
DiffEqNoiseProcess.VBT_BRIDGEFunction
VBT_BRIDGE(dW, W, W0, Wh, q, h, u, p, t, rng)

Generate noise for Virtual Brownian Tree (VBT) bridge interpolation.

The VBT bridge is a memory-efficient method for generating Brownian paths that can be evaluated at arbitrary time points without storing the entire path.

Arguments

  • dW: Noise increment container
  • W: Current noise value
  • W0: Starting noise value
  • Wh: Target noise value at end of interval
  • q: Interpolation parameter (0 to 1)
  • h: Total time interval
  • u, p, t: State, parameters, and time (for compatibility)
  • rng: Random number generator

Returns

Interpolated noise value using the VBT bridge formula

source
DiffEqNoiseProcess.INPLACE_WHITE_NOISE_DISTFunction
INPLACE_WHITE_NOISE_DIST(rand_vec, W, dt, u, p, t, rng)

Generate white noise distributed according to N(0, dt) in-place.

This is the in-place version of WHITENOISEDIST, modifying the provided array rather than allocating a new one.

Arguments

  • rand_vec: Array to fill with noise values
  • W: Current noise value (unused for white noise)
  • dt: Time step
  • u: Current state (for state-dependent noise)
  • p: Parameters
  • t: Current time
  • rng: Random number generator

Effects

Modifies rand_vec to contain random values distributed as N(0, dt)

source
DiffEqNoiseProcess.INPLACE_WHITE_NOISE_BRIDGEFunction
INPLACE_WHITE_NOISE_BRIDGE(rand_vec, W, W0, Wh, q, h, u, p, t, rng)

Generate white noise for Brownian bridge interpolation in-place.

This is the in-place version of WHITENOISEBRIDGE, modifying the provided array rather than allocating a new one.

Arguments

  • rand_vec: Array to fill with interpolated noise values
  • W: Current noise value
  • W0: Starting noise value
  • Wh: Target noise value at end of interval
  • q: Interpolation parameter (0 to 1)
  • h: Total time interval
  • u, p, t: State, parameters, and time (for compatibility)
  • rng: Random number generator

Effects

Modifies rand_vec to contain interpolated noise values

source
DiffEqNoiseProcess.INPLACE_VBT_BRIDGEFunction
INPLACE_VBT_BRIDGE(rand_vec, W, W0, Wh, q, h, u, p, t, rng)

Generate noise for Virtual Brownian Tree (VBT) bridge interpolation in-place.

This is the in-place version of VBT_BRIDGE, modifying the provided array rather than allocating a new one.

Arguments

  • rand_vec: Array to fill with interpolated noise values
  • W: Current noise value
  • W0: Starting noise value
  • Wh: Target noise value at end of interval
  • q: Interpolation parameter (0 to 1)
  • h: Total time interval
  • u, p, t: State, parameters, and time (for compatibility)
  • rng: Random number generator

Effects

Modifies rand_vec to contain VBT-interpolated noise values

source
DiffEqNoiseProcess.REAL_WHITE_NOISE_DISTFunction
REAL_WHITE_NOISE_DIST(dW, W, dt, u, p, t, rng)

Generate real-valued white noise distributed according to N(0, dt).

Unlike WHITENOISEDIST, this function always generates real-valued noise, even if the input type would normally support complex values.

Arguments

  • dW: Noise increment container
  • W: Current noise value (unused for white noise)
  • dt: Time step
  • u: Current state (for state-dependent noise)
  • p: Parameters
  • t: Current time
  • rng: Random number generator

Returns

Real-valued random numbers distributed as N(0, dt), scaled by √dt

source
DiffEqNoiseProcess.REAL_WHITE_NOISE_BRIDGEFunction
REAL_WHITE_NOISE_BRIDGE(dW, W, W0, Wh, q, h, u, p, t, rng)

Generate real-valued white noise for Brownian bridge interpolation.

This function ensures the generated noise is always real-valued, even for complex-valued endpoints.

Arguments

  • dW: Noise increment container
  • W: Current noise value
  • W0: Starting noise value
  • Wh: Target noise value at end of interval
  • q: Interpolation parameter (0 to 1)
  • h: Total time interval
  • u, p, t: State, parameters, and time (for compatibility)
  • rng: Random number generator

Returns

Real-valued interpolated noise value

source
DiffEqNoiseProcess.REAL_INPLACE_WHITE_NOISE_DISTFunction
REAL_INPLACE_WHITE_NOISE_DIST(rand_vec, W, dt, u, p, t, rng)

Generate real-valued white noise distributed according to N(0, dt) in-place.

This is the in-place version of REALWHITENOISE_DIST.

Arguments

  • rand_vec: Array to fill with noise values
  • W: Current noise value (unused for white noise)
  • dt: Time step
  • u: Current state (for state-dependent noise)
  • p: Parameters
  • t: Current time
  • rng: Random number generator

Effects

Modifies rand_vec to contain real-valued random numbers distributed as N(0, dt)

source
DiffEqNoiseProcess.REAL_INPLACE_WHITE_NOISE_BRIDGEFunction
REAL_INPLACE_WHITE_NOISE_BRIDGE(rand_vec, W, W0, Wh, q, h, u, p, t, rng)

Generate real-valued white noise for Brownian bridge interpolation in-place.

This is the in-place version of REALWHITENOISE_BRIDGE.

Arguments

  • rand_vec: Array to fill with interpolated noise values
  • W: Current noise value
  • W0: Starting noise value
  • Wh: Target noise value at end of interval
  • q: Interpolation parameter (0 to 1)
  • h: Total time interval
  • u, p, t: State, parameters, and time (for compatibility)
  • rng: Random number generator

Effects

Modifies rand_vec to contain real-valued interpolated noise values

source

Random Number Generation

DiffEqNoiseProcess.wiener_randnFunction
wiener_randn(rng::AbstractRNG, ::Type{T}) where {T}

Generate a random number from the standard normal distribution for type T.

Arguments

  • rng: Random number generator
  • T: Type of the random number to generate

Returns

A random number of type T from the standard normal distribution

source
wiener_randn(rng::AbstractRNG, proto::AbstractArray{T}) where {T <: Number}

Generate an array of random numbers from the standard normal distribution, matching the size of the prototype array.

Arguments

  • rng: Random number generator
  • proto: Prototype array whose size determines the output size

Returns

An array of random numbers from the standard normal distribution with the same size as proto

source
wiener_randn(rng::AbstractRNG, proto::T) where {T <: StaticArraysCore.SArray}

Generate a static array of random numbers from the standard normal distribution.

Arguments

  • rng: Random number generator
  • proto: Prototype static array

Returns

A static array of the same type as proto filled with standard normal random numbers

source
wiener_randn(rng::AbstractRNG, proto)

Generate random numbers from the standard normal distribution for arbitrary types.

Arguments

  • rng: Random number generator
  • proto: Prototype object whose type and size determine the output

Returns

Random values converted to the same type as proto

source
DiffEqNoiseProcess.wiener_randn!Function
wiener_randn!(rng::AbstractRNG, rand_vec::AbstractArray)

Fill an array with random numbers from the standard normal distribution in-place.

Arguments

  • rng: Random number generator
  • rand_vec: Array to fill with random values

Returns

The modified rand_vec filled with standard normal random numbers

source
wiener_randn!(rng::AbstractRNG, rand_vec)

Fill an arbitrary container with random numbers from the standard normal distribution in-place using broadcasting.

Arguments

  • rng: Random number generator (not used in this fallback)
  • rand_vec: Container to fill with random values

Returns

The modified rand_vec filled with standard normal random numbers

source
wiener_randn!(rng::AbstractRNG, rand_vec::GPUArraysCore.AbstractGPUArray)

Fill a GPU array with random numbers from the standard normal distribution in-place.

This specialized method works for GPUs because it doesn't pass the RNG to the GPU kernel, which may not be supported on all GPU backends.

Arguments

  • rng: Random number generator (not passed to GPU)
  • rand_vec: GPU array to fill with random values

Returns

The modified rand_vec filled with standard normal random numbers

source
wiener_randn!(y::AbstractRNG, x::AbstractArray{<:Complex{T}}) where {T <: Number}

Fill an array of complex numbers with random values from the standard complex normal distribution in-place.

Each complex number is generated as (a + bi)/√2 where a and b are independent standard normal random variables.

Arguments

  • y: Random number generator
  • x: Array of complex numbers to fill

Returns

The modified array filled with complex normal random numbers

source

Ornstein-Uhlenbeck Specific

DiffEqNoiseProcess.OrnsteinUhlenbeckType
OrnsteinUhlenbeck{T1, T2, T3}

Parameters for the Ornstein-Uhlenbeck process.

Fields

  • Θ: Mean reversion rate (higher values mean faster reversion)
  • μ: Long-term mean (the value the process reverts to)
  • σ: Volatility/diffusion coefficient

The process follows the SDE: dXt = Θ(μ - Xt)dt + σ dW_t

source
DiffEqNoiseProcess.OrnsteinUhlenbeck!Type
OrnsteinUhlenbeck!{T1, T2, T3}

In-place version of OrnsteinUhlenbeck parameters.

Fields

  • Θ: Mean reversion rate (higher values mean faster reversion)
  • μ: Long-term mean (the value the process reverts to)
  • σ: Volatility/diffusion coefficient

The process follows the SDE: dXt = Θ(μ - Xt)dt + σ dW_t

source
DiffEqNoiseProcess.ou_bridgeFunction
ou_bridge(dW, ou, W, W0, Wh, q, h, u, p, t, rng)

Generate Ornstein-Uhlenbeck bridge interpolation between two points.

Provides exact sampling from an OU process conditioned on both endpoints, useful for adaptive time-stepping and interpolation.

Arguments

  • dW: Noise increment container
  • ou: OrnsteinUhlenbeck parameters
  • W: Current noise process state
  • W0: Starting value
  • Wh: Target value at end of interval
  • q: Interpolation parameter (0 to 1)
  • h: Total time interval
  • u, p, t: State, parameters, and time (for compatibility)
  • rng: Random number generator

Returns

Interpolated OU process value that maintains the correct distribution

References

  • http://www.tandfonline.com/doi/pdf/10.1080/14697688.2014.941913
  • https://arxiv.org/pdf/1011.0067.pdf (page 18)
source
DiffEqNoiseProcess.ou_bridge!Function
ou_bridge!(rand_vec, ou, W, W0, Wh, q, h, u, p, t, rng)

Generate Ornstein-Uhlenbeck bridge interpolation in-place.

This is the in-place version of ou_bridge, modifying the provided array rather than allocating a new one.

Arguments

  • rand_vec: Array to fill with interpolated values
  • ou: OrnsteinUhlenbeck parameters
  • W: Current noise process state
  • W0: Starting value
  • Wh: Target value at end of interval
  • q: Interpolation parameter (0 to 1)
  • h: Total time interval
  • u, p, t: State, parameters, and time (for compatibility)
  • rng: Random number generator

Effects

Modifies rand_vec to contain the interpolated OU process values

source

Geometric Brownian Motion Specific

DiffEqNoiseProcess.GeometricBrownianMotionType
GeometricBrownianMotion{T1, T2}

Parameters for the geometric Brownian motion process.

Fields

  • μ: Drift parameter (expected return rate)
  • σ: Volatility parameter (standard deviation of returns)

The process follows the SDE: dXt = μXt dt + σXt dWt

This is commonly used in financial models, particularly the Black-Scholes model.

source
DiffEqNoiseProcess.GeometricBrownianMotion!Type
GeometricBrownianMotion!{T1, T2}

In-place version of GeometricBrownianMotion parameters.

Fields

  • μ: Drift parameter (expected return rate)
  • σ: Volatility parameter (standard deviation of returns)

The process follows the SDE: dXt = μXt dt + σXt dWt

source
DiffEqNoiseProcess.gbm_bridgeFunction
gbm_bridge(dW, gbm, W, W0, Wh, q, h, u, p, t, rng)

Generate geometric Brownian motion bridge interpolation between two points.

Provides exact sampling from a GBM process conditioned on both endpoints, useful for adaptive time-stepping and interpolation.

Arguments

  • dW: Noise increment container
  • gbm: GeometricBrownianMotion parameters
  • W: Current noise process state
  • W0: Starting value
  • Wh: Target value at end of interval
  • q: Interpolation parameter (0 to 1)
  • h: Total time interval
  • u, p, t: State, parameters, and time (for compatibility)
  • rng: Random number generator

Returns

Interpolated GBM process value that maintains the correct log-normal distribution

Reference

https://math.stackexchange.com/questions/412470/conditional-distribution-in-brownian-motion

source
DiffEqNoiseProcess.gbm_bridge!Function
gbm_bridge!(rand_vec, gbm, W, W0, Wh, q, h, u, p, t, rng)

Generate geometric Brownian motion bridge interpolation in-place.

This is the in-place version of gbm_bridge, modifying the provided array rather than allocating a new one.

Arguments

  • rand_vec: Array to fill with interpolated values
  • gbm: GeometricBrownianMotion parameters
  • W: Current noise process state
  • W0: Starting value
  • Wh: Target value at end of interval
  • q: Interpolation parameter (0 to 1)
  • h: Total time interval
  • u, p, t: State, parameters, and time (for compatibility)
  • rng: Random number generator

Effects

Modifies rand_vec to contain the interpolated GBM process values

source

Compound Poisson Specific

DiffEqNoiseProcess.CompoundPoissonProcessType
CompoundPoissonProcess{R, CR}

A compound Poisson process for modeling jump processes.

The process has jumps that occur according to a Poisson process with given rate, and jump sizes determined by a specified distribution.

Fields

  • rate: Jump rate function or constant (λ parameter)
  • currate: Current rate value (cached for efficiency)
  • computerates: Whether to recompute rates at each step

Constructor

CompoundPoissonProcess(rate, t0, W0; computerates = true, kwargs...)

Examples

# Constant rate
proc = CompoundPoissonProcess(2.0, 0.0, 0.0)

# Time-dependent rate
rate_func(u, p, t) = 1.0 + 0.5*sin(t)
proc = CompoundPoissonProcess(rate_func, 0.0, 0.0)

References

https://www.math.wisc.edu/~anderson/papers/AndPostleap.pdf Incorporating postleap checks in tau-leaping J. Chem. Phys. 128, 054103 (2008); https://doi.org/10.1063/1.2819665

source
DiffEqNoiseProcess.CompoundPoissonProcess!Type
CompoundPoissonProcess!{R, CR}

In-place version of CompoundPoissonProcess.

See CompoundPoissonProcess for details.

Constructor

CompoundPoissonProcess!(rate, t0, W0; computerates = true, kwargs...)
source
DiffEqNoiseProcess.cpp_bridgeFunction
cpp_bridge(dW, cpp, W, W0, Wh, q, h, u, p, t, rng)

Generate compound Poisson process bridge interpolation between two points.

Uses binomial thinning to distribute jumps appropriately between endpoints.

Arguments

  • dW: Noise increment container
  • cpp: CompoundPoissonProcess parameters
  • W: Current noise process state
  • W0: Starting value
  • Wh: Jump count difference (must be integer)
  • q: Interpolation parameter (0 to 1)
  • h: Total time interval
  • u, p, t: State, parameters, and time (for compatibility)
  • rng: Random number generator

Returns

Number of jumps distributed according to binomial thinning

source
DiffEqNoiseProcess.cpp_bridge!Function
cpp_bridge!(rand_vec, cpp, W, W0, Wh, q, h, u, p, t, rng)

In-place version of cpp_bridge.

Effects

Modifies rand_vec to contain binomially distributed jump counts

source