svcTMsPGOcc | R Documentation |
The function svcTMsPGOcc
fits multi-species multi-season spatially-varying coefficient 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. Models are implemented using a Nearest Neighbor Gaussian Process.
svcTMsPGOcc(occ.formula, det.formula, data, inits, priors, tuning,
svc.cols = 1, cov.model = 'exponential', NNGP = TRUE,
n.neighbors = 15, search.type = 'cb', std.by.sp = FALSE,
n.factors, svc.by.sp, n.batch, batch.length,
accept.rate = 0.43, n.omp.threads = 1,
verbose = TRUE, ar1 = FALSE, n.report = 100,
n.burn = round(.10 * n.batch * batch.length), n.thin = 1,
n.chains = 1, ...)
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 |
svc.cols |
a vector indicating the variables whose effects will be
estimated as spatially-varying coefficients. |
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: |
std.by.sp |
a logical value indicating whether the covariates are standardized
separately for each species within the corresponding range for each species ( |
n.factors |
the number of factors to use in the spatial factor model approach. Note this corresponds to the number of factors used for each spatially-varying coefficient that is estimated in the model. 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). |
svc.by.sp |
an optional list with length equal to |
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 within chains. The package must
be compiled for OpenMP support. For most Intel-based machines, we
recommend setting |
verbose |
if |
ar1 |
logical value indicating whether to include an AR(1) zero-mean
temporal random effect in the model. 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. |
... |
currently no additional arguments |
An object of class svcTMsPGOcc
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 four-dimensional array of posterior samples for the latent occurrence values for each species. Dimensions corresopnd to MCMC sample, species, site, and primary time period. |
psi.samples |
a four-dimensional array of posterior samples for the latent occupancy probability values for each species. Dimensions correspond to MCMC sample, species, site, and primary time period. |
w.samples |
a four-dimensional array of posterior samples for the latent spatial random effects for each spatial factor within each spatially-varying coefficient. Dimensions correspond to MCMC sample, factor, site, and spatially-varying coefficient. |
sigma.sq.psi.samples |
a |
sigma.sq.p.samples |
a |
beta.star.samples |
a |
alpha.star.samples |
a |
like.samples |
a four-dimensional array of posterior samples for the likelihood value used for calculating WAIC. Dimensions correspond to MCMC sample, species, site, and time period. |
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 |
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()
.
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.
Jeffrey W. Doser doserjef@msu.edu,
Andrew O. Finley finleya@msu.edu
Doser, J. W., Finley, A. O., Saunders, S. P., Kery, M., Weed, A. S., & Zipkin, E. F. (2024A). Modeling complex species-environment relationships through spatially-varying coefficient occupancy models. Journal of Agricultural, Biological and Environmental Statistics. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1007/s13253-023-00595-6")}.
Doser, J. W., Kery, M., Saunders, S. P., Finley, A. O., Bateman, B. L., Grand, J., Reault, S., Weed, A. S., & Zipkin, E. F. (2024B). Guidelines for the use of spatially varying coefficients in species distribution models. Global Ecology and Biogeography, 33(4), e13814. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1111/geb.13814")}.
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, \Sexpr[results=rd]{tools:::Rd_expr_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, \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1080/10618600.2018.1537924")}.
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. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.18637/jss.v067.i01")}.
# Simulate Data -----------------------------------------------------------
set.seed(500)
J.x <- 8
J.y <- 8
J <- J.x * J.y
# Years sampled
n.time <- sample(3:10, J, replace = TRUE)
# n.time <- rep(10, J)
n.time.max <- max(n.time)
# Replicates
n.rep <- matrix(NA, J, max(n.time))
for (j in 1:J) {
n.rep[j, 1:n.time[j]] <- sample(2:4, n.time[j], replace = TRUE)
}
N <- 7
# Community-level covariate effects
# Occurrence
beta.mean <- c(-3, -0.2, 0.5)
trend <- FALSE
sp.only <- 0
p.occ <- length(beta.mean)
tau.sq.beta <- c(0.6, 1.5, 1.4)
# Detection
alpha.mean <- c(0, 1.2, -1.5)
tau.sq.alpha <- c(1, 0.5, 2.3)
p.det <- length(alpha.mean)
# Random effects
psi.RE <- list()
p.RE <- list()
# 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]))
}
sp <- TRUE
svc.cols <- c(1, 2)
p.svc <- length(svc.cols)
n.factors <- 3
phi <- runif(p.svc * n.factors, 3 / .9, 3 / .3)
factor.model <- TRUE
cov.model <- 'exponential'
dat <- simTMsOcc(J.x = J.x, J.y = J.y, n.time = n.time, n.rep = n.rep, N = N,
beta = beta, alpha = alpha, sp.only = sp.only, trend = trend,
psi.RE = psi.RE, p.RE = p.RE, factor.model = factor.model,
svc.cols = svc.cols, n.factors = n.factors, phi = phi, sp = sp,
cov.model = cov.model)
y <- dat$y
X <- dat$X
X.p <- dat$X.p
coords <- dat$coords
X.re <- dat$X.re
X.p.re <- dat$X.p.re
occ.covs <- list(occ.cov.1 = X[, , 2],
occ.cov.2 = X[, , 3])
det.covs <- list(det.cov.1 = X.p[, , , 2],
det.cov.2 = X.p[, , , 3])
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 / .9, b = 3 / .1))
z.init <- apply(y, c(1, 2, 3), function(a) as.numeric(sum(a, na.rm = TRUE) > 0))
inits.list <- list(alpha.comm = 0, beta.comm = 0, beta = 0,
alpha = 0, tau.sq.beta = 1, tau.sq.alpha = 1,
phi = 3 / .5, z = z.init)
# Tuning
tuning.list <- list(phi = 1)
# Number of batches
n.batch <- 5
# Batch length
batch.length <- 25
n.burn <- 25
n.thin <- 1
n.samples <- n.batch * batch.length
# Note that this is just a test case and more iterations/chains may need to
# be run to ensure convergence.
out <- svcTMsPGOcc(occ.formula = ~ occ.cov.1 + occ.cov.2,
det.formula = ~ det.cov.1 + det.cov.2,
data = data.list,
inits = inits.list,
n.batch = n.batch,
batch.length = batch.length,
accept.rate = 0.43,
NNGP = TRUE,
n.neighbors = 5,
n.factors = n.factors,
svc.cols = svc.cols,
cov.model = 'exponential',
priors = prior.list,
tuning = tuning.list,
n.omp.threads = 1,
verbose = TRUE,
n.report = 1,
n.burn = n.burn,
n.thin = n.thin,
n.chains = 1)
summary(out)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.