
Function for prediction at new locations for single-species spatial N-mixture models
predict.spNMix.Rd
The function predict
collects posterior predictive samples for a set of new locations given an object of class `spNMix`. Prediction is possible for both the latent abundance state as well as detection.
Usage
# S3 method for spNMix
predict(object, X.0, coords.0, n.omp.threads = 1,
verbose = TRUE, n.report = 100, ignore.RE = FALSE,
type = 'abundance', include.sp = TRUE, ...)
Arguments
- object
an object of class spNMix
- X.0
the design matrix of covariates at the prediction locations. This should include a column of 1s for the intercept if an intercept is included in the model. If random effects are included in the abundance (or detection if
type = 'detection'
) portion of the model, the levels of the random effects at the new locations should be included as a column in the design matrix. The ordering of the levels should match the ordering used to fit the data inspNMix
. Columns should correspond to the order of how covariates were specified in the corresponding formula argument ofspNMix
. Column names of all variables must match the names of variables used when fitting the model (for the intercept, use'(Intercept)'
).- coords.0
the spatial coordinates corresponding to
X.0
. Note thatspAbundance
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 random abundance (or detection if
type = 'detection'
) effects from the subsequent predictions. IfTRUE
, random effects will be included. IfFALSE
, random effects will be set to 0 and predictions will only be generated from the fixed effects.- type
a quoted keyword indicating what type of prediction to produce. Valid keywords are 'abundance' to predict latent abundance and expected abundance values (this is the default), or 'detection' to predict detection probability given new values of detection covariates.
- 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 toFALSE
, predictions are given using the covariates and any unstructured random effects in the model. IfFALSE
, thecoords
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 intercept 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.spNMix
. When type = 'abundance'
, the list consists of:
- mu.0.samples
a
coda
object of posterior predictive samples for the expected abundance values. Note these will be per unit area if an offset was used when fitting the model withNMix()
- N.0.samples
a
coda
object of posterior predictive samples for the latent abundance values. These will be in the same units asmu.0.samples
.- w.0.samples
a
coda
object of posterior predictive samples for the latent spatial random effects.
When type = 'detection'
, the list consists of:
- p.0.samples
a
coda
object of posterior predictive samples for the detection probability values.
The return object will include additional objects used for standard extractor functions.
Examples
set.seed(200)
# Simulate Data -----------------------------------------------------------
J.x <- 15
J.y <- 15
J <- J.x * J.y
n.rep <- sample(3, J, replace = TRUE)
beta <- c(0.5, 1.5)
p.abund <- length(beta)
alpha <- c(0.5, 1.2, -0.5)
p.det <- length(alpha)
mu.RE <- list()
p.RE <- list()
phi <- runif(1, 3 / 1, 3 / .1)
sigma.sq <- runif(1, 0.2, 1.5)
kappa <- 0.5
sp <- TRUE
cov.model <- 'exponential'
dat <- simNMix(J.x = J.x, J.y = J.y, n.rep = n.rep, beta = beta, alpha = alpha,
kappa = kappa, mu.RE = mu.RE, p.RE = p.RE, sp = sp,
phi = phi, sigma.sq = sigma.sq, cov.model = cov.model,
family = 'NB')
# Split into fitting and prediction data set
pred.indx <- sample(1:J, round(J * .5), replace = FALSE)
y <- dat$y[-pred.indx, ]
# Abundance covariates
X <- dat$X[-pred.indx, ]
# Prediction covariates
X.0 <- dat$X[pred.indx, ]
# Detection covariates
X.p <- dat$X.p[-pred.indx, , ]
coords <- as.matrix(dat$coords[-pred.indx, ])
coords.0 <- as.matrix(dat$coords[pred.indx, ])
mu.0 <- dat$mu[pred.indx]
w.0 <- dat$w[pred.indx]
abund.covs <- X
colnames(abund.covs) <- c('int', 'abund.cov.1')
det.covs <- list(det.cov.1 = X.p[, , 2], det.cov.2 = X.p[, , 3])
data.list <- list(y = y,
abund.covs = abund.covs,
det.covs = det.covs,
coords = coords)
# Priors
prior.list <- list(beta.normal = list(mean = rep(0, p.abund),
var = rep(100, p.abund)),
alpha.normal = list(mean = rep(0, p.det),
var = rep(2.72, p.det)),
kappa.unif = c(0, 10))
# Starting values
inits.list <- list(alpha = alpha,
beta = beta,
kappa = kappa,
phi = 3 / 0.5,
sigma.sq = 1,
N = apply(y, 1, max, na.rm = TRUE))
# Tuning values
tuning.list <- list(phi = 0.5, kappa = 0.5, beta = 0.1, alpha = 0.1, w = 0.1)
n.batch <- 4
batch.length <- 25
n.burn <- 0
n.thin <- 1
n.chains <- 1
out <- spNMix(abund.formula = ~ abund.cov.1,
det.formula = ~ det.cov.1 + det.cov.2,
data = data.list,
n.batch = n.batch,
batch.length = batch.length,
inits = inits.list,
priors = prior.list,
NNGP = TRUE,
cov.model = 'spherical',
n.neighbors = 10,
accept.rate = 0.43,
n.omp.threads = 1,
verbose = TRUE,
n.report = 1,
n.burn = n.burn,
n.thin = n.thin,
n.chains = n.chains)
#> ----------------------------------------
#> Preparing to run the model
#> ----------------------------------------
#> No prior specified for phi.unif.
#> Setting uniform bounds based on the range of observed spatial coordinates.
#> No prior specified for sigma.sq.
#> Using an inverse-Gamma prior with the shape parameter set to 2 and scale parameter to 1.
#> 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 Poisson N-mixture model with 113 sites.
#>
#> Samples per Chain: 100 (4 batches of length 25)
#> Burn-in: 0
#> Thinning Rate: 1
#> Number of Chains: 1
#> Total Posterior Samples: 100
#>
#> Using the spherical spatial correlation model.
#>
#> Using 10 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 4, 25.00%
#> Parameter Acceptance Tuning
#> beta[1] 0.0 0.98020
#> beta[2] 0.0 0.98020
#> alpha[1] 8.0 0.98020
#> alpha[2] 12.0 0.98020
#> alpha[3] 8.0 0.98020
#> phi 24.0 0.98020
#> -------------------------------------------------
#> Batch: 2 of 4, 50.00%
#> Parameter Acceptance Tuning
#> beta[1] 8.0 0.97045
#> beta[2] 0.0 0.97045
#> alpha[1] 8.0 0.97045
#> alpha[2] 0.0 0.97045
#> alpha[3] 4.0 0.97045
#> phi 36.0 0.97045
#> -------------------------------------------------
#> Batch: 3 of 4, 75.00%
#> Parameter Acceptance Tuning
#> beta[1] 8.0 0.96079
#> beta[2] 0.0 0.96079
#> alpha[1] 24.0 0.96079
#> alpha[2] 20.0 0.96079
#> alpha[3] 8.0 0.96079
#> phi 44.0 0.98020
#> -------------------------------------------------
#> Batch: 4 of 4, 100.00%
summary(out)
#>
#> Call:
#> spNMix(abund.formula = ~abund.cov.1, det.formula = ~det.cov.1 +
#> det.cov.2, data = data.list, inits = inits.list, priors = prior.list,
#> cov.model = "spherical", NNGP = TRUE, n.neighbors = 10, n.batch = n.batch,
#> batch.length = batch.length, accept.rate = 0.43, n.omp.threads = 1,
#> verbose = TRUE, n.report = 1, n.burn = n.burn, n.thin = n.thin,
#> n.chains = n.chains)
#>
#> Samples per Chain: 100
#> Burn-in: 0
#> Thinning Rate: 1
#> Number of Chains: 1
#> Total Posterior Samples: 100
#> Run Time (min): 0.0033
#>
#> Abundance (log scale):
#> Mean SD 2.5% 50% 97.5% Rhat ESS
#> (Intercept) 0.5728 0.0925 0.4506 0.5528 0.6664 NA 2
#> abund.cov.1 1.5000 0.0000 1.5000 1.5000 1.5000 NA 0
#>
#> Detection (logit scale):
#> Mean SD 2.5% 50% 97.5% Rhat ESS
#> (Intercept) 0.6864 0.1798 0.4327 0.6910 0.9305 NA 3
#> det.cov.1 1.0259 0.0976 0.7995 1.0778 1.1366 NA 12
#> det.cov.2 -0.8295 0.1833 -1.2333 -0.8470 -0.5000 NA 5
#>
#> Spatial Covariance:
#> Mean SD 2.5% 50% 97.5% Rhat ESS
#> sigma.sq 3.5759 1.1783 0.5918 3.7347 5.3406 NA 4
#> phi 6.1687 0.9833 4.9977 5.8364 8.2912 NA 15
# Predict at new locations ------------------------------------------------
colnames(X.0) <- c('intercept', 'abund.cov')
out.pred <- predict(out, X.0, coords.0)
#> ----------------------------------------
#> Prediction description
#> ----------------------------------------
#> NNGP spatial abundance model with 113 observations.
#>
#> Number of covariates 2 (including intercept if specified).
#>
#> Using the spherical spatial correlation model.
#>
#> Using 10 nearest neighbors.
#>
#> Number of MCMC samples 100.
#>
#> Predicting at 112 locations.
#>
#> Source compiled with OpenMP support and model fit using 1 threads.
#> -------------------------------------------------
#> Predicting
#> -------------------------------------------------
#> Location: 100 of 112, 89.29%
#> Location: 112 of 112, 100.00%
#> Generating latent abundance state
str(out.pred)
#> List of 4
#> $ N.0.samples : 'mcmc' num [1:100, 1:112] 2 0 0 2 0 0 4 1 0 2 ...
#> ..- attr(*, "mcpar")= num [1:3] 1 100 1
#> $ w.0.samples : 'mcmc' num [1:100, 1:112] -0.312 0.512 0.109 1.356 0.281 ...
#> ..- attr(*, "mcpar")= num [1:3] 1 100 1
#> $ mu.0.samples: 'mcmc' num [1:100, 1:112] 0.396 0.98 0.655 2.28 0.778 ...
#> ..- attr(*, "mcpar")= num [1:3] 1 100 1
#> $ call : language predict.spNMix(object = out, X.0 = X.0, coords.0 = coords.0)
#> - attr(*, "class")= chr "predict.spNMix"