max_qEI | R Documentation |
Maximization of the qEI
criterion. This
is a refactoring of the function from the DiceOptim
package. However only the "genoud"
method can be used
for now. This function calls the function genoud
using analytical formulae of qEI
and its
gradient qEI.grad
.
max_qEI( model, npoints, lower, upper, crit = c("exact", "CL"), minimization = TRUE, optimcontrol = NULL, trace = 1 )
model |
An object in |
npoints |
Integer representing the desired number of iterations. |
lower, upper |
Numeric vectors of lower and upper bounds. |
crit |
Character |
minimization |
Logical specifying if the qEI to be maximized is used in minimization or in maximization. |
optimcontrol |
an optional list of control parameters for optimization. See details. |
trace |
Integer. Level of verbosity. |
The parameters of list optimcontrol
include
optimcontrol$method
, with character value "BFGS"
(default) or "genoud"
, specifying the method used to
maximize the criterion. This is irrelevant when crit
is
"CL"
because this method always uses genoud
.
When crit == "CL"
The elements of optimcontrol
can have the following
names and values.
parinit
Optional numeric
matrix of initial values. Must have
model@d
columns, the number of rows is not
constrained.
L
The "liar". Either a
character in c("max", "min", "mean")
, or
a scalar value specifying the liar. For the
character values: "min"
sets the liar to
min(model@y)
, "max"
set the liar to
max(model@y)
, and "mean" set it to the
prediction of the model. When L
is
NULL
, it is replaced by "min"
if
minimization
is TRUE
, or by "max"
otherwise.
Formal arguments of genoud
.
These include pop.size
with
default : 3 * 2^d
when d <
6
and 32 * d
otherwise, where d
is the number of inputs. One can also use max.generations
(default: 12), wait.generations
(default:
2) and BFGSburnin
(default: 2).
When optimcontrol$method == "BFGS"
The elements of optimcontrol
can have the
following names and values.
nStarts
(default: 4).
fastCompute
Logical. If TRUE
(default), a fast
approximation method based on a semi-analytic
formula is used, see Marmin (2014) for details.
samplingFun
A function which
samples a batch of starting point
Default : sampleFromEI
.
parinit
Optional 3d-array of
initial (or candidate)
batches. For each k
, the slice parinit[ , , k]
is a matrix of size c(npoints, d)
representing one batch. The number of initial
batches dim(parinit)[3]
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
Logical. If
TRUE
(default), a fast approximation method
based on a semi-analytic formula is used, see
Marmin(2014) for details.
parinit
Optional numeric matrix of candidate starting
points (one row corresponds to one point).
The formal arguments of the genoud
function
Main parameters are pop.size
(default:
50 * d * npoints
),
max.generations
(default: 5),
wait.generations
(default: 2),
BFGSburnin
(default: 2).
A list with components:
par A matrix with its rows containing the npoints
input vectors found.
value A value giving the qEI computed in par
.
this function may well be renamed
max_qEI_genoud
in a future version.
Sébastien Marmin, Clément Chevalier and David Ginsbourger.
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.
qEI_with_grad
.
set.seed(0) ## 3-points EI maximization. ## 9-points factorial design, and the corresponding response d <- 2; n <- 9 design.fact <- expand.grid(x1 = seq(0, 1, length = 3), x2 = seq(0, 1, length = 3)) response.branin <- apply(design.fact, 1, branin) lower <- c(0, 0); upper <- c(1, 1) ## number of point in the batch batchSize <- 3 ## model fit 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 ## ================================== ## Not run: maxBFGS <- max_qEI(model = fitted.model, npoints = batchSize, lower = lower, upper = upper, crit = "exact", optimcontrol = list(nStarts = 3, method = "BFGS")) print(maxBFGS$value) ## End(Not run) ## With a genetic algorithm 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)) print(maxGen$value) ## With the constant liar heuristic ## ================================ ## Not run: maxCL <- max_qEI(model = fitted.model, npoints = batchSize, lower = lower, upper = upper, crit = "CL", optimcontrol = list(pop.size = 20)) print(maxCL$value) ## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.