BM {Sim.DiffProc}R Documentation

Brownian motion, Brownian bridge, geometric Brownian motion, and arithmetic Brownian motion simulators

Description

The (S3) generic function for simulation of brownian motion, brownian bridge, geometric brownian motion, and arithmetic brownian motion.

Usage

BM(N, ...)
BB(N, ...)
GBM(N, ...)
ABM(N, ...)

## Default S3 method:
BM(N =1000,M=1,x0=0,t0=0,T=1,Dt=NULL, ...)
## Default S3 method:
BB(N =1000,M=1,x0=0,y=0,t0=0,T=1,Dt=NULL, ...)
## Default S3 method:
GBM(N =1000,M=1,x0=1,t0=0,T=1,Dt=NULL,theta=1,sigma=1, ...)
## Default S3 method:
ABM(N =1000,M=1,x0=0,t0=0,T=1,Dt=NULL,theta=1,sigma=1, ...)

Arguments

N

number of simulation steps.

M

number of trajectories.

x0

initial value of the process at time t_{0}.

y

terminal value of the process at time T of the BB.

t0

initial time.

T

final time.

Dt

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

theta

the interest rate of the ABM and GBM.

sigma

the volatility of the ABM and GBM.

...

potentially further arguments for (non-default) methods.

Details

The function BM returns a trajectory of the standard Brownian motion (Wiener process) in the time interval [t_{0},T]. Indeed, for W(dt) it holds true that W(dt) \rightarrow W(dt) - W(0) \rightarrow \mathcal{N}(0,dt), where \mathcal{N}(0,1) is normal distribution Normal.

The function BB returns a trajectory of the Brownian bridge starting at x_{0} at time t_{0} and ending at y at time T; i.e., the diffusion process solution of stochastic differential equation:

dX_{t}= \frac{y-X_{t}}{T-t} dt + dW_{t}

The function GBM returns a trajectory of the geometric Brownian motion starting at x_{0} at time t_{0}; i.e., the diffusion process solution of stochastic differential equation:

dX_{t}= \theta X_{t} dt + \sigma X_{t} dW_{t}

The function ABM returns a trajectory of the arithmetic Brownian motion starting at x_{0} at time t_{0}; i.e.,; the diffusion process solution of stochastic differential equation:

dX_{t}= \theta dt + \sigma dW_{t}

Value

X

an visible ts object.

Author(s)

A.C. Guidoum, K. Boukhetala.

References

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

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

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

See Also

This functions BM, BBridge and GBM are available in other packages such as "sde".

Examples


op <- par(mfrow = c(2, 2))

## Brownian motion
set.seed(1234)
X <- BM(M = 100)
plot(X,plot.type="single")
lines(as.vector(time(X)),rowMeans(X),col="red")

## Brownian bridge
set.seed(1234)
X <- BB(M =100)
plot(X,plot.type="single")
lines(as.vector(time(X)),rowMeans(X),col="red")

## Geometric Brownian motion
set.seed(1234)
X <- GBM(M = 100)
plot(X,plot.type="single")
lines(as.vector(time(X)),rowMeans(X),col="red")

## Arithmetic Brownian motion
set.seed(1234)
X <- ABM(M = 100)
plot(X,plot.type="single")
lines(as.vector(time(X)),rowMeans(X),col="red")

par(op)

[Package Sim.DiffProc version 4.9 Index]