Consider now a \(d\)-dimensional stochastic process \(X_{t}\) defined on a probability space \((\Omega, \mathfrak{F},\mathbb{P})\). We say that the bridge associated to \(X_{t}\) conditioned to the event \(\{X_{T}= a\}\) is the process: \[ \{\tilde{X}_{t}, t_{0} \leq t \leq T \}=\{X_{t}, t_{0} \leq t \leq T: X_{T}= a \} \] where \(T\) is a deterministic fixed time and \(a \in \mathbb{R}^d\) is fixed too.
bridgesdekd()
functionsThe (S3) generic function bridgesdekd()
(where k=1,2,3
) for simulation of 1,2 and 3-dim bridge stochastic differential equations,Itô or Stratonovich type, with different methods. The main arguments consist:
drift
and diffusion
coefficients as R expressions that depend on the state variable x
(y
and z
) and time variable t
.N
.M
(default: M=1
).x0
at initial time t0
.y
final time T
Dt
(default: Dt=(T-t0)/N
).type="ito"
for Ito or type="str"
for Stratonovich (by default type="ito"
).method
(default method="euler"
).By Monte-Carlo simulations, the following statistical measures (S3 method
) for class bridgesdekd()
(where k=1,2,3
) can be approximated for the process at any time \(t \in [t_{0},T]\) (default: at=(T-t0)/2
):
mean
.moment
with order=2
and center=TRUE
.Median
.Mode
.quantile
.min
and max
.skewness
and kurtosis
.cv
.moment
.summary
.We can just make use of the rsdekd()
function (where k=1,2,3
) to build our random number for class bridgesdekd()
(where k=1,2,3
) at any time \(t \in [t_{0},T]\). the main arguments consist:
object
an object inheriting from class bridgesdekd()
(where k=1,2,3
).at
time between \(s=t0\) and \(t=T\).The function dsde()
(where k=1,2,3
) approximate transition density for class bridgesdekd()
(where k=1,2,3
), the main arguments consist:
object
an object inheriting from class bridgesdekd()
(where k=1,2,3
).at
time between \(s=t0\) and \(t=T\).pdf
probability density function Joint
or Marginal
.The following we explain how to use this functions.
bridgesde1d()
We simulate a flow of \(5000\) trajectories, with integration step size \(\Delta t = 0.001\), and \(x_0 = 3\) at time \(t_0 = 0\), \(y = 1\) at terminal time \(T=1\).
R> f <- expression((1-x)/(1-t))
R> g <- expression(x)
R> mod <- bridgesde1d(drift=f,diffusion=g,x0=3,y=1,M=5000,method="milstein")
R> mod
Itô Bridge Sde 1D:
| dX(t) = (1 - X(t))/(1 - t) * dt + X(t) * dW(t)
Method:
| First-order Milstein scheme
Summary:
| Size of process | N = 1001.
| Crossing realized | C = 4867 among 5000.
| Initial value | x0 = 3.
| Ending value | y = 1.
| Time of process | t in [0,1].
| Discretization | Dt = 0.001.
R> summary(mod) ## default: summary at time = (T-t0)/2
Monte-Carlo Statistics for X(t) at time t = 0.5
| Crossing realized 4867 among 5000
Mean 1.98885
Variance 1.77324
Median 1.64320
Mode 1.17678
First quartile 1.13710
Third quartile 2.42955
Minimum 0.27403
Maximum 18.37158
Skewness 3.00933
Kurtosis 21.78707
Coef-variation 0.66955
3th-order moment 7.10590
4th-order moment 68.50650
5th-order moment 807.96130
6th-order moment 11063.58311
In Figure 1, we present the flow of trajectories, the mean path (red lines) of solution of \(X_{t}|X_{0}=3,X_{T}=1\):
R> plot(mod,ylab=expression(X[t]))
R> lines(time(mod),apply(mod$X,1,mean),col=2,lwd=2)
R> legend("topleft","mean path",inset = .01,col=2,lwd=2,cex=0.8,bty="n")
Bridge sde 1D
Figure 2, show approximation results for \(m(t)=\text{E}(X_{t}|X_{0}=3,X_{T}=1)\) and \(S(t)=\text{V}(X_{t}|X_{0}=3,X_{T}=1)\):
R> m <- apply(mod$X,1,mean)
R> S <- apply(mod$X,1,var)
R> out <- data.frame(m,S)
R> matplot(time(mod), out, type = "l", xlab = "time", ylab = "", col=2:3,lwd=2,lty=2:3,las=1)
R> legend("topright",c(expression(m(t),S(t))),col=2:3,lty=2:3,lwd=2,bty="n")
The following statistical measures (S3 method
) for class bridgesde1d()
can be approximated for the \(X_{t}|X_{0}=3,X_{T}=1\) process at any time \(t\), for example at=0.55
:
R> s = 0.55
R> mean(mod, at = s)
[1] 1.8894
R> moment(mod, at = s , center = TRUE , order = 2) ## variance
[1] 1.637
R> Median(mod, at = s)
[1] 1.5516
R> Mode(mod, at = s)
[1] 1.1551
R> quantile(mod , at = s)
0% 25% 50% 75% 100%
0.26933 1.07407 1.55160 2.30822 14.93212
R> kurtosis(mod , at = s)
[1] 16.668
R> skewness(mod , at = s)
[1] 2.7531
R> cv(mod , at = s )
[1] 0.67724
R> min(mod , at = s)
[1] 0.26933
R> max(mod , at = s)
[1] 14.932
R> moment(mod, at = s , center= TRUE , order = 4)
[1] 44.684
R> moment(mod, at = s , center= FALSE , order = 4)
[1] 136.08
The result summaries of the \(X_{t}|X_{0}=3,X_{T}=1\) process at time \(t=0.55\):
R> summary(mod, at = 0.55)
Monte-Carlo Statistics for X(t) at time t = 0.55
| Crossing realized 4867 among 5000
Mean 1.88942
Variance 1.63733
Median 1.55160
Mode 1.15509
First quartile 1.07407
Third quartile 2.30822
Minimum 0.26933
Maximum 14.93212
Skewness 2.75308
Kurtosis 16.66763
Coef-variation 0.67724
3th-order moment 5.76799
4th-order moment 44.68364
5th-order moment 395.80345
6th-order moment 4029.66373
Hence we can just make use of the rsde1d()
function to build our random number generator for \(X_{t}|X_{0}=3,X_{T}=1\) at time \(t=0.55\):
R> x <- rsde1d(object = mod, at = s)
R> head(x, n = 10)
[1] 1.3899 1.9269 1.1417 2.5318 1.8477 1.9571 1.1669 1.1350 2.3831
[10] 2.3246
R> summary(x)
Min. 1st Qu. Median Mean 3rd Qu. Max.
0.269 1.074 1.552 1.889 2.308 14.932
Display the random number generator for \(X_{t}|X_{0}=3,X_{T}=1\), see Figure 3:
R> plot(time(mod),mod$X[,1],type="l",ylab="X(t)",xlab="time",axes=F,lty=3)
R> points(s,x[1],pch=19,col=2,cex=0.5)
R> lines(c(s,s),c(0,x[1]),lty=2,col=2)
R> lines(c(0,s),c(x[1],x[1]),lty=2,col=2)
R> axis(1, s, bquote(at==.(s)),col=2,col.ticks=2)
R> axis(2, x[1], bquote(X[t==.(s)]),col=2,col.ticks=2)
R> legend('topright',col=2,pch=19,legend=bquote(X[t==.(s)]==.(x[1])),bty = 'n')
R> box()
The function dsde1d()
can be used to show the kernel density estimation for \(X_{t}|X_{0}=3,X_{T}=1\) at time \(t=0.55\) (hist=TRUE
based on truehist()
function in MASS package):
R> dens <- dsde1d(mod, at = s)
R> dens
Density of X(t-t0)|X(t0) = 3, X(T) = 1 at time t = 0.55
Data: x (4867 obs.); Bandwidth 'bw' = 0.1517
x f(x)
Min. :-0.1858 Min. :0.00000
1st Qu.: 3.7074 1st Qu.:0.00050
Median : 7.6007 Median :0.00284
Mean : 7.6007 Mean :0.06415
3rd Qu.:11.4940 3rd Qu.:0.04416
Max. :15.3873 Max. :0.54378
R> plot(dens,hist=TRUE) ## histgramme
R> plot(dens,add=TRUE) ## kernel density
Approximate the transitional densitie of \(X_{t}|X_{0}=3,X_{T}=1\) at \(t-s = \{0.25,0.75\}\):
R> plot(dsde1d(mod,at=0.75))
R> plot(dsde1d(mod,at=0.25),add=TRUE)
R> legend('topright',col=c('#0000FF4B','#FF00004B'),pch=15,legend=c("t-s=0.25","t-s=0.75"),bty = 'n')
Transitional densitie at time \(t-s = 0.25,0.75\)
bridgesde2d()
Assume that we want to describe the following \(2\)-dimensional bridge SDE’s in Stratonovich form:
\[\begin{equation}\label{eq:09} \begin{cases} dX_t = -(1+Y_{t}) X_{t} dt + 0.2 (1-Y_{t})\circ dW_{1,t},\quad X_{t_{0}}=1 \quad\text{and}\quad X_{T}=1\\ dY_t = -(1+X_{t}) Y_{t} dt + 0.1 (1-X_{t}) \circ dW_{2,t},\quad Y_{t_{0}}=-0.5 \quad\text{and}\quad Y_{T}=0.5 \end{cases} \end{equation}\]We simulate a flow of \(5000\) trajectories, with integration step size \(\Delta t = 0.01\), and using Runge-Kutta method order 1:
R> fx <- expression(-(1+y)*x , -(1+x)*y)
R> gx <- expression(0.2*(1-y),0.1*(1-x))
R> mod2 <- bridgesde2d(drift=fx,diffusion=gx,x0=c(1,-0.5),y=c(1,0.5),Dt=0.01,M=5000,type="str",method="rk1")
R> mod2
Stratonovich Bridge Sde 2D:
| dX(t) = expression(-(1 + y) * x) * dt + expression(0.2 * (1 - y)) o dW1(t)
| dY(t) = expression(-(1 + x) * y) * dt + expression(0.1 * (1 - x)) o dW2(t)
Method:
| Runge-Kutta method with order 1
Summary:
| Size of process | N = 1001.
| Crossing realized | C = 4998 among 5000.
| Initial values | x0 = (1,-0.5).
| Ending values | y = (1,0.5).
| Time of process | t in [0,10].
| Discretization | Dt = 0.01.
R> summary(mod2) ## default: summary at time = (T-t0)/2
Monte-Carlo Statistics for (X(t),Y(t)) at time t = 5
| Crossing realized 4998 among 5000
X Y
Mean 0.00831 -0.00254
Variance 0.01958 0.00512
Median 0.00589 -0.00209
Mode 0.00908 -0.00237
First quartile -0.08401 -0.04804
Third quartile 0.10458 0.04361
Minimum -0.51235 -0.28787
Maximum 0.63896 0.28438
Skewness 0.02005 -0.03346
Kurtosis 3.05726 3.52267
Coef-variation 16.84698 -28.14040
3th-order moment 0.00005 -0.00001
4th-order moment 0.00117 0.00009
5th-order moment 0.00001 0.00000
6th-order moment 0.00012 0.00000
In Figure 6, we present the flow of trajectories of \(X_{t}|X_{0}=1,X_{T}=1\) and \(Y_{t}|Y_{0}=-0.5,Y_{T}=0.5\):
R> plot(mod2,col=c('#FF00004B','#0000FF82'))
Bridge sde 2D
Figure 7, show approximation results for \(m_{1}(t)=\text{E}(X_{t}|X_{0}=1,X_{T}=1)\), \(m_{2}(t)=\text{E}(Y_{t}|Y_{0}=-0.5,Y_{T}=0.5)\),and \(S_{1}(t)=\text{V}(X_{t}|X_{0}=1,X_{T}=1)\), \(S_{2}(t)=\text{V}(Y_{t}|Y_{0}=-0.5,Y_{T}=0.5)\), and \(C_{12}(t)=\text{COV}(X_{t},Y_{t}|X_{0}=1,Y_{0}=-0.5,X_{T}=1,Y_{T}=0.5)\):
R> m1 <- apply(mod2$X,1,mean)
R> m2 <- apply(mod2$Y,1,mean)
R> S1 <- apply(mod2$X,1,var)
R> S2 <- apply(mod2$Y,1,var)
R> C12 <- sapply(1:dim(mod2$X)[1],function(i) cov(mod2$X[i,],mod2$Y[i,]))
R> out2 <- data.frame(m1,m2,S1,S2,C12)
R> matplot(time(mod2), out2, type = "l", xlab = "time", ylab = "", col=2:6,lwd=2,lty=2:6,las=1)
R> legend("top",c(expression(m[1](t),m[2](t),S[1](t),S[2](t),C[12](t))),col=2:6,lty=2:6,lwd=2,bty="n")
The following statistical measures (S3 method
) for class bridgesde2d()
can be approximated for the \(X_{t}|X_{0}=1,X_{T}=1\) and \(Y_{t}|Y_{0}=-0.5,Y_{T}=0.5\) process at any time \(t\), for example at=6.75
:
R> s = 6.75
R> mean(mod2, at = s)
[1] 0.0262114 0.0081418
R> moment(mod2, at = s , center = TRUE , order = 2) ## variance
[1] 0.0189430 0.0044029
R> Median(mod2, at = s)
[1] 0.027607 0.007586
R> Mode(mod2, at = s)
[1] 0.029730 0.016605
R> quantile(mod2 , at = s)
$x
0% 25% 50% 75% 100%
-0.526017 -0.065815 0.027607 0.119805 0.514165
$y
0% 25% 50% 75% 100%
-0.226050 -0.035261 0.007586 0.050428 0.269088
R> kurtosis(mod2 , at = s)
[1] 3.0379 3.2958
R> skewness(mod2 , at = s)
[1] -0.029599 0.077942
R> cv(mod2 , at = s )
[1] 5.2514 8.1506
R> min(mod2 , at = s)
[1] -0.52602 -0.22605
R> max(mod2 , at = s)
[1] 0.51417 0.26909
R> moment(mod2 , at = s , center= TRUE , order = 4)
[1] 0.001090540 0.000063916
R> moment(mod2 , at = s , center= FALSE , order = 4)
[1] 0.001161006 0.000066414
The result summaries of the \(X_{t}|X_{0}=1,X_{T}=1\) and \(Y_{t}|Y_{0}=-0.5,Y_{T}=0.5\) process at time \(t=6.75\):
R> summary(mod2, at = 6.75)
Monte-Carlo Statistics for (X(t),Y(t)) at time t = 6.75
| Crossing realized 4998 among 5000
X Y
Mean 0.02621 0.00814
Variance 0.01895 0.00440
Median 0.02761 0.00759
Mode 0.02973 0.01660
First quartile -0.06581 -0.03526
Third quartile 0.11981 0.05043
Minimum -0.52602 -0.22605
Maximum 0.51417 0.26909
Skewness -0.02960 0.07794
Kurtosis 3.03789 3.29578
Coef-variation 5.25142 8.15062
3th-order moment -0.00008 0.00002
4th-order moment 0.00109 0.00006
5th-order moment -0.00002 0.00000
6th-order moment 0.00010 0.00000
Hence we can just make use of the rsde2d()
function to build our random number generator for the couple \(X_{t},Y_{t}|X_{0}=1,Y_{0}=-0.5,X_{T}=1,Y_{T}=0.5\) at time \(t=6.75\):
R> x2 <- rsde2d(object = mod2, at = s)
R> head(x2, n = 10)
x y
1 0.043042 0.120642
2 0.125878 -0.035650
3 -0.077395 -0.029731
4 0.029373 0.034140
5 0.234090 -0.030233
6 -0.058481 -0.038801
7 0.082880 0.062691
8 0.161389 -0.033351
9 0.161178 0.023149
10 0.186901 -0.084523
R> summary(x2)
x y
Min. :-0.5260 Min. :-0.22605
1st Qu.:-0.0658 1st Qu.:-0.03526
Median : 0.0276 Median : 0.00759
Mean : 0.0262 Mean : 0.00814
3rd Qu.: 0.1198 3rd Qu.: 0.05043
Max. : 0.5142 Max. : 0.26909
Display the random number generator for the couple \(X_{t},Y_{t}|X_{0}=1,Y_{0}=-0.5,X_{T}=1,Y_{T}=0.5\), see Figure 8:
R> plot(ts.union(mod2$X[,1],mod2$Y[,1]),col=1:2,lty=3,plot.type="single",type="l",ylab= "",xlab="time",axes=F)
R> points(s,x2$x[1],pch=19,col=3,cex=0.8)
R> points(s,x2$y[1],pch=19,col=4,cex=0.8)
R> lines(c(s,s),c(-10,x2$x[1]),lty=2,col=6)
R> lines(c(0,s),c(x2$x[1],x2$x[1]),lty=2,col=3)
R> lines(c(0,s),c(x2$y[1],x2$y[1]),lty=2,col=4)
R> axis(1, s, bquote(at==.(s)),col=6,col.ticks=6)
R> axis(2, x2$x[1], bquote(X[t==.(s)]),col=3,col.ticks=3)
R> axis(2, x2$y[1], bquote(Y[t==.(s)]),col=4,col.ticks=4)
R> legend('topright',legend=bquote(c(X[t==.(s)]==.(x2$x[1]),Y[t==.(s)]==.(x2$y[1]))),bty = 'n')
R> box()
For each SDE type and for each numerical scheme, the density of \(X_{t}|X_{0}=1,X_{T}=1\) and \(Y_{t}|Y_{0}=-0.5,Y_{T}=0.5\) at time \(t=6.75\) are reported using dsde2d()
function, see e.g. Figure 9:
R> denM <- dsde2d(mod2,pdf="M",at =s)
R> denM
Marginal density of X(t-t0)|X(t0) = 1, X(T) = 1 at time t = 6.75
Data: x (4998 obs.); Bandwidth 'bw' = 0.02256
x f(x)
Min. :-0.59368 Min. :0.00004
1st Qu.:-0.29980 1st Qu.:0.02190
Median :-0.00593 Median :0.32544
Mean :-0.00593 Mean :0.84986
3rd Qu.: 0.28795 3rd Qu.:1.66051
Max. : 0.58183 Max. :2.98021
Marginal density of Y(t-t0)|Y(t0) = -0.5, Y(T) = 0.5 at time t = 6.75
Data: y (4998 obs.); Bandwidth 'bw' = 0.01048
y f(y)
Min. :-0.257486 Min. :0.0001
1st Qu.:-0.117984 1st Qu.:0.0658
Median : 0.021519 Median :0.6257
Mean : 0.021519 Mean :1.7903
3rd Qu.: 0.161021 3rd Qu.:3.2947
Max. : 0.300524 Max. :6.2143
R> plot(denM, main="Marginal Density")
Created using dsde2d()
plotted in (x, y)-space with dim = 2
. A contour
and image
plot of density obtained from a realization of the couple \(X_{t},Y_{t}|X_{0}=1,Y_{0}=-0.5,X_{T}=1,Y_{T}=0.5\) at time t=6.75
.
R> denJ <- dsde2d(mod2, pdf="J", n=100,at =s)
R> denJ
Joint density of (X(t-t0),Y(t-t0)|X(t0)=1,Y(t0)=-0.5,X(T)=1,Y(T)=0.5) at time t = 6.75
Data: (x,y) (2 x 4998 obs.);
x y f(x,y)
Min. :-0.52602 Min. :-0.226050 Min. : 0.0000
1st Qu.:-0.26597 1st Qu.:-0.102266 1st Qu.: 0.0187
Median :-0.00593 Median : 0.021519 Median : 0.2126
Mean :-0.00593 Mean : 0.021519 Mean : 1.9019
3rd Qu.: 0.25412 3rd Qu.: 0.145304 3rd Qu.: 1.8044
Max. : 0.51417 Max. : 0.269088 Max. :18.1845
R> plot(denJ,display="contour",main="Bivariate Transition Density at time t=6.755")
R> plot(denJ,display="image",main="Bivariate Transition Density at time t=6.755")
A \(3\)D plot of the transition density at \(t=6.75\) obtained with:
R> plot(denJ,main="Bivariate Transition Density at time t=6.75")
We approximate the bivariate transition density over the set transition horizons \(t\in [1,9]\) with \(\Delta t = 0.005\) using the code:
R> for (i in seq(1,9,by=0.005)){
+ plot(dsde2d(mod2, at = i,n=100),display="contour",main=paste0('Transition Density \n t = ',i))
+ }
bridgesde3d()
Assume that we want to describe the following bridges SDE’s (3D) in Itô form:
\[\begin{equation} \begin{cases} dX_t = -4 (1+X_{t}) Y_{t} dt + 0.2 dW_{1,t},\quad X_{t_{0}}=0 \quad\text{and}\quad X_{T}=0\\ dY_t = 4 (1-Y_{t}) X_{t} dt + 0.2 dW_{2,t},\quad Y_{t_{0}}=-1 \quad\text{and}\quad Y_{T}=-2\\ dZ_t = 4 (1-Z_{t}) Y_{t} dt + 0.2 dW_{3,t},\quad Z_{t_{0}}=0.5 \quad\text{and}\quad Z_{T}=0.5 \end{cases} \end{equation}\]We simulate a flow of \(5000\) trajectories, with integration step size \(\Delta t = 0.001\).
R> fx <- expression(-4*(1+x)*y, 4*(1-y)*x, 4*(1-z)*y)
R> gx <- rep(expression(0.2),3)
R> mod3 <- bridgesde3d(x0=c(0,-1,0.5),y=c(0,-2,0.5),drift=fx,diffusion=gx,M=5000)
R> mod3
Itô Bridge Sde 3D:
| dX(t) = -4 * (1 + X(t)) * Y(t) * dt + 0.2 * dW1(t)
| dY(t) = 4 * (1 - Y(t)) * X(t) * dt + 0.2 * dW2(t)
| dZ(t) = 4 * (1 - Z(t)) * Y(t) * dt + 0.2 * dW3(t)
Method:
| Euler scheme with order 0.5
Summary:
| Size of process | N = 1001.
| Crossing realized | C = 4999 among 5000.
| Initial values | x0 = (0,-1,0.5).
| Ending values | y = (0,-2,0.5).
| Time of process | t in [0,1].
| Discretization | Dt = 0.001.
R> summary(mod3) ## default: summary at time = (T-t0)/2
Monte-Carlo Statistics for (X(t),Y(t),Z(t)) at time t = 0.5
| Crossing realized 4999 among 5000
X Y Z
Mean 0.68653 0.50785 0.12107
Variance 0.01029 0.00698 0.01816
Median 0.68864 0.50841 0.12488
Mode 0.69242 0.50334 0.12559
First quartile 0.61908 0.45142 0.03415
Third quartile 0.75518 0.56523 0.20786
Minimum 0.25672 0.16737 -0.52257
Maximum 1.04315 0.86211 0.62638
Skewness -0.15098 0.01036 -0.11041
Kurtosis 3.38517 3.11495 3.39254
Coef-variation 0.14778 0.16452 1.11314
3th-order moment -0.00016 0.00001 -0.00027
4th-order moment 0.00036 0.00015 0.00112
5th-order moment -0.00002 0.00000 -0.00006
6th-order moment 0.00002 0.00001 0.00013
For plotting (back in time) using the command plot
, and plot3D
in space the results of the simulation are shown in Figure 12:
R> plot(mod3) ## in time
R> plot3D(mod3,display = "persp",main="3D Bridge SDE's") ## in space
Bridge sde 3D
Figure 13, show approximation results for \(m_{1}(t)=\text{E}(X_{t}|X_{0}=0,X_{T}=0)\), \(m_{2}(t)=\text{E}(Y_{t}|Y_{0}=-1,Y_{T}=-2)\), \(m_{3}(t)=\text{E}(Z_{t}|Z_{0}=0.5,Z_{T}=0.5)\) and \(S_{1}(t)=\text{V}(X_{t}|X_{0}=0,X_{T}=0)\), \(S_{2}(t)=\text{V}(Y_{t}|Y_{0}=-1,Y_{T}=-2)\), \(S_{3}(t)=\text{V}(Z_{t}|Z_{0}=0.5,Z_{T}=0.5)\),
R> m1 <- apply(mod3$X,1,mean)
R> m2 <- apply(mod3$Y,1,mean)
R> m3 <- apply(mod3$Z,1,mean)
R> S1 <- apply(mod3$X,1,var)
R> S2 <- apply(mod3$Y,1,var)
R> S3 <- apply(mod3$Z,1,var)
R> out3 <- data.frame(m1,m2,m3,S1,S2,S3)
R> matplot(time(mod3), out3, type = "l", xlab = "time", ylab = "", col=2:7,lwd=2,lty=2:7,las=1)
R> legend("bottom",c(expression(m[1](t),m[2](t),m[3](t),S[1](t),S[2](t),S[3](t))),col=2:7,lty=2:7,lwd=2,bty="n")
The following statistical measures (S3 method
) for class bridgesde3d()
can be approximated for the \(X_{t}|X_{0}=0,X_{T}=0\), \(Y_{t}|Y_{0}=-1,Y_{T}=-2\) and \(Z_{t}|Z_{0}=0.5,Z_{T}=0.5\) process at any time \(t\), for example at=0.75
:
R> s = 0.75
R> mean(mod3, at = s)
[1] 1.99799 0.12365 -0.50026
R> moment(mod3, at = s , center = TRUE , order = 2) ## variance
[1] 0.0111796 0.0046333 0.0303640
R> Median(mod3, at = s)
[1] 1.99805 0.12401 -0.49570
R> Mode(mod3, at = s)
[1] 1.99994 0.11729 -0.48679
R> quantile(mod3 , at = s)
$x
0% 25% 50% 75% 100%
1.5794 1.9274 1.9981 2.0685 2.4690
$y
0% 25% 50% 75% 100%
-0.119346 0.077376 0.124011 0.170480 0.342551
$z
0% 25% 50% 75% 100%
-1.146458 -0.616306 -0.495701 -0.380309 0.023288
R> kurtosis(mod3 , at = s)
[1] 3.1877 2.9253 3.0427
R> skewness(mod3 , at = s)
[1] 0.024288 -0.067700 -0.159695
R> cv(mod3 , at = s )
[1] 0.052925 0.550527 -0.348362
R> min(mod3 , at = s)
[1] 1.57943 -0.11935 -1.14646
R> max(mod3 , at = s)
[1] 2.469047 0.342551 0.023288
R> moment(mod3 , at = s , center= TRUE , order = 4)
[1] 0.000398562 0.000062823 0.002806433
R> moment(mod3 , at = s , center= FALSE , order = 4)
[1] 16.20417730 0.00071111 0.11271846
The result summaries of the \(X_{t}|X_{0}=0,X_{T}=0\), \(Y_{t}|Y_{0}=-1,Y_{T}=-2\) and \(Z_{t}|Z_{0}=0.5,Z_{T}=0.5\) process at time \(t=0.75\):
R> summary(mod3, at = 0.75)
Monte-Carlo Statistics for (X(t),Y(t),Z(t)) at time t = 0.75
| Crossing realized 4999 among 5000
X Y Z
Mean 1.99799 0.12365 -0.50026
Variance 0.01118 0.00463 0.03037
Median 1.99805 0.12401 -0.49570
Mode 1.99994 0.11729 -0.48679
First quartile 1.92736 0.07738 -0.61631
Third quartile 2.06850 0.17048 -0.38031
Minimum 1.57943 -0.11935 -1.14646
Maximum 2.46905 0.34255 0.02329
Skewness 0.02429 -0.06770 -0.15969
Kurtosis 3.18767 2.92532 3.04273
Coef-variation 0.05293 0.55053 -0.34836
3th-order moment 0.00003 -0.00002 -0.00085
4th-order moment 0.00040 0.00006 0.00281
5th-order moment 0.00001 0.00000 -0.00028
6th-order moment 0.00003 0.00000 0.00044
Hence we can just make use of the rsde3d()
function to build our random number generator for the triplet \(X_{t},Y_{t},Z_{t}|X_{0}=0,Y_{0}=-1,Z_{0}=0.5,X_{T}=0,Y_{T}=-2,Z_{T}=0.5\) at time \(t=0.75\):
R> x3 <- rsde3d(object = mod3, at = s)
R> head(x3, n = 10)
x y z
1 1.8355 0.131688 -0.69490
2 2.0015 0.089879 -0.93915
3 1.8763 0.100016 -0.21363
4 1.9862 0.096988 -0.33648
5 1.9481 0.088273 -0.51417
6 1.9816 0.091262 -0.55365
7 1.8656 -0.015809 -0.19584
8 1.9345 0.061476 -0.21002
9 2.0464 0.180440 -0.47656
10 1.8703 0.150539 -0.68193
R> summary(x3)
x y z
Min. :1.58 Min. :-0.1193 Min. :-1.1465
1st Qu.:1.93 1st Qu.: 0.0774 1st Qu.:-0.6163
Median :2.00 Median : 0.1240 Median :-0.4957
Mean :2.00 Mean : 0.1236 Mean :-0.5003
3rd Qu.:2.07 3rd Qu.: 0.1705 3rd Qu.:-0.3803
Max. :2.47 Max. : 0.3426 Max. : 0.0233
Display the random number generator for triplet \(X_{t},Y_{t},Z_{t}|X_{0}=0,Y_{0}=-1,Z_{0}=0.5,X_{T}=0,Y_{T}=-2,Z_{T}=0.5\) at time \(t=0.75\): , see Figure 14:
R> plot(ts.union(mod3$X[,1],mod3$Y[,1],mod3$Z[,1]),col=1:3,lty=3,plot.type="single",type="l",ylab= "",xlab="time",axes=F)
R> points(s,x3$x[1],pch=19,col=4,cex=0.8)
R> points(s,x3$y[1],pch=19,col=5,cex=0.8)
R> points(s,x3$z[1],pch=19,col=6,cex=0.8)
R> lines(c(s,s),c(-10,x3$x[1]),lty=2,col=7)
R> lines(c(0,s),c(x3$x[1],x3$x[1]),lty=2,col=4)
R> lines(c(0,s),c(x3$y[1],x3$y[1]),lty=2,col=5)
R> lines(c(0,s),c(x3$z[1],x3$z[1]),lty=2,col=6)
R> axis(1, s, bquote(at==.(s)),col=7,col.ticks=7)
R> axis(2, x3$x[1], bquote(X[t==.(s)]),col=4,col.ticks=4)
R> axis(2, x3$y[1], bquote(Y[t==.(s)]),col=5,col.ticks=5)
R> axis(2, x3$z[1], bquote(Z[t==.(s)]),col=6,col.ticks=6)
R> legend("bottomleft",legend=bquote(c(X[t==.(s)]==.(x3$x[1]),Y[t==.(s)]==.(x3$y[1]),Z[t==.(s)]==.(x3$z[1]))),bty = 'n',cex=0.75)
R> box()
For each SDE type and for each numerical scheme, the density of \(X_{t}|X_{0}=0,X_{T}=0\), \(Y_{t}|Y_{0}=-1,Y_{T}=-2\) and \(Z_{t}|Z_{0}=0.5,Z_{T}=0.5\) process at time \(t=0.75\) are reported using dsde3d()
function, see e.g. Figure 15:
R> denM <- dsde3d(mod3,pdf="M",at =s)
R> denM
Marginal density of X(t-t0)|X(t0) = 0, X(T) = 0 at time t = 0.75
Data: x (4999 obs.); Bandwidth 'bw' = 0.01726
x f(x)
Min. :1.5276 Min. :0.0001
1st Qu.:1.7759 1st Qu.:0.0153
Median :2.0242 Median :0.2646
Mean :2.0242 Mean :1.0059
3rd Qu.:2.2725 3rd Qu.:1.8871
Max. :2.5208 Max. :3.8189
Marginal density of Y(t-t0)|Y(t0) = -1, Y(T) = -2 at time t = 0.75
Data: y (4999 obs.); Bandwidth 'bw' = 0.01115
y f(y)
Min. :-0.15281 Min. :0.0001
1st Qu.:-0.02060 1st Qu.:0.0795
Median : 0.11160 Median :0.8635
Mean : 0.11160 Mean :1.8891
3rd Qu.: 0.24381 3rd Qu.:3.7539
Max. : 0.37601 Max. :5.6220
Marginal density of Z(t-t0)|Z(t0) = 0.5, Z(T) = 0.5 at time t = 0.75
Data: z (4999 obs.); Bandwidth 'bw' = 0.02856
z f(z)
Min. :-1.23212 Min. :0.00005
1st Qu.:-0.89685 1st Qu.:0.04242
Median :-0.56158 Median :0.38913
Mean :-0.56158 Mean :0.74494
3rd Qu.:-0.22632 3rd Qu.:1.44102
Max. : 0.10895 Max. :2.26437
R> plot(denM, main="Marginal Density")
For an approximate joint density for triplet \(X_{t},Y_{t},Z_{t}|X_{0}=0,Y_{0}=-1,Z_{0}=0.5,X_{T}=0,Y_{T}=-2,Z_{T}=0.5\) at time \(t=0.75\) (for more details, see package sm or ks.)
R> denJ <- dsde3d(mod3,pdf="J",at=0.75)
R> plot(denJ,display="rgl")
snssdekd()
& dsdekd()
& rsdekd()
- Monte-Carlo Simulation and Analysis of Stochastic Differential Equations.bridgesdekd()
& dsdekd()
& rsdekd()
- Constructs and Analysis of Bridges Stochastic Differential Equations.fptsdekd()
& dfptsdekd()
- Monte-Carlo Simulation and Kernel Density Estimation of First passage time.MCM.sde()
& MEM.sde()
- Parallel Monte-Carlo and Moment Equations for SDEs.fitsde()
- Parametric Estimation of 1-D Stochastic Differential Equation.Bladt, M. and Sorensen, M. (2007). Simple simulation of diffusion bridges with application to likelihood inference for diffusions. Working Paper, University of Copenhagen.
Guidoum AC, Boukhetala K (2017). Sim.DiffProc: Simulation of Diffusion Processes. R package version 3.8, URL https://cran.r-project.org/package=Sim.DiffProc.
Department of Probabilities & Statistics, Faculty of Mathematics, University of Science and Technology Houari Boumediene, BP 32 El-Alia, U.S.T.H.B, Algeria, E-mail (acguidoum@usthb.dz)↩
Faculty of Mathematics, University of Science and Technology Houari Boumediene, BP 32 El-Alia, U.S.T.H.B, Algeria, E-mail (kboukhetala@usthb.dz)↩