Bridge.jl

Bridge.jl

Documentation for Bridge.jl

Important concepts

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.

source
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 ContinuousTimeProcesses 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.

source
Base.valtypeFunction.
valtype(::ContinuousTimeProcess) -> T

Returns statespace (type) of a ContinuousTimeProcess{T].

source

Ordinary differential equations

ODESolver

Abstract (super-)type for solving methods for ordinary differential equations.

source
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."

source
Bridge.R3Type.
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$.

source
Bridge.BS3Type.
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.

source

Brownian motion

Stochastic differential equations

StatsBase.sampleFunction.
sample(tt, P, x1=zero(T))

Sample the process P on the grid tt exactly from its transitionprob(-ability) starting in x1.

source
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.

source
Bridge.quvarFunction.
quvar(X)

Computes quadratic variation of X.

source
Bridge.bracketFunction.
bracket(X)
bracket(X,Y)

Computes quadratic variation process of X (of X and Y).

source
Bridge.itoFunction.
ito(Y, X)

Integrate a valued stochastic process with respect to a stochastic differential.

source
Bridge.girsanovFunction.
girsanov(X::SamplePath, P::ContinuousTimeProcess, Pt::ContinuousTimeProcess)

Girsanov log likelihood $\mathrm{d}P/\mathrm{d}Pt(X)$

source
Bridge.lpFunction.
lp(s, x, t, y, P)

Log-transition density, shorthand for logpdf(transitionprob(s,x,t,P),y).

source
Bridge.llikelihoodFunction.
llikelihood(X::SamplePath, P::ContinuousTimeProcess)

Log-likelihood of observations X using transition density lp.

source
llikelihood(X::SamplePath, Pº::LocalGammaProcess, P::LocalGammaProcess)

Log-likelihood dPº/dP. (Up to proportionality.)

source
llikelihood(X::SamplePath, P::LocalGammaProcess)

Bridge log-likelihood with respect to reference measure P.P. (Up to proportionality.)

source
Bridge.eulerFunction.
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.

source
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.

source
Bridge.thetamethodFunction.
thetamethod(u, W, P, theta=0.5)

Solve stochastic differential equation using the theta method and Newton-Raphson steps.

source

Levy processes

GammaProcess

A GammaProcess with jump rate γ and inverse jump size λ has increments Gamma(t*γ, 1/λ) and Levy measure

\[ν(x)=γ x^{-1}\exp(-λ x), \]

Here Gamma(α,θ) is the Gamma distribution in julia's parametrization with shape parameter α and scale θ.

source
Bridge.nuFunction.
 nu(k,P)

(Bin-wise) integral of the Levy measure $\nu(B_k)$.

source

Miscellaneous

Bridge.endpoint!Function.
endpoint!(X::SamplePath, v)

Convenience functions setting the endpoint of X tov`.

source
Bridge.innerFunction.
inner(x[, y])

Short-hand for quadratic form x'x (or x'y).

source
Bridge.cumsum0Function.
cumsum0

Cumulative sum starting at 0,

source
Bridge.matFunction.
mat(X::SamplePath{SVector}) 
mat(yy::Vector{SVector})

Reinterpret X or yy to an array without change in memory.

source
Bridge.outerFunction.
outer(x[, y])

Short-hand for quadratic form xx' (or xy').

source
Bridge.CSplineType.
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$.

source
Bridge.integrateFunction.
integrate(cs::CSpline, s, t)

Integrate the cubic spline from s to t.

source
Bridge.logpdfnormalFunction.
logpdfnormal(x, A)

logpdf of centered gaussian with covariance A

source

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.mcstartFunction.
mcstart(x) -> state

Create state for random chain online statitics. The entries/value of x are ignored

source
Bridge.mcnextFunction.
mcnext(state, x) -> state

Update random chain online statistics when new chain value x was observed. Return new state.

source
Bridge.mcbandFunction.
mcband(mc)

Compute marginal 95% coverage interval for the chain from normal approximation.

source
Bridge.mcbandmeanFunction.
mcmeanband(mc)

Compute marginal confidence interval for the chain mean using normal approximation

source

Linear Processes

Bridge.LinProType.
LinPro(B, μ::T, σ)

Linear diffusion $dX = B(X - μ)dt + σdW$.

source
Bridge.PtildeType.
Ptilde(cs::CSpline, σ)

Affine diffusion $dX = cs(t) dt + σdW$ with cs a cubic spline ::CSpline.

source

Bridges

GuidedProp

General bridge proposal process

source
Bridge.bridgeFunction.
bridge(W, P, scheme! = euler!) -> Y

Integrate with scheme! and set Y[end] = P.v1.

source
Bridge.VsFunction.
Vs(s, T1, T2, v, P)

Time changed V for generation of U.

source
Bridge.mdbFunction.
mdb(u, W, P)
mdb!(copy(W), u, W, P)

Euler scheme with the diffusion coefficient correction of the modified diffusion bridge.

source
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.

source
Bridge.rFunction.
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$.

source
Bridge.gpK!Function.
gpK!(K::SamplePath, P)

Precompute $K = H^{-1}$ from $(d/dt)K = BK + KB' + a$ for a guided proposal.

source

Unsorted

LocalGammaProcess
source
Bridge.compensator0Function.
compensator0(kstart, P::LocalGammaProcess)

Compensator of GammaProcess approximating the LocalGammaProcess. For kstart == 1 (only choice) this is $\nu_0([b_1,\infty)$.

source
Bridge.compensatorFunction.
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).

source
Bridge.θFunction.
θ(x, P::LocalGammaProcess)

Inverse jump size compared to gamma process with same alpha and beta.

source
Bridge.softFunction.
soft(t, T1, T2)

Time change mapping s in [T1, T2](U-time) tot`in[T1, T2](X`-time).

source
Bridge.tofsFunction.
tofs(s, T1, T2)

Time change mapping t in [T1, T2] ($X$-time) to s in [T1, T2] (U-time).

source
Bridge.dotVsFunction.
dotVs (s, T1, T2, v, P)

Time changed time derivative of V for generation of U.

source