max_qEI {DiceOptim}R Documentation

Maximization of multipoint expected improvement criterion (qEI)

Description

Maximization of the qEI criterion. Two options are available : Constant Liar (CL), and brute force qEI maximization with Broyden-Fletcher-Goldfarb-Shanno (BFGS) algorithm, or GENetic Optimization Using Derivative (genoud) algorithm.

Usage

max_qEI(
  model,
  npoints,
  lower,
  upper,
  crit = "exact",
  minimization = TRUE,
  optimcontrol = NULL
)

Arguments

model

an object of class km ,

npoints

an integer representing the desired number of iterations,

lower

vector of lower bounds,

upper

vector of upper bounds,

crit

"exact", "CL" : a string specifying the criterion used. "exact" triggers the maximization of the multipoint expected improvement at each iteration (see max_qEI), "CL" applies the Constant Liar heuristic,

minimization

logical specifying if the qEI to be maximized is used in minimiziation or in maximization,

optimcontrol

an optional list of control parameters for optimization. See details.

Details

- CL is a heuristic method. First, the regular Expected Improvement EI is maximized (max_EI). Then, for the next points, the Expected Improvement is maximized again, but with an artificially updated Kriging model. Since the response values corresponding to the last best point obtained are not available, the idea of CL is to replace them by an arbitrary constant value L (a "lie") set by the user (default is the minimum of all currently available observations).

- The BFGS algorithm is implemented in the standard function optim. Analytical formulae of qEI and its gradient qEI.grad are used. The nStarts starting points are by default sampled with respect to the regular EI (sampleFromEI) criterion.

- The "genoud" method calls the function genoud using analytical formulae of qEI and its gradient qEI.grad.

The parameters of list optimcontrol are :

- optimcontrol$method : "BFGS" (default), "genoud" ; a string specifying the method used to maximize the criterion (irrelevant when crit is "CL" because this method always uses genoud),

- when crit="CL" :

+ optimcontrol$parinit : optional matrix of initial values (must have model@d columns, the number of rows is not constrained),

+ optimcontrol$L : "max", "min", "mean" or a scalar value specifying the liar ; "min" takes model@min, "max" takes model@max, "mean" takes the prediction of the model ; When L is NULL, "min" is taken if minimization==TRUE, else it is "max".

+ The parameters of function genoud. Main parameters are : "pop.size" (default : [N=3*2^model@d for dim<6 and N=32*model@d otherwise]), "max.generations" (default : 12), "wait.generations" (default : 2) and "BFGSburnin" (default : 2).

- when optimcontrol$method = "BFGS" :

+ optimcontrol$nStarts (default : 4),

+ optimcontrol$fastCompute : if TRUE (default), a fast approximation method based on a semi-analytic formula is used, see [Marmin 2014] for details,

+ optimcontrol$samplingFun : a function which sample a batch of starting point (default : sampleFromEI),

+ optimcontrol$parinit : optional 3d-array of initial (or candidate) batches (for all k, parinit[,,k] is a matrix of size npoints*model@d representing one batch). The number of initial batches (length(parinit[1,1,])) is not contrained and does not have to be equal to nStarts. If there is too few initial batches for nStarts, missing batches are drawn with samplingFun (default : NULL),

- when optimcontrol$method = "genoud" :

+ optimcontrol$fastCompute : if TRUE (default), a fast approximation method based on a semi-analytic formula is used, see [Marmin 2014] for details,

+ optimcontrol$parinit : optional matrix of candidate starting points (one row corresponds to one point),

+ The parameters of the genoud function. Main parameters are "pop.size" (default : [50*(model@d)*(npoints)]), "max.generations" (default : 5), "wait.generations" (default : 2), "BFGSburnin" (default : 2).

Value

A list with components:

par

A matrix containing the npoints input vectors found.

value

A value giving the qEI computed in par.

Author(s)

Sebastien Marmin

Clement Chevalier

David Ginsbourger

References

C. Chevalier and D. Ginsbourger (2014) Learning and Intelligent Optimization - 7th International Conference, Lion 7, Catania, Italy, January 7-11, 2013, Revised Selected Papers, chapter Fast computation of the multipoint Expected Improvement with applications in batch selection, pages 59-69, Springer.

D. Ginsbourger, R. Le Riche, L. Carraro (2007), A Multipoint Criterion for Deterministic Parallel Global Optimization based on Kriging. The International Conference on Non Convex Programming, 2007.

D. Ginsbourger, R. Le Riche, and L. Carraro. Kriging is well-suited to parallelize optimization (2010), In Lim Meng Hiot, Yew Soon Ong, Yoel Tenne, and Chi-Keong Goh, editors, Computational Intelligence in Expensive Optimization Problems, Adaptation Learning and Optimization, pages 131-162. Springer Berlin Heidelberg.

J. Mockus (1988), Bayesian Approach to Global Optimization. Kluwer academic publishers.

M. Schonlau (1997), Computer experiments and global optimization, Ph.D. thesis, University of Waterloo.

See Also

qEI, qEI.grad

Examples




set.seed(000)
# 3-points EI maximization.
# 9-points factorial design, and the corresponding response
d <- 2
n <- 9
design.fact <- expand.grid(seq(0,1,length=3), seq(0,1,length=3)) 
names(design.fact)<-c("x1", "x2")
design.fact <- data.frame(design.fact) 
names(design.fact)<-c("x1", "x2")
response.branin <- apply(design.fact, 1, branin)
response.branin <- data.frame(response.branin) 
names(response.branin) <- "y" 
lower <- c(0,0)
upper <- c(1,1)

# number of point in the bacth
batchSize <- 3

# model identification
fitted.model <- km(~1, design=design.fact, response=response.branin, 
                   covtype="gauss", control=list(pop.size=50,trace=FALSE), parinit=c(0.5, 0.5))

# maximization of qEI

# With a multistarted BFGS algorithm
maxBFGS <- max_qEI(model = fitted.model, npoints = batchSize, lower = lower, upper = upper, 
crit = "exact",optimcontrol=list(nStarts=3,method = "BFGS"))

# comparison
print(maxBFGS$value)
## Not run: 
# With a genetic algorithme using derivatives
maxGen  <- max_qEI(model = fitted.model, npoints = batchSize, lower = lower, upper = upper, 
crit = "exact", optimcontrol=list(nStarts=3,method = "genoud",pop.size=100,max.generations = 15))
# With the constant liar heuristic
maxCL   <- max_qEI(model = fitted.model, npoints = batchSize, lower = lower, upper = upper, 
crit = "CL",optimcontrol=list(pop.size=20))
print(maxGen$value)
print(maxCL$value)

## End(Not run)



[Package DiceOptim version 2.1.1 Index]