bridgesde1d | R Documentation |
The (S3) generic function bridgesde1d
for simulation of 1-dim bridge stochastic differential equations,Itô or Stratonovich type, with different methods.
bridgesde1d(N, ...)
## Default S3 method:
bridgesde1d(N = 1000, M=1, x0 = 0, y = 0, t0 = 0, T = 1, Dt,
drift, diffusion, 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 'bridgesde1d'
summary(object, at ,digits=NULL, ...)
## S3 method for class 'bridgesde1d'
time(x, ...)
## S3 method for class 'bridgesde1d'
mean(x, at, ...)
## S3 method for class 'bridgesde1d'
Median(x, at, ...)
## S3 method for class 'bridgesde1d'
Mode(x, at, ...)
## S3 method for class 'bridgesde1d'
quantile(x, at, ...)
## S3 method for class 'bridgesde1d'
kurtosis(x, at, ...)
## S3 method for class 'bridgesde1d'
skewness(x, at, ...)
## S3 method for class 'bridgesde1d'
min(x, at, ...)
## S3 method for class 'bridgesde1d'
max(x, at, ...)
## S3 method for class 'bridgesde1d'
moment(x, at, ...)
## S3 method for class 'bridgesde1d'
cv(x, at, ...)
## S3 method for class 'bridgesde1d'
bconfint(x, at, ...)
## S3 method for class 'bridgesde1d'
plot(x, ...)
## S3 method for class 'bridgesde1d'
lines(x, ...)
## S3 method for class 'bridgesde1d'
points(x, ...)
N |
number of simulation steps. |
M |
number of trajectories. |
x0 |
initial value of the process at time |
y |
terminal value of the process at time |
t0 |
initial time. |
T |
final time. |
Dt |
time step of the simulation (discretization). If it is |
drift |
drift coefficient: an |
diffusion |
diffusion coefficient: an |
alpha , mu |
weight of the predictor-corrector scheme; the default |
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. |
... |
potentially further arguments for (non-default) methods. |
The function bridgesde1d
returns a trajectory of the diffusion bridge starting at x
at time t0
and ending at y
at time T
.
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.
bridgesde1d
returns an object inheriting from class
"bridgesde1d"
.
X |
an invisible |
drift |
drift coefficient. |
diffusion |
diffusion coefficient. |
C |
indices of crossing realized of X(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
bridgesde2d
and bridgesde3d
for 2 and 3-dim.
DBridge
in package "sde".
## Example 1: Ito bridge sde
## Ito Bridge sde
## dX(t) = 2*(1-X(t)) *dt + dW(t)
## x0 = 2 at time t0=0 , and y = 1 at time T=1
set.seed(1234)
f <- expression( 2*(1-x) )
g <- expression( 1 )
mod1 <- bridgesde1d(drift=f,diffusion=g,x0=2,y=1,M=1000)
mod1
summary(mod1) ## Monte-Carlo statistics at T/2=0.5
summary(mod1,at=0.75) ## Monte-Carlo statistics at 0.75
## Not run:
plot(mod1)
lines(time(mod1),apply(mod1$X,1,mean),col=2,lwd=2)
lines(time(mod1),apply(mod1$X,1,bconfint,level=0.95)[1,],col=4,lwd=2)
lines(time(mod1),apply(mod1$X,1,bconfint,level=0.95)[2,],col=4,lwd=2)
legend("topleft",c("mean path",paste("bound of", 95," percent confidence")),
inset = .01,col=c(2,4),lwd=2,cex=0.8)
## End(Not run)
## Example 2: Stratonovich sde
## dX(t) = ((2-X(t))/(2-t)) dt + X(t) o dW(t)
## x0 = 2 at time t0=0 , and y = 2 at time T=1
set.seed(1234)
f <- expression((2-x)/(2-t))
g <- expression(x)
mod2 <- bridgesde1d(type="str",drift=f,diffusion=g,M=1000,x0=2,y=2)
mod2
summary(mod2,at = 0.25) ## Monte-Carlo statistics at 0.25
summary(mod2,at = 0.5) ## Monte-Carlo statistics at 0.5
summary(mod2,at = 0.75 )## Monte-Carlo statistics at 0.75
## Not run:
plot(mod2)
lines(time(mod2),apply(mod2$X,1,mean),col=2,lwd=2)
lines(time(mod2),apply(mod2$X,1,bconfint,level=0.95)[1,],col=4,lwd=2)
lines(time(mod2),apply(mod2$X,1,bconfint,level=0.95)[2,],col=4,lwd=2)
legend("topright",c("mean path",paste("bound of", 95," percent confidence")),
inset = .01,col=c(2,4),lwd=2,cex=0.8)
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.