easyGParetoptim {GPareto} | R Documentation |
EGO algorithm for multiobjective optimization
Description
User-friendly wrapper of the function GParetoptim
.
Generates initial DOEs and kriging models (objects of class km
),
and executes nsteps
iterations of multiobjective EGO methods.
Usage
easyGParetoptim(
fn,
...,
cheapfn = NULL,
budget,
lower,
upper,
par = NULL,
value = NULL,
noise.var = NULL,
control = list(method = "SMS", trace = 1, inneroptim = "pso", maxit = 100, seed = 42),
ncores = 1
)
Arguments
fn |
the multi-objective function to be minimized (vectorial output), found by a call to |
... |
additional parameters to be given to the objective |
cheapfn |
optional additional fast-to-evaluate objective function (handled next with class |
budget |
total number of calls to the objective function, |
lower |
vector of lower bounds for the variables to be optimized over, |
upper |
vector of upper bounds for the variables to be optimized over, |
par |
initial design of experiments. If not provided, |
value |
initial set of objective observations |
noise.var |
optional noise variance, for noisy objectives |
control |
an optional list of control parameters. See "Details", |
ncores |
number of CPU available (> 1 makes mean parallel |
Details
Does not require specific knowledge on kriging models (objects of class km
).
The problem considered is of the form: min f(x) = f_1(x), ..., f_p(x)
.
The control
argument is a list that can supply any of the following optional components:
-
method
: choice of multiobjective improvement function: "SMS
", "EHI
", "EMI
" or "SUR
" (seecrit_SMS
,crit_EHI
,crit_EMI
,crit_SUR
), -
trace
: if positive, tracing information on the progress of the optimization is produced (1
(default) for general progress,>1
for more details, e.g., warnings fromgenoud
), -
inneroptim
: choice of the inner optimization algorithm: "genoud
", "pso
" or "random
" (seegenoud
andpsoptim
), -
maxit
: maximum number of iterations of the inner loop, -
seed
: to fix the random variable generator, -
refPoint
: reference point for hypervolume computations (for "SMS
" and "EHI
" methods), -
extendper
: if no reference pointrefPoint
is provided, for each objective it is fixed to the maximum over the Pareto front plus extendper times the range. Default value to0.2
, corresponding to1.1
for a scaled objective with a Pareto front in[0,1]^n.obj
.
If noise.var="given_by_fn"
, fn
must return a list of two vectors, the first being the objective functions and the second
the corresponding noise variances. See examples in GParetoptim
.
For additional details or other possible arguments, see GParetoptim
.
Display of results and various post-processings are available with plotGPareto
.
Value
A list with components:
-
par
: all the non-dominated points found, -
value
: the matrix of objective values at the points given inpar
, -
history
: a list containing all the points visited by the algorithm (X
) and their corresponding objectives (y
), -
model
: a list of objects of classkm
, corresponding to the last kriging models fitted.
Note that in the case of noisy problems, value
and history$y.denoised
are denoised values. The original observations are available in the slot
history$y
.
Author(s)
Victor Picheny (INRA, Castanet-Tolosan, France)
Mickael Binois (Mines Saint-Etienne/Renault, France)
References
M. T. Emmerich, A. H. Deutz, J. W. Klinkenberg (2011), Hypervolume-based expected improvement: Monotonicity properties and exact computation,
Evolutionary Computation (CEC), 2147-2154.
V. Picheny (2015), Multiobjective optimization using Gaussian process emulators via stepwise uncertainty reduction,
Statistics and Computing, 25(6), 1265-1280.
T. Wagner, M. Emmerich, A. Deutz, W. Ponweiser (2010), On expected-improvement criteria for model-based multi-objective optimization.
Parallel Problem Solving from Nature, 718-727, Springer, Berlin.
J. D. Svenson (2011), Computer Experiments: Multiobjective Optimization and Sensitivity Analysis, Ohio State university, PhD thesis.
M. Binois, V. Picheny (2019), GPareto: An R Package for Gaussian-Process-Based Multi-Objective Optimization and Analysis,
Journal of Statistical Software, 89(8), 1-30, doi:10.18637/jss.v089.i08.
Examples
#---------------------------------------------------------------------------
# 2D objective function, 4 cases
#---------------------------------------------------------------------------
## Not run:
set.seed(25468)
n_var <- 2
fname <- ZDT3
lower <- rep(0, n_var)
upper <- rep(1, n_var)
#---------------------------------------------------------------------------
# 1- Expected Hypervolume Improvement optimization, using pso
#---------------------------------------------------------------------------
res <- easyGParetoptim(fn=fname, lower=lower, upper=upper, budget=15,
control=list(method="EHI", inneroptim="pso", maxit=20))
par(mfrow=c(1,2))
plotGPareto(res)
title("Pareto Front")
plot(res$history$X, main="Pareto set", col = "red", pch = 20)
points(res$par, col="blue", pch = 17)
#---------------------------------------------------------------------------
# 2- SMS Improvement optimization using random search, with initial DOE given
#---------------------------------------------------------------------------
library(DiceDesign)
design.init <- maximinESE_LHS(lhsDesign(10, n_var, seed = 42)$design)$design
response.init <- t(apply(design.init, 1, fname))
res <- easyGParetoptim(fn=fname, par=design.init, value=response.init, lower=lower, upper=upper,
budget=15, control=list(method="SMS", inneroptim="random", maxit=100))
par(mfrow=c(1,2))
plotGPareto(res)
title("Pareto Front")
plot(res$history$X, main="Pareto set", col = "red", pch = 20)
points(res$par, col="blue", pch = 17)
#---------------------------------------------------------------------------
# 3- Stepwise Uncertainty Reduction optimization, with one fast objective function
#---------------------------------------------------------------------------
fname <- camelback
cheapfn <- function(x) {
if (is.null(dim(x))) return(-sum(x))
else return(-rowSums(x))
}
res <- easyGParetoptim(fn=fname, cheapfn=cheapfn, lower=lower, upper=upper, budget=15,
control=list(method="SUR", inneroptim="pso", maxit=20))
par(mfrow=c(1,2))
plotGPareto(res)
title("Pareto Front")
plot(res$history$X, main="Pareto set", col = "red", pch = 20)
points(res$par, col="blue", pch = 17)
#---------------------------------------------------------------------------
# 4- Expected Hypervolume Improvement optimization, using pso, noisy fn
#---------------------------------------------------------------------------
noise.var <- c(0.1, 0.2)
funnoise <- function(x) {ZDT3(x) + sqrt(noise.var)*rnorm(n=2)}
res <- easyGParetoptim(fn=funnoise, lower=lower, upper=upper, budget=30, noise.var=noise.var,
control=list(method="EHI", inneroptim="pso", maxit=20))
par(mfrow=c(1,2))
plotGPareto(res)
title("Pareto Front")
plot(res$history$X, main="Pareto set", col = "red", pch = 20)
points(res$par, col="blue", pch = 17)
#---------------------------------------------------------------------------
# 5- Stepwise Uncertainty Reduction optimization, functional noise
#---------------------------------------------------------------------------
funnoise <- function(x) {ZDT3(x) + sqrt(abs(0.1*x))*rnorm(n=2)}
noise.var <- function(x) return(abs(0.1*x))
res <- easyGParetoptim(fn=funnoise, lower=lower, upper=upper, budget=30, noise.var=noise.var,
control=list(method="SUR", inneroptim="pso", maxit=20))
par(mfrow=c(1,2))
plotGPareto(res)
title("Pareto Front")
plot(res$history$X, main="Pareto set", col = "red", pch = 20)
points(res$par, col="blue", pch = 17)
## End(Not run)