optim_pn_MC: Maximizing excursion probabilities by sampling uniformly...

Description Usage Arguments Value Author(s) Examples

View source: R/optim_pn_MC.R

Description

Computes or maximizes excursion probability of points in Dinv based on either a fixed matrix of simulation points in D or on a Monte-Carlo optimization procedure.

Usage

1
2
optim_pn_MC(inv.integration.points, model, T, opt.index, inv.index, lower,
  upper, n.optpoints, pop.size, MC = TRUE, deterministicMat = NULL)

Arguments

inv.integration.points

Matrix of dimension n x dinv where n is the number of integration points and dinv the number of controlled parameters.

model

The current kriging model. km object.

T

Target threshold.

opt.index

Array with integers corresponding to the indices of the nuisance parameters.

inv.index

Array with integers corresponding to the indices of the controlled parameters.

lower

Array of size d. Lower bound of the input domain.

upper

Array of size d. Upper bound of the input domain.

n.optpoints

The number of simulation points in Dopt. Multivariate integrals in dimension n.optpoints are performed, so no value larger than 10 should be used.

pop.size

The number of Monte-Carlo samples used to maximize the excursion probability. Does not apply if MC = FALSE.

MC

Boolean. Is Monte-Carlo optimization required ?

deterministicMat

Applies only if MC = FALSE. Matrix of dimension (n*p) x d containing the simulation points used to compute (only once) pn through a multivariate integral.

Value

A list with the following fields. (i) integration.points: (n*p) x d matrix containing the simulation points optimizing the excursion probability. First p rows are associated to integration point 1, and so on. (ii) cov: (n*p) x p matrix containing the n different non-conditional covariance matrix in each batch of p points selected by the optimization procedure. (iii) mean: Array of size n*p with the kriging means of the simulation points. (iv) pn: Array of size n with the (underestimated) excursion probabilities for each integration point.

Author(s)

Clement Chevalier clement.chevalier@unine.ch

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
library(KrigInv)
library(mnormt)
library(randtoolbox)
myfun <- branin_robinv
d <- 3

set.seed(8)

n0 <- 30
T <- 10
opt.index <- c(3)
inv.index <- c(1,2)
lower <- rep(0,times=d)
upper <- rep(1,times=d)
d.inv <- length(inv.index);d.opt <- length(opt.index)
lower.inv <- lower[inv.index];upper.inv <- upper[inv.index]
lower.opt <- lower[opt.index];upper.opt <- upper[opt.index]

design <- matrix(runif(d*n0),nrow=n0)
response <- myfun(design)
model <- km(formula = ~1,design = design,response = response,covtype = "matern3_2")

n <- n.points <- 50 # number of integration points
inv.integration.points <- t(lower.inv + t(sobol(n=n.points,dim=d.inv))*(upper.inv-lower.inv))

p <- n.optpoints <- 2
pop.size <- 20
## Not run: 
result <- optim_pn_MC(inv.integration.points=inv.integration.points,model=model,T=T,
                      opt.index=opt.index,inv.index=inv.index,lower=lower,upper=upper,
                      n.optpoints=n.optpoints,pop.size=pop.size,MC=TRUE)

## End(Not run)
# an example with no MC
opt.simulation.points <- t(lower.opt + t(sobol(n=n.optpoints,dim=d.opt))*(upper.opt-lower.opt))
fullpoints <- matrix(c(0),nrow=n*p , ncol=d)
indices <- expand.grid(c(1:p),c(1:n)) 
index1 <- indices[,1]
index2 <- indices[,2]
fullpoints[,opt.index] <- opt.simulation.points[index1,]
fullpoints[,inv.index] <- inv.integration.points[index2,]
## Not run: 
result2 <- optim_pn_MC(inv.integration.points=inv.integration.points,model=model,T=T,
                       opt.index=opt.index,inv.index=inv.index,lower=lower,upper=upper,
                       n.optpoints=n.optpoints,pop.size=1,MC=FALSE,
                       deterministicMat = fullpoints)

## End(Not run)                        

IRSN/RobustInv documentation built on Nov. 20, 2019, 10:46 p.m.