Processing math: 100%

Monte-Carlo Simulations and Analysis of Stochastic Differential Equations

A.C. Guidoum1 and K. Boukhetala2

2020-05-04

snssde1d()

Assume that we want to describe the following SDE:

Ito form3:

dXt=12θ2Xtdt+θXtdWt,X0=x0>0 Stratonovich form: dXt=12θ2Xtdt+θXtdWt,X0=x0>0

In the above f(t,x)=12θ2x and g(t,x)=θx (θ>0), Wt 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=1000,type="ito") # Using Ito
R> mod2 <- snssde1d(drift=f,diffusion=g,x0=10,M=1000,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  = 1000.
    | 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  = 1000.
    | 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 Xt process at any time t:

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.19043
Variance              34.77438
Median                 9.80795
Mode                   7.81562
First quartile         7.05936
Third quartile        13.70165
Minimum                2.30727
Maximum               48.28898
Skewness               1.50480
Kurtosis               6.28741
Coef-variation         0.52697
3th-order moment     308.57990
4th-order moment    7603.10190
5th-order moment  167122.80676
6th-order moment 4610438.55545
R> summary(mod2, at = 1)

    Monte-Carlo Statistics for X(t) at time t = 1
                               
Mean                   12.94687
Variance               44.96980
Median                 11.53663
Mode                    9.47033
First quartile          8.21275
Third quartile         16.07308
Minimum                 2.34762
Maximum                56.04992
Skewness                1.47591
Kurtosis                6.41717
Coef-variation          0.51796
3th-order moment      445.08292
4th-order moment    12977.32870
5th-order moment   332826.52899
6th-order moment 10773943.92221

Hence we can just make use of the rsde1d() function to build our random number generator for the conditional density of the Xt|X0 (Xmod1t|X0 and Xmod2t|X0) at time t=1.

R> x1 <- rsde1d(object = mod1, at = 1)  # X(t=1) | X(0)=x0 (Ito SDE)
R> x2 <- rsde1d(object = mod2, at = 1)  # X(t=1) | X(0)=x0 (Stratonovich SDE)
R> head(data.frame(x1,x2),n=5)
       x1      x2
1 18.8246 14.6275
2 13.2760  7.1711
3  6.5486  7.8493
4 14.7511  5.3086
5 23.5182 35.2471

The function dsde1d() can be used to show the Approximate transitional density for Xt|X0 at time ts=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))
Approximate transitional density for $X_{t}|X_{0}$ at time $t-s=1$ with log-normal curves. mod1: Ito and mod2: Stratonovich  Approximate transitional density for $X_{t}|X_{0}$ at time $t-s=1$ with log-normal curves. mod1: Ito and mod2: Stratonovich

Figure 1: Approximate transitional density for Xt|X0 at time ts=1 with log-normal curves. mod1: Ito and mod2: Stratonovich

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> ## Ito
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: Ito and mod2: Stratonovich  mod1: Ito and mod2: Stratonovich

Figure 2: mod1: Ito and mod2: Stratonovich

Return to snssde1d()

snssde2d()

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

Ito form: {dXt=fx(t,Xt,Yt)dt+gx(t,Xt,Yt)dW1,tdYt=fy(t,Xt,Yt)dt+gy(t,Xt,Yt)dW2,t Stratonovich form: {dXt=fx(t,Xt,Yt)dt+gx(t,Xt,Yt)dW1,tdYt=fy(t,Xt,Yt)dt+gy(t,Xt,Yt)dW2,t

where (W1,t,W2,t) are a two independent standard Wiener process if corr = NULL. 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 Xt. In mathematical terms, the equation is written as an Ito equation: dXt=1μXtdt+σdWt,X0=x0 In these equations, μ and σ are positive constants called, respectively, the relaxation time and the diffusion constant. The time integral of the OU process Xt (or indeed of any process Xt) is defined to be the process Yt that satisfies: Yt=Y0+XtdtdYt=Xtdt,Y0=y0 Yt is not itself a Markov process; however, Xt and Yt together comprise a bivariate continuous Markov process. We wish to find the solutions Xt and Yt to the coupled time-evolution equations: {dXt=1μXtdt+σdWtdYt=Xtdt

We simulate a flow of 1000 trajectories of (Xt,Yt), with integration step size Δ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=1000,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  = 1000.
    | Initial values    | (x0,y0) = (5,0).
    | Time of process   | t in [0,10].
    | Discretization    | Dt = 0.01.

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

R> summary(mod2d, at = 10)

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

R> ## in time
R> plot(mod2d)
R> ## in plane (O,X,Y)
R> plot2d(mod2d,type="n") 
R> points2d(mod2d,col=rgb(0,100,0,50,maxColorValue=255), pch=16)
 Ornstein-Uhlenbeck process and its integral  Ornstein-Uhlenbeck process and its integral

Figure 3: Ornstein-Uhlenbeck process and its integral

Hence we can just make use of the rsde2d() function to build our random number for (Xt,Yt) at time t=10.

R> out <- rsde2d(object = mod2d, at = 10) 
R> head(out,n=3)
          x       y
1  0.835404 25.1479
2 -0.033169 19.9462
3 -1.169357  9.2241

The density of Xt and Yt at time t=10 are reported using dsde2d() function, see e.g. Figure 4: the marginal density of Xt and Yt at time t=10. For plotted in (x, y)-space with dim = 2. A contour and image plot of density obtained from a realization of system (Xt,Yt) at time t=10, see:

R> ## the marginal density
R> denM <- dsde2d(mod2d,pdf="M",at =10)
R> plot(denM, main="Marginal Density")
R> ## the Joint density
R> denJ <- dsde2d(mod2d, pdf="J", n=100,at =10)
R> plot(denJ,display="contour",main="Bivariate Transition Density at time t=10")
Marginal and Joint density at time t=10 Marginal and Joint density at time t=10

Figure 4: Marginal and Joint density at time t=10

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

R> plot(denJ,display="persp",main="Bivariate Transition Density at time t=10")
Marginal and Joint density at time t=10

Figure 5: Marginal and Joint density at time t=10

We approximate the bivariate transition density over the set transition horizons t[1,10] by Δ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 ˙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: ¨Xμ(1X2)˙X+X=0 where x is the position coordinate (which is a function of the time t), and μ 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 ξt, with delta-type correlation function E(ξtξt+h)=2σδ(h) ¨Xμ(1X2)˙X+X=ξt, where μ>0 . It’s solution cannot be obtained in terms of elementary functions, even in the phase plane. The White noise ξt is formally derivative of the Wiener process Wt. The representation of a system of two first order equations follows the same idea as in the deterministic case by letting ˙x=y, from physical equation we get the above system: {˙X=Y˙Y=μ(1X2)YX+ξt The system can mathematically explain by a Stratonovitch equations: {dXt=YtdtdYt=(μ(1X2t)YtXt)dt+2σdW2,t

Implemente in R as follows, with integration step size Δ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")

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

R> plot(mod2d,ylim=c(-8,8))   ## back in time
R> plot2d(mod2d)              ## in plane (O,X,Y)
The stochastic Van-der-Pol equationThe stochastic Van-der-Pol equation

Figure 6: The stochastic Van-der-Pol equation

Return to snssde2d()

The Heston Model

Consider a system of stochastic differential equations:

{dXt=μXtdt+XtYtdB1,tdYt=ν(θYt)dt+σYtdB2,t

Conditions to ensure positiveness of the volatility process are that 2νθ>σ2, and the two Brownian motions (B1,t,B2,t) are correlated. Σ to describe the correlation structure, for example: Σ=(10.30.32)

R> mu = 1.2; sigma=0.1;nu=2;theta=0.5
R> fx <- expression( mu*x ,nu*(theta-y)) 
R> gx <- expression( x*sqrt(y) ,sigma*sqrt(y))
R> Sigma <- matrix(c(1,0.3,0.3,2),nrow=2,ncol=2) # correlation matrix
R> HM <- snssde2d(drift=fx,diffusion=gx,Dt=0.001,x0=c(100,1),corr=Sigma,M=1000)
R> HM
Itô Sde 2D:
    | dX(t) = mu * X(t) * dt + X(t) * sqrt(Y(t)) * dB1(t)
    | dY(t) = nu * (theta - Y(t)) * dt + sigma * sqrt(Y(t)) * dB2(t)
    | Correlation structure:                    
             1.0 0.3
             0.3 2.0
Method:
    | Euler scheme with order 0.5
Summary:
    | Size of process   | N  = 1001.
    | Number of simulation  | M  = 1000.
    | Initial values    | (x0,y0) = (100,1).
    | Time of process   | t in [0,1].
    | Discretization    | Dt = 0.001.

Hence we can just make use of the rsde2d() function to build our random number for (Xt,Yt) at time t=1.

R> out <- rsde2d(object = HM, at = 1) 
R> head(out,n=3)
       x       y
1 667.47 0.59944
2 152.76 0.59491
3 414.38 0.58505

The density of Xt and Yt at time t=1 are reported using dsde2d() function. See:

R> denJ <- dsde2d(HM,pdf="J",at =1,lims=c(-100,900,0.4,0.75))
R> plot(denJ,display="contour",main="Bivariate Transition Density at time t=10")
R> plot(denJ,display="persp",main="Bivariate Transition Density at time t=10")

Return to snssde2d()

snssde3d()

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

Ito form: {dXt=fx(t,Xt,Yt,Zt)dt+gx(t,Xt,Yt,Zt)dW1,tdYt=fy(t,Xt,Yt,Zt)dt+gy(t,Xt,Yt,Zt)dW2,tdZt=fz(t,Xt,Yt,Zt)dt+gz(t,Xt,Yt,Zt)dW3,t Stratonovich form: {dXt=fx(t,Xt,Yt,Zt)dt+gx(t,Xt,Yt,Zt)dW1,tdYt=fy(t,Xt,Yt,Zt)dt+gy(t,Xt,Yt,Zt)dW2,tdZt=fz(t,Xt,Yt,Zt)dt+gz(t,Xt,Yt,Zt)dW3,t

(W1,t,W2,t,W3,t) are three independents standard Wiener process if corr = NULL. 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 Ito form: {dXt=4(1Xt)Ytdt+0.2dW1,tdYt=4(1Yt)Xtdt+0.2dW2,tdZt=4(1Zt)Ytdt+0.2dW3,t

with (W1,t,W2,t,W3,t) are three indpendant standard Wiener process.

We simulate a flow of 1000 trajectories, with integration step size Δ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=1000)
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  = 1000.
    | 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 (Xt,Yt,Zt) process at any time t, for example at=1:

R> s = 1
R> mean(mod3d, at = s)
R> moment(mod3d, at = s , center = TRUE , order = 2) ## variance
R> Median(mod3d, at = s)
R> Mode(mod3d, at = s)
R> quantile(mod3d , at = s)
R> kurtosis(mod3d , at = s)
R> skewness(mod3d , at = s)
R> cv(mod3d , at = s )
R> min(mod3d , at = s)
R> max(mod3d , at = s)
R> moment(mod3d, at = s , center= TRUE , order = 4)
R> moment(mod3d, at = s , center= FALSE , order = 4)

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

R> summary(mod3d, at = t)

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

R> plot(mod3d,union = TRUE)         ## back in time
R> plot3D(mod3d,display="persp")    ## in space (O,X,Y,Z)
 Flow of $1000$ trajectories of $(X_t ,Y_t ,Z_t)$  Flow of $1000$ trajectories of $(X_t ,Y_t ,Z_t)$

Figure 7: Flow of 1000 trajectories of (Xt,Yt,Zt)

Hence we can just make use of the rsde3d() function to build our random number for (Xt,Yt,Zt) at time t=1.

R> out <- rsde3d(object = mod3d, at = 1) 
R> head(out,n=3)
         x        y       z
1 -0.83626 1.005921 0.85364
2 -0.76947 0.768522 0.77156
3 -0.59477 0.018633 0.58222

For each SDE type and for each numerical scheme, the marginal density of Xt, Yt and Zt at time t=1 are reported using dsde3d() function, see e.g. Figure 8.

R> den <- dsde3d(mod3d,pdf="M",at =1)
R> plot(den, main="Marginal Density") 
 Marginal density of $X_t$, $Y_t$ and $Z_t$ at time $t=1$

Figure 8: Marginal density of Xt, Yt and Zt at time t=1

For an approximate joint transition density for (Xt,Yt,Zt) (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 Uw(x,y,z,t), Vw(x,y,z,t) and Sw(x,y,z,t) are neglected and the dispersion coefficient D(x,y,z) is constant. A system becomes (see Boukhetala,1996): {dXt=(KXtX2t+Y2t+Z2t)dt+σdW1,tdYt=(KYtX2t+Y2t+Z2t)dt+σdW2,tdZt=(KZtX2t+Y2t+Z2t)dt+σdW3,t

with initial conditions (X0,Y0,Z0)=(1,1,1), by specifying the drift and diffusion coefficients of three processes Xt, Yt and Zt 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))

The results of simulation (3D) are shown in Figure 9:

R> plot3D(mod3d,display="persp",col="blue")
 Attractive model for 3D diffusion processes

Figure 9: Attractive model for 3D diffusion processes

Return to snssde3d()

Transformation of an SDE one-dimensional

Next is an example of one-dimensional SDE driven by three correlated Wiener process (B1,t,B2,t,B3,t), as follows: dXt=B1,tdt+B2,tdB3,t with: Σ=(10.20.50.210.70.50.71) To simulate the solution of the process Xt, we make a transformation to a system of three equations as follows: {dXt=Ytdt+ZtdB3,tdYt=dB1,tdZt=dB2,t

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

R> fx <- expression(y,0,0) 
R> gx <- expression(z,1,1)
R> Sigma <-matrix(c(1,0.2,0.5,0.2,1,-0.7,0.5,-0.7,1),nrow=3,ncol=3)
R> modtra <- snssde3d(drift=fx,diffusion=gx,M=1000,corr=Sigma)
R> modtra
Itô Sde 3D:
    | dX(t) = Y(t) * dt + Z(t) * dB1(t)
    | dY(t) = 0 * dt + 1 * dB2(t)
    | dZ(t) = 0 * dt + 1 * dB3(t)
    | Correlation structure:                          
             1.0  0.2  0.5
             0.2  1.0 -0.7
             0.5 -0.7  1.0
Method:
    | Euler scheme with order 0.5
Summary:
    | Size of process   | N  = 1001.
    | Number of simulation  | M  = 1000.
    | Initial values    | (x0,y0,z0) = (0,0,0).
    | Time of process   | t in [0,1].
    | Discretization    | Dt = 0.001.

The histogram and kernel density of Xt at time t=1 are reported using rsde3d() function, and we calculate emprical variance-covariance matrix C(s,t)=Cov(Xs,Xt), see e.g. Figure 10.

R> X <- rsde3d(modtra,at=1)$x
R> MASS::truehist(X,xlab = expression(X[t==1]));box()
R> lines(density(X),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)
R> ## Cov-Matrix
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 = ""))
The histogram and kernel density of $X_t$ at time $t=1$. Emprical variance-covariance matrixThe histogram and kernel density of $X_t$ at time $t=1$. Emprical variance-covariance matrix

Figure 10: The histogram and kernel density of Xt at time t=1. Emprical variance-covariance matrix

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 (2020). Performing Parallel Monte Carlo and Moment Equations Methods for Ito and Stratonovich Stochastic Differential Systems: R Package Sim.DiffProc. Accept Submission to Journal of Statistical Software.

  3. Guidoum AC, Boukhetala K (2020). Sim.DiffProc: Simulation of Diffusion Processes. R package version 4.6, 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 Xmod1t the following Stratonovich SDE: dXt=θXtdWt.