Monte-Carlo Simulation and Analysis of Stochastic Differential Equations

A.C. Guidoum1 and K. Boukhetala2

2017-11-28

snssde1d()

Assume that we want to describe the following SDE:

Itô form3:

\[\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:

R> theta = 0.5
R> f <- expression( (0.5*theta^2*x) )
R> g <- expression( theta*x )
R> mod1 <- snssde1d(drift=f,diffusion=g,x0=10,M=10000,type="ito") # Using Itô
R> mod2 <- snssde1d(drift=f,diffusion=g,x0=10,M=10000,type="str") # Using Stratonovich 
R> mod1
Itô Sde 1D:
    | dX(t) = (0.5 * theta^2 * X(t)) * dt + theta * X(t) * dW(t)
Method:
    | Euler scheme with order 0.5
Summary:
    | Size of process   | N  = 1001.
    | Number of simulation  | M  = 10000.
    | Initial value     | x0 = 10.
    | Time of process   | t in [0,1].
    | Discretization    | Dt = 0.001.
R> mod2
Stratonovich Sde 1D:
    | dX(t) = (0.5 * theta^2 * X(t)) * dt + theta * X(t) o dW(t)
Method:
    | Euler scheme with order 0.5
Summary:
    | Size of process   | N  = 1001.
    | Number of simulation  | M  = 10000.
    | 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\):

The following statistical measures (S3 method) for class snssde1d() can be approximated for the \(X_{t}\) process at any time \(t\), for example at=1:

R> s = 1
R> mean(mod1, at = s)
[1] 11.363
R> moment(mod1, at = s , center = TRUE , order = 2) ## variance
[1] 37.589
R> Median(mod1, at = s)
[1] 9.9825
R> Mode(mod1, at =s)
[1] 7.9492
R> quantile(mod1 , at = s)
     0%     25%     50%     75%    100% 
 1.6273  7.0831  9.9825 14.0993 75.2783 
R> kurtosis(mod1 , at = s)
[1] 8.6708
R> skewness(mod1 , at = s)
[1] 1.7559
R> cv(mod1 , at = s )
[1] 0.5396
R> min(mod1 , at = s)
[1] 1.6273
R> max(mod1 , at = s)
[1] 75.278
R> moment(mod1, at = s , center= TRUE , order = 4)
[1] 12254
R> moment(mod1, at = s , center= FALSE , order = 4)
[1] 76435

The summary of the results of mod1 and mod2 at time \(t=1\) of class snssde1d() is given by:

R> summary(mod1, at = 1)

    Monte-Carlo Statistics for X(t) at time t = 1
                              
Mean                   11.3626
Variance               37.5926
Median                  9.9825
Mode                    7.9492
First quartile          7.0831
Third quartile         14.0993
Minimum                 1.6273
Maximum                75.2783
Skewness                1.7559
Kurtosis                8.6708
Coef-variation          0.5396
3th-order moment      404.7246
4th-order moment    12253.5885
5th-order moment   383645.0793
6th-order moment 15583336.9981
R> summary(mod2, at = 1)

    Monte-Carlo Statistics for X(t) at time t = 1
                              
Mean                  10.03733
Variance              27.63061
Median                 8.89745
Mode                   7.46607
First quartile         6.39703
Third quartile        12.44638
Minimum                1.42120
Maximum               49.15364
Skewness               1.57450
Kurtosis               7.03141
Coef-variation         0.52369
3th-order moment     228.68076
4th-order moment    5368.13144
5th-order moment  113997.58809
6th-order moment 2980394.32500

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\).

R> x1 <- rsde1d(object = mod1, at = 1)  # X(t=1) | X(0)=x0 (Itô SDE)
R> x2 <- rsde1d(object = mod2, at = 1)  # X(t=1) | X(0)=x0 (Stratonovich SDE)
R> head(x1,n=10)
 [1] 17.1014 20.0566 10.2885 14.1245 19.7357  9.8151 11.0635  8.9138
 [9] 17.8666  9.1985
R> head(x2,n=10)
 [1]  8.5531  6.7120  4.2750  4.9903  6.3169 15.3767 14.0683  7.8177
 [9]  5.9116  5.6135
R> summary(data.frame(x1,x2))
       x1              x2       
 Min.   : 1.63   Min.   : 1.42  
 1st Qu.: 7.08   1st Qu.: 6.40  
 Median : 9.98   Median : 8.90  
 Mean   :11.36   Mean   :10.04  
 3rd Qu.:14.10   3rd Qu.:12.45  
 Max.   :75.28   Max.   :49.15  

The function dsde1d() can be used to show the Approximate transitional density for \(X_{t}|X_{0}\) at time \(t-s=1\) with log-normal curves:

R> mu1 = log(10); sigma1= sqrt(theta^2)  # log mean and log variance for mod1 
R> mu2 = log(10)-0.5*theta^2 ; sigma2 = sqrt(theta^2) # log mean and log variance for mod2
R> AppdensI <- dsde1d(mod1, at = 1)
R> AppdensS <- dsde1d(mod2, at = 1)
R> plot(AppdensI , dens = function(x) dlnorm(x,meanlog=mu1,sdlog = sigma1))
R> 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):

R> ## Itô
R> plot(mod1,ylab=expression(X^mod1))
R> lines(time(mod1),apply(mod1$X,1,mean),col=2,lwd=2)
R> lines(time(mod1),apply(mod1$X,1,bconfint,level=0.95)[1,],col=4,lwd=2)
R> lines(time(mod1),apply(mod1$X,1,bconfint,level=0.95)[2,],col=4,lwd=2)
R> legend("topleft",c("mean path",paste("bound of", 95,"% confidence")),inset = .01,col=c(2,4),lwd=2,cex=0.8)
R> ## Stratonovich
R> plot(mod2,ylab=expression(X^mod2))
R> lines(time(mod2),apply(mod2$X,1,mean),col=2,lwd=2)
R> lines(time(mod2),apply(mod2$X,1,bconfint,level=0.95)[1,],col=4,lwd=2)
R> lines(time(mod2),apply(mod2$X,1,bconfint,level=0.95)[2,],col=4,lwd=2)
R> legend("topleft",c("mean path",paste("bound of",95,"% confidence")),col=c(2,4),inset =.01,lwd=2,cex=0.8)
mod1: Itô and mod2: Stratonovich

mod1: Itô and mod2: Stratonovich

Return to snssde1d()

snssde2d()

The following \(2\)-dimensional SDE’s with a vector of drift and a diagonal matrix of diffusion coefficients:

Itô 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:

Ornstein-Uhlenbeck process and its integral

The Ornstein-Uhlenbeck (OU) process has a long history in physics. Introduced in essence by Langevin in his famous 1908 paper on Brownian motion, the process received a more thorough mathematical examination several decades later by Uhlenbeck and Ornstein (1930). The OU process is understood here to be the univariate continuous Markov process \(X_t\). In mathematical terms, the equation is written as an Ito equation: \[\begin{equation}\label{eq016} dX_t = -\frac{1}{\mu} X_t dt + \sqrt{\sigma} dW_t,\quad X_{0}=x_{0} \end{equation}\] In these equations, \(\mu\) and \(\sigma\) are positive constants called, respectively, the relaxation time and the diffusion constant. The time integral of the OU process \(X_t\) (or indeed of any process \(X_t\)) is defined to be the process \(Y_t\) that satisfies: \[\begin{equation}\label{eq017} Y_{t} = Y_{0}+\int X_{t} dt \Leftrightarrow dY_t = X_{t} dt ,\quad Y_{0}=y_{0} \end{equation}\] \(Y_t\) is not itself a Markov process; however, \(X_t\) and \(Y_t\) together comprise a bivariate continuous Markov process. We wish to find the solutions \(X_t\) and \(Y_t\) to the coupled time-evolution equations: \[\begin{equation}\label{eq018} \begin{cases} dX_t = -\frac{1}{\mu} X_t dt + \sqrt{\sigma} dW_t\\ dY_t = X_{t} dt \end{cases} \end{equation}\]

We simulate a flow of \(10000\) trajectories of \((X_{t},Y_{t})\), with integration step size \(\Delta t = 0.01\), and using second Milstein method.

R> x0=5;y0=0
R> mu=3;sigma=0.5
R> fx <- expression(-(x/mu),x)  
R> gx <- expression(sqrt(sigma),0)
R> mod2d <- snssde2d(drift=fx,diffusion=gx,Dt=0.01,M=10000,x0=c(x0,y0),method="smilstein")
R> mod2d
Itô Sde 2D:
    | dX(t) = -(X(t)/mu) * dt + sqrt(sigma) * dW1(t)
    | dY(t) = X(t) * dt + 0 * dW2(t)
Method:
    | Second-order Milstein scheme
Summary:
    | Size of process   | N  = 1001.
    | Number of simulation  | M  = 10000.
    | Initial values    | (x0,y0) = (5,0).
    | Time of process   | t in [0,10].
    | Discretization    | Dt = 0.01.

The following statistical measures (S3 method) for class snssde2d() can be approximated for the \((X_{t},Y_{t})\) process at any time \(t\), for example at=5:

R> s = 5
R> mean(mod2d, at = s)
[1]  0.93766 12.14376
R> moment(mod2d, at = s , center = TRUE , order = 2) ## variance
[1] 0.72106 7.16658
R> Median(mod2d, at = s)
[1]  0.93424 12.17686
R> Mode(mod2d, at = s)
[1]  0.9398 12.3973
R> quantile(mod2d , at = s)
$x
      0%      25%      50%      75%     100% 
-2.34076  0.36344  0.93424  1.51820  4.23466 

$y
     0%     25%     50%     75%    100% 
 2.0998 10.3602 12.1769 13.9129 21.6283 
R> kurtosis(mod2d , at = s)
[1] 3.0216 3.0434
R> skewness(mod2d , at = s)
[1]  0.008700 -0.013332
R> cv(mod2d , at = s )
[1] 0.90565 0.22046
R> min(mod2d , at = s)
[1] -2.3408  2.0998
R> max(mod2d , at = s)
[1]  4.2347 21.6283
R> moment(mod2d, at = s , center= TRUE , order = 4)
[1]   1.5713 156.3417
R> moment(mod2d, at = s , center= FALSE , order = 4)
[1]     6.168 28232.729

The summary of the results of mod2d at time \(t=5\) of class snssde2d() is given by:

R> summary(mod2d, at = s)

    Monte-Carlo Statistics for (X(t),Y(t)) at time t = 5
                        X          Y
Mean              0.93766   12.14376
Variance          0.72113    7.16730
Median            0.93424   12.17686
Mode              0.93980   12.39728
First quartile    0.36344   10.36019
Third quartile    1.51820   13.91290
Minimum          -2.34076    2.09978
Maximum           4.23466   21.62835
Skewness          0.00870   -0.01333
Kurtosis          3.02159    3.04343
Coef-variation    0.90565    0.22046
3th-order moment  0.00533   -0.25582
4th-order moment  1.57131  156.34170
5th-order moment  0.05265  -25.90355
6th-order moment  5.65976 5559.77008

For plotting (back in time) using the command plot, the results of the simulation are shown in Figure 3.

R> plot(mod2d)
Ornstein-Uhlenbeck process and its integral

Ornstein-Uhlenbeck process and its integral

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\).

R> out <- rsde2d(object = mod2d, at = 10) 
R> head(out,n=10)
           x      y
1  -0.125743 20.720
2   0.198664 13.106
3  -1.249980 17.829
4  -0.360778 18.707
5   0.278290 11.122
6  -0.692971 22.814
7   0.752602 23.364
8  -0.087329 14.902
9   0.724739 12.780
10  0.398091 11.419
R> summary(out)
       x                y        
 Min.   :-3.524   Min.   :-5.69  
 1st Qu.:-0.375   1st Qu.:11.11  
 Median : 0.197   Median :14.49  
 Mean   : 0.199   Mean   :14.48  
 3rd Qu.: 0.783   3rd Qu.:17.90  
 Max.   : 3.680   Max.   :31.74  
R> cov(out)
        x       y
x 0.75136  2.0579
y 2.05786 25.5791

Figure 4, show simulation results for moments of system :

R> mx <- apply(mod2d$X,1,mean)
R> my <- apply(mod2d$Y,1,mean)
R> Sx <- apply(mod2d$X,1,var)
R> Sy <- apply(mod2d$Y,1,var)
R> Cxy <- sapply(1:1001,function(i) cov(mod2d$X[i,],mod2d$Y[i,]))
R> out_b <- data.frame(mx,my,Sx,Sy,Cxy)
R> matplot(time(mod2d), out_b, type = "l", xlab = "time", ylab = "",col=2:6,lwd=2,lty=2:6,las=1)
R> legend("topleft",c(expression(hat(E)(X[t]),hat(E)(Y[t]),hat(Var)(X[t]),hat(Var)(Y[t]),hat(Cov)(X[t],Y[t]))),inset = .05,col=2:6,lty=2:6,lwd=2,cex=0.9)

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 5: the marginal density of \(X_t\) and \(Y_t\) at time \(t=10\).

R> denM <- dsde2d(mod2d,pdf="M",at =10)
R> denM

 Marginal density of X(t-t0)|X(t0)=5 at time t = 10

Data: x (10000 obs.);   Bandwidth 'bw' = 0.1232

       x                f(x)        
 Min.   :-3.8940   Min.   :0.00000  
 1st Qu.:-1.9081   1st Qu.:0.00153  
 Median : 0.0778   Median :0.03565  
 Mean   : 0.0778   Mean   :0.12576  
 3rd Qu.: 2.0637   3rd Qu.:0.24133  
 Max.   : 4.0496   Max.   :0.47596  

 Marginal density of Y(t-t0)|Y(t0)=0 at time t = 10

Data: y (10000 obs.);   Bandwidth 'bw' = 0.7214

       y               f(y)         
 Min.   :-7.853   Min.   :0.000001  
 1st Qu.: 2.586   1st Qu.:0.000739  
 Median :13.026   Median :0.009363  
 Mean   :13.026   Mean   :0.023925  
 3rd Qu.:23.465   3rd Qu.:0.047221  
 Max.   :33.904   Max.   :0.078524  
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 system \((X_{t},Y_{t})\) at time t=10.

R> denJ <- dsde2d(mod2d, pdf="J", n=100,at =10)
R> denJ

 Joint density of (X(t-t0),Y(t-t0)|X(t0)=5,Y(t0)=0) at time t = 10

Data: (x,y) (2 x 10000 obs.);

       x                 y              f(x,y)        
 Min.   :-3.5243   Min.   :-5.689   Min.   :0.000000  
 1st Qu.:-1.7233   1st Qu.: 3.668   1st Qu.:0.000000  
 Median : 0.0778   Median :13.026   Median :0.000194  
 Mean   : 0.0778   Mean   :13.026   Mean   :0.003634  
 3rd Qu.: 1.8788   3rd Qu.:22.383   3rd Qu.:0.002543  
 Max.   : 3.6799   Max.   :31.740   Max.   :0.042365  
R> plot(denJ,display="contour",main="Bivariate Transition Density at time t=10")
R> plot(denJ,display="image",main="Bivariate Transition Density at time t=10")

A \(3\)D plot of the transition density at \(t=10\) obtained with:

R> plot(denJ,main="Bivariate Transition Density at time t=10")

We approximate the bivariate transition density over the set transition horizons \(t\in [1,10]\) by \(\Delta t = 0.005\) using the code:

R> for (i in seq(1,10,by=0.005)){ 
+ plot(dsde2d(mod2d, at = i,n=100),display="contour",main=paste0('Transition Density \n t = ',i))
+ }

Return to snssde2d()

The stochastic Van-der-Pol equation

The Van der Pol (1922) equation is an ordinary differential equation that can be derived from the Rayleigh differential equation by differentiating and setting \(\dot{x}=y\), see Naess and Hegstad (1994); Leung (1995) and for more complex dynamics in Van-der-Pol equation see Jing et al. (2006). It is an equation describing self-sustaining oscillations in which energy is fed into small oscillations and removed from large oscillations. This equation arises in the study of circuits containing vacuum tubes and is given by: \[\begin{equation}\label{eq:12} \ddot{X}-\mu (1-X^{2}) \dot{X} + X = 0 \end{equation}\] where \(x\) is the position coordinate (which is a function of the time \(t\)), and \(\mu\) is a scalar parameter indicating the nonlinearity and the strength of the damping, to simulate the deterministic equation see Grayling (2014) for more details. Consider stochastic perturbations of the Van-der-Pol equation, and random excitation force of such systems by White noise \(\xi_{t}\), with delta-type correlation function \(\text{E}(\xi_{t}\xi_{t+h})=2\sigma \delta (h)\) \[\begin{equation}\label{eq:13} \ddot{X}-\mu (1-X^{2}) \dot{X} + X = \xi_{t}, \end{equation}\] where \(\mu > 0\) . It’s solution cannot be obtained in terms of elementary functions, even in the phase plane. The White noise \(\xi_{t}\) is formally derivative of the Wiener process \(W_{t}\). The representation of a system of two first order equations follows the same idea as in the deterministic case by letting \(\dot{x}=y\), from physical equation we get the above system: \[\begin{equation}\label{eq:14} \begin{cases} \dot{X} = Y \\ \dot{Y} = \mu \left(1-X^{2}\right) Y - X + \xi_{t} \end{cases} \end{equation}\] The system can mathematically explain by a Stratonovitch equations: \[\begin{equation}\label{eq:15} \begin{cases} dX_{t} = Y_{t} dt \\ dY_{t} = \left(\mu (1-X^{2}_{t}) Y_{t} - X_{t}\right) dt + 2 \sigma \circ dW_{2,t} \end{cases} \end{equation}\]

Implemente in R as follows, with integration step size \(\Delta t = 0.01\) and using stochastic Runge-Kutta methods 1-stage.

R> mu = 4; sigma=0.1
R> fx <- expression( y ,  (mu*( 1-x^2 )* y - x)) 
R> gx <- expression( 0 ,2*sigma)
R> mod2d <- snssde2d(drift=fx,diffusion=gx,N=10000,Dt=0.01,type="str",method="rk1")
R> 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 with 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.

For plotting (back in time) using the command plot, and plot2d in plane the results of the simulation are shown in Figure 8.

R> plot2d(mod2d) ## in plane (O,X,Y)
R> plot(mod2d)   ## back in time

Return to snssde2d()

snssde3d()

The following \(3\)-dimensional SDE’s with a vector of drift and a diagonal matrix of diffusion coefficients:

Itô 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:

Basic example

Assume that we want to describe the following SDE’s (3D) in Itô form: \[\begin{equation}\label{eq0166} \begin{cases} dX_t = 4 (-1-X_{t}) Y_{t} dt + 0.2 dW_{1,t}\\ dY_t = 4 (1-Y_{t}) X_{t} dt + 0.2 dW_{2,t}\\ dZ_t = 4 (1-Z_{t}) Y_{t} dt + 0.2 dW_{3,t} \end{cases} \end{equation}\]

We simulate a flow of \(10000\) 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> mod3d <- snssde3d(x0=c(x=2,y=-2,z=-2),drift=fx,diffusion=gx,M=10000)
R> mod3d
Itô 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.
    | Number of simulation  | M  = 10000.
    | Initial values    | (x0,y0,z0) = (2,-2,-2).
    | Time of process   | t in [0,1].
    | Discretization    | Dt = 0.001.

The following statistical measures (S3 method) for class snssde3d() can be approximated for the \((X_{t},Y_{t},Z_{t})\) process at any time \(t\), for example at=1:

R> s = 1
R> mean(mod3d, at = s)
[1] -0.79605  0.87793  0.79497
R> moment(mod3d, at = s , center = TRUE , order = 2) ## variance
[1] 0.009681 0.103446 0.010391
R> Median(mod3d, at = s)
[1] -0.80262  0.85968  0.80243
R> Mode(mod3d, at = s)
[1] -0.81559  0.83478  0.81109
R> quantile(mod3d , at = s)
$x
      0%      25%      50%      75%     100% 
-1.10478 -0.86437 -0.80262 -0.73359 -0.35522 

$y
      0%      25%      50%      75%     100% 
-0.11828  0.65401  0.85968  1.07801  2.31093 

$z
     0%     25%     50%     75%    100% 
0.25421 0.73398 0.80243 0.86626 1.08362 
R> kurtosis(mod3d , at = s)
[1] 3.1122 3.2932 3.5399
R> skewness(mod3d , at = s)
[1]  0.34678  0.38076 -0.53128
R> cv(mod3d , at = s )
[1] -0.12361  0.36637  0.12823
R> min(mod3d , at = s)
[1] -1.10478 -0.11828  0.25421
R> max(mod3d , at = s)
[1] -0.35522  2.31093  1.08362
R> moment(mod3d, at = s , center= TRUE , order = 4)
[1] 0.00029174 0.03524778 0.00038230
R> moment(mod3d, at = s , center= FALSE , order = 4)
[1] 0.43762 1.15222 0.43740

The summary of the results of mod3d at time \(t=1\) of class snssde3d() is given by:

R> summary(mod3d, at = s)

  Monte-Carlo Statistics for (X(t),Y(t),Z(t)) at time t = 1
                        X        Y        Z
Mean             -0.79605  0.87793  0.79497
Variance          0.00968  0.10346  0.01039
Median           -0.80262  0.85968  0.80243
Mode             -0.81559  0.83478  0.81109
First quartile   -0.86437  0.65401  0.73398
Third quartile   -0.73359  1.07801  0.86626
Minimum          -1.10478 -0.11828  0.25421
Maximum          -0.35522  2.31093  1.08362
Skewness          0.34678  0.38076 -0.53128
Kurtosis          3.11222  3.29316  3.53988
Coef-variation   -0.12361  0.36637  0.12823
3th-order moment  0.00033  0.01267 -0.00056
4th-order moment  0.00029  0.03525  0.00038
5th-order moment  0.00003  0.01353 -0.00006
6th-order moment  0.00002  0.02214  0.00003

For plotting (back in time) using the command plot, and plot3D in space the results of the simulation are shown in Figure 9.

R> plot(mod3d,union = TRUE)         ## back in time
R> plot3D(mod3d,display="persp")    ## in space (O,X,Y,Z)
3D SDEs

3D SDEs

Hence we can just make use of the rsde3d() function to build our random number for \((X_{t},Y_{t},Z_{t})\) at time \(t = 1\).

R> out <- rsde3d(object = mod3d, at = s) 
R> head(out,n=10)
          x       y       z
1  -0.86795 1.03414 0.77641
2  -0.80683 0.78026 0.78566
3  -0.86290 1.49416 0.82761
4  -0.70053 0.45556 0.74760
5  -0.74101 0.82970 0.68647
6  -0.98061 0.42199 0.66585
7  -0.79099 0.56426 0.76548
8  -0.84643 0.78707 0.71232
9  -0.69390 0.44960 0.64746
10 -0.61427 0.68495 0.56993
R> summary(out)
       x                y                z        
 Min.   :-1.105   Min.   :-0.118   Min.   :0.254  
 1st Qu.:-0.864   1st Qu.: 0.654   1st Qu.:0.734  
 Median :-0.803   Median : 0.860   Median :0.802  
 Mean   :-0.796   Mean   : 0.878   Mean   :0.795  
 3rd Qu.:-0.734   3rd Qu.: 1.078   3rd Qu.:0.866  
 Max.   :-0.355   Max.   : 2.311   Max.   :1.084  
R> cov(out)
          x         y         z
x  0.009682 -0.018004 -0.004178
y -0.018004  0.103457  0.019585
z -0.004178  0.019585  0.010392

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 10.

R> den <- dsde3d(mod3d,pdf="M",at =1)
R> den

 Marginal density of X(t-t0)|X(t0)=2 at time t = 1

Data: x (10000 obs.);   Bandwidth 'bw' = 0.01392

       x                 f(x)       
 Min.   :-1.14654   Min.   :0.0000  
 1st Qu.:-0.93827   1st Qu.:0.0173  
 Median :-0.73000   Median :0.4388  
 Mean   :-0.73000   Mean   :1.1992  
 3rd Qu.:-0.52173   3rd Qu.:2.2751  
 Max.   :-0.31346   Max.   :4.1212  

 Marginal density of Y(t-t0)|Y(t0)=-2 at time t = 1

Data: y (10000 obs.);   Bandwidth 'bw' = 0.04513

       y                 f(y)        
 Min.   :-0.25368   Min.   :0.00001  
 1st Qu.: 0.42132   1st Qu.:0.01324  
 Median : 1.09633   Median :0.13243  
 Mean   : 1.09633   Mean   :0.37001  
 3rd Qu.: 1.77133   3rd Qu.:0.71624  
 Max.   : 2.44633   Max.   :1.28773  

 Marginal density of Z(t-t0)|Z(t0)=-2 at time t = 1

Data: z (10000 obs.);   Bandwidth 'bw' = 0.01408

       z                f(z)       
 Min.   :0.21196   Min.   :0.0000  
 1st Qu.:0.44044   1st Qu.:0.0162  
 Median :0.66891   Median :0.2836  
 Mean   :0.66891   Mean   :1.0931  
 3rd Qu.:0.89739   3rd Qu.:1.9952  
 Max.   :1.12586   Max.   :4.0831  
R> plot(den, main="Marginal Density") 

For an approximate joint transition density for \((X_t,Y_t,Z_t)\) (for more details, see package sm or ks.)

R> denJ <- dsde3d(mod3d,pdf="J")
R> plot(denJ,display="rgl")

Return to snssde3d()

Attractive model for 3D diffusion processes

If we assume that \(U_w( x , y , z , t )\), \(V_w( x , y , z , t )\) and \(S_w( x , y , z , t )\) are neglected and the dispersion coefficient \(D( x , y , z )\) is constant. A system becomes (see Boukhetala,1996): \[\begin{eqnarray}\label{eq19} % \nonumber to remove numbering (before each equation) \begin{cases} dX_t = \left(\frac{-K X_{t}}{\sqrt{X^{2}_{t} + Y^{2}_{t} + Z^{2}_{t}}}\right) dt + \sigma dW_{1,t} \nonumber\\ dY_t = \left(\frac{-K Y_{t}}{\sqrt{X^{2}_{t} + Y^{2}_{t} + Z^{2}_{t}}}\right) dt + \sigma dW_{2,t} \\ dZ_t = \left(\frac{-K Z_{t}}{\sqrt{X^{2}_{t} + Y^{2}_{t} + Z^{2}_{t}}}\right) dt + \sigma dW_{3,t} \nonumber \end{cases} \end{eqnarray}\]

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.

R> K = 4; s = 1; sigma = 0.2
R> 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)) ) 
R> gx <- rep(expression(sigma),3)
R> mod3d <- snssde3d(drift=fx,diffusion=gx,N=10000,x0=c(x=1,y=1,z=1))
R> mod3d
Itô 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 with 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 = 0.0001.

The results of simulation are shown:

R> plot3D(mod3d,display="persp",col="blue")

Return to snssde3d()

Transformation of an SDE one-dimensional

Next is an example of one-dimensional SDE driven by three independent Brownian motions (\(W_{1,t}\),\(W_{2,t}\),\(W_{3,t}\)), as follows: \[\begin{equation}\label{eq20} dX_{t} = \mu W_{1,t} dt + \sigma W_{2,t} \circ dW_{3,t} \end{equation}\] To simulate the solution of the process \(X_t\), we make a transformation to a system of three equations as follows: \[\begin{eqnarray}\label{eq21} \begin{cases} % \nonumber to remove numbering (before each equation) dX_t = \mu Y_{t} dt + \sigma Z_{t} \circ dW_{3,t} \nonumber\\ dY_t = dW_{1,t} \\ dZ_t = dW_{2,t} \nonumber \end{cases} \end{eqnarray}\]

run by calling the function snssde3d() to produce a simulation of the solution, with \(\mu = 1\) and \(\sigma = 1\).

R> fx <- expression(y,0,0) 
R> gx <- expression(z,1,1)
R> modtra <- snssde3d(drift=fx,diffusion=gx,M=10000,type="str")
R> modtra
Stratonovich Sde 3D:
    | dX(t) = Y(t) * dt + Z(t) o dW1(t)
    | dY(t) = 0 * dt + 1 o dW2(t)
    | dZ(t) = 0 * dt + 1 o dW3(t)
Method:
    | Euler scheme with order 0.5
Summary:
    | Size of process   | N  = 1001.
    | Number of simulation  | M  = 10000.
    | Initial values    | (x0,y0,z0) = (0,0,0).
    | Time of process   | t in [0,1].
    | Discretization    | Dt = 0.001.
R> summary(modtra)

  Monte-Carlo Statistics for (X(t),Y(t),Z(t)) at time t = 1
                         X         Y         Z
Mean              -0.01044  -0.01647  -0.01930
Variance           0.78983   1.01403   0.96885
Median            -0.01837  -0.01114  -0.01814
Mode              -0.13435   0.02920  -0.01126
First quartile    -0.55623  -0.68685  -0.68549
Third quartile     0.53909   0.65842   0.65707
Minimum           -4.98033  -4.31390  -3.40609
Maximum            4.33912   3.99552   3.45899
Skewness           0.00603  -0.02683  -0.02529
Kurtosis           4.07151   3.07348   2.96449
Coef-variation   -85.13068 -61.15766 -51.00351
3th-order moment   0.00424  -0.02740  -0.02412
4th-order moment   2.53991   3.16031   2.78267
5th-order moment  -0.01685  -0.37694  -0.15523
6th-order moment  17.68542  16.73578  13.12140

the following code produces the result in Figure 12.

R> plot(modtra$X,plot.type="single",ylab=expression(X[t]))
R> lines(time(modtra),apply(modtra$X,1,mean),col=2,lwd=2)
R> legend("topleft",c("mean path"),col=2,lwd=2,cex=0.8)
Simulation of X_t

Simulation of \(X_t\)

The histogram and kernel density of \(X_t\) at time \(t=1\) are reported using dsde3d() function, see e.g. Figure 13.

R> den <- dsde3d(modtra,pdf="Marginal",at=1)
R> den$resx

Call:
    density.default(x = x, na.rm = TRUE)

Data: x (10000 obs.);   Bandwidth 'bw' = 0.1166

       x                y          
 Min.   :-5.330   Min.   :0.00000  
 1st Qu.:-2.825   1st Qu.:0.00122  
 Median :-0.321   Median :0.01272  
 Mean   :-0.321   Mean   :0.09971  
 3rd Qu.: 2.184   3rd Qu.:0.14684  
 Max.   : 4.689   Max.   :0.49431  
R> MASS::truehist(den$ech$x,xlab = expression(X[t==1]));box()
R> lines(den$resx,col="red",lwd=2)
R> 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)

Figure 14 and 15, show approximation results for \(m_{1}(t)= \text{E}(X_{t})\), \(S_{1}(t)=\text{V}(X_{t})\) and \(C(s,t)=\text{Cov}(X_{s},X_{t})\):

R> m1  <- apply(modtra$X,1,mean) ## m1(t)
R> S1  <- apply(modtra$X,1,var)  ## s1(t)
R> out_a <- data.frame(m1,S1)
R> matplot(time(modtra), out_a, type = "l", xlab = "time", ylab = "", col=2:3,lwd=2,lty=2:3,las=1)
R> legend("topleft",c(expression(m[1](t),S[1](t))),inset = .09,col=2:3,lty=2:3,lwd=2,cex=0.9)

R> color.palette=colorRampPalette(c('white','green','blue','red'))
R> filled.contour(time(modtra), time(modtra), cov(t(modtra$X)), color.palette=color.palette,plot.title = title(main = expression(paste("Covariance empirique:",cov(X[s],X[t]))),xlab = "time", ylab = "time"),key.title = title(main = ""))

Return to snssde3d()

Further reading

  1. snssdekd() & dsdekd() & rsdekd()- Monte-Carlo Simulation and Analysis of Stochastic Differential Equations.
  2. bridgesdekd() & dsdekd() & rsdekd() - Constructs and Analysis of Bridges Stochastic Differential Equations.
  3. fptsdekd() & dfptsdekd() - Monte-Carlo Simulation and Kernel Density Estimation of First passage time.
  4. MCM.sde() & MEM.sde() - Parallel Monte-Carlo and Moment Equations for SDEs.
  5. fitsde() - Parametric Estimation of 1-D Stochastic Differential Equation.

References

  1. Boukhetala K (1996). Modelling and Simulation of a Dispersion Pollutant with Attractive Centre, volume 3, pp. 245-252. Computer Methods and Water Resources, Computational Mechanics Publications, Boston, USA.

  2. Guidoum AC, Boukhetala K (2017). Performing Parallel Monte Carlo and Moment Equations Methods for Itô and Stratonovich Stochastic Differential Systems: R Package Sim.DiffProc. Preprint submitted to Journal of Statistical Software.

  3. Guidoum AC, Boukhetala K (2017). Sim.DiffProc: Simulation of Diffusion Processes. R package version 4.0, URL https://cran.r-project.org/package=Sim.DiffProc.


  1. 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)

  2. 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)

  3. The equivalently of \(X_{t}^{\text{mod1}}\) the following Stratonovich SDE: \(dX_{t} = \theta X_{t} \circ dW_{t}\).