Stochastic calculus and univariate and multivariate stochastic processes/Markov processes in continuous time. See ./example/tutorial.jl for an introduction. I am personally interested in simulating diffusion bridges and doing Bayesian inference on discretely observed diffusion processes, but this package is written to be of general use and contributions are welcome.
- Define and simulate diffusion processes in one or more dimension
- Continuous and discrete likelihood using Girsanovs theorem and transition densities
- Monte Carlo sample diffusion bridges, diffusion processes conditioned to hit a point v at a prescribed time T
- Brownian motion in one and more dimensions
- Ornstein-Uhlenbeck processes
- Geometric Brownian motion
- Bessel processes
- Gamma processes
- Basic stochastic calculus functionality (Ito integral, quadratic variation)
- Euler-Scheme and implicit methods (Rungekutta)
The layout/api was originally written to be compatible with Simon Danisch's package FixedSizeArrays.jl. It was refactored to be compatible with StaticArrays.jl by Dan Getz.
The example programs in the example/ directory have additional dependencies: ConjugatePriors and a plotting library.
The key objects introduced are the abstract type ContinuousTimeProcess{T}
parametrised by the state space of the path, for example T == Float64
and various structs
suptyping it, for example Wiener{Float64}
for a real Brownian motion. These play roughly a similar role as types subtyping Distribution
in the Distributions.jl package.
Secondly, the struct
struct SamplePath{T}
tt::Vector{Float64}
yy::Vector{T}
SamplePath{T}(tt, yy) where {T} = new(tt, yy)
end
serves as container for sample path returned by direct and approximate samplers (sample
, euler
, ...).
tt
is the vector of the grid points of the simulation and yy
the corresponding vector of states.
Help is available at the REPL:
help?> euler
search: euler euler! eulergamma default_worker_pool schedule @schedule
euler(u, W, P) -> X
Solve stochastic differential equation ``dX_t = b(t, X_t)dt + σ(t, X_t)dW_t, X_0 = u``
using the Euler scheme.
Pre-defined processes defined are
Wiener
, WienerBridge
, Gamma
, LinPro
(linear diffusion/generalized Ornstein-Uhlenbeck) and others.
It is also quite transparent how to add a new process:
# Define a diffusion process
struct OrnsteinUhlenbeck <: ContinuousTimeProcess{Float64}
β::Float64 # drift parameter (also known as inverse relaxation time)
σ::Float64 # diffusion parameter
function OrnsteinUhlenbeck(β::Float64, σ::Float64)
isnan(β) || β > 0. || error("Parameter λ must be positive.")
isnan(σ) || σ > 0. || error("Parameter σ must be positive.")
new(β, σ)
end
end
# define drift and diffusion coefficient of OrnsteinUhlenbeck
import Bridge: b, σ, a, transitionprob
Bridge.b(t,x, P::OrnsteinUhlenbeck) = -P.β*x
Bridge.σ(t, x, P::OrnsteinUhlenbeck) = P.σ
Bridge.a(t, x, P::OrnsteinUhlenbeck) = P.σ^2
# simulate OrnsteinUhlenbeck using Euler scheme
W = sample(0:0.01:10, Wiener())
X = euler(0.1, W, OrnsteinUhlenbeck(20.0, 1.0))
See the documentation for more functionality and issue #12 (Feedback and Contribution) for coordination of the development. Bridge is free software under the MIT licence. If you use Bridge.jl in a closed environment I’d be happy to hear about your use case in a mail to [email protected] and able to give some support.
F. v. d. Meulen, M. Schauer: Bayesian estimation of discretely observed multi-dimensional diffusion processes using guided proposals. Electronic Journal of Statistics 11 (1), 2017, doi:10.1214/17-EJS1290.
M. Schauer, F. v. d. Meulen, H. v. Zanten: Guided proposals for simulating multi-dimensional diffusion bridges. Bernoulli 23 (4A), 2017, doi:10.3150/16-BEJ833.