snssde1d: Simulation of 1-D Stochastic Differential Equation

View source: R/snssde.R

snssde1dR Documentation

Simulation of 1-D Stochastic Differential Equation

Description

The (S3) generic function snssde1d of simulation of solution to 1-dim stochastic differential equation of Itô or Stratonovich type, with different methods.

Usage

snssde1d(N, ...)
## Default S3 method:
snssde1d(N = 1000, M = 1, x0 = 0, t0 = 0, T = 1, Dt, 
   drift, diffusion, alpha = 0.5, mu = 0.5, type = c("ito", "str"), 
   method = c("euler", "milstein", "predcorr", "smilstein", "taylor", 
   "heun", "rk1", "rk2", "rk3"), ...)	
   
## S3 method for class 'snssde1d'
summary(object, at ,digits=NULL, ...)
## S3 method for class 'snssde1d'
time(x, ...)
## S3 method for class 'snssde1d'
mean(x, at, ...)
## S3 method for class 'snssde1d'
Median(x, at, ...)
## S3 method for class 'snssde1d'
Mode(x, at, ...)
## S3 method for class 'snssde1d'
quantile(x, at, ...)
## S3 method for class 'snssde1d'
kurtosis(x, at, ...)
## S3 method for class 'snssde1d'
min(x, at, ...)
## S3 method for class 'snssde1d'
max(x, at, ...)
## S3 method for class 'snssde1d'
skewness(x, at, ...)
## S3 method for class 'snssde1d'
moment(x, at,  ...)
## S3 method for class 'snssde1d'
cv(x, at,  ...)
## S3 method for class 'snssde1d'
bconfint(x, at,  ...)

## S3 method for class 'snssde1d'
plot(x, ...)
## S3 method for class 'snssde1d'
lines(x, ...)
## S3 method for class 'snssde1d'
points(x, ...)		

Arguments

N

number of simulation steps.

M

number of trajectories (Monte-Carlo).

x0

initial value of the process at time t0.

t0

initial time.

T

ending time.

Dt

time step of the simulation (discretization). If it is missing a default \Delta t = \frac{T-t_{0}}{N}.

drift

drift coefficient: an expression of two variables t and x.

diffusion

diffusion coefficient: an expression of two variables t and x.

alpha, mu

weight of the predictor-corrector scheme; the default alpha = 0.5 and mu = 0.5.

type

if type="ito" simulation sde of Itô type, else type="str" simulation sde of Stratonovich type; the default type="ito".

method

numerical methods of simulation, the default method = "euler".

x, object

an object inheriting from class "snssde1d".

at

time between t0 and T. Monte-Carlo statistics of the solution X_{t} at time at. The default at = T.

digits

integer, used for number formatting.

...

potentially further arguments for (non-default) methods.

Details

The function snssde1d returns a ts x of length N+1; i.e. solution of the sde of Ito or Stratonovich types; If Dt is not specified, then the best discretization \Delta t = \frac{T-t_{0}}{N}.

The Ito stochastic differential equation is:

dX(t) = a(t,X(t)) dt + b(t,X(t)) dW(t)

Stratonovich sde :

dX(t) = a(t,X(t)) dt + b(t,X(t)) \circ dW(t)

The methods of approximation are classified according to their different properties. Mainly two criteria of optimality are used in the literature: the strong and the weak (orders of) convergence. The method of simulation can be one among: Euler-Maruyama Order 0.5, Milstein Order 1, Milstein Second-Order, Predictor-Corrector method, Itô-Taylor Order 1.5, Heun Order 2 and Runge-Kutta Order 1, 2 and 3.

An overview of this package, see browseVignettes('Sim.DiffProc') for more informations.

Value

snssde1d returns an object inheriting from class "snssde1d".

X

an invisible ts object.

drift

drift coefficient.

diffusion

diffusion coefficient.

type

type of sde.

method

the numerical method used.

Author(s)

A.C. Guidoum, K. Boukhetala.

References

Guidoum AC, Boukhetala K (2020). "Performing Parallel Monte Carlo and Moment Equations Methods for Itô and Stratonovich Stochastic Differential Systems: R Package Sim.DiffProc". Journal of Statistical Software, 96(2), 1–82. doi:10.18637/jss.v096.i02

Friedman, A. (1975). Stochastic differential equations and applications. Volume 1, ACADEMIC PRESS.

Henderson, D. and Plaschko,P. (2006). Stochastic differential equations in science and engineering. World Scientific.

Allen, E. (2007). Modeling with Ito stochastic differential equations. Springer-Verlag.

Jedrzejewski, F. (2009). Modeles aleatoires et physique probabiliste. Springer-Verlag.

Iacus, S.M. (2008). Simulation and inference for stochastic differential equations: with R examples. Springer-Verlag, New York.

Kloeden, P.E, and Platen, E. (1989). A survey of numerical methods for stochastic differential equations. Stochastic Hydrology and Hydraulics, 3, 155–178.

Kloeden, P.E, and Platen, E. (1991a). Relations between multiple ito and stratonovich integrals. Stochastic Analysis and Applications, 9(3), 311–321.

Kloeden, P.E, and Platen, E. (1991b). Stratonovich and ito stochastic taylor expansions. Mathematische Nachrichten, 151, 33–50.

Kloeden, P.E, and Platen, E. (1995). Numerical Solution of Stochastic Differential Equations. Springer-Verlag, New York.

Oksendal, B. (2000). Stochastic Differential Equations: An Introduction with Applications. 5th edn. Springer-Verlag, Berlin.

Platen, E. (1980). Weak convergence of approximations of ito integral equations. Z Angew Math Mech. 60, 609–614.

Platen, E. and Bruti-Liberati, N. (2010). Numerical Solution of Stochastic Differential Equations with Jumps in Finance. Springer-Verlag, New York

Saito, Y, and Mitsui, T. (1993). Simulation of Stochastic Differential Equations. The Annals of the Institute of Statistical Mathematics, 3, 419–432.

See Also

snssde2d and snssde3d for 2 and 3-dim sde.

sde.sim in package "sde".

simulate in package "yuima".

Examples


## Example 1: Ito sde
## dX(t) = 2*(3-X(t)) dt + 2*X(t) dW(t)
set.seed(1234)

f <- expression(2*(3-x) )
g <- expression(1)
mod1 <- snssde1d(drift=f,diffusion=g,M=4000,x0=10,Dt=0.01)
mod1
summary(mod1)
## Not run: 
plot(mod1)
lines(time(mod1),apply(mod1$X,1,mean),col=2,lwd=2)
lines(time(mod1),apply(mod1$X,1,bconfint,level=0.95)[1,],col=4,lwd=2)
lines(time(mod1),apply(mod1$X,1,bconfint,level=0.95)[2,],col=4,lwd=2)
legend("topright",c("mean path",paste("bound of", 95," percent confidence")),
       inset = .01,col=c(2,4),lwd=2,cex=0.8)

## End(Not run)
## Example 2: Stratonovich sde
## dX(t) = ((2-X(t))/(2-t)) dt + X(t) o dW(t)
set.seed(1234)

f <- expression((2-x)/(2-t))
g <- expression(x)
mod2 <- snssde1d(type="str",drift=f,diffusion=g,M=4000,x0=1, method="milstein")
mod2
summary(mod2,at = 0.25)
summary(mod2,at = 1)
## Not run: 
plot(mod2)
lines(time(mod2),apply(mod2$X,1,mean),col=2,lwd=2)
lines(time(mod2),apply(mod2$X,1,bconfint,level=0.95)[1,],col=4,lwd=2)
lines(time(mod2),apply(mod2$X,1,bconfint,level=0.95)[2,],col=4,lwd=2)
legend("topleft",c("mean path",paste("bound of", 95," percent confidence")),
       inset = .01,col=c(2,4),lwd=2,cex=0.8)

## End(Not run)

Sim.DiffProc documentation built on May 29, 2024, 8:09 a.m.