easyEGO.cst {DiceOptim}R Documentation

EGO algorithm with constraints

Description

User-friendly wrapper of the function EGO.cst Generates initial DOEs and kriging models (objects of class km), and executes nsteps iterations of EGO methods integrating constraints.

Usage

easyEGO.cst(
  fun,
  constraint,
  n.cst = 1,
  budget,
  lower,
  upper,
  cheapfun = FALSE,
  equality = FALSE,
  X = NULL,
  y = NULL,
  C = NULL,
  control = list(method = "EFI", trace = 1, inneroptim = "genoud", maxit = 100, seed =
    42),
  ...
)

Arguments

fun

scalar function to be minimized,

constraint

vectorial function corresponding to the constraints, see details below,

n.cst

number of constraints,

budget

total number of calls to the objective and constraint functions,

lower

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

upper

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

cheapfun

optional boolean, TRUE if the objective is a fast-to-evaluate function that does not need a kriging model

equality

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

X

initial design of experiments. If not provided, X is taken as a maximin LHD with budget/3 points

y

initial set of objective observations f(X). Computed if not provided.

C

initial set of constraint observations g(X). Computed if not provided.

control

an optional list of control parameters. See "Details".

...

additional parameters to be given to BOTH the objective fun and constraints.

Details

Does not require knowledge on kriging models (objects of class km)

The problem considered is of the form: min f(x) s.t. g(x) \le 0, g having a vectorial output. By default all its components are supposed to be inequalities, but one can use a Boolean vector in equality to specify which are equality constraints, hence of the type g(x) = 0. The control argument is a list that can supply any of the following components:

For additional details, see EGO.cst.

Value

A list with components:

Author(s)

Victor Picheny

Mickael Binois

References

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.

M. Schonlau, W.J. Welch, and D.R. Jones (1998), Global versus local search in constrained optimization of computer models, Lecture Notes-Monograph Series, 11-25.

M.J. Sasena, P. Papalambros, and P.Goovaerts (2002), Exploration of metamodeling sampling criteria for constrained global optimization, Engineering optimization, 34, 263-278.

R.B. Gramacy, G.A. Gray, S. Le Digabel, H.K.H Lee, P. Ranjan, G. Wells, Garth, and S.M. Wild (2014+), Modeling an augmented Lagrangian for improved blackbox constrained optimization, arXiv preprint arXiv:1403.4890.

J.M. Parr (2012), Improvement criteria for constraint handling and multiobjective optimization, University of Southampton.

V. Picheny (2014), A stepwise uncertainty reduction approach to constrained global optimization, Proceedings of the 17th International Conference on Artificial Intelligence and Statistics, JMLR W&CP 33, 787-795.

Examples

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

set.seed(25468)
library(DiceDesign)

n_var <- 2 
fun <- goldsteinprice
fun1.cst <- function(x){return(-branin(x) + 25)}
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
constraint <- function(x){return(c(fun1.cst(x), fun2.cst(x)))}

# For illustration purposes
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)

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

#---------------------------------------------------------------------------
# 1- Expected Feasible Improvement criterion, expensive objective function,
# two inequality constraints, 15 observations budget, using genoud
#---------------------------------------------------------------------------
res <- easyEGO.cst(fun=fun, constraint=constraint, n.cst=2, lower=lower, upper=upper, budget=15, 
                   control=list(method="EFI", inneroptim="genoud", maxit=20))

cat("best design found:", res$par, "\n")
cat("corresponding objective and constraints:", res$value, "\n")

# Objective function in colour, constraint boundaries in red
# Initial DoE: white circles, added points: blue crosses, best solution: red cross

filled.contour(seq(0, 1, length.out = n.grid), seq(0, 1, length.out = n.grid), nlevels = 50,
               matrix(obj.grid, n.grid), main = "Two inequality constraints",
               xlab = expression(x[1]), ylab = expression(x[2]), color = terrain.colors, 
               plot.axes = {axis(1); axis(2);
                            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(res$history$X, col = "blue", pch = 4, lwd = 2)       
                            points(res$par[1], res$par[2], col = "red", pch = 4, lwd = 2, cex=2) 
               }
)

#---------------------------------------------------------------------------
# 2- Augmented Lagrangian Improvement criterion, expensive objective function,
# one inequality and one equality constraint, 25 observations budget, using random search
#---------------------------------------------------------------------------
res2 <- easyEGO.cst(fun=fun, constraint=constraint, n.cst=2, lower=lower, upper=upper, budget=25,
                   equality = c(TRUE, FALSE),
                   control=list(method="AL", inneroptim="random", maxit=100))

cat("best design found:", res2$par, "\n")
cat("corresponding objective and constraints:", res2$value, "\n")

# Objective function in colour, inequality constraint boundary in red, equality
# constraint in orange
# Initial DoE: white circles, added points: blue crosses, best solution: red cross

filled.contour(seq(0, 1, length.out = n.grid), seq(0, 1, length.out = n.grid), nlevels = 50,
               matrix(obj.grid, n.grid), xlab = expression(x[1]), ylab = expression(x[2]),
               main = "Inequality (red) and equality (orange) constraints", color = terrain.colors, 
               plot.axes = {axis(1); axis(2);
                            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 = "orange")
                            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(res2$history$X, col = "blue", pch = 4, lwd = 2)
                            points(res2$par[1], res2$par[2], col = "red", pch = 4, lwd = 2, cex=2) 
               }
)

#---------------------------------------------------------------------------
# 3- Stepwise Uncertainty Reduction criterion, fast objective function,
# single inequality constraint, with initial DOE given + 10 observations,
# using genoud
#---------------------------------------------------------------------------
n.init <- 12 
design.grid <- round(maximinESE_LHS(lhsDesign(n.init, n_var, seed = 42)$design)$design, 1)
cst2.init <- apply(design.grid, 1, fun2.cst)

res3 <- easyEGO.cst(fun=fun, constraint=fun2.cst, n.cst=1, lower=lower, upper=upper, budget=10,
                    X=design.grid, C=cst2.init,
                   cheapfun=TRUE, control=list(method="SUR", inneroptim="genoud", maxit=20))
        
cat("best design found:", res3$par, "\n")
cat("corresponding objective and constraint:", res3$value, "\n")

# Objective function in colour, inequality constraint boundary in red
# Initial DoE: white circles, added points: blue crosses, best solution: red cross
                            
filled.contour(seq(0, 1, length.out = n.grid), seq(0, 1, length.out = n.grid), nlevels = 50,
               matrix(obj.grid, n.grid), main = "Single constraint, fast objective",
               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)
                            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(res3$history$X, col = "blue", pch = 4, lwd = 2)
                            points(res3$par[1], res3$par[2], col = "red", pch = 4, lwd = 2, cex=2) 
                                           }
                )



[Package DiceOptim version 2.1.1 Index]