| coordProf_UQ {profExtrema} | R Documentation | 
Coordinate profiles UQ from a kriging model
Description
The function coordProf_UQ computes the profile extrema functions for posterior realizations of a Gaussian process and its confidence bounds
Usage
coordProf_UQ(object, threshold, allResMean = NULL,
  quantiles_uq = c(0.05, 0.95), options_approx = NULL,
  options_full_sims = NULL, options_sims = NULL,
  options_bound = NULL, plot_level = 0, plot_options = NULL,
  return_level = 1)
Arguments
| object | either a km model or a list containing partial results. If  | 
| threshold | the threshold of interest | 
| allResMean | a list resulting from  | 
| quantiles_uq | a vector containing the quantiles to be computed | 
| options_approx | an optional list of options for approxMaxMin, see approxMaxMin for details. | 
| options_full_sims | an optional list of options for getAllMaxMin, see getAllMaxMin for details. If NULL the full computations are not excuted. NOTE: this computations might be very expensive! | 
| options_sims | an optional list of options for the posterior simulations. 
 | 
| options_bound | an optional list containing  | 
| plot_level | an integer to select the plots to return (0=no plots, 1=basic plots, 2= all plots) | 
| plot_options | an optional list of parameters for plots. See setPlotOptions for currently available options. | 
| return_level | an integer to select the amount of details returned | 
Value
If return_level=1 a list containing
- profSups:an array- dxfullDesignSizexnsimscontaining the profile sup for each coordinate for each realization.
- profInfs:an array- dxfullDesignSizexnsimscontaining the profile inf for each coordinate for each realization.
- prof_quantiles_approx:a list containing the quantiles (levels set by- quantiles_uq) of the profile extrema functions.
 if return_level=2 the same list as above but also including more: a list containing 
- times:a list containing- tSpts:computational time for selecting pilot points.
- tApprox1ord:vector containing the computational time required for profile extrema computation for each realization
 
- simuls:a matrix containing the value of the field simulated at the pilot points
- sPts:the pilot points
Author(s)
Dario Azzimonti
Examples
if (!requireNamespace("DiceKriging", quietly = TRUE)) {
stop("DiceKriging needed for this example to work. Please install it.",
     call. = FALSE)
}
# Compute a kriging model from 50 evaluations of the Branin function
# Define the function
g<-function(x){
  return(-branin(x))
}
gp_des<-lhs::maximinLHS(20,2)
reals<-apply(gp_des,1,g)
kmModel<-km(design = gp_des,response = reals,covtype = "matern3_2")
threshold=-10
d<-2
# Compute coordinate profiles UQ starting from GP model
# define simulation options
options_sims<-list(algorithm="B", lower=rep(0,d), upper=rep(1,d),
                   batchsize=80, optimcontrol = list(method="genoud",pop.size=100,print.level=0),
                   integcontrol = list(distrib="sobol",n.points=1000), nsim=150)
# define 1 order approximation options
init_des<-lhs::maximinLHS(15,d)
options_approx<- list(multistart=4,heavyReturn=TRUE,
                      initDesign=init_des,fullDesignSize=100,
                      smoother="1order")
# define plot options
options_plots<-list(save=FALSE, titleProf = "Coordinate profiles",
                    title2d = "Posterior mean",qq_fill=TRUE)
## Not run: 
# profile UQ on approximate coordinate profiles
cProfiles_UQ<-coordProf_UQ(object = kmModel,threshold = threshold,allResMean = NULL,
                            quantiles_uq = c(0.05,0.95),options_approx = options_approx,
                            options_full_sims = NULL,options_sims = options_sims,
                            options_bound = NULL,plot_level = 3,
                            plot_options = options_plots,return_level = 3)
# profile UQ on full optim coordinate profiles
options_full_sims<-list(multistart=4,heavyReturn=TRUE)
cProfiles_UQ_full<-coordProf_UQ(object = cProfiles_UQ,threshold = threshold,allResMean = NULL,
                            quantiles_uq = c(0.05,0.95),options_approx = options_approx,
                            options_full_sims = options_full_sims,options_sims = options_sims,
                            options_bound = NULL,plot_level = 3,
                            plot_options = options_plots,return_level = 3)
# profile UQ on full optim coordinate profiles with bound
cProfiles_UQ_full_bound<-coordProf_UQ(object = cProfiles_UQ_full,threshold = threshold,
                                      allResMean = NULL, quantiles_uq = c(0.05,0.95),
                                      options_approx = options_approx,
                                      options_full_sims = options_full_sims,
                                      options_sims = options_sims,
                                      options_bound = list(beta=0.024,alpha=0.05),
                                      plot_level = 3, plot_options = options_plots,
                                      return_level = 3)
## End(Not run)