opre_l: Financial options using an Euler-Maruyama discretisation

View source: R/opre_l.R

opre_lR Documentation

Financial options using an Euler-Maruyama discretisation

Description

Financial options based on scalar geometric Brownian motion and Heston models, similar to Mike Giles' original 2008 Operations Research paper, Giles (2008), using an Euler-Maruyama discretisation

Usage

opre_l(l, N, option)

Arguments

l

the level to be simulated.

N

the number of samples to be computed.

option

the option type, between 1 and 5. The options are:

1 = European call;
2 = Asian call;
3 = lookback call;
4 = digital call;
5 = Heston model.

Details

This function is based on GPL-2 'Matlab' code by Mike Giles.

Value

A named list containing:

sums

is a vector of length six \left(\sum Y_i, \sum Y_i^2, \sum Y_i^3, \sum Y_i^4, \sum X_i, \sum X_i^2\right) where Y_i are iid simulations with expectation E[P_0] when l=0 and expectation E[P_l-P_{l-1}] when l>0, and X_i are iid simulations with expectation E[P_l]. Note that only the first two components of this are used by the main mlmc() driver, the full vector is used by mlmc.test() for convergence tests etc;

cost

is a scalar with the total cost of the paths simulated, computed as N \times 4^l for level l.

Author(s)

Louis Aslett <louis.aslett@durham.ac.uk>

Mike Giles <Mike.Giles@maths.ox.ac.uk>

Tigran Nagapetyan <nagapetyan@stats.ox.ac.uk>

References

Giles, M.B. (2008) 'Multilevel Monte Carlo Path Simulation', Operations Research, 56(3), pp. 607–617. Available at: \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1287/opre.1070.0496")}.

Examples


# These are similar to the MLMC tests for the original
# 2008 Operations Research paper, using an Euler-Maruyama
# discretisation with 4^l timesteps on level l.
#
# The differences are:
# -- the plots do not have the extrapolation results
# -- two plots are log_2 rather than log_4
# -- the new MLMC driver is a little different
# -- switch to X_0=100 instead of X_0=1
#
# Note the following takes quite a while to run, for a toy example see after
# this block.

N0   <- 1000 # initial samples on coarse levels
Lmin <- 2 # minimum refinement level
Lmax <- 6 # maximum refinement level

test.res <- list()
for(option in 1:5) {
  if(option == 1) {
    cat("\n ---- Computing European call ---- \n")
    N      <- 1000000 # samples for convergence tests
    L      <- 5 # levels for convergence tests
    Eps    <- c(0.005, 0.01, 0.02, 0.05, 0.1)
  } else if(option == 2) {
    cat("\n ---- Computing Asian call ---- \n")
    N      <- 1000000 # samples for convergence tests
    L      <- 5 # levels for convergence tests
    Eps    <- c(0.005, 0.01, 0.02, 0.05, 0.1)
  } else if(option == 3) {
    cat("\n ---- Computing lookback call ---- \n")
    N      <- 1000000 # samples for convergence tests
    L      <- 5 # levels for convergence tests
    Eps    <- c(0.01, 0.02, 0.05, 0.1, 0.2)
  } else if(option == 4) {
    cat("\n ---- Computing digital call ---- \n")
    N      <- 4000000 # samples for convergence tests
    L      <- 5 # levels for convergence tests
    Eps    <- c(0.02, 0.05, 0.1, 0.2, 0.5)
  } else if(option == 5) {
    cat("\n ---- Computing Heston model ---- \n")
    N      <- 2000000 # samples for convergence tests
    L      <- 5 # levels for convergence tests
    Eps    <- c(0.005, 0.01, 0.02, 0.05, 0.1)
  }

  test.res[[option]] <- mlmc.test(opre_l, N, L, N0, Eps, Lmin, Lmax, option = option)

  # print exact analytic value, based on S0=K
  T   <- 1
  r   <- 0.05
  sig <- 0.2
  K   <- 100

  k   <- 0.5*sig^2/r;
  d1  <- (r+0.5*sig^2)*T / (sig*sqrt(T))
  d2  <- (r-0.5*sig^2)*T / (sig*sqrt(T))

  if(option == 1) {
    val <- K*( pnorm(d1) - exp(-r*T)*pnorm(d2) )
  } else if(option == 2) {
    val <- NA
  } else if(option == 3) {
    val <- K*( pnorm(d1) - pnorm(-d1)*k - exp(-r*T)*(pnorm(d2) - pnorm(d2)*k) )
  } else if(option == 4) {
    val <- K*exp(-r*T)*pnorm(d2)
  } else if(option == 5) {
    val <- NA
  }

  if(is.na(val)) {
    cat(sprintf("\n Exact value unknown, MLMC value: %f \n", test.res[[option]]$P[1]))
  } else {
    cat(sprintf("\n Exact value: %f, MLMC value: %f \n", val, test.res[[option]]$P[1]))
  }

  # plot results
  plot(test.res[[option]])
}


# The level sampler can be called directly to retrieve the relevant level sums:
opre_l(l = 7, N = 10, option = 1)


louisaslett/mlmc documentation built on Sept. 3, 2024, 10:36 p.m.