mcqmc06_l | R Documentation |
Financial options based on scalar geometric Brownian motion, similar to Mike Giles' MCQMC06 paper, Giles (2008), using a Milstein discretisation.
mcqmc06_l(l, N, option)
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:
|
This function is based on GPL-2 C++ code by Mike Giles.
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 2^l
for level l
.
Louis Aslett <louis.aslett@durham.ac.uk>
Mike Giles <Mike.Giles@maths.ox.ac.uk>
Giles, M. (2008) 'Improved Multilevel Monte Carlo Convergence using the Milstein Scheme', in A. Keller, S. Heinrich, and H. Niederreiter (eds) Monte Carlo and Quasi-Monte Carlo Methods 2006. Berlin, Heidelberg: Springer, pp. 343–358. Available at: \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1007/978-3-540-74496-2_20")}.
# 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
#
# Note the following takes quite a while to run, for a toy example see after
# this block.
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, 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
B <- 0.85*K
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))
d3 <- (2*log(B/K) + (r+0.5*sig^2)*T) / (sig*sqrt(T))
d4 <- (2*log(B/K) + (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 <- K*( pnorm(d1) - exp(-r*T)*pnorm(d2) -
((K/B)^(1-1/k))*((B^2)/(K^2)*pnorm(d3) - exp(-r*T)*pnorm(d4)) )
}
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:
mcqmc06_l(l = 7, N = 10, option = 1)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.