spMvGLM {spBayes}R Documentation

Function for fitting multivariate Bayesian generalized linear spatial regression models

Description

The function spMvGLM fits multivariate stationary Bayesian generalized linear spatial regression models. Given a set of knots, spMvGLM fits a predictive process model (see references below).

Usage

spMvGLM(formula, family="binomial", data = parent.frame(), coords, knots,
      starting, tuning, priors, cov.model, n.samples,
      verbose=TRUE, n.report=100, ...)

Arguments

formula for a multivariate model with q response variables, this is a list of univariate formulas. See example below.
family currently only supports binomial and poisson data using the logit and log link functions, respectively.
data an optional data frame containing the variables in the model. If not found in data, the variables are taken from environment(formula), typically the environment from which spMvGLM is called.
coords an n x 2 matrix of the observation coordinates in R^2 (e.g., easting and northing).
knots either a m x 2 matrix of the predictive process knot coordinates in R^2 (e.g., easting and northing) or a vector of length two or three with the first and second elements recording the number of columns and rows in the desired knot grid. The third, optional, element sets the offset of the outermost knots from the extent of the coords extent.
starting a list with each tag corresponding to a parameter name. Valid list tags are beta, A, phi, and nu. The value portion of each tag is a vector of parameter's starting value. For A the vector is of length q(q-q)/2+q and phi and nu are of length q. Here, A holds the the lower-triangle elements in column major ordering of the Cholesky square root of the spatial cross-covariance matrix. If the predictive process is used then w must be of length qm; otherwise, it must be of length qn. Alternatively, w can be set as a scalar, in which case the value is repeated.
tuning a list with each tag corresponding to a parameter name. Valid list tags are beta, A, phi, nu, and w. The value portion of each tag defines the variance of the Metropolis normal proposal distribution. The tuning value for beta can be a vector of length p or the lower-triangle of the pxp Cholesky square-root of the desired proposal variance matrix. For A, the vector is of length q(q-q)/2+q and phi and nu are of length q. If the predictive process is used then w must be of length m; otherwise, it must be of length n. Alternatively, w can be set as a scalar, in which the value is repeated.
priors a list with each tag corresponding to a parameter name. Valid list tags are beta.flat, beta.normal, K.IW, Psi.IW, phi.unif, and nu.unif. If beta.normal then covariate specific mean and variance hyperparameters are passed as the first and second list elements, respectively. K is assumed to follow an inverse-Wishart distribution, whereas the spatial range phi and smoothness nu parameters are assumed to follow Uniform distributions. The hyperparameters of the inverse-Wishart are passed as a list of length two, with the first and second elements corresponding to the df and qxq scale matrix, respectively. The hyperparameters of the Uniform are also passed as a vector of length 2xq with consecutive elements representing the first and second elements corresponding to the lower and upper support in the order of the univariate models given in formula.
cov.model a quoted key word 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". See below for details.
n.samples the number of MCMC iterations.
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 Metropolis acceptance and MCMC progress.
... currently no additional arguments.

Value

An object of class spMvGLM, which is a list with the following tags:

coords the n x 2 matrix specified by coords.
knot.coords the m x 2 matrix as specified by knots.
p.samples a coda object of posterior samples for the defined parameters.
acceptance the Metropolis sampling acceptance rate.
sp.effects a matrix that holds samples from the posterior distribution of the spatial random effects. The rows of this matrix correspond to the n point observations and the columns are the posterior samples.

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

Author(s)

Andrew O. Finley finleya@msu.edu,
Sudipto Banerjee baner009@umn.edu

References

Finley, A.O., S. Banerjee, and R.E. McRoberts. (2008) A Bayesian approach to quantifying uncertainty in multi-source forest area estimates. Environmental and Ecological Statistics, 15:241–258.

Banerjee, S., A.E. Gelfand, A.O. Finley, and H. Sang. (2008) Gaussian Predictive Process Models for Large Spatial Datasets. Journal of the Royal Statistical Society Series B, 70:825–848.

Finley, A.O., S. Banerjee, P. Waldmann, and T. Ericsson. (2008). Hierarchical spatial modeling of additive and dominance genetic variance for large spatial trial datasets. Biometrics. DOI: 10.1111/j.1541-0420.2008.01115.x

Finley, A.O,. H. Sang, S. Banerjee, and A.E. Gelfand. (2008). Improving the performance of predictive process modeling for large datasets. Computational Statistics and Data Analysis, DOI: 10.1016/j.csda.2008.09.008

Banerjee, S., Carlin, B.P., and Gelfand, A.E. (2004). Hierarchical modeling and analysis for spatial data. Chapman and Hall/CRC Press, Boca Raton, Fla.

See Also

spGGT,spMvLM, spMvGLM

Examples

## Not run: 
################################
##Spatial multivariate poisson
################################

##Generate some data
n <- 100
q <- 3   
nltr <- q*(q-1)/2+q

coords <- cbind(runif(n,1,100),runif(n,1,100))

theta <- rep(3/50,q)

A <- matrix(0,q,q)
A[lower.tri(A,TRUE)] <- rnorm(nltr,1,1)
K <- A%*%t(A)

Psi <- diag(0,q)

c1 <- mvCovInvLogDet(coords=coords, knots=coords,
                     cov.model="exponential",
                     V=K, Psi=Psi, theta=theta,
                     modified.pp=TRUE, SWM=FALSE)

w <- mvrnorm(1,rep(0,nrow(c1$C)),c1$C)

X <- mkMvX(list(matrix(1,n,1), matrix(1,n,1), matrix(1,n,1)))
beta <- c(-1,0,1)
y <- rpois(n*q, exp(X%*%beta+w))

y.1 <- y[seq(1,length(y),q)]
y.2 <- y[seq(2,length(y),q)]
y.3 <- y[seq(3,length(y),q)]

##Specify starting values and collect samples. For
##a true analysis, several longer chains should be
##run.
A.starting <- diag(1,q)[lower.tri(diag(1,q), TRUE)]

beta.starting <- coefficients(glm(y~X-1, family="poisson"))
beta.tuning <- t(chol(vcov(glm(y~X-1, family="poisson"))))

n.samples <- 15000

m.1 <- spMvGLM(list(y.1~1,y.2~1,y.3~1), family="poisson", coords=coords,
               knots=c(8,8,0),
               starting=
               list("beta"=beta.starting, "phi"=rep(0.06,q),
                    "A"=A.starting, "w"=0), 
               tuning=
               list("beta"=beta.tuning, "phi"=rep(0.01,q),
                    "A"=rep(0.005,nltr), "w"=0.001),
               priors=
               list("beta.Flat", "phi.Unif"=rep(c(0.03, 0.3),q),
                    "K.IW"=list(q+1, diag(0.1,q))),
               cov.model="exponential",
               n.samples=n.samples, verbose=TRUE, n.report=500)

m.1$p.samples[,paste("phi_",1:q,sep="")] <-
  3/m.1$p.samples[,paste("phi_",1:q,sep="")]

burn.in <- 0.75*n.samples
print(summary(mcmc(m.1$p.samples[burn.in:n.samples,])))

beta.hat <- colMeans(m.1$p.samples[burn.in:n.samples,1:q])
w.hat <- rowMeans(m.1$sp.effects[,burn.in:n.samples])

y.hat <- exp(X%*%beta.hat+w.hat)

y.hat.1 <- y.hat[seq(1,length(y.hat),q)]
y.hat.2 <- y.hat[seq(2,length(y.hat),q)]
y.hat.3 <- y.hat[seq(3,length(y.hat),q)]

##Take a look
par(mfrow=c(3,2))
surf <- mba.surf(cbind(coords,y.1),no.X=100, no.Y=100, extend=TRUE)$xyz.est
image(surf, main="Observed counts")
contour(surf, drawlabels=FALSE, add=TRUE)
text(coords, labels=y.1, cex=1, col="blue")

surf <- mba.surf(cbind(coords,y.hat.1),no.X=100, no.Y=100, extend=TRUE)$xyz.est
image(surf, main="Fitted counts")
contour(surf, add=TRUE)
contour(surf, drawlabels=FALSE, add=TRUE)
text(coords, labels=round(y.hat.1,0), cex=1, col="blue")

surf <- mba.surf(cbind(coords,y.2),no.X=100, no.Y=100, extend=TRUE)$xyz.est
image(surf)
contour(surf, drawlabels=FALSE, add=TRUE)
text(coords, labels=y.2, cex=1, col="blue")

surf <- mba.surf(cbind(coords,y.hat.2),no.X=100, no.Y=100, extend=TRUE)$xyz.est
image(surf)
contour(surf, drawlabels=FALSE, add=TRUE)
text(coords, labels=round(y.hat.2,0), cex=1, col="blue")

surf <- mba.surf(cbind(coords,y.3),no.X=100, no.Y=100, extend=TRUE)$xyz.est
image(surf)
contour(surf, drawlabels=FALSE, add=TRUE)
text(coords, labels=y.3, cex=1, col="blue")

surf <- mba.surf(cbind(coords,y.hat.3),no.X=100, no.Y=100, extend=TRUE)$xyz.est
image(surf)
contour(surf, drawlabels=FALSE, add=TRUE)
text(coords, labels=round(y.hat.3,0), cex=1, col="blue")



## End(Not run)

[Package spBayes version 0.1-2 Index]