| control.georob {georob} | R Documentation |
Control Parameters for georob
Description
This page documents parameters used to control georob. It
describes the arguments of the functions control.georob,
param.transf, fwd.transf, dfwd.transf,
bwd.transf, control.rq, control.nleqslv,
control.nlminb and control.optim, which all serve to
control the behaviour of georob.
Usage
control.georob(ml.method = c("REML", "ML"), reparam = TRUE,
maximizer = c("nlminb", "optim"), initial.param = TRUE,
initial.fixef = c("lmrob", "rq", "lm"), bhat = NULL,
min.rweight = 0.25,
param.tf = param.transf(), fwd.tf = fwd.transf(),
deriv.fwd.tf = dfwd.transf(), bwd.tf = bwd.transf(),
psi.func = c("logistic", "t.dist", "huber"),
irwls.maxiter = 50,
irwls.ftol = 1.e-5, force.gradient = FALSE,
min.condnum = 1.e-12, zero.dist = sqrt(.Machine[["double.eps"]]),
error.family.estimation = c("gaussian", "long.tailed"),
error.family.cov.effects = c("gaussian", "long.tailed"),
error.family.cov.residuals = c("gaussian", "long.tailed"),
cov.bhat = TRUE, full.cov.bhat = FALSE, cov.betahat = TRUE,
cov.delta.bhat = TRUE, full.cov.delta.bhat = TRUE,
cov.delta.bhat.betahat = TRUE,
cov.ehat = TRUE, full.cov.ehat = FALSE,
cov.ehat.p.bhat = FALSE, full.cov.ehat.p.bhat = FALSE,
hessian = TRUE,
rq = control.rq(), lmrob = lmrob.control(),
nleqslv = control.nleqslv(),
optim = control.optim(), nlminb = control.nlminb(),
pcmp = control.pcmp(), ...)
param.transf(variance = "log", snugget = "log", nugget = "log", scale = "log",
alpha = c(
RMaskey = "log", RMdewijsian = "logit2", RMfbm = "logit2", RMgencauchy = "logit2",
RMgenfbm = "logit2", RMlgd = "identity", RMqexp = "logit1", RMstable = "logit2"
),
beta = c(RMdagum = "logit1", RMgencauchy = "log", RMlgd = "log"),
delta = "logit1", gamma = c(RMcauchy = "log", RMdagum = "logit1"),
kappa = "logit3", lambda = "log", mu = "log", nu = "log",
f1 = "log", f2 ="log", omega = "identity", phi = "identity", zeta = "identity")
fwd.transf(...)
dfwd.transf(...)
bwd.transf(...)
control.rq(tau = 0.5, rq.method = c("br", "fnb", "pfn"),
rq.alpha = 0.1, ci = FALSE, iid = TRUE,
interp = TRUE, tcrit = TRUE, rq.beta = 0.99995, eps = 1e-06,
Mm.factor = 0.8, max.bad.fixup = 3, ...)
control.nleqslv(method = c("Broyden", "Newton"),
global = c("dbldog", "pwldog", "qline", "gline", "none"),
xscalm = c("fixed", "auto"), control = list(ftol = 1e-04), ...)
control.optim(method = c("BFGS", "Nelder-Mead", "CG",
"L-BFGS-B", "SANN", "Brent"), lower = -Inf, upper = Inf,
control = list(reltol = 1e-05), ...)
control.nlminb(control = list(rel.tol = 1.e-5), lower = -Inf,
upper = Inf, ...)
Arguments
ml.method |
a character keyword defining whether Gaussian maximum
likelihood ( |
reparam |
a logical scalar. If |
maximizer |
a character keyword defining whether the Gaussian
(restricted) log-likelihood is maximized by |
initial.param |
a logical scalar, controlling whether initial values
of variogram parameters are computed for solving the robustified
estimating equations of the variogram and anisotropy parameters. If
|
initial.fixef |
a character keyword defining whether the function
|
bhat |
a numeric vector with initial values for the spatial random
effects |
min.rweight |
a positive numeric with the “robustness
weight” of the initial |
param.tf |
a function such as |
fwd.tf |
a function such as |
deriv.fwd.tf |
a function such as |
bwd.tf |
a function such as |
psi.func |
a character keyword defining what |
irwls.maxiter |
a positive integer equal to the maximum number of
IRWLS iterations to solve the estimating equations of
|
irwls.ftol |
a positive numeric with the convergence criterion for
IRWLS. Convergence is assumed if the objective function change of a IRWLS
iteration does not exceed |
force.gradient |
a logical scalar controlling whether the estimating
equations or the gradient of the Gaussian restricted log-likelihood are
evaluated even if all variogram parameters are fixed (default
|
min.condnum |
a positive numeric with the minimum acceptable
ratio of smallest to largest singular value of the model matrix
|
zero.dist |
a positive numeric equal to the maximum distance, separating two sampling locations that are still considered as being coincident. |
error.family.estimation |
a character keyword, defining the
probability distribution for |
error.family.cov.effects |
a character keyword, defining the
probability distribution for |
error.family.cov.residuals |
a character keyword, defining the
probability distribution for |
cov.bhat |
a logical scalar controlling whether the covariances of
|
full.cov.bhat |
a logical scalar controlling whether the full
covariance matrix ( |
cov.betahat |
a logical scalar controlling whether the covariance
matrix of |
cov.delta.bhat |
a logical scalar controlling whether the covariances of
|
full.cov.delta.bhat |
a logical scalar controlling whether the full covariance
matrix ( |
cov.delta.bhat.betahat |
a logical scalar controlling whether the covariance
matrix of |
cov.ehat |
a logical scalar controlling whether the covariances of
|
full.cov.ehat |
a logical scalar controlling whether the full covariance
matrix ( |
cov.ehat.p.bhat |
a logical scalar controlling whether the covariances of
|
full.cov.ehat.p.bhat |
a logical scalar controlling whether the full
covariance matrix ( |
hessian |
a logical scalar controlling whether for Gaussian (RE)ML the Hessian should be computed at the MLEs. |
rq |
a list of arguments passed to |
lmrob |
a list of arguments passed to the |
nleqslv |
a list of arguments passed to
|
nlminb |
a list of arguments passed to |
optim |
a list of arguments passed to |
pcmp |
a list of arguments, passed e.g. to |
... |
for |
variance, snugget, nugget, scale, alpha, beta, delta, gamma, kappa, lambda, mu, nu |
character strings with names of transformation functions of the variogram parameters. |
f1, f2, omega, phi, zeta |
character strings with names of transformation functions of the anisotropy variogram parameters. |
tau, rq.method, rq.alpha, ci, iid, interp, tcrit |
arguments passed
as |
rq.beta, eps, Mm.factor, max.bad.fixup |
arguments passed as
|
method, global, xscalm, control, lower, upper, reltol, rel.tol |
arguments passed to related arguments of
|
Details
Parameter transformations
The arguments param.tf, fwd.tf, deriv.fwd.tf,
bwd.tf define the transformations of the variogram parameters for
RE(ML) estimation. Implemented are currently "log",
"logit1", "logit2", "logit3" (various variants of
logit-transformation, see code of function fwd.transf) and "identity" (= no)
transformations. These are the possible values that the many arguments
of the function param.transf accept (as quoted character strings)
and these are the names of the list components returned by
fwd.transf, dfwd.transf and bwd.transf. Additional
transformations can be implemented by:
Extending the function definitions by arguments like
fwd.tf = fwd.transf(my.fun = function(x) your transformation),
deriv.fwd.tf = dfwd.transf(my.fun = function(x) your derivative),
bwd.tf = bwd.transf(my.fun = function(x) your back-transformation),Assigning to a given argument of
param.transfthe name of the new function, e.g.
variance = "my.fun".
Note that the values given for the arguments of param.transf
must match the names of the functions returned by fwd.transf,
dfwd.transf and bwd.transf.
Approximation of covariances of fixed and random effects and residuals
The robustified estimating equations of robust REML depend on the
covariances of \widehat{\boldsymbol{B}}.
These covariances (and the covariances of
\boldsymbol{B}-\widehat{\boldsymbol{B}},
\widehat{\boldsymbol{\beta}},
\widehat{\boldsymbol{\varepsilon}},
\widehat{\boldsymbol{\varepsilon}} +
\widehat{\boldsymbol{B}}) are
approximated by expressions that in turn depend on the variances of
\varepsilon,
\psi(\varepsilon/\tau) and the expectation
of \psi'(\varepsilon/\tau) (= \partial / \partial \varepsilon \,
\psi(\varepsilon/\tau)). The arguments
error.family.estimation, error.family.cov.effects and
error.family.cov.residuals control what parametric distribution
for \varepsilon is used to compute the variance of
\varepsilon,
\psi(\varepsilon/\tau) and the expectation
of \psi'(\varepsilon/\tau) when
solving the estimating equations (
error.family.estimation),computing the covariances of
\widehat{\boldsymbol{\beta}},\widehat{\boldsymbol{B}}and\boldsymbol{B}-\widehat{\boldsymbol{B}}(error.family.cov.effects) andcomputing the covariances of
\widehat{\boldsymbol{\varepsilon}}=\boldsymbol{Y} - \boldsymbol{X} \widehat{\boldsymbol{\beta}} - \widehat{\boldsymbol{B}}and\widehat{\boldsymbol{\varepsilon}} + \widehat{\boldsymbol{B}} =\boldsymbol{Y} - \boldsymbol{X} \widehat{\boldsymbol{\beta}}
(error.family.cov.residuals).
Possible options are: "gaussian" or "long.tailed". In
the latter case the probability density function of
\varepsilon is assumed to be proportional to
1/\tau \, \exp(-\rho_c(\varepsilon/\tau)), where
\psi_c(x)=\rho_c^\prime(x).
Value
control.georob, control.rq, control.nleqslv,
control.optim and control.nlminb all create lists with
control parameters passed to georob,
rq, nleqslv,
optim or nlminb, see arguments
above and the help pages of the respective functions for information
about the components of these lists. Note that the list returned by
control.georob contains some components (irwls.initial,
tuning.psi.nr, cov.bhat.betahat,
aux.cov.pred.target) that cannot be changed by the user.
param.transf generates a list with character strings that
define what transformations are used for estimating the variogram
parameters, and fwd.transf, bwd.transf and
dfwd.transf return lists of functions with forward and backward
transformations and the first derivatives of the forward
transformations, see section Parameter transformations above.
Author(s)
Andreas Papritz papritz@retired.ethz.ch.
See Also
georobPackage for a description of the model and a brief summary of the algorithms;
georob for (robust) fitting of spatial linear models;
georobObject for a description of the class georob;
profilelogLik for computing profiles of Gaussian likelihoods;
plot.georob for display of RE(ML) variogram estimates;
georobModelBuilding for stepwise building models of class georob;
cv.georob for assessing the goodness of a fit by georob;
georobMethods for further methods for the class georob;
predict.georob for computing robust Kriging predictions;
lgnpp for unbiased back-transformation of Kriging prediction
of log-transformed data;
georobSimulation for simulating realizations of a Gaussian process
from model fitted by georob; and finally
sample.variogram and fit.variogram.model
for robust estimation and modelling of sample variograms.
Examples
data(meuse)
r.logzn.rob <- georob(log(zinc) ~ sqrt(dist), data = meuse, locations = ~ x + y,
variogram.model = "RMexp",
param = c(variance = 0.15, nugget = 0.05, scale = 200),
tuning.psi = 1, control = control.georob(cov.bhat = TRUE,
cov.ehat.p.bhat = TRUE, initial.fixef = "rq"), verbose = 2)
qqnorm(rstandard(r.logzn.rob, level = 0)); abline(0, 1)
qqnorm(ranef(r.logzn.rob, standard = TRUE)); abline(0, 1)