Skip to contents

The function predict collects posterior predictive samples for a set of new locations given an object of class `spDS`. Prediction is possible for both the latent abundance state as well as detection.

Usage

# S3 method for spDS
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 spDS

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 in spDS. Columns should correspond to the order of how covariates were specified in the corresponding formula argument of spDS. 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 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 random abundance (or detection if type = 'detection') effects from the subsequent predictions. If TRUE, random effects will be included. If FALSE, 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 to FALSE, predictions are given using the covariates and any unstructured random effects in the model. If FALSE, the coords 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.spDS. When type = 'abundance', the list consists of:

mu.0.samples

a coda object of posterior predictive samples for the expected abundance values, or expected abundance per unit area (i.e., density) values when an offset was used when fitting the model with spDS().

N.0.samples

a coda object of posterior predictive samples for the latent abundance values. These will be in the same units as mu.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:

sigma.0.samples

a coda object of posterior predictive samples for sigma (the parameter controlling detection probability).

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

Examples

set.seed(123)
J.x <- 10
J.y <- 10 
J <- J.x * J.y
# Number of distance bins from which to simulate data. 
n.bins <- 5
# Length of each bin. This should be of length n.bins
bin.width <- c(.10, .10, .20, .3, .1)
# Abundance coefficients
beta <- c(1.0, 0.2, 0.3, -0.2)
p.abund <- length(beta)
# Detection coefficients
alpha <- c(-1.0, -0.3)
p.det <- length(alpha)
# Detection decay function
det.func <- 'halfnormal'
mu.RE <- list()
p.RE <- list()
sp <- TRUE
phi <- 3 / .5
sigma.sq <- 0.8
cov.model <- 'exponential'
family <- 'NB'
kappa <- 0.1
offset <- 1.8
transect <- 'point'

dat <- simDS(J.x = J.x, J.y = J.y, n.bins = n.bins, bin.width = bin.width,
             beta = beta, alpha = alpha, det.func = det.func, kappa = kappa, 
             mu.RE = mu.RE, p.RE = p.RE, sp = sp, family = family, 
             offset = offset, transect = transect, phi = phi, sigma.sq = sigma.sq,
             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, ]
# Abundance covariates
X <- dat$X[-pred.indx, ]
# Prediction covariates
X.0 <- dat$X[pred.indx, ]
# Detection covariates
X.p <- dat$X.p[-pred.indx, ]
dist.breaks <- dat$dist.breaks
coords <- dat$coords[-pred.indx, ]
coords.0 <- dat$coords[pred.indx, ]

covs <- cbind(X, X.p)
colnames(covs) <- c('int.abund', 'abund.cov.1', 'abund.cov.2', 'abund.cov.3', 
                    'int.det', 'det.cov.1')

data.list <- list(y = y, 
                  covs = covs,
                  dist.breaks = dist.breaks, 
                  coords = coords,
                  offset = offset)

# Priors
prior.list <- list(beta.normal = list(mean = 0, var = 10),
                   alpha.normal = list(mean = 0,
                                       var = 10), 
                   kappa.unif = c(0, 100), 
                   phi.unif = c(3 / 1, 3 / .1),
                   sigma.sq.ig = c(2, 1)) 
# Starting values
inits.list <- list(alpha = 0,
                   beta = 0,
                   kappa = 1, 
                   phi = 3 / .5, 
                   sigma.sq = 1)
# Tuning values
tuning <- list(beta = 0.1, alpha = 0.1, beta.star = 0.3, alpha.star = 0.1, 
               kappa = 0.2, phi = 1, w = 1) 

out <- spDS(abund.formula = ~ abund.cov.1 + abund.cov.2 + abund.cov.3,
            det.formula = ~ det.cov.1,
            data = data.list, 
            n.batch = 10, 
            batch.length = 25, 
            inits = inits.list, 
            family = 'NB',
            det.func = 'halfnormal', 
            transect = 'point', 
            cov.model = 'exponential', 
            NNGP = TRUE,
            n.neighbors = 5,
            tuning = tuning,
            priors = prior.list, 
            accept.rate = 0.43, 
            n.omp.threads = 1, 
            verbose = TRUE, 
            n.report = 100,
            n.burn = 100,
            n.thin = 1,
            n.chains = 1) 
#> ----------------------------------------
#> 	Preparing to run the model
#> ----------------------------------------
#> **NOTE**: spatial negative binomial models can be difficult to
#> estimate as they contain two forms of overdispersion. If experiencing
#> very poor mixing/convergence of MCMC chains (particularly kappa and phi),
#> consider using a spatial Poisson model or more informative
#> priors on kappa or phi.
#> N is not specified in initial values.
#> Setting initial values based on observed data
#> 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 Negative Binomial HDS model with 75 sites.
#> 
#> Samples per Chain: 250 (10 batches of length 25)
#> Burn-in: 100 
#> Thinning Rate: 1 
#> Number of Chains: 1 
#> Total Posterior Samples: 150 
#> 
#> 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%
summary(out)
#> 
#> Call:
#> spDS(abund.formula = ~abund.cov.1 + abund.cov.2 + abund.cov.3, 
#>     det.formula = ~det.cov.1, data = data.list, inits = inits.list, 
#>     priors = prior.list, tuning = tuning, cov.model = "exponential", 
#>     NNGP = TRUE, n.neighbors = 5, n.batch = 10, batch.length = 25, 
#>     accept.rate = 0.43, family = "NB", transect = "point", det.func = "halfnormal", 
#>     n.omp.threads = 1, verbose = TRUE, n.report = 100, n.burn = 100, 
#>     n.thin = 1, n.chains = 1)
#> 
#> Samples per Chain: 250
#> Burn-in: 100
#> Thinning Rate: 1
#> Number of Chains: 1
#> Total Posterior Samples: 150
#> Run Time (min): 0.0034
#> 
#> Abundance (log scale): 
#>               Mean     SD    2.5%    50%  97.5% Rhat ESS
#> (Intercept) 1.3766 0.2580  0.9201 1.4176 1.8294   NA   6
#> abund.cov.1 0.0625 0.5817 -1.0732 0.0797 1.0241   NA   2
#> abund.cov.2 1.2440 0.6769  0.2924 0.9531 2.3067   NA   2
#> abund.cov.3 0.6211 0.3685 -0.1991 0.7685 1.0279   NA   3
#> 
#> Detection (log scale): 
#>                Mean     SD    2.5%     50%   97.5% Rhat ESS
#> (Intercept) -1.0568 0.0675 -1.1547 -1.0573 -0.9311   NA   3
#> det.cov.1   -0.6874 0.1247 -0.9101 -0.7022 -0.4674   NA   9
#> 
#> Spatial Covariance: 
#>            Mean     SD   2.5%    50%   97.5% Rhat ESS
#> sigma.sq 0.6580 0.3474 0.1737 0.6621  1.3966   NA   4
#> phi      9.7562 5.7922 3.3895 7.9988 22.7218   NA   8
#> 
#> NB overdispersion: 
#>         Mean     SD   2.5%    50%  97.5% Rhat ESS
#> kappa 0.0598 0.0158 0.0354 0.0578 0.0915   NA  10

# Predict at new locations ------------------------------------------------
colnames(X.0) <- c('intercept', 'abund.cov.1', 'abund.cov.2', 'abund.cov.3')
out.pred <- predict(out, X.0, coords.0)
#> ----------------------------------------
#> 	Prediction description
#> ----------------------------------------
#> NNGP spatial abundance model with 75 observations.
#> 
#> Number of covariates 4 (including intercept if specified).
#> 
#> Using the exponential spatial correlation model.
#> 
#> Using 5 nearest neighbors.
#> 
#> Number of MCMC samples 150.
#> 
#> Predicting at 25 locations.
#> 
#> Source compiled with OpenMP support and model fit using 1 threads.
#> -------------------------------------------------
#> 		Predicting
#> -------------------------------------------------
#> Location: 25 of 25, 100.00%
#> Generating latent abundance state
mu.0.quants <- apply(out.pred$mu.0.samples, 2, quantile, c(0.025, 0.5, 0.975))
plot(dat$mu[pred.indx], mu.0.quants[2, ], pch = 19, xlab = 'True', 
     ylab = 'Fitted', ylim = c(min(mu.0.quants), max(mu.0.quants)))
segments(dat$mu[pred.indx], mu.0.quants[1, ], dat$mu[pred.indx], mu.0.quants[3, ])
lines(dat$mu[pred.indx], dat$mu[pred.indx])