Assume that we want to describe the following SDE:
Ito form1:
\[\begin{equation}\label{eq:05} dX_{t} = \frac{1}{2}\theta^{2} X_{t} dt + \theta X_{t} dW_{t},\qquad X_{0}=x_{0} > 0 \end{equation}\] Stratonovich form: \[\begin{equation}\label{eq:06} dX_{t} = \frac{1}{2}\theta^{2} X_{t} dt +\theta X_{t} \circ dW_{t},\qquad X_{0}=x_{0} > 0 \end{equation}\]In the above \(f(t,x)=\frac{1}{2}\theta^{2} x\) and \(g(t,x)= \theta x\) (\(\theta > 0\)), \(W_{t}\) is a standard Wiener process. To simulate this models using snssde1d() function we need to specify:
drift and diffusion coefficients as R expressions that depend on the state variable x and time variable t.N=1000 (by default: N=1000).M=500 (by default: M=1).t0=0, x0=10 and end time T=1 (by default: t0=0, x0=0 and T=1).Dt=0.001 (by default: Dt=(T-t0)/N).type="ito" for Ito or type="str" for Stratonovich (by default type="ito").method (by default method="euler").theta = 0.5
f <- expression( (0.5*theta^2*x) )
g <- expression( theta*x )
mod1 <- snssde1d(drift=f,diffusion=g,x0=10,M=500,type="ito") # Using Ito
mod2 <- snssde1d(drift=f,diffusion=g,x0=10,M=500,type="str") # Using Stratonovich
mod1## Ito Sde 1D:
## | dX(t) = (0.5 * theta^2 * X(t)) * dt + theta * X(t) * dW(t)
## Method:
## | Euler scheme of order 0.5
## Summary:
## | Size of process | N = 1001.
## | Number of simulation | M = 500.
## | Initial value | x0 = 10.
## | Time of process | t in [0,1].
## | Discretization | Dt = 0.001.
mod2## Stratonovich Sde 1D:
## | dX(t) = (0.5 * theta^2 * X(t)) * dt + theta * X(t) o dW(t)
## Method:
## | Euler scheme of order 0.5
## Summary:
## | Size of process | N = 1001.
## | Number of simulation | M = 500.
## | Initial value | x0 = 10.
## | Time of process | t in [0,1].
## | Discretization | Dt = 0.001.
Using Monte-Carlo simulations, the following statistical measures (S3 method) for class snssde1d() can be approximated for the \(X_{t}\) process at any time \(t\):
mean.median.quantile.skewness and kurtosis.moment.bconfint.The summary of the results of mod1 and mod2 at time \(t=1\) of class snssde1d() is given by:
summary(mod1, at = 1)##
## Monte-Carlo Statistics for X(t) at time t = 1
##
## Mean 11.04319
## Variance 34.51348
## Median 9.67614
## First quartile 6.88003
## Third quartile 13.35324
## Skewness 1.41619
## Kurtosis 5.60120
## Moment of order 3 287.14814
## Moment of order 4 6672.03563
## Moment of order 5 125490.75202
## Int.conf Inf (95%) 3.62732
## Int.conf Sup (95%) 25.94828
summary(mod2, at = 1)##
## Monte-Carlo Statistics for X(t) at time t = 1
##
## Mean 9.61779
## Variance 28.14247
## Median 8.48501
## First quartile 6.11644
## Third quartile 12.16451
## Skewness 2.61720
## Kurtosis 19.52529
## Moment of order 3 390.73323
## Moment of order 4 15464.00809
## Moment of order 5 668137.79437
## Int.conf Inf (95%) 2.98779
## Int.conf Sup (95%) 20.63778
Hence we can just make use of the rsde1d() function to build our random number generator for the conditional density of the \(X_{t}|X_{0}\) (\(X_{t}^{\text{mod1}}| X_{0}\) and \(X_{t}^{\text{mod2}}|X_{0}\)) at time \(t = 1\).
x1 <- rsde1d(object = mod1, at = 1) # X(t=1) | X(0)=x0 (Itô SDE)
x2 <- rsde1d(object = mod2, at = 1) # X(t=1) | X(0)=x0 (Stratonovich SDE)
summary(data.frame(x1,x2))## x1 x2
## Min. : 2.204 Min. : 1.500
## 1st Qu.: 6.880 1st Qu.: 6.116
## Median : 9.676 Median : 8.485
## Mean :11.043 Mean : 9.618
## 3rd Qu.:13.353 3rd Qu.:12.165
## Max. :38.898 Max. :59.373
The function dsde1d() can be used to show the kernel density estimation for \(X_{t}|X_{0}\) at time \(t=1\) with log-normal curves:
mu1 = log(10); sigma1= sqrt(theta^2) # log mean and log variance for mod1
mu2 = log(10)-0.5*theta^2 ; sigma2 = sqrt(theta^2) # log mean and log variance for mod2
AppdensI <- dsde1d(mod1, at = 1)
AppdensS <- dsde1d(mod2, at = 1)
plot(AppdensI , dens = function(x) dlnorm(x,meanlog=mu1,sdlog = sigma1))
plot(AppdensS , dens = function(x) dlnorm(x,meanlog=mu2,sdlog = sigma2))
In Figure 2, we present the flow of trajectories, the mean path (red lines) of solution of and , with their empirical \(95\%\) confidence bands, that is to say from the \(2.5th\) to the \(97.5th\) percentile for each observation at time \(t\) (blue lines):
plot(mod1,plot.type="single",ylab=expression(X^mod1))
lines(time(mod1),mean(mod1),col=2,lwd=2)
lines(time(mod1),bconfint(mod1,level=0.95)[,1],col=4,lwd=2)
lines(time(mod1),bconfint(mod1,level=0.95)[,2],col=4,lwd=2)
legend("topleft",c("mean path",paste("bound of",95,"% confidence")),col=c(2,4),lwd=2,cex=0.8)
plot(mod2,plot.type="single",ylab=expression(X^mod2))
lines(time(mod2),mean(mod2),col=2,lwd=2)
lines(time(mod2),bconfint(mod2,level=0.95)[,1],col=4,lwd=2)
lines(time(mod2),bconfint(mod2,level=0.95)[,2],col=4,lwd=2)
legend("topleft",c("mean path",paste("bound of",95,"% confidence")),col=c(2,4),lwd=2,cex=0.8)
The following \(2\)-dimensional SDE’s with a vector of drift and a diagonal matrix of diffusion coefficients:
Ito form: \[\begin{equation}\label{eq:09} \begin{cases} dX_t = f_{x}(t,X_{t},Y_{t}) dt + g_{x}(t,X_{t},Y_{t}) dW_{1,t}\\ dY_t = f_{y}(t,X_{t},Y_{t}) dt + g_{y}(t,X_{t},Y_{t}) dW_{2,t} \end{cases} \end{equation}\] Stratonovich form: \[\begin{equation}\label{eq:10} \begin{cases} dX_t = f_{x}(t,X_{t},Y_{t}) dt + g_{x}(t,X_{t},Y_{t}) \circ dW_{1,t}\\ dY_t = f_{y}(t,X_{t},Y_{t}) dt + g_{y}(t,X_{t},Y_{t}) \circ dW_{2,t} \end{cases} \end{equation}\]\(W_{1,t}\) and \(W_{2,t}\) is a two independent standard Wiener process. To simulate \(2d\) models using snssde2d() function we need to specify:
drift (2d) and diffusion (2d) coefficients as R expressions that depend on the state variable x, y and time variable t.N (default: N=1000).M (default: M=1).t0, x0 and end time T (default: t0=0, x0=c(0,0) and T=1).Dt (default: Dt=(T-t0)/N).type="ito" for Ito or type="str" for Stratonovich (default type="ito").method (default method="euler").We simulate a flow of \(500\) trajectories of \((X_{t},Y_{t})\), with integration step size \(\Delta t = 0.01\), and using second Milstein method.
x=5;y=0
mu=3;sigma=0.5
fx <- expression(-(x/mu),x)
gx <- expression(sqrt(sigma),0)
mod2d <- snssde2d(drift=fx,diffusion=gx,Dt=0.01,M=500,x0=c(x,y),method="smilstein")
mod2d## Ito Sde 2D:
## | dX(t) = -(X(t)/mu) * dt + sqrt(sigma) * dW1(t)
## | dY(t) = X(t) * dt + 0 * dW2(t)
## Method:
## | Second Milstein scheme of order 1.5
## Summary:
## | Size of process | N = 1001.
## | Number of simulation | M = 500.
## | Initial values | (x0,y0) = (5,0).
## | Time of process | t in [0,10].
## | Discretization | Dt = 0.01.
summary(mod2d)##
## Monte-Carlo Statistics for (X(t),Y(t)) at time t = 10
## X Y
## Mean 0.16435 14.47619
## Variance 0.82224 24.76162
## Median 0.12762 14.46409
## First quartile -0.50351 10.83284
## Third quartile 0.76750 17.93925
## Skewness 0.18123 0.02360
## Kurtosis 2.89768 2.58953
## Moment of order 3 0.13512 2.90818
## Moment of order 4 1.95904 1587.73848
## Moment of order 5 0.79038 108.72301
## Int.conf Inf (95%) -1.46128 4.66876
## Int.conf Sup (95%) 2.06418 24.03663
For plotting (back in time) using the command plot, the results of the simulation are shown in Figure 3.
plot(mod2d)
Take note of the well known result, which can be derived from either this equations. That for any \(t > 0\) the OU process \(X_t\) and its integral \(Y_t\) will be the normal distribution with mean and variance given by: \[ \begin{cases} \text{E}(X_{t}) =x_{0} e^{-t/\mu} &\text{and}\quad\text{Var}(X_{t})=\frac{\sigma \mu}{2} \left (1-e^{-2t/\mu}\right )\\ \text{E}(Y_{t}) = y_{0}+x_{0}\mu \left (1-e^{-t/\mu}\right ) &\text{and}\quad\text{Var}(Y_{t})=\sigma\mu^{3}\left (\frac{t}{\mu}-2\left (1-e^{-t/\mu}\right )+\frac{1}{2}\left (1-e^{-2t/\mu}\right )\right ) \end{cases} \]
Hence we can just make use of the rsde2d() function to build our random number for \((X_{t},Y_{t})\) at time \(t = 10\).
out <- rsde2d(object = mod2d, at = 10)
summary(out)## x y
## Min. :-2.6374 Min. : 1.749
## 1st Qu.:-0.5035 1st Qu.:10.833
## Median : 0.1276 Median :14.464
## Mean : 0.1643 Mean :14.476
## 3rd Qu.: 0.7675 3rd Qu.:17.939
## Max. : 3.1062 Max. :27.377
For each SDE type and for each numerical scheme, the density of \(X_t\) and \(Y_t\) at time \(t=10\) are reported using dsde2d() function, see e.g. Figure 4: the marginal density of \(X_t\) and \(Y_t\) at time \(t=10\).
denM <- dsde2d(mod2d,pdf="M",at =10)
denM##
## Marginal density of X(t)|X(0)=x0 at time t = 10
##
## Data: x (500 obs.); Bandwidth 'bw' = 0.2355
##
## x f(x)
## Min. :-3.343802 Min. :0.0000383
## 1st Qu.:-1.554705 1st Qu.:0.0059189
## Median : 0.234393 Median :0.0645287
## Mean : 0.234393 Mean :0.1395981
## 3rd Qu.: 2.023491 3rd Qu.:0.2782026
## Max. : 3.812589 Max. :0.4003288
##
## Marginal density of Y(t)|Y(0)=y0 at time t = 10
##
## Data: y (500 obs.); Bandwidth 'bw' = 1.292
##
## y f(y)
## Min. :-2.127983 Min. :0.00000732
## 1st Qu.: 6.217316 1st Qu.:0.00361625
## Median :14.562615 Median :0.02345091
## Mean :14.562615 Mean :0.02992753
## 3rd Qu.:22.907915 3rd Qu.:0.05506848
## Max. :31.253214 Max. :0.07816690
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 system \((X_{t},Y_{t})\) at time t=10.
denJ <- dsde2d(mod2d,pdf="J",at =10)
denJ##
## Joint density of (X(t),Y(t)|X(0)=x0,Y(0)=y0) at time t = 10
##
## Data: (x,y) (2 x 500 obs.);
##
## x y f(x,y)
## Min. :-2.6373740 Min. : 1.748691 Min. :0.00000000
## 1st Qu.:-1.2014904 1st Qu.: 8.155653 1st Qu.:0.00030789
## Median : 0.2343932 Median :14.562615 Median :0.00267498
## Mean : 0.2343932 Mean :14.562615 Mean :0.00660300
## 3rd Qu.: 1.6702768 3rd Qu.:20.969578 3rd Qu.:0.00991682
## Max. : 3.1061604 Max. :27.376540 Max. :0.03481699
plot(denJ,display="contour",main="Bivariate Density")
plot(denJ,display="image",drawpoints=TRUE,col.pt="green",cex=0.25,pch=19,main="Bivariate Density")
A \(3\)D plot of the density obtained with:
plot(denJ,main="Bivariate Density")
Implemente in R as follows, with integration step size \(\Delta t = 0.01\) and using stochastic Runge-Kutta methods 1-stage.
mu = 4; sigma=0.1
fx <- expression( y , (mu*( 1-x^2 )* y - x))
gx <- expression( 0 ,2*sigma)
mod2d <- snssde2d(drift=fx,diffusion=gx,N=10000,Dt=0.01,type="str",method="rk1")
mod2d## Stratonovich Sde 2D:
## | dX(t) = Y(t) * dt + 0 o dW1(t)
## | dY(t) = (mu * (1 - X(t)^2) * Y(t) - X(t)) * dt + 2 * sigma o dW2(t)
## Method:
## | Runge-Kutta method of order 1
## Summary:
## | Size of process | N = 10001.
## | Number of simulation | M = 1.
## | Initial values | (x0,y0) = (0,0).
## | Time of process | t in [0,100].
## | Discretization | Dt = 0.01.
plot2d(mod2d) ## in plane (O,X,Y)
plot(mod2d) ## back in time
The following \(3\)-dimensional SDE’s with a vector of drift and a diagonal matrix of diffusion coefficients:
Ito form: \[\begin{equation}\label{eq17} \begin{cases} dX_t = f_{x}(t,X_{t},Y_{t},Z_{t}) dt + g_{x}(t,X_{t},Y_{t},Z_{t}) dW_{1,t}\\ dY_t = f_{y}(t,X_{t},Y_{t},Z_{t}) dt + g_{y}(t,X_{t},Y_{t},Z_{t}) dW_{2,t}\\ dZ_t = f_{z}(t,X_{t},Y_{t},Z_{t}) dt + g_{z}(t,X_{t},Y_{t},Z_{t}) dW_{3,t} \end{cases} \end{equation}\] Stratonovich form: \[\begin{equation}\label{eq18} \begin{cases} dX_t = f_{x}(t,X_{t},Y_{t},Z_{t}) dt + g_{x}(t,X_{t},Y_{t},Z_{t}) \circ dW_{1,t}\\ dY_t = f_{y}(t,X_{t},Y_{t},Z_{t}) dt + g_{y}(t,X_{t},Y_{t},Z_{t}) \circ dW_{2,t}\\ dZ_t = f_{z}(t,X_{t},Y_{t},Z_{t}) dt + g_{z}(t,X_{t},Y_{t},Z_{t}) \circ dW_{3,t} \end{cases} \end{equation}\]\(W_{1,t}\), \(W_{2,t}\) and \(W_{3,t}\) is a 3 independent standard Wiener process. To simulate this system using snssde3d() function we need to specify:
drift (3d) and diffusion (3d) coefficients as R expressions that depend on the state variables x, y , z and time variable t.N (default: N=1000).M (default: M=1).t0, x0 and end time T (default: t0=0, x0=c(0,0,0) and T=1).Dt (default: Dt=(T-t0)/N).type="ito" for Ito or type="str" for Stratonovich (default type="ito").method (default method="euler").We simulate a flow of 500 trajectories, with integration step size \(\Delta t = 0.001\).
fx <- expression(4*(-1-x)*y , 4*(1-y)*x , 4*(1-z)*y)
gx <- rep(expression(0.2),3)
mod3d <- snssde3d(x0=c(x=2,y=-2,z=-2),drift=fx,diffusion=gx,N=1000,M=500)
mod3d## Ito 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 of order 0.5
## Summary:
## | Size of process | N = 1001.
## | Number of simulation | M = 500.
## | Initial values | (x0,y0,z0) = (2,-2,-2).
## | Time of process | t in [0,1].
## | Discretization | Dt = 0.001.
summary(mod3d)##
## Monte-Carlo Statistics for (X(t),Y(t),Z(t)) at time t = 1
## X Y Z
## Mean -0.79248 0.86862 0.79436
## Variance 0.00870 0.09680 0.00967
## Median -0.79519 0.81592 0.80148
## First quartile -0.86393 0.63552 0.72762
## Third quartile -0.73416 1.05511 0.86827
## Skewness 0.33787 0.58001 -0.22169
## Kurtosis 3.19287 2.96294 2.77273
## Moment of order 3 0.00027 0.01747 -0.00021
## Moment of order 4 0.00024 0.02776 0.00026
## Moment of order 5 0.00003 0.01234 -0.00002
## Int.conf Inf (95%) -0.95983 0.38219 0.59257
## Int.conf Sup (95%) -0.59059 1.58563 0.97425
plot(mod3d,union = TRUE) ## back in time
plot3D(mod3d,display="persp") ## in space (O,X,Y,Z)
For each SDE type and for each numerical scheme, the marginal density of \(X_t\), \(Y_t\) and \(Z_t\) at time \(t=1\) are reported using dsde3d() function, see e.g. Figure 8.
den <- dsde3d(mod3d,pdf="M",at =1)
den##
## Marginal density for X(t)|X(0)=x0 at time t = 1
##
## Data: x (500 obs.); Bandwidth 'bw' = 0.02422
##
## x f(x)
## Min. :-1.0999820 Min. :0.000374
## 1st Qu.:-0.9165744 1st Qu.:0.049528
## Median :-0.7331667 Median :0.581942
## Mean :-0.7331667 Mean :1.361743
## 3rd Qu.:-0.5497591 3rd Qu.:2.762533
## Max. :-0.3663515 Max. :4.147440
##
## Marginal density for Y(t)|Y(0)=y0 at time t = 1
##
## Data: y (500 obs.); Bandwidth 'bw' = 0.08079
##
## y f(y)
## Min. :-0.0270619 Min. :0.0001683
## 1st Qu.: 0.4918822 1st Qu.:0.0685150
## Median : 1.0108264 Median :0.3049756
## Mean : 1.0108264 Mean :0.4812734
## 3rd Qu.: 1.5297705 3rd Qu.:0.9105846
## Max. : 2.0487147 Max. :1.3148288
##
## Marginal density for Z(t)|Z(0)=z0 at time t = 1
##
## Data: z (500 obs.); Bandwidth 'bw' = 0.02553
##
## z f(z)
## Min. :0.4249639 Min. :0.000355
## 1st Qu.:0.6016392 1st Qu.:0.135101
## Median :0.7783146 Median :0.799373
## Mean :0.7783146 Mean :1.413634
## 3rd Qu.:0.9549899 3rd Qu.:2.807504
## Max. :1.1316652 Max. :3.909531
plot(den, main="Marginal Density")
For Joint density for \((X_t,Y_t,Z_t)\) (for more details, see package sm or ks.)
denJ <- dsde3d(mod3d,pdf="J")
denJ
plot(denJ,display="persp")
mtext("Contour surface of kernel density estimate")with initial conditions \((X_{0},Y_{0},Z_{0})=(1,1,1)\), by specifying the drift and diffusion coefficients of three processes \(X_{t}\), \(Y_{t}\) and \(Z_{t}\) as R expressions which depends on the three state variables (x,y,z) and time variable t, with integration step size Dt=0.0001.
K = 4; s = 1; sigma = 0.2
fx <- expression( (-K*x/sqrt(x^2+y^2+z^2)) , (-K*y/sqrt(x^2+y^2+z^2)) , (-K*z/sqrt(x^2+y^2+z^2)) )
gx <- rep(expression(sigma),3)
mod3d <- snssde3d(drift=fx,diffusion=gx,N=10000,x0=c(x=1,y=1,z=1))
mod3d## Ito Sde 3D:
## | dX(t) = (-K * X(t)/sqrt(X(t)^2 + Y(t)^2 + Z(t)^2)) * dt + sigma * dW1(t)
## | dY(t) = (-K * Y(t)/sqrt(X(t)^2 + Y(t)^2 + Z(t)^2)) * dt + sigma * dW2(t)
## | dZ(t) = (-K * Z(t)/sqrt(X(t)^2 + Y(t)^2 + Z(t)^2)) * dt + sigma * dW3(t)
## Method:
## | Euler scheme of order 0.5
## Summary:
## | Size of process | N = 10001.
## | Number of simulation | M = 1.
## | Initial values | (x0,y0,z0) = (1,1,1).
## | Time of process | t in [0,1].
## | Discretization | Dt = 1e-04.
The results of simulation are shown:
plot3D(mod3d,display="persp",col="blue")
run by calling the function snssde3d() to produce a simulation of the solution, with \(\mu = 1\) and \(\sigma = 1\).
fx <- expression(y,0,0)
gx <- expression(z,1,1)
modtra <- snssde3d(drift=fx,diffusion=gx,M=500)
modtra## Ito Sde 3D:
## | dX(t) = Y(t) * dt + Z(t) * dW1(t)
## | dY(t) = 0 * dt + 1 * dW2(t)
## | dZ(t) = 0 * dt + 1 * dW3(t)
## Method:
## | Euler scheme of order 0.5
## Summary:
## | Size of process | N = 1001.
## | Number of simulation | M = 500.
## | Initial values | (x0,y0,z0) = (0,0,0).
## | Time of process | t in [0,1].
## | Discretization | Dt = 0.001.
summary(modtra)##
## Monte-Carlo Statistics for (X(t),Y(t),Z(t)) at time t = 1
## X Y Z
## Mean -0.00637 0.01584 0.01956
## Variance 0.82162 1.06164 0.98480
## Median 0.03692 0.08289 0.05760
## First quartile -0.54098 -0.66921 -0.60428
## Third quartile 0.56493 0.67993 0.64193
## Skewness -0.35442 -0.22477 -0.10933
## Kurtosis 4.02014 3.02166 3.15964
## Moment of order 3 -0.26395 -0.24587 -0.10684
## Moment of order 4 2.71381 3.40564 3.06433
## Moment of order 5 -3.13266 -2.24807 -1.09098
## Int.conf Inf (95%) -1.96767 -2.10310 -2.06506
## Int.conf Sup (95%) 1.72106 1.94134 2.03285
the following code produces the result in Figure 9.
plot(modtra$X,plot.type="single",ylab="X")
lines(time(modtra),mean(modtra)$X,col=2,lwd=2)
lines(time(modtra),bconfint(modtra,level=0.95)$X[,1],col=4,lwd=2)
lines(time(modtra),bconfint(modtra,level=0.95)$X[,2],col=4,lwd=2)
legend("topleft",c("mean path",paste("bound of",95,"% confidence")),col=c(2,4),lwd=2,cex=0.8)
The histogram and kernel density of \(X_t\) at time \(t=1\) are reported using dsde3d() function, see e.g. Figure 10.
den <- dsde3d(modtra,at=1)
den$resx## NULL
MASS::truehist(den$ech$x,xlab = expression(X[t==1]));box()
lines(den$resx,col="red",lwd=2)
legend("topleft",c("Distribution histogram","Kernel Density"),inset =.01,pch=c(15,NA),lty=c(NA,1),col=c("cyan","red"),lwd=2,cex=0.8)
The equivalently of \(X_{t}^{\text{mod1}}\) the following Stratonovich SDE: \(dX_{t} = \theta X_{t} \circ dW_{t}\).↩