Bridge.jl
Documentation for Bridge.jl
Important concepts
Bridge.ContinuousTimeProcess
— Type.ContinuousTimeProcess{T}
Types inheriting from the abstract type ContinuousTimeProcess{T}
characterize the properties of a T
-valued stochastic process, play a similar role as distribution types like Exponential
in the package Distributions
.
Bridge.SamplePath
— Type.SamplePath{T} <: AbstractPath{T}
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 discretely observed ContinuousTimeProcess
es and for the sample path returned by direct and approximate samplers. tt
is the vector of the grid points of the observation/simulation and yy
is the corresponding vector of states.
It supports getindex, setindex!, length, copy, vcat, start, next, done, endof
.
Base.valtype
— Function.valtype(::ContinuousTimeProcess) -> T
Returns statespace (type) of a ContinuousTimeProcess{T]
.
Ordinary differential equations
Bridge.ODESolver
— Type.ODESolver
Abstract (super-)type for solving methods for ordinary differential equations.
Bridge.solve!
— Function.solve!(method, F, X::SamplePath, x0, P) -> X, [err]
Solve ordinary differential equation $(d/dx) x(t) = F(t, x(t), P)$ on the fixed grid X.tt
writing into X.yy
method::R3
- using a non-adaptive Ralston (1965) update (order 3).
method::BS3
use non-adaptive Bogacki–Shampine method to give error estimate.
Call _solve!
to inline. "Pretty fast if x
is a bitstype or a StaticArray."
Bridge.R3
— Type.R3
Ralston (1965) update (order 3 step of the Bogacki–Shampine 1989 method) to solve $y(t + dt) - y(t) = \int_t^{t+dt} F(s, y(s)) ds$.
Bridge.BS3
— Type.BS3
Ralston (1965) update (order 3 step of the Bogacki–Shampine 1989 method) to solve $y(t + dt) - y(t) = \int_t^{t+dt} F(s, y(s)) ds$. Uses Bogacki–Shampine method to give error estimate.
Brownian motion
Stochastic differential equations
StatsBase.sample
— Function.sample(tt, P, x1=zero(T))
Sample the process P
on the grid tt
exactly from its transitionprob
(-ability) starting in x1
.
StatsBase.sample!
— Function.sample!(X, P, x1=zero(T))
Sample the process P
on the grid X.tt
exactly from its transitionprob
(-ability) starting in x1
writing into X.yy
.
Bridge.quvar
— Function.quvar(X)
Computes quadratic variation of X
.
Bridge.bracket
— Function.bracket(X)
bracket(X,Y)
Computes quadratic variation process of X
(of X
and Y
).
Bridge.ito
— Function.ito(Y, X)
Integrate a valued stochastic process with respect to a stochastic differential.
Bridge.girsanov
— Function.girsanov(X::SamplePath, P::ContinuousTimeProcess, Pt::ContinuousTimeProcess)
Girsanov log likelihood $\mathrm{d}P/\mathrm{d}Pt(X)$
Bridge.lp
— Function.lp(s, x, t, y, P)
Log-transition density, shorthand for logpdf(transitionprob(s,x,t,P),y)
.
Bridge.llikelihood
— Function.llikelihood(X::SamplePath, P::ContinuousTimeProcess)
Log-likelihood of observations X
using transition density lp
.
llikelihood(X::SamplePath, Pº::LocalGammaProcess, P::LocalGammaProcess)
Log-likelihood dPº/dP
. (Up to proportionality.)
llikelihood(X::SamplePath, P::LocalGammaProcess)
Bridge log-likelihood with respect to reference measure P.P
. (Up to proportionality.)
Bridge.euler
— Function.euler(u, W, P) -> X
Solve stochastic differential equation $dX_t = b(t,X_t)dt + σ(t,X_t)dW_t$ using the Euler scheme.
Bridge.euler!
— Function.euler!(Y, u, W, P) -> X
Solve stochastic differential equation $dX_t = b(t,X_t)dt + σ(t,X_t)dW_t$ using the Euler scheme in place.
Bridge.thetamethod
— Function.thetamethod(u, W, P, theta=0.5)
Solve stochastic differential equation using the theta method and Newton-Raphson steps.
Levy processes
Bridge.GammaProcess
— Type.GammaProcess
A GammaProcess with jump rate γ
and inverse jump size λ
has increments Gamma(t*γ, 1/λ)
and Levy measure
Here Gamma(α,θ)
is the Gamma distribution in julia's parametrization with shape parameter α
and scale θ
.
Bridge.nu
— Function. nu(k,P)
(Bin-wise) integral of the Levy measure $\nu(B_k)$.
Miscellaneous
Bridge.endpoint!
— Function.endpoint!(X::SamplePath, v)
Convenience functions setting the endpoint of X to
v`.
Bridge.inner
— Function.inner(x[, y])
Short-hand for quadratic form x'x (or x'y).
Bridge.cumsum0
— Function.cumsum0
Cumulative sum starting at 0,
Bridge.mat
— Function.mat(X::SamplePath{SVector})
mat(yy::Vector{SVector})
Reinterpret X
or yy
to an array without change in memory.
Bridge.outer
— Function.outer(x[, y])
Short-hand for quadratic form xx' (or xy').
Bridge.CSpline
— Type.CSpline(s, t, x, y = x, m0 = (y-x)/(t-s), m1 = (y-x)/(t-s))
Cubic spline parametrized by $f(s) = x$ and $f(t) = y$, $f'(s) = m_0$, $f'(t) = m_1$.
Bridge.integrate
— Function.integrate(cs::CSpline, s, t)
Integrate the cubic spline from s
to t
.
Bridge.logpdfnormal
— Function.logpdfnormal(x, A)
logpdf of centered gaussian with covariance A
Online statistics
Online updating of the tuple state = (m, m2, n)
where
m
- mean(x[1:n])
m2
- sum of squares of differences from the current mean, $\textstyle\sum_{i=1}^n (x_i - \bar x_n)^2$
n
- number of iterations
Bridge.mcstart
— Function.mcstart(x) -> state
Create state for random chain online statitics. The entries/value of x
are ignored
Bridge.mcnext
— Function.mcnext(state, x) -> state
Update random chain online statistics when new chain value x
was observed. Return new state
.
Bridge.mcband
— Function.mcband(mc)
Compute marginal 95% coverage interval for the chain from normal approximation.
Bridge.mcbandmean
— Function.mcmeanband(mc)
Compute marginal confidence interval for the chain mean using normal approximation
Linear Processes
Bridge.LinPro
— Type.LinPro(B, μ::T, σ)
Linear diffusion $dX = B(X - μ)dt + σdW$.
Bridge.Ptilde
— Type.Ptilde(cs::CSpline, σ)
Affine diffusion $dX = cs(t) dt + σdW$ with cs a cubic spline ::CSpline
.
Bridges
Bridge.GuidedProp
— Type.GuidedProp
General bridge proposal process
Bridge.bridge
— Function.bridge(W, P, scheme! = euler!) -> Y
Integrate with scheme!
and set Y[end] = P.v1
.
Bridge.Vs
— Function.Vs(s, T1, T2, v, P)
Time changed V
for generation of U
.
Bridge.mdb
— Function.mdb(u, W, P)
mdb!(copy(W), u, W, P)
Euler scheme with the diffusion coefficient correction of the modified diffusion bridge.
Bridge.mdb!
— Function.mdb(u, W, P)
mdb!(copy(W), u, W, P)
Euler scheme with the diffusion coefficient correction of the modified diffusion bridge.
Bridge.r
— Function.r(t, x, T, v, P)
Returns $r(t,x) = \operatorname{grad}_x \log p(t,x; T, v)$ where $p$ is the transition density of the process $P$.
Bridge.gpK!
— Function.gpK!(K::SamplePath, P)
Precompute $K = H^{-1}$ from $(d/dt)K = BK + KB' + a$ for a guided proposal.
Unsorted
Bridge.LocalGammaProcess
— Type.LocalGammaProcess
Bridge.compensator0
— Function.compensator0(kstart, P::LocalGammaProcess)
Compensator of GammaProcess approximating the LocalGammaProcess. For kstart == 1
(only choice) this is $\nu_0([b_1,\infty)$.
Bridge.compensator
— Function.compensator(kstart, P::LocalGammaProcess)
Compensator of LocalGammaProcess For kstart = 1
, this is $\sum_{k=1}^N \nu(B_k)$, for kstart = 0
, this is $\sum_{k=0}^N \nu(B_k) - C$ (where $C$ is a constant).
Bridge.θ
— Function.θ(x, P::LocalGammaProcess)
Inverse jump size compared to gamma process with same alpha and beta.
Bridge.soft
— Function.soft(t, T1, T2)
Time change mapping s
in [T1, T2]
(
U
-time) to
t`in
[T1, T2](
X`-time).
Bridge.tofs
— Function.tofs(s, T1, T2)
Time change mapping t
in [T1, T2]
($X$-time) to s
in [T1, T2]
(U
-time).
Bridge.dotVs
— Function.dotVs (s, T1, T2, v, P)
Time changed time derivative of V
for generation of U
.