solve_game {GPGame} | R Documentation |
Main solver
Description
Main function to solve games.
Usage
solve_game(
fun,
...,
equilibrium = "NE",
crit = "sur",
model = NULL,
n.init = NULL,
n.ite,
d,
nobj,
x.to.obj = NULL,
noise.var = NULL,
Nadir = NULL,
Shadow = NULL,
integcontrol = NULL,
simucontrol = NULL,
filtercontrol = NULL,
kmcontrol = NULL,
returncontrol = NULL,
ncores = 1,
trace = 1,
seed = NULL,
calibcontrol = NULL,
freq.exploit = 1000
)
Arguments
fun |
fonction with vectorial output |
... |
additional parameter to be passed to |
equilibrium |
either ' |
crit |
' |
model |
list of |
n.init |
number of points of the initial design of experiments if no model is given |
n.ite |
number of iterations of sequential optimization |
d |
variable dimension |
nobj |
number of objectives (players) |
x.to.obj |
for |
noise.var |
noise variance. Either a scalar (same noise for all objectives), a vector (constant noise, different for each objective),
a function (type closure) with vectorial output (variable noise, different for each objective) or |
Nadir , Shadow |
optional vectors of size |
integcontrol |
optional list for handling integration points. See Details. |
simucontrol |
optional list for handling conditional simulations. See Details. |
filtercontrol |
optional list for handling filters. See Details. |
kmcontrol |
optional list for handling |
returncontrol |
optional list for choosing return options. See Details. |
ncores |
number of CPU available (> 1 makes mean parallel |
trace |
controls the level of printing: |
seed |
to fix the random variable generator |
calibcontrol |
an optional list for calibration problems, containing |
freq.exploit |
an optional integer to force exploitation (i.e. evaluation of the predicted equilibrium) every |
Details
If noise.var="given_by_fn"
, fn
returns a list of two vectors, the first being the objective functions and the second
the corresponding noise variances.
integcontrol
controls the way the design space is discretized. One can directly provide a set of points integ.pts
with
corresponding indices expanded.indices
(for NE
). Otherwise, the points are generated according to the number of strategies n.s
.
If n.s
is a scalar, it corresponds to the total number of strategies (to be divided equally among players),
otherwise it corresponds to the nb of strategies per player. In addition, one may choose the type of discretization with gridtype
.
Options are 'lhs
' or 'cartesian
'. Finally, lb
and ub
are vectors specifying the bounds for the design variables.
By default the design space is [0,1]^d
.
A renew
slot is available, if TRUE
, then integ.pts
are changed at each iteration. Available only for KSE
and CKSE
.
For CKSE
, setting the slot kweights=TRUE
allows to increase the number of integration points, with nsamp
(default to 1e4
) virtual simulation points.
simucontrol
controls options on conditional GP simulations. Options are IS
: if TRUE
, importance sampling is used for ynew
;
n.ynew
number of samples of Y(x_{n+1})
and n.sim
number of sample path generated.
filtercontrol
controls filtering options. filter
sets how to select a subset of simulation and candidate points,
either either a single value or a vector of two to use different filters for simulation and candidate points.
Possible values are 'window
', 'Pnash
' (for NE
), 'PND
' (probability of non domination), 'none
'.
nsimPoints
and ncandPoints
set the maximum number of simulation/candidate points wanted
(use with filter 'Pnash
' for now). Default values are 800
and 200
, resp.
randomFilter
(TRUE
by default except for filter window
) sets whereas the filter acts randomly or deterministically.
For more than 3 objectives, PND
is estimated by sampling; the number of samples is controled by nsamp
(default to max(20, 5 * nobj)
).
kmcontrol
Options for handling nobj
km
models.
cov.reestim
(Boolean, TRUE
by default) specifies if the kriging hyperparameters
should be re-estimated at each iteration,
returncontrol
sets options for the last iterations and what is returned by the algorithm.
track.Eq
allows to estimate the equilibrium at each iteration; options are 'none
' to do nothing,
"mean
" (default) to compute the equilibrium of the prediction mean (all candidates),
"empirical
" (for KSE
) and "pex
"/"psim
" (NE
only)
for using Pnash
estimate (along with mean estimate, on integ.pts only, NOT reestimated if filter.simu
or crit
is Pnash
).
The boolean force.exploit.last
(default to TRUE
) allows to evaluate the equilibrium on the predictive mean - if not already evaluated -
instead of using crit
(i.e., sur
) for KSE
and CKSE
.
Value
A list with components:
model
: a list of objects of classkm
corresponding to the last kriging models fitted.Jplus
: recorded values of the acquisition function maximizerinteg.pts
andexpanded.indices
: the discrete space used,predEq
: a list containing the recorded values of the estimated best solution,Eq.design, Eq.poff
: estimated equilibrium and corresponding pay-off
Note: with CKSE, kweights are not used when the mean on integ.pts is used. Also, CKSE does not support non-constant mean at this stage.
References
V. Picheny, M. Binois, A. Habbal (2016+), A Bayesian optimization approach to find Nash equilibria, https://arxiv.org/abs/1611.02440.
Examples
################################################################
# Example 1: Nash equilibrium, 2 variables, 2 players, no filter
################################################################
# Define objective function (R^2 -> R^2)
fun1 <- function (x)
{
if (is.null(dim(x))) x <- matrix(x, nrow = 1)
b1 <- 15 * x[, 1] - 5
b2 <- 15 * x[, 2]
return(cbind((b2 - 5.1*(b1/(2*pi))^2 + 5/pi*b1 - 6)^2 + 10*((1 - 1/(8*pi)) * cos(b1) + 1),
-sqrt((10.5 - b1)*(b1 + 5.5)*(b2 + 0.5)) - 1/30*(b2 - 5.1*(b1/(2*pi))^2 - 6)^2-
1/3 * ((1 - 1/(8 * pi)) * cos(b1) + 1)))
}
# To use parallel computation (turn off on Windows)
library(parallel)
parallel <- FALSE #TRUE #
if(parallel) ncores <- detectCores() else ncores <- 1
# Simple configuration: no filter, discretization is a 21x21 grid
# Grid definition
n.s <- rep(21, 2)
x.to.obj <- c(1,2)
gridtype <- 'cartesian'
# Run solver with 6 initial points, 4 iterations
# Increase n.ite to at least 10 for better results
res <- solve_game(fun1, equilibrium = "NE", crit = "sur", n.init=6, n.ite=4,
d = 2, nobj=2, x.to.obj = x.to.obj,
integcontrol=list(n.s=n.s, gridtype=gridtype),
ncores = ncores, trace=1, seed=1)
# Get estimated equilibrium and corresponding pay-off
NE <- res$Eq.design
Poff <- res$Eq.poff
# Draw results
plotGame(res)
################################################################
# Example 2: same example, KS equilibrium with given Nadir
################################################################
# Run solver with 6 initial points, 4 iterations
# Increase n.ite to at least 10 for better results
res <- solve_game(fun1, equilibrium = "KSE", crit = "sur", n.init=6, n.ite=4,
d = 2, nobj=2, x.to.obj = x.to.obj,
integcontrol=list(n.s=400, gridtype="lhs"),
ncores = ncores, trace=1, seed=1, Nadir=c(Inf, -20))
# Get estimated equilibrium and corresponding pay-off
NE <- res$Eq.design
Poff <- res$Eq.poff
# Draw results
plotGame(res, equilibrium = "KSE", Nadir=c(Inf, -20))
################################################################
# Example 3: Nash equilibrium, 4 variables, 2 players, filtering
################################################################
fun2 <- function(x, nobj = 2){
if (is.null(dim(x))) x <- matrix(x, 1)
y <- matrix(x[, 1:(nobj - 1)], nrow(x))
z <- matrix(x[, nobj:ncol(x)], nrow(x))
g <- rowSums((z - 0.5)^2)
tmp <- t(apply(cos(y * pi/2), 1, cumprod))
tmp <- cbind(t(apply(tmp, 1, rev)), 1)
tmp2 <- cbind(1, t(apply(sin(y * pi/2), 1, rev)))
return(tmp * tmp2 * (1 + g))
}
# Grid definition: player 1 plays x1 and x2, player 2 x3 and x4
# The grid is a lattice made of two LHS designs of different sizes
n.s <- c(44, 43)
x.to.obj <- c(1,1,2,2)
gridtype <- 'lhs'
# Set filtercontrol: window filter applied for integration and candidate points
# 500 simulation and 200 candidate points are retained.
filtercontrol <- list(nsimPoints=500, ncandPoints=200,
filter=c("window", "window"))
# Set km control: lower bound is specified for the covariance range
# Covariance type and model trend are specified
kmcontrol <- list(lb=rep(.2,4), model.trend=~1, covtype="matern3_2")
# Run solver with 20 initial points, 4 iterations
# Increase n.ite to at least 20 for better results
res <- solve_game(fun2, equilibrium = "NE", crit = "psim", n.init=20, n.ite=2,
d = 4, nobj=2, x.to.obj = x.to.obj,
integcontrol=list(n.s=n.s, gridtype=gridtype),
filtercontrol=filtercontrol,
kmcontrol=kmcontrol,
ncores = 1, trace=1, seed=1)
# Get estimated equilibrium and corresponding pay-off
NE <- res$Eq.design
Poff <- res$Eq.poff
# Draw results
plotGame(res)
################################################################
# Example 4: same example, KS equilibrium
################################################################
# Grid definition: simple lhs
integcontrol=list(n.s=1e4, gridtype='lhs')
# Run solver with 20 initial points, 4 iterations
# Increase n.ite to at least 20 for better results
res <- solve_game(fun2, equilibrium = "KSE", crit = "sur", n.init=20, n.ite=2,
d = 4, nobj=2,
integcontrol=integcontrol,
filtercontrol=filtercontrol,
kmcontrol=kmcontrol,
ncores = 1, trace=1, seed=1)
# Get estimated equilibrium and corresponding pay-off
NE <- res$Eq.design
Poff <- res$Eq.poff
# Draw results
plotGame(res, equilibrium = "KSE")