Basic Usage Tutorial
This tutorial covers the fundamental usage patterns of DiffEqNoiseProcess.jl, including creating noise processes, using them in problems, and direct simulation.
Creating and Using Noise Processes
Basic Wiener Process
A Wiener process (Brownian motion) is the most fundamental noise process:
using DiffEqNoiseProcess, SciMLBase
# Create a standard Wiener process
# Parameters: t0 (start time), W0 (initial value), Z0 (end value)
W = WienerProcess(0.0, 0.0, 1.0)
# You can specify additional options
W_custom = WienerProcess(0.0, 0.0, 1.0;
save_everystep = true,
rswm = RSWM())t: 1-element Vector{Float64}:
0.0
u: 1-element Vector{Float64}:
0.0Real-Valued Wiener Process
For scalar problems, use a real-valued Wiener process:
# Real-valued Wiener process
W_real = RealWienerProcess(0.0, 0.0, 1.0)t: 1-element Vector{Float64}:
0.0
u: 1-element Vector{Float64}:
0.0Geometric Brownian Motion
A geometric Brownian motion process is useful for financial modeling:
# Parameters: μ (drift), σ (volatility), t0, W0, Z0
μ = 0.05 # 5% drift
σ = 0.2 # 20% volatility
GBM = GeometricBrownianMotionProcess(μ, σ, 0.0, 1.0, 1.0)t: 1-element Vector{Float64}:
0.0
u: 1-element Vector{Float64}:
1.0Using Noise in SDE Problems
Noise processes can be passed directly to SDE problems:
# Define a simple SDE: du = u*dt + 0.1*u*dW
function f(u, p, t)
return u
end
function g(u, p, t)
return 0.1 * u
end
u0 = 1.0
tspan = (0.0, 1.0)
# Create the problem with custom noise
prob = SDEProblem(f, g, u0, tspan, noise = W)SDEProblem with uType Float64 and tType Float64. In-place: false
Non-trivial mass matrix: false
timespan: (0.0, 1.0)
u0: 1.0Direct Simulation
You can simulate noise processes directly without solving an SDE:
# Create a noise problem
prob = NoiseProblem(W, (0.0, 1.0))
# Solve it with a specified timestep
sol = solve(prob; dt = 0.01)
# Access the noise values
println("Final noise value: ", sol.u[end])
println("Number of steps: ", length(sol.t))Final noise value: -1.5102604952421415
Number of steps: 101Ensemble Simulations
Generate multiple realizations of the same noise process:
# Create ensemble problem
enprob = EnsembleProblem(prob)
# Generate 100 trajectories
ensol = solve(enprob; dt = 0.01, trajectories = 100)
println("Generated $(length(ensol.u)) noise trajectories")Generated 100 noise trajectoriesDirect Interface for Advanced Usage
For advanced users, noise processes support a direct stepping interface:
# Create a fresh Wiener process
W_direct = WienerProcess(0.0, 0.0, 1.0)
# Set the timestep
dt = 0.1
W_direct.dt = dt
# For state-dependent noise (not applicable here, but required by interface)
u = nothing
p = nothing
# Calculate and accept steps
calculate_step!(W_direct, dt, u, p)
for i in 1:5
accept_step!(W_direct, dt, u, p)
println("Step $i: W = $(W_direct.curW)")
endStep 1: W = -0.10671619714290562
Step 2: W = 0.14156516417076823
Step 3: W = -0.03729584725419072
Step 4: W = -0.13657936130589254
Step 5: W = -0.24017182215070354This direct interface is primarily used internally by the SDE solvers but can be useful for custom implementations.