sfMsPGOcc {spOccupancy} | R Documentation |
Function for Fitting Spatial Factor Multi-Species Occupancy Models
Description
The function sfMsPGOcc
fits multi-species spatial occupancy models with species correlations (i.e., a spatially-explicit joint species distribution model with imperfect detection). We use Polya-Gamma latent variables and a spatial factor modeling approach. Currently, models are implemented using a Nearest Neighbor Gaussian Process.
Usage
sfMsPGOcc(occ.formula, det.formula, data, inits, priors, tuning,
cov.model = 'exponential', NNGP = TRUE,
n.neighbors = 15, search.type = 'cb', n.factors, n.batch,
batch.length, accept.rate = 0.43, n.omp.threads = 1,
verbose = TRUE, n.report = 100,
n.burn = round(.10 * n.batch * batch.length), n.thin = 1,
n.chains = 1, k.fold, k.fold.threads = 1, k.fold.seed,
k.fold.only = FALSE, ...)
Arguments
occ.formula |
a symbolic description of the model to be fit for the occurrence portion of the model using R's model syntax. Random intercepts are allowed using lme4 syntax (Bates et al. 2015). Only right-hand side of formula is specified. See example below. |
det.formula |
a symbolic description of the model to be fit for the detection portion of the model using R's model syntax. Only right-hand side of formula is specified. See example below. Random intercepts are allowed using lme4 syntax (Bates et al. 2015). |
data |
a list containing data necessary for model fitting.
Valid tags are |
inits |
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 |
tuning |
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:
|
NNGP |
if |
n.neighbors |
number of neighbors used in the NNGP. Only used if
|
search.type |
a quoted keyword that specifies the type of nearest
neighbor search algorithm. Supported method key words are: |
n.factors |
the number of factors to use in the spatial factor model approach. Typically, the number of factors is set to be small (e.g., 4-5) relative to the total number of species in the community, which will lead to substantial decreases in computation time. However, the value can be anywhere between 1 and N (the number of species in the community). |
n.batch |
the number of MCMC batches in each chain to run for the Adaptive MCMC sampler. See Roberts and Rosenthal (2009) for details. |
batch.length |
the length of each MCMC batch to run for the Adaptive MCMC sampler. See Roberts and Rosenthal (2009) for details. |
accept.rate |
target acceptance rate for Adaptive MCMC. Defaul is 0.43. See Roberts and Rosenthal (2009) for details. |
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 |
verbose |
if |
n.report |
the interval to report Metropolis sampler acceptance and MCMC progress. Note this is specified in terms of batches and not overall samples for spatial models. |
n.burn |
the number of samples out of the total |
n.thin |
the thinning interval for collection of MCMC samples. The
thinning occurs after the |
n.chains |
the number of chains to run in sequence. |
k.fold |
specifies the number of k folds for cross-validation.
If not specified as an argument, then cross-validation is not performed
and |
k.fold.threads |
number of threads to use for cross-validation. If
|
k.fold.seed |
seed used to split data set into |
k.fold.only |
a logical value indicating whether to only perform
cross-validation ( |
... |
currently no additional arguments |
Value
An object of class sfMsPGOcc
that is a list comprised of:
beta.comm.samples |
a |
alpha.comm.samples |
a |
tau.sq.beta.samples |
a |
tau.sq.alpha.samples |
a |
beta.samples |
a |
alpha.samples |
a |
theta.samples |
a |
lambda.samples |
a |
z.samples |
a three-dimensional array of posterior samples for the latent occurrence values for each species. |
psi.samples |
a three-dimensional array of posterior samples for the latent occupancy probability values for each species. |
w.samples |
a three-dimensional array of posterior samples for the latent spatial random effects for each latent factor. |
sigma.sq.psi.samples |
a |
sigma.sq.p.samples |
a |
beta.star.samples |
a |
alpha.star.samples |
a |
like.samples |
a three-dimensional array of posterior samples for the likelihood value associated with each site and species. Used for calculating WAIC. |
rhat |
a list of Gelman-Rubin diagnostic values for some of the model parameters. |
ESS |
a list of effective sample sizes for some of the model parameters. |
run.time |
MCMC sampler execution time reported using |
k.fold.deviance |
vector of scoring rules (deviance) from k-fold cross-validation.
A separate value is reported for each species.
Only included if |
The return object will include additional objects used for
subsequent prediction and/or model fit evaluation. Note that detection
probability estimated values are not included in the model object, but can
be extracted using fitted()
.
Note
Some of the underlying code used for generating random numbers from the Polya-Gamma distribution is taken from the pgdraw package written by Daniel F. Schmidt and Enes Makalic. Their code implements Algorithm 6 in PhD thesis of Jesse Bennett Windle (2013) https://repositories.lib.utexas.edu/handle/2152/21842.
Author(s)
Jeffrey W. Doser doserjef@msu.edu,
Andrew O. Finley finleya@msu.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, 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.
Finley, A. O., Datta, A., and Banerjee, S. (2020). spNNGP R package for nearest neighbor Gaussian process models. arXiv preprint arXiv:2001.09111.
Polson, N.G., J.G. Scott, and J. Windle. (2013) Bayesian Inference for Logistic Models Using Polya-Gamma Latent Variables. Journal of the American Statistical Association, 108:1339-1349.
Roberts, G.O. and Rosenthal J.S. (2009) Examples of adaptive MCMC. Journal of Computational and Graphical Statistics, 18(2):349-367.
Bates, Douglas, Martin Maechler, Ben Bolker, Steve Walker (2015). Fitting Linear Mixed-Effects Models Using lme4. Journal of Statistical Software, 67(1), 1-48. doi:10.18637/jss.v067.i01.
Hooten, M. B., and Hobbs, N. T. (2015). A guide to Bayesian model selection for ecologists. Ecological Monographs, 85(1), 3-28.
Christensen, W. F., and Amemiya, Y. (2002). Latent variable analysis of multivariate spatial data. Journal of the American Statistical Association, 97(457), 302-317.
Examples
set.seed(400)
# Simulate Data -----------------------------------------------------------
J.x <- 7
J.y <- 7
J <- J.x * J.y
n.rep <- sample(2:4, size = J, replace = TRUE)
N <- 8
# Community-level covariate effects
# Occurrence
beta.mean <- c(0.2, -0.15)
p.occ <- length(beta.mean)
tau.sq.beta <- c(0.6, 0.3)
# Detection
alpha.mean <- c(0.5, 0.2, -.2)
tau.sq.alpha <- c(0.2, 0.3, 0.8)
p.det <- length(alpha.mean)
# Random effects
psi.RE <- list()
# Include a non-spatial random effect on occurrence
psi.RE <- list(levels = c(20),
sigma.sq.psi = c(0.5))
p.RE <- list()
# Include a random effect on detection
p.RE <- list(levels = c(40),
sigma.sq.p = c(2))
# Draw species-level effects from community means.
beta <- matrix(NA, nrow = N, ncol = p.occ)
alpha <- matrix(NA, nrow = N, ncol = p.det)
for (i in 1:p.occ) {
beta[, i] <- rnorm(N, beta.mean[i], sqrt(tau.sq.beta[i]))
}
for (i in 1:p.det) {
alpha[, i] <- rnorm(N, alpha.mean[i], sqrt(tau.sq.alpha[i]))
}
n.factors <- 4
phi <- runif(n.factors, 3/1, 3/.4)
dat <- simMsOcc(J.x = J.x, J.y = J.y, n.rep = n.rep, N = N, beta = beta, alpha = alpha,
phi = phi, sp = TRUE, cov.model = 'exponential',
factor.model = TRUE, n.factors = n.factors, psi.RE = psi.RE,
p.RE = p.RE)
# Number of batches
n.batch <- 10
# Batch length
batch.length <- 25
n.samples <- n.batch * batch.length
y <- dat$y
X <- dat$X
X.p <- dat$X.p
X.p.re <- dat$X.p.re
X.re <- dat$X.re
coords <- as.matrix(dat$coords)
# Package all data into a list
occ.covs <- cbind(X, X.re)
colnames(occ.covs) <- c('int', 'occ.cov', 'occ.re')
det.covs <- list(det.cov.1 = X.p[, , 2],
det.cov.2 = X.p[, , 3],
det.re = X.p.re[, , 1])
data.list <- list(y = y,
occ.covs = occ.covs,
det.covs = det.covs,
coords = coords)
# Priors
prior.list <- list(beta.comm.normal = list(mean = 0, var = 2.72),
alpha.comm.normal = list(mean = 0, var = 2.72),
tau.sq.beta.ig = list(a = 0.1, b = 0.1),
tau.sq.alpha.ig = list(a = 0.1, b = 0.1),
phi.unif = list(a = 3/1, b = 3/.1))
# Initial values
lambda.inits <- matrix(0, N, n.factors)
diag(lambda.inits) <- 1
lambda.inits[lower.tri(lambda.inits)] <- rnorm(sum(lower.tri(lambda.inits)))
inits.list <- list(alpha.comm = 0,
beta.comm = 0,
beta = 0,
alpha = 0,
tau.sq.beta = 1,
tau.sq.alpha = 1,
phi = 3 / .5,
lambda = lambda.inits,
z = apply(y, c(1, 2), max, na.rm = TRUE))
# Tuning
tuning.list <- list(phi = 1)
out <- sfMsPGOcc(occ.formula = ~ occ.cov + (1 | occ.re),
det.formula = ~ det.cov.1 + det.cov.2 + (1 | det.re),
data = data.list,
inits = inits.list,
n.batch = n.batch,
batch.length = batch.length,
accept.rate = 0.43,
priors = prior.list,
cov.model = "exponential",
tuning = tuning.list,
n.omp.threads = 1,
verbose = TRUE,
NNGP = TRUE,
n.neighbors = 5,
n.factors = n.factors,
search.type = 'cb',
n.report = 10,
n.burn = 50,
n.thin = 1,
n.chains = 1)
summary(out)