spNNGP {spNNGP} | R Documentation |
Function for Fitting Univariate Bayesian Spatial Regression Models
Description
The function spNNGP
fits Gaussian and non-Gaussian univariate Bayesian spatial
regression models using Nearest Neighbor Gaussian Processes (NNGP).
Usage
spNNGP(formula, data = parent.frame(), coords, method = "response",
family="gaussian", weights, n.neighbors = 15,
starting, tuning, priors, cov.model = "exponential",
n.samples, n.omp.threads = 1, search.type = "cb", ord,
return.neighbor.info = FALSE, neighbor.info,
fit.rep = FALSE, sub.sample, verbose = TRUE, n.report = 100, ...)
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 |
method |
a quoted keyword that specifies the NNGP sampling
algorithm. Supported method keywords are: |
family |
a quoted keyword that specifies the data likelihood. Choices are "gaussian" for continuous outcome and "binomial" for discrete outcome which assumes a logistic link modeled using Polya-Gamma latent variables. |
weights |
specifies the number of trials for each observation when
|
n.neighbors |
number of neighbors used in the NNGP. |
starting |
a list with each tag corresponding to a parameter name. Valid tags are |
tuning |
a list with each tag corresponding to a parameter
name. Valid tags are |
priors |
a list with each tag corresponding to a
parameter name. Valid tags are |
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:
|
n.samples |
the number of posterior samples to collect. |
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 |
fit.rep |
if |
sub.sample |
an optional list that specifies the samples used
for |
verbose |
if |
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 |
n.report |
the interval to report Metropolis sampler acceptance and MCMC progress. |
... |
currently no additional arguments. |
Details
Model parameters can be fixed at their starting
values by setting their
tuning
values to zero.
The no nugget model is specified by setting tau.sq
to zero
in the starting
and tuning
lists.
Value
An object of class NNGP
with additional class designations for
method
and family
. The return object is a list comprising:
p.beta.samples |
a |
p.theta.samples |
a |
p.w.samples |
is a matrix of posterior samples for the spatial
random effects, where rows correspond to locations in |
y.hat.samples |
if |
y.hat.quants |
if |
y.rep.samples |
if |
y.rep.quants |
if |
s.indx |
if |
neighbor.info |
returned if |
run.time |
execution time for parameter estimation reported using
|
The return object will include additional objects used for subsequent prediction and/or model fit evaluation.
Author(s)
Andrew O. Finley finleya@msu.edu,
Abhirup Datta abhidatta@jhu.edu,
Sudipto Banerjee sudipto@ucla.edu
References
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.
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, 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.
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 <- 100
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))
##Fit a Response and Latent NNGP model
n.samples <- 500
starting <- list("phi"=phi, "sigma.sq"=5, "tau.sq"=1)
tuning <- list("phi"=0.5, "sigma.sq"=0.5, "tau.sq"=0.5)
priors <- list("phi.Unif"=c(3/1, 3/0.01), "sigma.sq.IG"=c(2, 5), "tau.sq.IG"=c(2, 1))
cov.model <- "exponential"
m.s <- spNNGP(y~x-1, coords=coords, starting=starting, method="latent", n.neighbors=10,
tuning=tuning, priors=priors, cov.model=cov.model,
n.samples=n.samples, n.omp.threads=1)
summary(m.s)
plot(apply(m.s$p.w.samples, 1, median), w)
m.r <- spNNGP(y~x-1, coords=coords, starting=starting, method="response", n.neighbors=10,
tuning=tuning, priors=priors, cov.model=cov.model,
n.samples=n.samples, n.omp.threads=1)
summary(m.r)
##Fit with some user defined neighbor ordering
##ord <- order(coords[,2]) ##y-axis
ord <- order(coords[,1]+coords[,2]) ##x+y-axis
##ord <- sample(1:n, n) ##random
m.r.xy <- spNNGP(y~x-1, coords=coords, starting=starting, method="response", n.neighbors=10,
tuning=tuning, priors=priors, cov.model=cov.model,
ord=ord, return.neighbor.info=TRUE,
n.samples=n.samples, n.omp.threads=1)
summary(m.r.xy)
## Not run:
##Visualize the neighbor sets and ordering constraint
n.indx <- m.r.xy$neighbor.info$n.indx
ord <- m.r.xy$neighbor.info$ord
##This is how the data are ordered internally for model fitting
coords.ord <- coords[ord,]
for(i in 1:n){
plot(coords.ord, cex=1, xlab="Easting", ylab="Northing")
points(coords.ord[i,,drop=FALSE], col="blue", pch=19, cex=1)
points(coords.ord[n.indx[[i]],,drop=FALSE], col="red", pch=19, cex=1)
readline(prompt = "Pause. Press <Enter> to continue...")
}
## End(Not run)