Skip to contents

The function svcAbund fits univariate spatially-varying coefficient GLMMs.

Usage

svcAbund(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, family = 'Gaussian',
         n.omp.threads = 1, verbose = TRUE, n.report = 100, 
         n.burn = round(.10 * n.batch * batch.length), n.thin = 1, 
         n.chains = 1, ...)

Arguments

formula

a symbolic description of the model to be fit for the model using R's model syntax. Only right-hand side of formula is specified. See example below. Random intercepts and slopes are allowed using lme4 syntax (Bates et al. 2015).

data

a list containing data necessary for model fitting. Valid tags are y, covs, coords, and z (for family = 'zi-Gaussian' only. y is a vector of the observed count values, where the values represent the observed values at each site. covs is a list, matrix, or data frame of covariates used in the model, where each column (or list element) represents a different covariate. coords is a \(J \times 2\) matrix of the observation coordinates. Note that svcAbundance assumes coordinates are specified in a projected coordinate system. z is used for fitting a zero-inflated Gaussian model. It is a vector where each value indicates the binary component of the model. In the context of abundance models, this can be thought of as the component of the model that indicates whether the species is present at each location, and then the supplied values in y are the observed abundance values at those locations where z = 1.

inits

a list with each tag corresponding to a parameter name. Valid tags are beta, sigma.sq, phi, w, nu, tau.sq, sigma.sq.mu. nu is only specified if cov.model = "matern", sigma.sq.mu is only specified if there are random effects in formula, and The value portion of each tag is the parameter's initial value. See priors description for definition of each parameter name. Additionally, the tag fix can be set to TRUE to fix the starting values across all chains. If fix is not specified (the default), starting values are varied randomly across chains.

priors

a list with each tag corresponding to a parameter name. Valid tags are beta.normal, phi.unif, sigma.sq.ig, nu.unif, tau.sq.ig, sigma.sq.mu.ig. Abundance (beta) regression coefficients are assumed to follow a normal distribution. The hyperparameters of the normal distribution are passed as a list of length two with the first and second elements corresponding to the mean and variance of the normal distribution, which are each specified as vectors of length equal to the number of coefficients to be estimated or of length one if priors are the same for all coefficients. If not specified, prior means are set to 0 and prior variances are set to 100. The spatial variance parameter, sigma.sq, and the Gaussian residual variance parameter, tau.sq, are assumed to follow an inverse-Gamma distribution. The spatial decay phi and spatial smoothness nu, parameters are assumed to follow Uniform distributions. The hyperparameters of the inverse-Gamma for sigma.sq is passed as a list of length two with the first and second elements corresponding to the shape and scale parameters of the inverse-Gamma distribution either for each spatially-varying coefficient, or a single value if assumign the same values for all spatially-varying coefficients. The hyperparameters of the inverse-Gamma for tau.sq is passed as a vector of length two, with the first and second elements corresponding to the shape and scale, respectively. The hyperparameters of the Uniform are also passed as a list of length two with the first and second elements corresponding to the lower and upper support, respectively, for each SVC or a single value if giving the same prior for each SVC. sigma.sq.mu are the random effect variances for any random effects, and are assumed to follow an inverse-Gamma distribution. The hyperparameters of the inverse-Gamma distribution are passed as a list of length two with the first and second elements corresponding to the shape and scale parameters, respectively, which are each specified as vectors of length equal to the number of random effects or of length one if priors are the same for all random effect variances.

svc.cols

a vector indicating the variables whose effects will be estimated as spatially-varying coefficients. svc.cols can be an integer vector with values indicating the order of covariates specified in the model formula (with 1 being the intercept if specified), or it can be specified as a character vector with names corresponding to variable names in occ.covs (for the intercept, use '(Intercept)'). svc.cols default argument of 1 results in a univariate spatial GLMM analogous to spAbund (assuming an intercept is included in the model).

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: "exponential", "matern", "spherical", and "gaussian".

tuning

a single numeric value representing the initial variance of the adaptive sampler for phi and nu. See Roberts and Rosenthal (2009) for details.

NNGP

if TRUE, model is fit with an NNGP. See Datta et al. (2016) and Finley et al. (2019) for more information. Currently only NNGP is supported, functionality for a full GP may be addded in future package development.

n.neighbors

number of neighbors used in the NNGP. Only used if NNGP = TRUE. Datta et al. (2016) showed that 15 neighbors is usually sufficient, but that as few as 5 neighbors can be adequate for certain data sets, which can lead to even greater decreases in run time. We recommend starting with 15 neighbors (the default) and if additional gains in computation time are desired, subsequently compare the results with a smaller number of neighbors using WAIC.

search.type

a quoted keyword that specifies the type of nearest neighbor search algorithm. Supported method key words are: "cb" and "brute". The "cb" should generally be much faster. If locations do not have identical coordinate values on the axis used for the nearest neighbor ordering then "cb" and "brute" should produce identical neighbor sets. However, if there are identical coordinate values on the axis used for nearest neighbor ordering, then "cb" and "brute" might produce different, but equally valid, neighbor sets, e.g., if data are on a grid.

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.

family

the distribution to use for abundance. Currently, spatially-varying coefficient models are available for family = 'Gaussian' and family = 'zi-Gaussian'.

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 n.omp.threads up to the number of hyperthreaded cores. Note, n.omp.threads > 1 might not work on some systems.

verbose

if TRUE, messages about data preparation, model specification, and progress of the sampler are printed to the screen. Otherwise, no messages are printed.

n.report

the interval to report Metropolis sampler acceptance and MCMC progress.

n.burn

the number of samples out of the total n.batch * batch.length samples in each chain to discard as burn-in. By default, the first 10% of samples is discarded.

n.thin

the thinning interval for collection of MCMC samples. The thinning occurs after the n.burn samples are discarded. Default value is set to 1.

n.chains

the number of MCMC chains to run in sequence.

...

currently no additional arguments

References

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 .

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 .

Roberts, G.O. and Rosenthal J.S. (2009) Examples of adaptive MCMC. Journal of Computational and Graphical Statistics, 18(2):349-367.

Author

Jeffrey W. Doser doserjef@msu.edu,
Andrew O. Finley finleya@msu.edu

Value

An object of class svcAbund that is a list comprised of:

beta.samples

a coda object of posterior samples for the abundance regression coefficients.

tau.sq.samples

a coda object of posterior samples for the residual variance parameter.

y.rep.samples

a coda object of posterior samples for the abundance replicate (fitted) values with dimensions corresponding to MCMC samples and site.

mu.samples

a coda object of posterior samples for the expected abundance samples with dimensions corresponding to MCMC samples and site.

theta.samples

a coda object of posterior samples for spatial covariance parameters.

w.samples

a three-dimensional array of posterior samples for the spatially-varying coefficients with dimensions corresponding to MCMC sample, SVC, and site.

sigma.sq.mu.samples

a coda object of posterior samples for variances of random effects included in the model. Only included if random effects are specified in formula.

beta.star.samples

a coda object of posterior samples for the abundance random effects. Only included if random effects are specified in formula.

like.samples

a coda object of posterior samples for the likelihood value associated with each site. 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 proc.time().

The return object will include additional objects used for subsequent prediction and/or model fit evaluation.

Examples

set.seed(1000)
# Sites
J.x <- 10
J.y <- 10
J <- J.x * J.y
# Abundance ---------------------------
beta <- c(5, 0.5, -0.2, 0.75)
p <- length(beta)
mu.RE <- list()
mu.RE <- list(levels = c(35, 40),
              sigma.sq.mu = c(0.7, 1.5),
              beta.indx = list(1, 1))
# Spatial parameters ------------------
sp <- TRUE
svc.cols <- c(1, 2)
p.svc <- length(svc.cols)
cov.model <- "exponential"
sigma.sq <- runif(p.svc, 0.4, 4)
phi <- runif(p.svc, 3/1, 3/0.6)
tau.sq <- 2
z <- rbinom(J, 1, 0.5)

# Get all the data
dat <- simAbund(J.x = J.x, J.y = J.y, beta = beta, tau.sq = tau.sq,
                mu.RE = mu.RE, sp = sp, svc.cols = svc.cols, 
                family = 'zi-Gaussian', cov.model = cov.model, 
                sigma.sq = sigma.sq, phi = phi, z = z)
# Get data in format for spAbundance --------------------------------------
y <- dat$y
X <- dat$X
X.re <- dat$X.re
coords <- dat$coords

# Package all data into a list
covs <- cbind(X, X.re)
colnames(covs) <- c('int', 'cov.1', 'cov.2', 'cov.3', 'factor.1', 'factor.2')

# Data list bundle
data.list <- list(y = y, covs = covs, coords = coords, z = z)
# Priors
prior.list <- list(beta.normal = list(mean = 0, var = 1000),
                   sigma.sq.ig = list(a = 2, b = 1), tau.sq = c(2, 1),
                   sigma.sq.mu.ig = list(a = 2, b = 1),
                   phi.unif = list(a = 3 / 1, b = 3 / 0.1))

# Starting values
inits.list <- list(beta = 0, alpha = 0,
                   sigma.sq = 1, phi = 3 / 0.5, 
                   tau.sq = 2, sigma.sq.mu = 0.5)
# Tuning
tuning.list <- list(phi = 1)

n.batch <- 10
batch.length <- 25
n.burn <- 100
n.thin <- 1

out <- svcAbund(formula = ~ cov.1 + cov.2 + cov.3 + 
                            (1 | factor.1) + (1 | factor.2),
                svc.cols = c(1, 2),
                data = data.list,
                n.batch = n.batch,
                batch.length = batch.length,
                inits = inits.list,
                priors = prior.list,
                accept.rate = 0.43,
                family = 'zi-Gaussian',
                cov.model = "exponential",
                tuning = tuning.list,
                n.omp.threads = 1,
                verbose = TRUE,
                NNGP = TRUE,
                n.neighbors = 5,
                n.report = 25,
                n.burn = n.burn,
                n.thin = n.thin,
                n.chains = 3)
#> ----------------------------------------
#> 	Preparing to run the model
#> ----------------------------------------
#> No prior specified for tau.sq.
#> Using an inverse-Gamma prior with the shape parameter set to 2 and scale parameter to 0.5.
#> w is not specified in initial values.
#> Setting initial value to 0
#> ----------------------------------------
#> 	Building the neighbor list
#> ----------------------------------------
#> ----------------------------------------
#> Building the neighbors of neighbors list
#> ----------------------------------------
#> ----------------------------------------
#> 	Model description
#> ----------------------------------------
#> Spatial NNGP model with 53 sites.
#> 
#> Samples per chain: 250 (10 batches of length 25)
#> Burn-in: 100 
#> Thinning Rate: 1 
#> Number of Chains: 3 
#> Total Posterior Samples: 450 
#> 
#> Number of spatially-varying coefficients: 2 
#> Using the exponential spatial correlation model.
#> 
#> Using 5 nearest neighbors.
#> 
#> Source compiled with OpenMP support and model fit using 1 thread(s).
#> 
#> Adaptive Metropolis with target acceptance rate: 43.0
#> ----------------------------------------
#> 	Chain 1
#> ----------------------------------------
#> Sampling ... 
#> Batch: 10 of 10, 100.00%
#> ----------------------------------------
#> 	Chain 2
#> ----------------------------------------
#> Sampling ... 
#> Batch: 10 of 10, 100.00%
#> ----------------------------------------
#> 	Chain 3
#> ----------------------------------------
#> Sampling ... 
#> Batch: 10 of 10, 100.00%