krige.bayes | R Documentation |
The function krige.bayes
performs Bayesian analysis of
geostatistical data allowing specifications of
different levels of uncertainty in the model parameters.
It returns results on the posterior distributions for the model
parameters and on the predictive distributions for prediction
locations (if provided).
krige.bayes(geodata, coords = geodata$coords, data = geodata$data,
locations = "no", borders, model, prior, output)
model.control(trend.d = "cte", trend.l = "cte", cov.model = "matern",
kappa = 0.5, aniso.pars = NULL, lambda = 1)
prior.control(beta.prior = c("flat", "normal", "fixed"),
beta = NULL, beta.var.std = NULL,
sigmasq.prior = c("reciprocal", "uniform",
"sc.inv.chisq", "fixed"),
sigmasq = NULL, df.sigmasq = NULL,
phi.prior = c("uniform", "exponential","fixed",
"squared.reciprocal", "reciprocal"),
phi = NULL, phi.discrete = NULL,
tausq.rel.prior = c("fixed", "uniform", "reciprocal"),
tausq.rel, tausq.rel.discrete = NULL)
post2prior(obj)
geodata |
a list containing elements |
coords |
an |
data |
a vector with n data values. By default it takes the
component |
locations |
an |
borders |
optional. If missing,
by default reads the element |
model |
a list defining the fixed components of the model.
It can take an output to a call to |
prior |
a list with the specification of priors for the model
parameters.
It can take an output to a call to |
output |
a list specifying output options.
It can take an output to a call to |
trend.d |
specifies the trend (covariates) values at the data
locations. See documentation
of |
trend.l |
specifies the trend (covariates) at the prediction
locations. Must be of the same type as defined for |
cov.model |
string indicating the name of the model for the correlation function. Further details in the
documentation for |
kappa |
additional smoothness parameter. Only used if the
correlation function is one of: |
aniso.pars |
fixed parameters for geometric anisotropy
correction. If |
lambda |
numerical value of the Box-Cox transformation parameter.
The value |
beta.prior |
prior distribution for the mean (vector)
parameter |
beta |
mean hyperparameter for the distribution of the mean (vector) parameter |
beta.var.std |
standardised (co)variance hyperparameter(s)
for the prior for the mean
(vector) parameter |
sigmasq.prior |
specifies the prior for the parameter
|
sigmasq |
fixed value of the sill parameter
|
df.sigmasq |
numerical. Number of degrees of freedom for the
prior for the parameter |
phi.prior |
prior distribution for the range parameter
|
phi |
fixed value of the range parameter |
phi.discrete |
support points of the discrete prior
for the range parameter |
tausq.rel.prior |
specifies a prior distribution for the
relative nugget parameter
|
tausq.rel |
fixed value for the relative nugget parameter.
Only used if
|
tausq.rel.discrete |
support points of the discrete prior
for the relative nugget parameter |
obj |
an object of the class |
krige.bayes
is a generic function for Bayesian geostatistical
analysis of (transformed) Gaussian where predictions take into account the parameter
uncertainty.
It can be set to run conventional kriging methods which
use known parameters or plug-in estimates. However, the
functions krige.conv
and ksline
are preferable for
prediction with fixed parameters.
PRIOR SPECIFICATION
The basis of the Bayesian algorithm is the discretisation of the prior
distribution for the parameters \phi
and \tau^2_{rel}
=\frac{\tau^2}{\sigma^2}
.
The Tech. Report (see References
below)
provides details on the results used in the current implementation.
The expressions of the implemented priors for the parameter \phi
are:
p(\phi) \propto 1
.
p(\phi) = \frac{1}{\nu}
\exp(-\frac{1}{\nu} * \phi)
.
p(\phi) \propto \frac{1}{\phi}
.
p(\phi) \propto
\frac{1}{\phi^2}
.
fixed known or estimated value of \phi
.
The expressions of the implemented priors for the parameter \tau^2_{rel}
are:
fixed known or estimated value of
\tau^2_{rel}
. Defaults to zero.
p(\tau^2_{rel}) \propto 1
.
p(\tau^2_{rel}) \propto \frac{1}{\tau^2_{rel}}
.
Apart from those a user defined prior can be specifyed by
entering a vector of probabilities for a discrete distribution
with suport points given by the argument phi.discrete
and/or
tausq.rel.discrete
.
CONTROL FUNCTIONS
The function call includes auxiliary control functions which allows
the user to specify and/or change the specification of model
components
(using model.control
), prior
distributions (using prior.control
) and
output options (using output.control
).
Default options are available in most of the cases.
An object with class
"krige.bayes"
and
"kriging"
.
The attribute prediction.locations
containing the name of the
object with the coordinates of the prediction locations (argument
locations
) is assigned to the object.
Returns a list with the following components:
posterior |
results on on the posterior distribution of the
model parameters. A list with the following possible components: |
summary information on the posterior distribution
of the mean parameter \beta
.
summary information on the posterior distribution
of the variance parameter \sigma^2
(partial sill).
summary information on the posterior distribution
of the correlation parameter \phi
(range parameter) .
summary information on the posterior distribution
of the relative nugget variance parameter
\tau^2_{rel}
.
information on discrete the joint distribution of these parameters.
a data.frame with a sample from the posterior distribution. Each column corresponds to one of the basic model parameters.
predictive |
results on the predictive distribution at the prediction locations, if provided. A list with the following possible components: |
expected values.
expected variance.
type of posterior distribution.
mean of the simulations at each locations.
variance of the simulations at each locations.
quantiles computed from the the simulations.
probabilities computed from the simulations.
simulations from the predictive distribution.
prior |
a list with information on the prior distribution
and hyper-parameters of the model parameters ( |
model |
model specification as defined by |
.Random.seed |
system random seed before running the function.
Allows reproduction of results. If
the |
max.dist |
maximum distance found between two data locations. |
call |
the function call. |
Paulo J. Ribeiro Jr. paulojus@leg.ufpr.br,
Peter J. Diggle p.diggle@lancaster.ac.uk.
Diggle, P.J. & Ribeiro Jr, P.J. (2002) Bayesian inference in Gaussian model-based geostatistics. Geographical and Environmental Modelling, Vol. 6, No. 2, 129-146.
The technical details about the implementation of krige.bayes
can be
found at:
Ribeiro, P.J. Jr. and Diggle, P.J. (1999) Bayesian inference in
Gaussian model-based geostatistics. Tech. Report ST-99-08, Dept
Maths and Stats, Lancaster University.
Available at:
http://www.leg.ufpr.br/geoR/geoRdoc/bayeskrige.pdf
Further information about geoR can be found at:
http://www.leg.ufpr.br/geoR/.
For a extended list of examples of the usage see http://www.leg.ufpr.br/geoR/tutorials/examples.krige.bayes.R and/or the geoR tutorials page at http://www.leg.ufpr.br/geoR/tutorials/.
lines.variomodel.krige.bayes
,
plot.krige.bayes
for outputs related to the
parameters in the model,
image.krige.bayes
and
persp.krige.bayes
for graphical output of
prediction results.
krige.conv
and
ksline
for conventional kriging methods.
## Not run:
# generating a simulated data-set
ex.data <- grf(70, cov.pars=c(10, .15), cov.model="matern", kappa=2)
#
# defining the grid of prediction locations:
ex.grid <- as.matrix(expand.grid(seq(0,1,l=21), seq(0,1,l=21)))
#
# computing posterior and predictive distributions
# (warning: the next command can be time demanding)
ex.bayes <- krige.bayes(ex.data, loc=ex.grid,
model = model.control(cov.m="matern", kappa=2),
prior = prior.control(phi.discrete=seq(0, 0.7, l=51),
phi.prior="reciprocal"))
#
# Prior and posterior for the parameter phi
plot(ex.bayes, type="h", tausq.rel = FALSE, col=c("red", "blue"))
#
# Plot histograms with samples from the posterior
par(mfrow=c(3,1))
hist(ex.bayes)
par(mfrow=c(1,1))
# Plotting empirical variograms and some Bayesian estimates:
# Empirical variogram
plot(variog(ex.data, max.dist = 1), ylim=c(0, 15))
# Since ex.data is a simulated data we can plot the line with the "true" model
lines.variomodel(ex.data, lwd=2)
# adding lines with summaries of the posterior of the binned variogram
lines(ex.bayes, summ = mean, lwd=1, lty=2)
lines(ex.bayes, summ = median, lwd=2, lty=2)
# adding line with summary of the posterior of the parameters
lines(ex.bayes, summary = "mode", post = "parameters")
# Plotting again the empirical variogram
plot(variog(ex.data, max.dist=1), ylim=c(0, 15))
# and adding lines with median and quantiles estimates
my.summary <- function(x){quantile(x, prob = c(0.05, 0.5, 0.95))}
lines(ex.bayes, summ = my.summary, ty="l", lty=c(2,1,2), col=1)
# Plotting some prediction results
op <- par(no.readonly = TRUE)
par(mfrow=c(2,2), mar=c(4,4,2.5,0.5), mgp = c(2,1,0))
image(ex.bayes, main="predicted values")
image(ex.bayes, val="variance", main="prediction variance")
image(ex.bayes, val= "simulation", number.col=1,
main="a simulation from the \npredictive distribution")
image(ex.bayes, val= "simulation", number.col=2,
main="another simulation from \nthe predictive distribution")
#
par(op)
## End(Not run)
##
## For a extended list of exemples of the usage of krige.bayes()
## see http://www.leg.ufpr.br/geoR/tutorials/examples.krige.bayes.R
##
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.