bridgesde3d | R Documentation |
The (S3) generic function bridgesde3d
for simulation of 3-dim bridge stochastic differential equations,Itô or Stratonovich type, with different methods.
bridgesde3d(N, ...)
## Default S3 method:
bridgesde3d(N=1000,M=1, x0=c(0,0,0),
y=c(0,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 'bridgesde3d'
summary(object, at,
digits=NULL, ...)
## S3 method for class 'bridgesde3d'
time(x, ...)
## S3 method for class 'bridgesde3d'
mean(x, at, ...)
## S3 method for class 'bridgesde3d'
Median(x, at, ...)
## S3 method for class 'bridgesde3d'
Mode(x, at, ...)
## S3 method for class 'bridgesde3d'
quantile(x, at, ...)
## S3 method for class 'bridgesde3d'
kurtosis(x, at, ...)
## S3 method for class 'bridgesde3d'
skewness(x, at, ...)
## S3 method for class 'bridgesde3d'
min(x, at, ...)
## S3 method for class 'bridgesde3d'
max(x, at, ...)
## S3 method for class 'bridgesde3d'
moment(x, at, ...)
## S3 method for class 'bridgesde3d'
cv(x, at, ...)
## S3 method for class 'bridgesde3d'
bconfint(x, at, ...)
## S3 method for class 'bridgesde3d'
plot(x, ...)
## S3 method for class 'bridgesde3d'
lines(x, ...)
## S3 method for class 'bridgesde3d'
points(x, ...)
## S3 method for class 'bridgesde3d'
plot3D(x, display = c("persp","rgl"), ...)
N |
number of simulation steps. |
M |
number of trajectories. |
x0 |
initial value (numeric vector of length 3) of the process |
y |
terminal value (numeric vector of length 3) of the process |
t0 |
initial time. |
T |
final time. |
Dt |
time step of the simulation (discretization). If it is |
drift |
drift coefficient: an |
diffusion |
diffusion coefficient: an |
corr |
the correlation structure of three Brownian motions W1(t), W2(t) and W3(t); must be a real symmetric positive-definite square matrix of dimension 3. |
alpha |
weight |
mu |
weight |
type |
if |
method |
numerical methods of simulation, the default |
x , object |
an object inheriting from class |
at |
time between |
digits |
integer, used for number formatting. |
display |
|
... |
potentially further arguments for (non-default) methods. |
The function bridgesde3d
returns a mts
of the diffusion bridge starting at x
at time t0
and ending at y
at time T
. W1(t), W2(t) and W3(t) three 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.
bridgesde3d
returns an object inheriting from class
"bridgesde3d"
.
X , Y , Z |
an invisible |
driftx , drifty , driftz |
drift coefficient of X(t), Y(t) and Z(t). |
diffx , diffy , diffz |
diffusion coefficient of X(t), Y(t) and Z(t). |
Cx , Cy , Cz |
indices of crossing realized of X(t), Y(t)) and Z(t). |
type |
type of sde. |
method |
the numerical method used. |
A.C. Guidoum, K. Boukhetala.
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
bridgesde1d
for simulation of 1-dim SDE. DBridge
in package "sde".
bridgesde2d
for simulation of 2-dim SDE.
## dX(t) = 4*(-1-X(t))*Y(t) dt + 0.2 * dW1(t) ; x01 = 0 and y01 = 0
## dY(t) = 4*(1-Y(t)) *X(t) dt + 0.2 * dW2(t) ; x02 = -1 and y02 = -2
## dZ(t) = 4*(1-Z(t)) *Y(t) dt + 0.2 * dW3(t) ; x03 = 0.5 and y03 = 0.5
## W1(t), W2(t) and W3(t) are three correlated Brownian motions with Sigma
set.seed(1234)
fx <- expression(4*(-1-x)*y, 4*(1-y)*x, 4*(1-z)*y)
gx <- rep(expression(0.2),3)
# correlation matrix
Sigma <-matrix(c(1,0.3,0.5,0.3,1,0.2,0.5,0.2,1),nrow=3,ncol=3)
res <- bridgesde3d(x0=c(0,-1,0.5),y=c(0,-2,0.5),drift=fx,diffusion=gx,corr=Sigma,M=200)
res
summary(res) ## Monte-Carlo statistics at time T/2=0.5
summary(res,at=0.25) ## Monte-Carlo statistics at time 0.25
summary(res,at=0.75) ## Monte-Carlo statistics at time 0.75
##
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)
lines(time(res),apply(res$Z,1,mean),col=5,lwd=2)
legend("topleft",c(expression(E(X[t])),expression(E(Y[t])),
expression(E(Z[t]))),lty=1,inset = .01,col=c(3,4,5))
##
plot3D(res,display = "persp",main="3-dim bridge sde")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.