Skip to contents

The function predict collects posterior predictive samples for a set of new locations given an object of class `sfMsAbund`.

Usage

# S3 method for sfMsAbund
predict(object, X.0, coords.0, n.omp.threads = 1, 
        verbose = TRUE, n.report = 100, ignore.RE = FALSE, 
        z.0.samples, include.sp = TRUE, ...)

Arguments

object

an object of class sfMsAbund

X.0

the design matrix of covariates at the prediction locations. This can be either a two-dimensional matrix with rows corresponding to sites and columns corresponding to covariates, or can be a three-dimensional array, with dimensions corresponding to site, replicate, and covariate, respectively. Note that the first covariate should consist of all 1s for the intercept if an intercept is included in the model. If random effects are included in the the model, the levels of the random effects at the new locations/time periods should be included as an element of the three-dimensional array. The ordering of the levels should match the ordering used to fit the data in sfMsAbund. The covariates should be organized in the same order as they were specified in the corresponding formula argument of sfMsAbund. Names of the third dimension (covariates) of any random effects in X.0 must match the name of the random effects used to fit the model, if specified in the corresponding formula argument of sfMsAbund. See example below.

coords.0

the spatial coordinates corresponding to X.0. Note that spAbundance assumes coordinates are specified in a projected coordinate system.

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, model specification and progress of the sampler is printed to the screen. Otherwise, nothing is printed to the screen.

n.report

the interval to report sampling progress.

ignore.RE

logical value that specifies whether or not to remove unstructured random effects from the subsequent predictions. If TRUE, unstructured random effects will be included. If FALSE, unstructured random effects will be set to 0 and predictions will only be generated from the fixed effects.

z.0.samples

a three-dimensional array with dimensions corresponding to MCMC samples, species, and prediction locations. The array contains the full posterior samples of the predicted binary portion of a zero-inflated Gaussian model. In the context of abundance models, this typically corresponds to estimates of the presence or absence of each species at the location. When using spOccupancy to generate the first stage samples of the zero-inflated Gaussian model, this is the object contained in the z.0.samples object of the predition function for the spOccupancy object. Ignored for all model types other than zero-inflated Gaussian.

include.sp

a logical value used to indicate whether spatial random effects should be included in the predictions. By default, this is set to TRUE. If set to FALSE, predictions are given using the covariates and any unstructured random effects in the model. If FALSE, the coords.0 argument is not required.

...

currently no additional arguments

Note

When ignore.RE = FALSE, both sampled levels and non-sampled levels of random effects are supported for prediction. For sampled levels, the posterior distribution for the random effect corresponding to that level of the random effect will be used in the prediction. For non-sampled levels, random values are drawn from a normal distribution using the posterior samples of the random effect variance, which results in fully propagated uncertainty in predictions with models that incorporate random effects.

Author

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

Value

A list object of class predict.sfMsAbund. The list consists of:

mu.0.samples

a three or four-dimensional object of posterior predictive samples for the expected abundance values with dimensions corresponding to posterior predictive sample, species, site, and replicate. Note if an offset was used when fitting the model with sfMsAbund, the abundance values are reported per unit of the offset.

y.0.samples

a three or four-dimensional object of posterior predictive samples for the abundance values with dimensions corresponding to posterior predictive sample, species, site, and replicate. These will be in the same units as mu.0.samples.

w.0.samples

a three-dimensional array of posterior predictive samples for the latent factors. Array dimensions correspond to MCMC sample, latent factor, and site.

The return object will include additional objects used for standard extractor functions.

Examples

set.seed(408)
J.x <- 8
J.y <- 8
J <- J.x * J.y
n.rep <- sample(3, size = J, replace = TRUE)
n.sp <- 6
# Community-level covariate effects
beta.mean <- c(-2, 0.5)
p.abund <- length(beta.mean)
tau.sq.beta <- c(0.2, 1.2)
mu.RE <- list()
# Draw species-level effects from community means.
beta <- matrix(NA, nrow = n.sp, ncol = p.abund)
for (i in 1:p.abund) {
  beta[, i] <- rnorm(n.sp, beta.mean[i], sqrt(tau.sq.beta[i]))
}
sp <- TRUE 
kappa <- runif(n.sp, 0.1, 1)
factor.model <- TRUE
n.factors <- 3
cov.model <- 'spherical'
phi <- runif(n.factors, 3 / 1, 3 / .1)

dat <- simMsAbund(J.x = J.x, J.y = J.y, n.rep = n.rep, n.sp = n.sp, beta = beta,
                  mu.RE = mu.RE, sp = sp, kappa = kappa, family = 'NB', 
                  factor.model = factor.model, n.factors = n.factors, 
                  phi = phi, cov.model = cov.model)

# Split into fitting and prediction data set
pred.indx <- sample(1:J, round(J * .25), replace = FALSE)
y <- dat$y[, -pred.indx, , drop = FALSE]
# Occupancy covariates
X <- dat$X[-pred.indx, , , drop = FALSE]
# Coordinates
coords <- dat$coords[-pred.indx, ]
# Prediction values
y.0 <- dat$y[, pred.indx, , drop = FALSE]
X.0 <- dat$X[pred.indx, , , drop = FALSE]
coords.0 <- dat$coords[pred.indx, ]

# Package all data into a list
covs <- list(int = X[, , 1], abund.cov.1 = X[, , 2])
data.list <- list(y = y, covs = covs, coords = coords)
prior.list <- list(beta.comm.normal = list(mean = 0, var = 100),
                   kappa.unif = list(a = 0, b = 10),
                   phi.unif = list(a = 3 / 1, b = 3 / .1),
                   tau.sq.beta.ig = list(a = .1, b = .1))
inits.list <- list(beta.comm = 0,
                   beta = 0,
                   kappa = 0.5,
                   phi = 3 / .5,
                   tau.sq.beta = 1)
tuning.list <- list(kappa = 0.3, beta = 0.1, lambda = 0.5, w = 0.5, 
                    phi = 1)

# Small
n.batch <- 2
batch.length <- 25
n.burn <- 20
n.thin <- 1
n.chains <- 1

out <- sfMsAbund(formula = ~ abund.cov.1,
                 data = data.list,
                 n.batch = n.batch,
                 inits = inits.list,
                 priors = prior.list,
                 tuning = tuning.list,
                 batch.length = batch.length,
                 n.omp.threads = 1,
                 n.factors = 1,
                 verbose = TRUE,
                 n.neighbors = 5, 
                 n.report = 1,
                 n.burn = n.burn,
                 n.thin = n.thin,
                 n.chains = n.chains)
#> ----------------------------------------
#> 	Preparing to run the model
#> ----------------------------------------
#> lambda is not specified in initial values.
#> Setting initial values of the lower triangle to 0
#> 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 Factor NNGP Multi-species Poisson Abundance
#> model fit with 48 sites and 6 species.
#> 
#> Samples per Chain: 50 
#> Burn-in: 20 
#> Thinning Rate: 1 
#> Number of Chains: 1 
#> Total Posterior Samples: 30 
#> 
#> Using the exponential spatial correlation model.
#> 
#> Using 1 latent spatial factors.
#> 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: 1 of 2, 50.00%
#> 	Number	Parameter	Acceptance	Tuning
#> 	1	beta[1]		84.0		0.10202
#> 	2	beta[1]		88.0		0.10202
#> 	3	beta[1]		68.0		0.10202
#> 	4	beta[1]		60.0		0.10202
#> 	5	beta[1]		76.0		0.10202
#> 	6	beta[1]		68.0		0.10202
#> 	1	phi		40.0		1.00000
#> -------------------------------------------------
#> Batch: 2 of 2, 100.00%

# Predict at new locations
out.pred <- predict(out, X.0, coords.0)
#> ----------------------------------------
#> 	Prediction description
#> ----------------------------------------
#> Spatial Factor NNGP GLMM with 48 observations.
#> 
#> Number of covariates 2 (including intercept if specified).
#> 
#> Using the exponential spatial correlation model.
#> 
#> Using 5 nearest neighbors.
#> Using 1 latent spatial factors.
#> 
#> Number of MCMC samples 30.
#> 
#> Predicting at 16 non-sampled locations.
#> 
#> 
#> Source compiled with OpenMP support and model fit using 1 threads.
#> -------------------------------------------------
#> 		Predicting
#> -------------------------------------------------
#> Location: 16 of 16, 100.00%
#> Generating abundance predictions
str(out.pred)
#> List of 6
#>  $ y.0.samples : num [1:30, 1:6, 1:16, 1:3] NA NA NA NA NA NA NA NA NA NA ...
#>  $ w.0.samples : num [1:30, 1, 1:16] 0.631 -1.318 0.163 0.437 -0.592 ...
#>  $ mu.0.samples: num [1:30, 1:6, 1:16, 1:3] NA NA NA NA NA NA NA NA NA NA ...
#>  $ run.time    : 'proc_time' Named num [1:5] 0.002 0.001 0.003 0 0
#>   ..- attr(*, "names")= chr [1:5] "user.self" "sys.self" "elapsed" "user.child" ...
#>  $ call        : language predict.sfMsAbund(object = out, X.0 = X.0, coords.0 = coords.0)
#>  $ object.class: chr "sfMsAbund"
#>  - attr(*, "class")= chr "predict.sfMsAbund"