critcst_optimizer {DiceOptim}R Documentation

Maximization of constrained Expected Improvement criteria

Description

Given objects of class km for the objective and constraints, and a set of tuning parameters (lower, upper and critcontrol), critcst_optimizer performs the maximization of a constrained Expected Improvement or SUR criterion and delivers the next point to be visited in an EGO-like procedure.

The latter maximization relies either on a genetic algorithm using derivatives, genoud or exhaustive search at pre-specified points. It is important to remark that the information needed about the objective and constraint functions reduces here to the vector of response values embedded in the models (no call to the objective/constraint functions or simulators (except possibly for the objective)).

Usage

critcst_optimizer(
  crit = "EFI",
  model.fun,
  model.constraint,
  equality = FALSE,
  lower,
  upper,
  type = "UK",
  critcontrol = NULL,
  optimcontrol = NULL
)

Arguments

crit

sampling criterion. Three choices are available : "EFI", "AL" and "SUR",

model.fun

object of class km corresponding to the objective function, or, if the objective function is fast-to-evaluate, either the objective function to be minimized or a fastfun object, see details and examples below,

model.constraint

either one or a list of models of class km, one for each constraint,

equality

either FALSE if all constraints are inequalities, or a Boolean vector indicating which are equalities,

lower

vector of lower bounds for the variables to be optimized over,

upper

vector of upper bounds for the variables to be optimized over,

type

"SK" or "UK" (default), depending whether uncertainty related to trend estimation has to be taken into account.

critcontrol

optional list of control parameters for criterion crit, see details.
Options for the checkPredict function: threshold (1e-4) and distance (covdist) are used to avoid numerical issues occuring when adding points too close to the existing ones.

optimcontrol

optional list of control parameters for optimization of the selected infill criterion. "method" set the optimization method; one can choose between "discrete" and "genoud". For each method, further parameters can be set.
For "discrete", one has to provide the argument "candidate.points".
For "genoud", one can control, among others, "pop.size" (default : [N = 3*2^dim for dim < 6 and N = 32*dim otherwise]), "max.generations" (12), "wait.generations" (2)), see genoud. Numbers into brackets are the default values. @return A list with components:

  • par: The best set of parameters found,

  • value: The value of expected improvement at par.

Details

Extension of the function max_EI for constrained optimization.

Available infill criteria with crit are :

Depending on the selected criterion, parameters can be given with critcontrol. Also options for checkPredict are available. More precisions are given in the corresponding help pages.

If the objective function to minimize is inexpensive, i.e. no need of a kriging model, then one can provide it in model.obj, which is handled next with class fastfun (or directly as a fastfun object). See example below.

In the case of equality constraints, it is possible to define them with equality. Additionally, one can modify the tolerance on the constraints using the tolConstraints component of critcontrol: an optional vector giving a tolerance for each of the constraints (equality or inequality). It is highly recommended to use it when there are equality constraints since the default tolerance of 0.05 (resp. 0 for inequality constraints) in such case might not be suited.

Author(s)

Victor Picheny

Mickael Binois

References

W.R. Jr. Mebane and J.S. Sekhon (2011), Genetic optimization using derivatives: The rgenoud package for R, Journal of Statistical Software, 42(11), 1-26

D.R. Jones, M. Schonlau, and W.J. Welch (1998), Efficient global optimization of expensive black-box functions, Journal of Global Optimization, 13, 455-492.

See Also

critcst_optimizer, crit_EFI, crit_AL, crit_SUR_cst

Examples

#---------------------------------------------------------------------------
# 2D objective function, 2 cases
#---------------------------------------------------------------------------

set.seed(2546)
library(DiceDesign)

n_var <- 2
fun <- branin

fun1.cst <- function(x){return(goldsteinprice(x)+.5)}
fun2.cst <- function(x){return(3/2 - x[1] - 2*x[2] - .5*sin(2*pi*(x[1]^2 - 2*x[2])))}

# Constraint function with vectorial output
cstfun <- function(x){return(c(fun1.cst(x), fun2.cst(x)))}

n.grid <- 31
test.grid <- expand.grid(X1 = seq(0, 1, length.out = n.grid), X2 = seq(0, 1, length.out = n.grid))
obj.grid <- apply(test.grid, 1, fun)
cst1.grid <- apply(test.grid, 1, fun1.cst)
cst2.grid <- apply(test.grid, 1, fun2.cst)

n_appr <- 12 
design.grid <- round(maximinESE_LHS(lhsDesign(n_appr, n_var, seed = 2)$design)$design, 1)
obj.init <- apply(design.grid, 1, fun)
cst1.init <- apply(design.grid, 1, fun1.cst)
cst2.init <- apply(design.grid, 1, fun2.cst)
model.fun <- km(~., design = design.grid, response = obj.init)
model.constraint1 <- km(~., design = design.grid, response = cst1.init, lower=c(.2,.2))
model.constraint2 <- km(~., design = design.grid, response = cst2.init, lower=c(.2,.2))
models.cst <- list(model.constraint1, model.constraint2)

lower <- rep(0, n_var)
upper <- rep(1, n_var)

#---------------------------------------------------------------------------
# Augmented Lagrangian Improvement, fast objective function, two ineq constraints,
# optimization with genoud
#---------------------------------------------------------------------------
critcontrol <- list(lambda=c(.5,2), rho=.5)
optimcontrol <- list(method = "genoud", max.generations=10, pop.size=20)

AL_grid <- apply(test.grid, 1, crit_AL, model.fun = fastfun(fun, design.grid),
                 model.constraint = models.cst, critcontrol=critcontrol)

cstEGO1 <- critcst_optimizer(crit = "AL", model.fun = fun,
                             model.constraint = models.cst, equality = FALSE,
                             lower = lower, upper = upper, 
                             optimcontrol = optimcontrol, critcontrol=critcontrol)

filled.contour(seq(0, 1, length.out = n.grid), seq(0, 1, length.out = n.grid), nlevels = 50,
               matrix(AL_grid, n.grid), main = "AL map and its maximizer (blue)",
               xlab = expression(x[1]), ylab = expression(x[2]), color = terrain.colors, 
               plot.axes = {axis(1); axis(2);
                            points(design.grid[,1], design.grid[,2], pch = 21, bg = "white")
                            contour(seq(0, 1, length.out = n.grid), seq(0, 1, length.out = n.grid), 
                            matrix(obj.grid, n.grid), nlevels = 10, add=TRUE,drawlabels=TRUE,
                                   col = "black")
                            contour(seq(0, 1, length.out = n.grid), seq(0, 1, length.out = n.grid), 
                            matrix(cst1.grid, n.grid), level = 0, add=TRUE,drawlabels=FALSE,
                                   lwd=1.5, col = "red")
                            contour(seq(0, 1, length.out = n.grid), seq(0, 1, length.out = n.grid), 
                            matrix(cst2.grid, n.grid), level = 0, add=TRUE,drawlabels=FALSE,
                                   lwd=1.5, col = "red")
                            points(cstEGO1$par, col = "blue", pch = 4, lwd = 2)
                            }
              )
#---------------------------------------------------------------------------
# SUR, expensive objective function, one equality constraint,
# optimization with genoud, integration on a regular grid
#---------------------------------------------------------------------------
optimcontrol <- list(method = "genoud", s = 40, maxit = 40)
critcontrol  <- list(tolConstraints = .15, integration.points=as.matrix(test.grid))

SUR_grid <- apply(test.grid, 1, crit_SUR_cst, model.fun = model.fun,
                  model.constraint = model.constraint1, equality = TRUE, critcontrol = critcontrol)

cstEGO2 <- critcst_optimizer(crit = "SUR", model.fun = model.fun,
                             model.constraint = model.constraint1, equality = TRUE,
                             lower = lower, upper = upper, 
                             optimcontrol = optimcontrol, critcontrol = critcontrol)

filled.contour(seq(0, 1, length.out = n.grid), seq(0, 1, length.out = n.grid), nlevels = 50,
               matrix(SUR_grid, n.grid), main = "SUR map and its maximizer (blue)",
               xlab = expression(x[1]), ylab = expression(x[2]), color = terrain.colors, 
               plot.axes = {axis(1); axis(2);
                            points(design.grid[,1], design.grid[,2], pch = 21, bg = "white")
                            contour(seq(0, 1, length.out = n.grid), seq(0, 1, length.out = n.grid), 
                            matrix(obj.grid, n.grid), nlevels = 10, add=TRUE,
                            drawlabels=TRUE, col = "black")
                            contour(seq(0, 1, length.out = n.grid), seq(0, 1, length.out = n.grid), 
                            matrix(cst1.grid, n.grid), level = c(-critcontrol$tolConstraints,
                            critcontrol$tolConstraints), 
                            add=TRUE, drawlabels=FALSE,lwd=1.5, col = "orange")
                            points(cstEGO2$par, col = "blue", pch = 4, lwd = 2)
                            }
              )


[Package DiceOptim version 2.1.1 Index]