svcTPGOcc {spOccupancy} | R Documentation |
Function for Fitting Multi-Season Single-Species Spatially-Varying Coefficient Occupancy Models Using Polya-Gamma Latent Variables
Description
Function for fitting multi-season single-species spatially-varying coefficient occupancy models using Polya-Gamma latent variables. Models are fit using Nearest Neighbor Gaussian Processes.
Usage
svcTPGOcc(occ.formula, det.formula, data, inits, priors,
tuning, svc.cols = 1, cov.model = 'exponential', NNGP = TRUE,
n.neighbors = 15, search.type = 'cb', 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, k.fold, k.fold.threads = 1,
k.fold.seed = 100, 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. Only right-hand side of formula is specified. See example below. Random intercepts are allowed using lme4 syntax (Bates et al. 2015). |
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: |
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 in each chain to run for the Adaptive MCMC sampler. See Roberts and Rosenthal (2009) for details. |
accept.rate |
target acceptance rate for Adaptive MCMC. Default 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 |
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 MCMC progress. |
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 svcTPGOcc
that is a list comprised of:
beta.samples |
a |
alpha.samples |
a |
z.samples |
a three-dimensional array of posterior samples for the latent occupancy values, with dimensions corresponding to posterior sample, site, and primary time period. |
psi.samples |
a three-dimensional array of posterior samples for the latent occupancy probability values, with dimensions corresponding to posterior sample, site, and primary time period. |
theta.samples |
a |
w.samples |
a three-dimensional array of posterior samples for the latent spatial random effects for all spatially-varying coefficients. Dimensions correspond to MCMC sample, coefficient, and sites. |
sigma.sq.psi.samples |
a |
sigma.sq.p.samples |
a |
beta.star.samples |
a |
alpha.star.samples |
a |
eta.samples |
a |
like.samples |
a three-dimensional array of posterior samples for the likelihood values associated with each site and primary time period. 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 |
execution time reported using |
k.fold.deviance |
scoring rule (deviance) from k-fold cross-validation.
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 that if k.fold.only = TRUE
, the
return list object will only contain run.time
and k.fold.deviance
.
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
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.
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.
Kery, M., & Royle, J. A. (2021). Applied hierarchical modeling in ecology: Analysis of distribution, abundance and species richness in R and BUGS: Volume 2: Dynamic and advanced models. Academic Press. Section 4.6.
Finley, A. O., and Banerjee, S. (2020). Bayesian spatially varying coefficient models in the spBayes R package. Environmental Modelling and Software, 125, 104608.
Hooten, M. B., and Hobbs, N. T. (2015). A guide to Bayesian model selection for ecologists. Ecological monographs, 85(1), 3-28.
MacKenzie, D. I., J. D. Nichols, G. B. Lachman, S. Droege, J. Andrew Royle, and C. A. Langtimm. 2002. Estimating Site Occupancy Rates When Detection Probabilities Are Less Than One. Ecology 83: 2248-2255.
Examples
set.seed(1000)
# Sites
J.x <- 15
J.y <- 15
J <- J.x * J.y
# Years sampled
n.time <- sample(10, J, replace = TRUE)
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(4, n.time[j], replace = TRUE)
}
# Occurrence --------------------------
beta <- c(-2, -0.5, -0.2, 0.75)
trend <- TRUE
sp.only <- 0
psi.RE <- list()
# Detection ---------------------------
alpha <- c(1, 0.7, -0.5)
p.RE <- list()
# Spatial parameters ------------------
sp <- TRUE
svc.cols <- c(1, 2, 3)
p.svc <- length(svc.cols)
cov.model <- "exponential"
sigma.sq <- runif(p.svc, 0.1, 1)
phi <- runif(p.svc, 3 / 1, 3 / 0.2)
rho <- 0.8
sigma.sq.t <- 1
ar1 <- TRUE
x.positive <- FALSE
# Get all the data
dat <- simTOcc(J.x = J.x, J.y = J.y, n.time = n.time, n.rep = n.rep,
beta = beta, alpha = alpha, sp.only = sp.only, trend = trend,
psi.RE = psi.RE, p.RE = p.RE,
sp = sp, cov.model = cov.model, sigma.sq = sigma.sq, phi = phi,
svc.cols = svc.cols, ar1 = ar1, rho = rho, sigma.sq.t = sigma.sq.t,
x.positive = x.positive)
# Prep the data for svcTPGOcc ---------------------------------------------
# Full data set
y <- dat$y
X <- dat$X
X.re <- dat$X.re
X.p <- dat$X.p
X.p.re <- dat$X.p.re
coords <- dat$coords
# Package all data into a list
occ.covs <- list(int = X[, , 1],
trend = X[, , 2],
occ.cov.1 = X[, , 3],
occ.cov.2 = X[, , 4])
# Detection
det.covs <- list(det.cov.1 = X.p[, , , 2],
det.cov.2 = X.p[, , , 3])
# Data list bundle
data.list <- list(y = y,
occ.covs = occ.covs,
det.covs = det.covs,
coords = coords)
# Priors
prior.list <- list(beta.normal = list(mean = 0, var = 2.72),
alpha.normal = list(mean = 0, var = 2.72),
phi.unif = list(a = 3/1, b = 3/.1))
# Starting values
z.init <- apply(y, c(1, 2), function(a) as.numeric(sum(a, na.rm = TRUE) > 0))
inits.list <- list(beta = 0, alpha = 0,
sigma.sq = 1, phi = 3 / 0.5,
z = z.init, nu = 1)
# Tuning
tuning.list <- list(phi = 0.4, nu = 0.3, rho = 0.5, sigma.sq = 0.5)
# MCMC settings
n.batch <- 2
n.burn <- 0
n.thin <- 1
# Run the model
out <- svcTPGOcc(occ.formula = ~ trend + occ.cov.1 + occ.cov.2,
det.formula = ~ det.cov.1 + det.cov.2,
data = data.list,
inits = inits.list,
tuning = tuning.list,
priors = prior.list,
cov.model = "exponential",
svc.cols = svc.cols,
NNGP = TRUE,
ar1 = TRUE,
n.neighbors = 5,
n.batch = n.batch,
batch.length = 25,
verbose = TRUE,
n.report = 25,
n.burn = n.burn,
n.thin = n.thin,
n.chains = 1)