spConjNNGP {spNNGP} | R Documentation |
Function for Fitting Univariate Bayesian Conjugate Spatial Regression Models
Description
The function spConjNNGP
fits Gaussian univariate Bayesian conjugate spatial
regression models using Nearest Neighbor Gaussian Processes (NNGP).
Usage
spConjNNGP(formula, data = parent.frame(), coords, knots, n.neighbors = 15,
theta.alpha, sigma.sq.IG, cov.model = "exponential",
k.fold = 5, score.rule = "crps",
X.0, coords.0, n.omp.threads = 1, search.type = "cb",
ord, return.neighbor.info = TRUE,
neighbor.info, fit.rep = FALSE, n.samples, verbose = TRUE, ...)
Arguments
formula |
a symbolic description of the regression model to be fit. See example below. |
data |
an optional data frame containing the variables in the
model. If not found in data, the variables are taken from
|
coords |
an |
knots |
an |
n.neighbors |
number of neighbors used in the NNGP. |
theta.alpha |
a vector or matrix of parameter values for
|
sigma.sq.IG |
a vector of length two that holds the
hyperparameters, shape and scale respectively, for the
inverse-Gamma prior on |
cov.model |
a quoted keyword that specifies the covariance
function used to model the spatial dependence structure among the
observations. Supported covariance model key words are:
|
k.fold |
specifies the number of
k folds for cross-validation. If |
score.rule |
a quoted keyword |
X.0 |
the design matrix for prediction locations. An
intercept should be provided in the first column if one is specified
in |
coords.0 |
the spatial coordinates corresponding to
|
n.omp.threads |
a positive integer indicating
the number of threads to use for SMP parallel processing. The package must
be compiled for OpenMP support. For most Intel-based machines, we recommend setting
|
search.type |
a quoted keyword that specifies type of nearest
neighbor search algorithm. Supported method key words are: |
ord |
an index vector of length |
return.neighbor.info |
if |
neighbor.info |
see the |
fit.rep |
if |
n.samples |
gives the number of posterior samples
returned. Note, point and associated variance estimates for model
parameters are not based on posterior samples. Only specify
|
verbose |
if |
... |
currently no additional arguments. |
Value
An object of class NNGP
and conjugate
, and, if knots
is provided, SLGP
. Among other elements, the object contains:
theta.alpha |
the input |
beta.hat |
a matrix of regression coefficient estimates
corresponding to the returned |
beta.var |
|
sigma.sq.hat |
estimate of |
sigma.sq.var |
|
k.fold.scores |
results from the k-fold cross-validation if
|
y.0.hat |
prediction if |
y.0.var.hat |
|
n.neighbors |
number of neighbors used in the NNGP. |
neighbor.info |
returned if |
run.time |
execution time for parameter estimation reported using
|
Author(s)
Andrew O. Finley finleya@msu.edu,
Abhirup Datta abhidatta@jhu.edu,
Sudipto Banerjee sudipto@ucla.edu
References
Datta, A., S. Banerjee, A.O. Finley, and A.E. Gelfand. (2016) Hierarchical Nearest-Neighbor Gaussian process models for large geostatistical datasets. Journal of the American Statistical Association, doi: 10.1080/01621459.2015.1044091.
Finley, A.O., A. Datta, S. Banerjee (2022) spNNGP R Package for Nearest Neighbor Gaussian Process Models. Journal of Statistical Software, doi: 10.18637/jss.v103.i05.
Finley, A.O., A. Datta, B.D. Cook, D.C. Morton, H.E. Andersen, and S. Banerjee. (2019) Efficient algorithms for Bayesian Nearest Neighbor Gaussian Processes. Journal of Computational and Graphical Statistics, doi: 10.1080/10618600.2018.1537924.
Gneiting, T and A.E. Raftery. (2007) Strictly proper scoring rules, prediction, and estimation. Journal of the American Statistical Association, doi: 10.1198/016214506000001437.
Shirota, S., A.O. Finley, B.D. Cook, and S. Banerjee (2019) Conjugate Nearest Neighbor Gaussian Process models for efficient statistical interpolation of large spatial data. https://arxiv.org/abs/1907.10109.
Examples
rmvn <- function(n, mu=0, V = matrix(1)){
p <- length(mu)
if(any(is.na(match(dim(V),p))))
stop("Dimension problem!")
D <- chol(V)
t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p)))
}
##Make some data
set.seed(1)
n <- 2000
coords <- cbind(runif(n,0,1), runif(n,0,1))
x <- cbind(1, rnorm(n))
B <- as.matrix(c(1,5))
sigma.sq <- 5
tau.sq <- 1
phi <- 3/0.5
D <- as.matrix(dist(coords))
R <- exp(-phi*D)
w <- rmvn(1, rep(0,n), sigma.sq*R)
y <- rnorm(n, x%*%B + w, sqrt(tau.sq))
ho <- sample(1:n, 1000)
y.ho <- y[ho]
x.ho <- x[ho,,drop=FALSE]
w.ho <- w[ho]
coords.ho <- coords[ho,]
y <- y[-ho]
x <- x[-ho,,drop=FALSE]
w <- w[-ho,,drop=FALSE]
coords <- coords[-ho,]
##Fit a Conjugate NNGP model and predict for the holdout
sigma.sq.IG <- c(2, sigma.sq)
cov.model <- "exponential"
g <- 10
theta.alpha <- cbind(seq(phi,30,length.out=g), seq(tau.sq/sigma.sq,5,length.out=g))
colnames(theta.alpha) <- c("phi", "alpha")
m.c <- spConjNNGP(y~x-1, coords=coords, n.neighbors = 10,
X.0 = x.ho, coords.0 = coords.ho,
k.fold = 5, score.rule = "crps",
n.omp.threads = 1,
theta.alpha = theta.alpha, sigma.sq.IG = sigma.sq.IG, cov.model = cov.model)
m.c