spSVC {spBayes} | R Documentation |
Function for fitting univariate Bayesian spatially-varying coefficient regression models
Description
The function spSVC
fits Gaussian univariate Bayesian spatially-varying
coefficient regression models.
Usage
spSVC(formula, data = parent.frame(), svc.cols=1, coords,
priors, starting, tuning, cov.model, center.scale=FALSE,
amcmc, n.samples, n.omp.threads = 1,
verbose=TRUE, n.report=100, ...)
Arguments
formula |
a symbolic description of the regression model to be fit. See example below. |
data |
an optional data frame containing the variables in the
model. If not found in data, the variables are taken from
|
svc.cols |
a vector indicating which columns of the regression
design matrix |
coords |
an |
priors |
a list with each tag corresponding to a
parameter name. Valid tags are There are two specification for the Gaussian Process (GP) on the If the regression coefficients, i.e., |
starting |
a list with each tag corresponding to a
parameter name. Valid tags are |
tuning |
a list with each tag corresponding to a
parameter name. Valid tags are |
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:
|
center.scale |
if |
amcmc |
a list with tags |
n.samples |
the number of MCMC iterations. This argument is
ignored if |
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 |
verbose |
if |
n.report |
the interval to report Metropolis sampler acceptance and MCMC progress. |
... |
currently no additional arguments. |
Details
Model parameters can be fixed at their starting
values by setting their
tuning
values to zero.
The no nugget model is specified by removing tau.sq
from the starting
list.
Value
An object of class spSVC
, which is a list comprising:
coords |
the |
p.theta.samples |
a |
acceptance |
the Metropolis sampling
acceptance percent. Reported at |
The return object will include additional objects used for subsequent
parameter recovery, prediction, and model fit evaluation using
spRecover
, spPredict
,
spDiag
, respectively.
Author(s)
Andrew O. Finley finleya@msu.edu,
Sudipto Banerjee sudipto@ucla.edu
References
Finley, A.O., S. Banerjee, and A.E. Gelfand. (2015) spBayes for large univariate and multivariate point-referenced spatio-temporal data models. Journal of Statistical Software, 63:1–28. https://www.jstatsoft.org/article/view/v063i13.
Roberts G.O. and Rosenthal J.S. (2006). Examples of Adaptive MCMC. http://probability.ca/jeff/ftpdir/adaptex.pdf.
Finley, A.O. and S. Banerjee (2019) Bayesian spatially varying coefficient models in the spBayes R package. https://arxiv.org/abs/1903.03028.
See Also
Examples
## Not run:
library(Matrix)
rmvn <- function(n, mu=0, V = matrix(1)){
p <- length(mu)
if(any(is.na(match(dim(V),p))))
stop("Dimension problem!")
D <- chol(V)
t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p)))
}
##Assume both columns of X are space-varying and the two GPs don't covary
set.seed(1)
n <- 200
coords <- cbind(runif(n,0,1), runif(n,0,1))
X <- as.matrix(cbind(1, rnorm(n)))
colnames(X) <- c("x.1", "x.2")
Z <- t(bdiag(as.list(as.data.frame(t(X)))))
B <- as.matrix(c(1,5))
p <- length(B)
sigma.sq <- c(1,5)
tau.sq <- 1
phi <- 3/0.5
D <- as.matrix(dist(coords))
C <- exp(-phi*D)%x%diag(sigma.sq)
w <- rmvn(1, rep(0,p*n), C)
mu <- as.vector(X%*%B + Z%*%w)
y <- rnorm(n, mu, sqrt(tau.sq))
##fit a model to the simulated dat
starting <- list("phi"=rep(3/0.5, p), "sigma.sq"=rep(1, p), "tau.sq"=1)
tuning <- list("phi"=rep(0.1, p), "sigma.sq"=rep(0.1, p), "tau.sq"=0.1)
cov.model <- "exponential"
priors <- list("phi.Unif"=list(rep(3/2, p), rep(3/0.0001, p)),
"sigma.sq.IG"=list(rep(2, p), rep(2, p)),
"tau.sq.IG"=c(2, 1))
##fit model
n.samples <- 2000
m.1 <- spSVC(y~X-1, coords=coords, starting=starting, svc.cols=c(1,2),
tuning=tuning, priors=priors, cov.model=cov.model,
n.samples=n.samples, n.omp.threads=4)
plot(m.1$p.theta.samples, density=FALSE)
##recover posterior samples
m.1 <- spRecover(m.1, start=floor(0.75*n.samples), thin=2, n.omp.threads=4)
summary(m.1$p.beta.recover.samples)
summary(m.1$p.theta.recover.samples)
##check fitted values
quant <- function(x){quantile(x, prob=c(0.025, 0.5, 0.975))}
##fitted y
y.hat <- apply(m.1$p.y.samples, 1, quant)
rng <- c(-15, 20)
plot(y, y.hat[2,], pch=19, cex=0.5, xlab="Fitted y", ylab="Observed y",
xlim=rng, ylim=rng)
arrows(y, y.hat[2,], y, y.hat[1,], angle=90, length=0.05)
arrows(y, y.hat[2,], y, y.hat[3,], angle=90, length=0.05)
lines(rng, rng, col="blue")
##recovered w
w.hat <- apply(m.1$p.w.recover.samples, 1, quant)
w.1.indx <- seq(1, p*n, p)
w.2.indx <- seq(2, p*n, p)
par(mfrow=c(1,2))
rng <- c(-5,5)
plot(w[w.1.indx], w.hat[2,w.1.indx], pch=19, cex=0.5, xlab="Fitted w.1", ylab="Observed w.1",
xlim=rng, ylim=rng)
arrows(w[w.1.indx], w.hat[2,w.1.indx], w[w.1.indx], w.hat[1,w.1.indx], angle=90, length=0.05)
arrows(w[w.1.indx], w.hat[2,w.1.indx], w[w.1.indx], w.hat[3,w.1.indx], angle=90, length=0.05)
lines(rng, rng, col="blue")
rng <- c(-10,10)
plot(w[w.2.indx], w.hat[2,w.2.indx], pch=19, cex=0.5, xlab="Fitted w.2", ylab="Observed w.2",
xlim=rng, ylim=rng)
arrows(w[w.2.indx], w.hat[2,w.2.indx], w[w.2.indx], w.hat[1,w.2.indx], angle=90, length=0.05)
arrows(w[w.2.indx], w.hat[2,w.2.indx], w[w.2.indx], w.hat[3,w.2.indx], angle=90, length=0.05)
lines(rng, rng, col="blue")
## End(Not run)