bridgesde2d: Simulation of 2-D Bridge SDE's

View source: R/bridgesde.R

bridgesde2dR Documentation

Simulation of 2-D Bridge SDE's

Description

The (S3) generic function bridgesde2d for simulation of 2-dim bridge stochastic differential equations,Itô or Stratonovich type, with different methods.

Usage

bridgesde2d(N, ...)
## Default S3 method:
bridgesde2d(N = 1000, M = 1, x0 = c(0, 0), 
   y = c(0, 0),t0 = 0, T = 1, Dt,drift, diffusion, corr = NULL,
   alpha = 0.5, mu = 0.5,type = c("ito", "str"),method = 
   c("euler", "milstein","predcorr", "smilstein", "taylor",
   "heun", "rk1", "rk2", "rk3"), ...)

## S3 method for class 'bridgesde2d'
summary(object, at,
     digits=NULL, ...)								  
## S3 method for class 'bridgesde2d'
time(x, ...)
## S3 method for class 'bridgesde2d'
mean(x, at, ...)
## S3 method for class 'bridgesde2d'
Median(x, at, ...)
## S3 method for class 'bridgesde2d'
Mode(x, at, ...)
## S3 method for class 'bridgesde2d'
quantile(x, at, ...)
## S3 method for class 'bridgesde2d'
kurtosis(x, at, ...)
## S3 method for class 'bridgesde2d'
skewness(x, at, ...)
## S3 method for class 'bridgesde2d'
min(x, at, ...)
## S3 method for class 'bridgesde2d'
max(x, at, ...)
## S3 method for class 'bridgesde2d'
moment(x, at, ...)
## S3 method for class 'bridgesde2d'
cv(x, at, ...)
## S3 method for class 'bridgesde2d'
bconfint(x, at, ...)

## S3 method for class 'bridgesde2d'
plot(x, ...)
## S3 method for class 'bridgesde2d'
lines(x, ...)
## S3 method for class 'bridgesde2d'
points(x, ...)	
## S3 method for class 'bridgesde2d'
plot2d(x, ...)
## S3 method for class 'bridgesde2d'
lines2d(x, ...)
## S3 method for class 'bridgesde2d'
points2d(x, ...)								  

Arguments

N

number of simulation steps.

M

number of trajectories.

x0

initial value (numeric vector of length 2) of the process X_t and Y_t at time t_0.

y

terminal value (numeric vector of length 2) of the process X_t and Y_t at time T.

t0

initial time.

T

final time.

Dt

time step of the simulation (discretization). If it is missing a default \Delta t = \frac{T-t_{0}}{N}.

drift

drift coefficient: an expression of three variables t, x and y for process X_t and Y_t.

diffusion

diffusion coefficient: an expression of three variables t, x and y for process X_t and Y_t.

corr

the correlation structure of two Brownian motions W1(t) and W2(t); must be a real symmetric positive-definite square matrix of dimension 2.

alpha, mu

weight of the predictor-corrector scheme; the default alpha = 0.5 and mu = 0.5.

type

if type="ito" simulation diffusion bridge of Itô type, else type="str" simulation diffusion bridge of Stratonovich type; the default type="ito".

method

numerical methods of simulation, the default method = "euler"; see snssde2d.

x, object

an object inheriting from class "bridgesde2d".

at

time between t0 and T. Monte-Carlo statistics of the solution (X_{t},Y_{t}) at time at. The default at = T/2.

digits

integer, used for number formatting.

...

potentially further arguments for (non-default) methods.

Details

The function bridgesde2d returns a mts of the diffusion bridge starting at x at time t0 and ending at y at time T. W1(t) and W2(t) are two standard Brownian motion independent if corr=NULL.

The methods of approximation are classified according to their different properties. Mainly two criteria of optimality are used in the literature: the strong and the weak (orders of) convergence. The method of simulation can be one among: Euler-Maruyama Order 0.5, Milstein Order 1, Milstein Second-Order, Predictor-Corrector method, Itô-Taylor Order 1.5, Heun Order 2 and Runge-Kutta Order 1, 2 and 3.

An overview of this package, see browseVignettes('Sim.DiffProc') for more informations.

Value

bridgesde2d returns an object inheriting from class "bridgesde2d".

X, Y

an invisible mts (2-dim) object (X(t),Y(t)).

driftx, drifty

drift coefficient of X(t) and Y(t).

diffx, diffy

diffusion coefficient of X(t) and Y(t).

Cx, Cy

indices of crossing realized of X(t) and Y(t).

type

type of sde.

method

the numerical method used.

Author(s)

A.C. Guidoum, K. Boukhetala.

References

Bladt, M. and Sorensen, M. (2007). Simple simulation of diffusion bridges with application to likelihood inference for diffusions. Working Paper, University of Copenhagen.

Iacus, S.M. (2008). Simulation and inference for stochastic differential equations: with R examples. Springer-Verlag, New York

See Also

bridgesde1d for simulation of 1-dim SDE.

DBridge in package "sde".

Examples

## dX(t) = 4*(-1-X(t)) dt + 0.2 dW1(t)
## dY(t) = X(t) dt + 0 dW2(t)
## x01 = 0 , y01 = 0
## x02 = 0, y02 = 0 
## W1(t) and W2(t) two correlated Brownian motion with matrix Sigma=matrix(c(1,0.7,0.7,1),nrow=2)
set.seed(1234)

fx <- expression(4*(-1-x) , x)
gx <- expression(0.2 , 0)
Sigma= matrix(c(1,0.7,0.7,1),nrow=2)
res <- bridgesde2d(drift=fx,diffusion=gx,Dt=0.005,M=500,corr=Sigma)
res
summary(res) ## Monte-Carlo statistics at time T/2=2.5
summary(res,at=1) ## Monte-Carlo statistics at time 1
summary(res,at=4) ## Monte-Carlo statistics at time 4
##
plot(res,type="n")
lines(time(res),apply(res$X,1,mean),col=3,lwd=2)
lines(time(res),apply(res$Y,1,mean),col=4,lwd=2)
legend("topright",c(expression(E(X[t])),expression(E(Y[t]))),lty=1,inset = .7,col=c(3,4))
##
plot2d(res)

Sim.DiffProc documentation built on May 29, 2024, 8:09 a.m.