get_active_in {funGp}R Documentation

Extraction of active inputs in a given model structure

Description

The fgpm_factory function returns an object of class "Xfgpm" with the function call of all the evaluated models stored in the @log.success@args and @log.crashes@args slots. The get_active_in function interprets the arguments linked to any structural configuration and returns a list with two elements: (i) a matrix of scalar input variables kept active; and (ii) a list of functional input variables kept active.

Usage

get_active_in(sIn = NULL, fIn = NULL, args)

Arguments

sIn

An optional matrix of scalar input coordinates with all the orignal scalar input variables.

fIn

An optional list of functional input coordinates with all the original functional input variables.

args

An object of class "modelCall", which specifies the model structure for which the active inputs should be extracted.

Value

An object of class "list", containing the following information extracted from the args parameter: (i) a matrix of scalar input variables kept active; and (ii) a list of functional input variables kept active.

Author(s)

José Betancourt, François Bachoc, Thierry Klein and Jérémy Rohmer

References

Betancourt, J., Bachoc, F., Klein, T., Idier, D., Rohmer, J., and Deville, Y. (2024), "funGp: An R Package for Gaussian Process Regression with Scalar and Functional Inputs". Journal of Statistical Software, 109, 5, 1–51. (doi:10.18637/jss.v109.i05)

Betancourt, J., Bachoc, F., and Klein, T. (2020), R Package Manual: "Gaussian Process Regression for Scalar and Functional Inputs with funGp - The in-depth tour". RISCOPE project. [HAL]

See Also

* which_on for details on how to obtain only the indices of the active inputs.

* modelCall for details on the args argument.

* fgpm_factory for funGp heuristic model selection.

* Xfgpm for details on object delivered by fgpm_factory.

Examples

# Use precalculated Xfgpm object named xm
# indices of active inputs in the best model
xm@log.success@args[[1]] # the full fgpm call
set.seed(100)
n.tr <- 32
sIn <- expand.grid(x1 = seq(0,1,length = n.tr^(1/5)), x2 = seq(0,1,length = n.tr^(1/5)),
x3 = seq(0,1,length = n.tr^(1/5)), x4 = seq(0,1,length = n.tr^(1/5)),
x5 = seq(0,1,length = n.tr^(1/5)))
fIn <- list(f1 = matrix(runif(n.tr*10), ncol = 10), f2 = matrix(runif(n.tr*22), ncol = 22))
which_on(sIn, fIn, xm@log.success@args[[1]]) # only the indices extracted by which_on

# data structures of active inputs
active <- get_active_in(sIn, fIn, xm@log.success@args[[1]])
active$sIn.on # scalar data structures
active$fIn.on # functional data structures
# identifying selected model and corresponding fgpm arguments
opt.model <- xm@model
opt.args <- xm@log.success@args[[1]]

# generating new input data for prediction
n.pr <- 243
sIn.pr <- expand.grid(x1 = seq(0,1,length = n.pr^(1/5)), x2 = seq(0,1,length = n.pr^(1/5)),
                      x3 = seq(0,1,length = n.pr^(1/5)), x4 = seq(0,1,length = n.pr^(1/5)),
                      x5 = seq(0,1,length = n.pr^(1/5)))
fIn.pr <- list(f1 = matrix(runif(n.pr*10), ncol = 10), f2 = matrix(runif(n.pr*22), ncol = 22))

# pruning data structures for prediction to keep only active inputs!!
active <- get_active_in(sIn.pr, fIn.pr, opt.args)

# making predictions
preds <- predict(opt.model, sIn.pr = active$sIn.on, fIn.pr = active$fIn.on)

# plotting predictions
plot(preds)


# preparing new data for simulation based on inputs kept active____________________________
opt.model <- xm@model
opt.args <- xm@log.success@args[[1]]

# generating new input data for simulation
n.sm <- 243
sIn.sm <- expand.grid(x1 = seq(0,1,length = n.pr^(1/5)), x2 = seq(0,1,length = n.pr^(1/5)),
                      x3 = seq(0,1,length = n.pr^(1/5)), x4 = seq(0,1,length = n.pr^(1/5)),
                      x5 = seq(0,1,length = n.pr^(1/5)))
fIn.sm <- list(f1 = matrix(runif(n.sm*10), ncol = 10), f2 = matrix(runif(n.sm*22), ncol = 22))

# pruning data structures for simulation to keep only active inputs!!
active <- get_active_in(sIn.sm, fIn.sm, opt.args)

# making light simulations
sims_l <- simulate(opt.model, nsim = 10, sIn.sm = active$sIn.on, fIn.sm = active$fIn.on)

# plotting light simulations
plot(sims_l)

## Not run: 
# rebuilding of 3 best models using new data_______________________________________________
# NOTE: this example is of higher complexity than the previous ones. We recomend you run
#       the previous examples and understand the @log.success and @log.crashes slots in
#       the Xfgpm object delivered by fgpm_factory.
#
#       In the second example above we showed how to use get_active_in to prune the input
#       data structures for prediction based on the fgpm arguments of the best model found
#       by fgpm_factory. In this new example we generalize that concept by: (i) rebuilding
#       the 3 best models found by fgpm_factory using new data, (ii) pruning the input
#       data structures used for prediction with each of the models, and (iii) plotting
#       the predictions made by the three models. The key ingredient here is that the
#       three best models might have different scalar and functional inputs active. The
#       get_active_in function will allow to process the data structures in order to
#       extract only the scalar inputs required to re-build the model and then to make
#       predictions with each model. Check also the funGp manual for further details
#
#       funGp manual: https://doi.org/10.18637/jss.v109.i05


# <<<<<<< PART 1: calling fgpm_factory to perform the structural optimization >>>>>>>
#         -------------------------------------------------------------------
# this part is precalculated and loaded via data("precalculated_Xfgpm_objects")
summary(xm)

# <<<<<<< PART 2: re-building the three best models found by fgpm_factory >>>>>>>
#         ---------------------------------------------------------------
# recovering the fgpm arguments of the three best models
argStack <- xm@log.success@args[1:3]

# new data arrived, now we have 243 observations
n.nw <- 243 # more points!
sIn.nw <- expand.grid(x1 = seq(0,1,length = n.nw^(1/5)), x2 = seq(0,1,length = n.nw^(1/5)),
                      x3 = seq(0,1,length = n.nw^(1/5)), x4 = seq(0,1,length = n.nw^(1/5)),
                      x5 = seq(0,1,length = n.nw^(1/5)))
fIn.nw <- list(f1 = matrix(runif(n.nw*10), ncol = 10), f2 = matrix(runif(n.nw*22), ncol = 22))
sOut.nw <- fgp_BB7(sIn.nw, fIn.nw, n.nw)

# the second best model
modelDef(xm,2)
# re-building the three best models based on the new data (compact code with all 3 calls)
newEnv <- list(sIn = sIn.nw, fIn = fIn.nw, sOut = sOut.nw)
modStack <- lapply(1:3, function(i) eval(parse(text =  modelDef(xm,i)), env = newEnv))


# <<<<<<< PART 3: making predictions from the three best models found by fgpm_factory >>>>>>>
#         ---------------------------------------------------------------------------
# generating input data for prediction
n.pr <- 32
sIn.pr <- expand.grid(x1 = seq(0,1,length = n.pr^(1/5)), x2 = seq(0,1,length = n.pr^(1/5)),
                      x3 = seq(0,1,length = n.pr^(1/5)), x4 = seq(0,1,length = n.pr^(1/5)),
                      x5 = seq(0,1,length = n.pr^(1/5)))
fIn.pr <- list(f1 = matrix(runif(n.pr*10), ncol = 10), matrix(runif(n.pr*22), ncol = 22))

# making predictions based on the three best models (compact code with all 3 calls)
preds <- do.call(cbind, Map(function(model, args) {
  active <- get_active_in(sIn.pr, fIn.pr, args)
  predict(model, sIn.pr = active$sIn.on, fIn.pr = active$fIn.on)$mean
}, modStack, argStack))


# <<<<<<< PART 4: plotting predictions from the three best models found by fgpm_factory >>>>>>>
#         -----------------------------------------------------------------------------
# plotting predictions made by the three models
plot(1, xlim = c(1,nrow(preds)), ylim = range(preds), xaxt = "n",
     xlab = "Prediction point index", ylab = "Output",
     main = "Predictions with best 3 structural configurations")
axis(1, 1:nrow(preds))
for (i in seq_len(n.pr)) {lines(rep(i,2), range(preds[i,1:3]), col = "grey35", lty = 3)}
points(preds[,1], pch = 21, bg = "black")
points(preds[,2], pch = 23, bg = "red")
points(preds[,3], pch = 24, bg = "green")
legend("bottomleft", legend = c("Model 1", "Model 2", "Model 3"),
       pch = c(21, 23, 24), pt.bg = c("black", "red", "green"), inset = c(.02,.08))

## End(Not run)


[Package funGp version 1.0.0 Index]