Monte-Carlo Simulation and Analysis of Stochastic Differential Equations

A.C. Guidoum1 and K. Boukhetala2

2018-09-10

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.294
R> moment(mod1, at = s , center = TRUE , order = 2) ## variance
[1] 35.496
R> Median(mod1, at = s)
[1] 9.9441
R> Mode(mod1, at =s)
[1] 7.8018
R> quantile(mod1 , at = s)
     0%     25%     50%     75%    100% 
 1.6445  7.1602  9.9441 13.8679 59.7271 
R> kurtosis(mod1 , at = s)
[1] 7.2934
R> skewness(mod1 , at = s)
[1] 1.621
R> cv(mod1 , at = s )
[1] 0.52756
R> min(mod1 , at = s)
[1] 1.6445
R> max(mod1 , at = s)
[1] 59.727
R> moment(mod1, at = s , center= TRUE , order = 4)
[1] 9191.3
R> moment(mod1, at = s , center= FALSE , order = 4)
[1] 68116

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.29388
Variance              35.49969
Median                 9.94412
Mode                   7.80175
First quartile         7.16021
Third quartile        13.86789
Minimum                1.64453
Maximum               59.72708
Skewness               1.62104
Kurtosis               7.29338
Coef-variation         0.52756
3th-order moment     342.87069
4th-order moment    9191.32518
5th-order moment  230302.31460
6th-order moment 7122226.23084
R> summary(mod2, at = 1)

    Monte-Carlo Statistics for X(t) at time t = 1
                              
Mean                   9.97015
Variance              27.21940
Median                 8.81041
Mode                   7.08137
First quartile         6.33736
Third quartile        12.32062
Minimum                1.41600
Maximum               61.05984
Skewness               1.69339
Kurtosis               8.48428
Coef-variation         0.52328
3th-order moment     240.47720
4th-order moment    6285.96643
5th-order moment  168367.21002
6th-order moment 5803105.06475

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]  4.5617 15.2902  5.7172  7.8927 21.2781 15.1225  6.7414 14.9453
 [9] 27.0423 14.1436
R> head(x2,n=10)
 [1]  5.4088 13.2887  4.6744  6.3361 11.2292  6.5593 11.6521 10.4583
 [9] 11.3059  6.9635
R> summary(data.frame(x1,x2))
       x1              x2       
 Min.   : 1.64   Min.   : 1.42  
 1st Qu.: 7.16   1st Qu.: 6.34  
 Median : 9.94   Median : 8.81  
 Mean   :11.29   Mean   : 9.97  
 3rd Qu.:13.87   3rd Qu.:12.32  
 Max.   :59.73   Max.   :61.06  

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.94922 12.20826
R> moment(mod2d, at = s , center = TRUE , order = 2) ## variance
[1] 0.72484 7.16594
R> Median(mod2d, at = s)
[1]  0.93207 12.19928
R> Mode(mod2d, at = s)
[1]  0.83543 12.28881
R> quantile(mod2d , at = s)
$x
      0%      25%      50%      75%     100% 
-2.52964  0.38431  0.93207  1.51615  4.24452 

$y
     0%     25%     50%     75%    100% 
 2.0354 10.3969 12.1993 14.0398 21.8081 
R> kurtosis(mod2d , at = s)
[1] 2.9848 2.9321
R> skewness(mod2d , at = s)
[1] 0.033821 0.014161
R> cv(mod2d , at = s )
[1] 0.89696 0.21928
R> min(mod2d , at = s)
[1] -2.5296  2.0354
R> max(mod2d , at = s)
[1]  4.2445 21.8081
R> moment(mod2d, at = s , center= TRUE , order = 4)
[1]   1.5685 150.5964
R> moment(mod2d, at = s , center= FALSE , order = 4)
[1]     6.3782 28785.4283

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.94922   12.20826
Variance          0.72491    7.16666
Median            0.93207   12.19928
Mode              0.83543   12.28881
First quartile    0.38431   10.39693
Third quartile    1.51615   14.03979
Minimum          -2.52964    2.03544
Maximum           4.24452   21.80815
Skewness          0.03382    0.01416
Kurtosis          2.98479    2.93212
Coef-variation    0.89696    0.21928
3th-order moment  0.02087    0.27168
4th-order moment  1.56848  150.59645
5th-order moment  0.10568    4.77462
6th-order moment  5.53409 5174.30078

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.025937 12.5317
2   1.089531 19.3829
3   0.790466  7.0203
4  -0.719261 19.4002
5   0.970252 27.9771
6   0.933042 16.8523
7   1.038908 18.4764
8   0.814775 21.3297
9  -1.107891 12.9421
10  0.591950 11.4679
R> summary(out)
       x                y        
 Min.   :-3.567   Min.   :-6.34  
 1st Qu.:-0.410   1st Qu.:11.04  
 Median : 0.190   Median :14.49  
 Mean   : 0.180   Mean   :14.53  
 3rd Qu.: 0.772   3rd Qu.:17.89  
 Max.   : 3.426   Max.   :30.81  
R> cov(out)
        x       y
x 0.75571  2.1272
y 2.12724 25.9638

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

       x                f(x)        
 Min.   :-3.9387   Min.   :0.00000  
 1st Qu.:-2.0046   1st Qu.:0.00172  
 Median :-0.0705   Median :0.03835  
 Mean   :-0.0705   Mean   :0.12913  
 3rd Qu.: 1.8636   3rd Qu.:0.23996  
 Max.   : 3.7977   Max.   :0.45557  

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

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

       y               f(y)         
 Min.   :-8.516   Min.   :0.000001  
 1st Qu.: 1.861   1st Qu.:0.000821  
 Median :12.238   Median :0.009732  
 Mean   :12.238   Mean   :0.024068  
 3rd Qu.:22.616   3rd Qu.:0.047036  
 Max.   :32.993   Max.   :0.078066  
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.5667   Min.   :-6.3358   Min.   :0.000000  
 1st Qu.:-1.8186   1st Qu.: 2.9512   1st Qu.:0.000000  
 Median :-0.0705   Median :12.2383   Median :0.000206  
 Mean   :-0.0705   Mean   :12.2383   Mean   :0.003772  
 3rd Qu.: 1.6776   3rd Qu.:21.5254   3rd Qu.:0.002818  
 Max.   : 3.4257   Max.   :30.8124   Max.   :0.039369  
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.79530  0.87844  0.79514
R> moment(mod3d, at = s , center = TRUE , order = 2) ## variance
[1] 0.0098268 0.1069379 0.0102464
R> Median(mod3d, at = s)
[1] -0.80159  0.85335  0.80060
R> Mode(mod3d, at = s)
[1] -0.81126  0.81668  0.80312
R> quantile(mod3d , at = s)
$x
      0%      25%      50%      75%     100% 
-1.07356 -0.86605 -0.80159 -0.73201 -0.33451 

$y
       0%       25%       50%       75%      100% 
-0.072332  0.651842  0.853354  1.080896  2.870194 

$z
     0%     25%     50%     75%    100% 
0.29802 0.73177 0.80060 0.86770 1.13628 
R> kurtosis(mod3d , at = s)
[1] 3.2611 3.4961 3.2981
R> skewness(mod3d , at = s)
[1]  0.38751  0.45523 -0.40903
R> cv(mod3d , at = s )
[1] -0.12465  0.37229  0.12731
R> min(mod3d , at = s)
[1] -1.073563 -0.072332  0.298022
R> max(mod3d , at = s)
[1] -0.33451  2.87019  1.13628
R> moment(mod3d, at = s , center= TRUE , order = 4)
[1] 0.00031497 0.03998848 0.00034634
R> moment(mod3d, at = s , center= FALSE , order = 4)
[1] 0.43646 1.18650 0.43761

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.79530  0.87844  0.79514
Variance          0.00983  0.10695  0.01025
Median           -0.80159  0.85335  0.80060
Mode             -0.81126  0.81668  0.80312
First quartile   -0.86605  0.65184  0.73177
Third quartile   -0.73201  1.08090  0.86770
Minimum          -1.07356 -0.07233  0.29802
Maximum          -0.33451  2.87019  1.13628
Skewness          0.38751  0.45523 -0.40903
Kurtosis          3.26113  3.49611  3.29813
Coef-variation   -0.12465  0.37229  0.12731
3th-order moment  0.00038  0.01592 -0.00042
4th-order moment  0.00031  0.03999  0.00035
5th-order moment  0.00004  0.02089 -0.00005
6th-order moment  0.00002  0.03489  0.00002

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.91235 1.11431 0.80904
2  -0.82417 0.78085 0.71303
3  -0.69025 0.87803 0.79209
4  -0.76321 0.49149 0.65043
5  -0.76783 0.48448 0.62280
6  -0.67763 0.39713 0.60449
7  -0.85963 1.05202 0.77641
8  -0.74741 1.04911 0.78306
9  -0.79319 0.88463 0.83531
10 -0.86931 0.87753 0.77403
R> summary(out)
       x                y                 z        
 Min.   :-1.074   Min.   :-0.0723   Min.   :0.298  
 1st Qu.:-0.866   1st Qu.: 0.6518   1st Qu.:0.732  
 Median :-0.802   Median : 0.8534   Median :0.801  
 Mean   :-0.795   Mean   : 0.8784   Mean   :0.795  
 3rd Qu.:-0.732   3rd Qu.: 1.0809   3rd Qu.:0.868  
 Max.   :-0.335   Max.   : 2.8702   Max.   :1.136  
R> cov(out)
           x         y          z
x  0.0098277 -0.018378 -0.0042096
y -0.0183785  0.106949  0.0197329
z -0.0042096  0.019733  0.0102475

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

       x                 f(x)       
 Min.   :-1.11598   Min.   :0.0000  
 1st Qu.:-0.91001   1st Qu.:0.0434  
 Median :-0.70404   Median :0.4680  
 Mean   :-0.70404   Mean   :1.2126  
 3rd Qu.:-0.49807   3rd Qu.:2.3853  
 Max.   :-0.29209   Max.   :4.0569  

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

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

       y                 f(y)        
 Min.   :-0.20935   Min.   :0.00001  
 1st Qu.: 0.59479   1st Qu.:0.00099  
 Median : 1.39893   Median :0.06041  
 Mean   : 1.39893   Mean   :0.31059  
 3rd Qu.: 2.20307   3rd Qu.:0.55081  
 Max.   : 3.00721   Max.   :1.27133  

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

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

       z                f(z)       
 Min.   :0.25470   Min.   :0.0000  
 1st Qu.:0.48593   1st Qu.:0.0121  
 Median :0.71715   Median :0.2991  
 Mean   :0.71715   Mean   :1.0802  
 3rd Qu.:0.94837   3rd Qu.:2.0534  
 Max.   :1.17959   Max.   :3.9860  
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.00675  0.01140   -0.00362
Variance            0.82947  1.02639    1.01124
Median              0.00356  0.01360   -0.00226
Mode                0.02129  0.11475    0.10404
First quartile     -0.57297 -0.67813   -0.68853
Third quartile      0.57052  0.69228    0.67940
Minimum            -5.37654 -3.82762   -3.53156
Maximum             4.80013  4.04002    4.19328
Skewness           -0.04712  0.00645   -0.00476
Kurtosis            4.25285  2.92843    2.98956
Coef-variation   -134.86881 88.86827 -277.70256
3th-order moment   -0.03560  0.00670   -0.00484
4th-order moment    2.92603  3.08501    3.05716
5th-order moment   -0.77761  0.04420    0.02058
6th-order moment   24.55746 15.26236   15.43598

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

       x                y          
 Min.   :-5.742   Min.   :0.00000  
 1st Qu.:-3.015   1st Qu.:0.00070  
 Median :-0.288   Median :0.00728  
 Mean   :-0.288   Mean   :0.09159  
 3rd Qu.: 2.438   3rd Qu.:0.11802  
 Max.   : 5.165   Max.   :0.48233  
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. TEX.sde() - Converting Sim.DiffProc Objects to LaTeX.
  6. 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 (2018). Sim.DiffProc: Simulation of Diffusion Processes. R package version 4.1, 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}\).