| crit_SUR_Eq {GPGame} | R Documentation | 
SUR criterion for equilibria
Description
Computes the SUR criterion associated to an equilibrium for a given xnew and a set of trajectories of objective functions
on a predefined grid.
Usage
crit_SUR_Eq(
  idx,
  model,
  integcontrol,
  Simu,
  precalc.data = NULL,
  equilibrium,
  n.ynew = NULL,
  cross = FALSE,
  IS = FALSE,
  plot = FALSE,
  kweights = NULL,
  Nadir = NULL,
  Shadow = NULL,
  calibcontrol = NULL
)
Arguments
| idx | is the index on the grid of the strategy evaluated | 
| model | is a list of  | 
| integcontrol | is a list containing:  | 
| Simu | is a matrix of size [ | 
| precalc.data | is a list of length  | 
| equilibrium | equilibrium type: either " | 
| n.ynew | is the number of  | 
| cross | if  | 
| IS | if  | 
| plot | if  | 
| kweights | kriging weights for  | 
| Nadir,Shadow | optional vectors of size  | 
| calibcontrol | an optional list for calibration problems, containing  | 
Value
Criterion value.
References
V. Picheny, M. Binois, A. Habbal (2016+), A Bayesian optimization approach to find Nash equilibria, https://arxiv.org/abs/1611.02440.
See Also
crit_PNash for an alternative infill criterion
Examples
##############################################
# 2 variables, 2 players
##############################################
library(DiceKriging)
set.seed(42)
# Objective function (R^2 -> R^2)
fun <- 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)))
}
# Grid definition
n.s <- rep(14, 2)
x.to.obj   <- c(1,2)
gridtype <- 'cartesian'
integcontrol <- generate_integ_pts(n.s=n.s, d=4, nobj=2, x.to.obj = x.to.obj, gridtype=gridtype)
integ.pts <- integcontrol$integ.pts
expanded.indices <- integcontrol$expanded.indices
# Kriging models
n.init <- 11
design <- integ.pts[sample.int(n=nrow(integ.pts), size=n.init, replace=FALSE),]
response <- t(apply(design, 1, fun))
mf1 <- km(~., design = design, response = response[,1], lower=c(.1,.1))
mf2 <- km(~., design = design, response = response[,2], lower=c(.1,.1))
model <- list(mf1, mf2)
# Conditional simulations
Simu <- t(Reduce(rbind, lapply(model, simulate, nsim=10, newdata=integ.pts, cond=TRUE,
                                   checkNames=FALSE, nugget.sim = 10^-8)))
# Useful precalculations with the package KrigInv can be reused for computational speed.
# library(KrigInv)
# precalc.data <- lapply(model, FUN=KrigInv:::precomputeUpdateData, integration.points=integ.pts)
# Compute criterion for all points on the grid
crit_grid <- lapply(X=1:prod(n.s), FUN=crit_SUR_Eq, model=model,
                    integcontrol=integcontrol, equilibrium = "NE",
                    # precalc.data=precalc.data, # Uncomment if precalc.data is computed
                    Simu=Simu, n.ynew=10, IS=FALSE, cross=FALSE)
crit_grid <- unlist(crit_grid)
# Draw contour of the criterion
filled.contour(seq(0, 1, length.out = n.s[1]), seq(0, 1, length.out = n.s[2]),
               matrix(pmax(0, crit_grid), n.s[1], n.s[2]), main = "SUR criterion",
               xlab = expression(x[1]), ylab = expression(x[2]), color = terrain.colors,
               plot.axes = {axis(1); axis(2);
                            points(design[,1], design[,2], pch = 21, bg = "white")
                           }
)