bhomGP {bhetGP} | R Documentation |
MCMC sampling for Homoskedastic GP
Description
Conducts MCMC sampling of hyperparameters for a homGP.
Separate length scale parameters theta_y
govern the correlation strength of the response. g
governs
noise for the noise. tau2_y
control the amplitude of the mean process.
In Matern covariance, v
governs smoothness.
Usage
bhomGP(
x = NULL,
y = NULL,
reps_list = NULL,
nmcmc = 1000,
sep = TRUE,
inits = NULL,
priors = NULL,
cov = c("exp2", "matern", "ARD matern"),
v = 2.5,
stratergy = c("default", "flat"),
vecchia = FALSE,
m = min(25, length(y) - 1),
ordering = NULL,
reps = TRUE,
verb = TRUE,
cores = 4
)
Arguments
x |
vector or matrix of input locations |
y |
vector of response values |
reps_list |
list object from hetGP::find_reps |
nmcmc |
number of MCMC iterations |
sep |
logical indicating whether to fit isotropic GP ( |
inits |
set initial values for hyparameters:
|
priors |
hyperparameters for priors and proposals (see details) |
cov |
covariance kernel, either Matern, ARD Matern
or squared exponential ( |
v |
Matern smoothness parameter (only used if |
stratergy |
choose initialization stratergy; "default" uses hetGP for
|
vecchia |
logical indicating whether to use Vecchia approximation |
m |
size of Vecchia conditioning sets (only used if
|
ordering |
optional ordering for Vecchia approximation, must correspond
to rows of |
reps |
logical; if |
verb |
logical indicating whether to print progress |
cores |
if |
Details
Maps inputs x
to mean response y
. Utilizes Metropolis Hastings
sampling of the length scale and nugget parameters with proposals and priors
controlled by priors
. g
is estimated by default but may be specified and fixed.
When vecchia = TRUE
, all calculations leverage the Vecchia approximation with
specified conditioning set size m
. tau2_y
is inferred by default but
may be pre-specified and fixed.
NOTE on OpenMP: The Vecchia implementation relies on OpenMP parallelization for efficient computation. This function will produce a warning message if the package was installed without OpenMP (this is the default for CRAN packages installed on Apple machines). To set up OpenMP parallelization, download the package source code and install using the gcc/g++ compiler.
Proposals for g
and theta
follow a uniform sliding window
scheme, e.g.
theta_star <- runif(1, l * theta_t / u, u * theta_t / l)
,
with defaults l = 1
and u = 2
provided in priors
.
To adjust these, set priors = list(l = new_l, u = new_u)
.
Priors on g
, theta_y
follow Gamma
distributions with shape parameters (alpha
) and rate parameters
(beta
) controlled within the priors
list object.
Defaults are
-
priors$alpha$theta_lam <- 1.5
-
priors$beta$theta_lam <- 4
-
priors$alpha$theta_y <- 1.5
-
priors$beta$theta_y <- 4
-
priors$alpha$g <- 1.5
-
priors$beta$g <- 4
tau2_y
is not sampled; rather directly inferred
under conjugate Inverse Gamma prior with shape (alpha
) and scale
parameters (beta
) within the priors
list object
-
priors$a$tau2_y <- 10
-
priors$a$tau2_y <- 4
These priors are designed for x
scaled to
[0, 1] and y
having mean mean_y
. These may be
adjusted using the priors
input.
Initial values for theta_y
, and g
may be
specified by the user. If no initial values are specified, stratergy
will determine the initialization method. stratergy = "default"
leverages mleHomGP for initial values of hyper-parameters if
vecchia = FALSE
and scaled vecchia approach if vecchia = TRUE
.
For scaled Vecchia code from https://github.com/katzfuss-group/scaledVecchia/blob/master/vecchia_scaled.R is used to fit a vecchia approximated GP to (x, y). A script is leveraged internally within this package that fits this method.
Optionally, choose stratergy = "flat" which will start at uninformative initial values or specify initial values.
The output object of class bhomgp
or bhomgp_vec
is designed for
use with trim
, predict
, and plot
.
Value
a list of the S3 class bhomgp
or bhomgp_vec
with elements:
-
x
: copy of input matrix -
y
: copy of response vector -
Ylist
: list of all responses observed per location (x) -
A
: number of replicates at each location -
nmcmc
: number of MCMC iterations -
priors
: copy of proposal/prior settings -
settings
: setting for predictions usingbhetgp
orbhetgp_vec
object -
g_y
: vector of MCMC samples forg_y
-
theta_y
: vector of MCMC samples fortheta_y
(length scale of mean process) -
tau2
: vector of MAP estimates fortau2
(scale parameter of outer layer) -
llik_y
: vector of MVN log likelihood of Y for reach Gibbs iteration -
time
: computation time in seconds -
x_approx
: conditioning set, NN and ordering forvecchia = TRUE
-
m
: copy of size of conditioning set forvecchia = TRUE
References
Binois, Mickael, Robert B. Gramacy, and Mike Ludkovski. "Practical heteroscedastic Gaussian process modeling for large simulation experiments." Journal of Computational and Graphical Statistics 27.4 (2018): 808-821.
Katzfuss, Matthias, Joseph Guinness, and Earl Lawrence. "Scaled Vecchia approximation for fast computer-model emulation." SIAM/ASA Journal on Uncertainty Quantification 10.2 (2022): 537-554.
Sauer, Annie Elizabeth. "Deep Gaussian process surrogates for computer experiments." (2023).
Examples
fx <- function(x){
result <- (6 * x - 2)^2* sin(12 * x - 4)
}
# Training data
r <- 10
xn <- seq(0, 1, length = 25)
x <- rep(xn, r)
f <- fx(x)
y <- f + rnorm(length(x))
# Testing data
xx <- seq(0, 1, length = 100)
yy <- fx(xx)
# Example 1: Full model, no Vecchia
fit <- bhomGP(x, y, nmcmc = 100)
fit <- trim(fit, 50, 10)
fit <- predict(fit, xx)
plot(fit) # can run with trace = TRUE to view trace plots
# Example 2: Vecchia approximated model
fit <- bhomGP(x, y, nmcmc = 100, vecchia = TRUE, m = 5)
fit <- trim(fit, 50, 10)
fit <- predict(fit, xx, vecchia = TRUE)
plot(fit)