mcqmc06_l: Financial options using a Milstein discretisation

Description Usage Arguments Details Author(s) References Examples

Description

Financial options based on scalar geometric Brownian motion, similar to Mike Giles' MCQMC06 paper, using a Milstein discretisation

Usage

1
mcqmc06_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 = barrier call.

Details

This function is based on GPL-2 C++ code by Mike Giles.

Author(s)

Louis Aslett <aslett@stats.ox.ac.uk>

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

References

M.B. Giles. 'Improved multilevel Monte Carlo convergence using the Milstein scheme', p.343-358 in Monte Carlo and Quasi-Monte Carlo Methods 2006, Springer, 2007.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
## Not run: 
# These are similar to the MLMC tests for the MCQMC06 paper
# using a Milstein discretisation with 2^l timesteps on level l
#
# The figures are slightly different due to:
# -- change in MSE split
# -- change in cost calculation
# -- different random number generation
# -- switch to S_0=100

M    <- 2 # refinement cost factor
N0   <- 200 # initial samples on coarse levels
Lmin <- 2 # minimum refinement level
Lmax <- 10 # maximum refinement level

test.res <- list()
for(option in 1:5) {
  if(option==1) {
    cat("\n ---- Computing European call ---- \n")
    N      <- 20000 # samples for convergence tests
    L      <- 8 # 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      <- 20000 # samples for convergence tests
    L      <- 8 # 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      <- 20000 # samples for convergence tests
    L      <- 10 # levels for convergence tests
    Eps    <- c(0.005, 0.01, 0.02, 0.05, 0.1)
  } else if(option==4) {
    cat("\n ---- Computing digital call ---- \n")
    N      <- 200000 # samples for convergence tests
    L      <- 8 # levels for convergence tests
    Eps    <- c(0.01, 0.02, 0.05, 0.1, 0.2)
  } else if(option==5) {
    cat("\n ---- Computing barrier call ---- \n")
    N      <- 200000 # samples for convergence tests
    L      <- 8 # levels for convergence tests
    Eps    <- c(0.005, 0.01, 0.02, 0.05, 0.1)
  }

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

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

## End(Not run)

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

mlmc documentation built on May 1, 2019, 9:45 p.m.

Related to mcqmc06_l in mlmc...