krigeTg function

TransGaussian kriging using Box-Cox transforms

TransGaussian kriging using Box-Cox transforms

TransGaussian (ordinary) kriging function using Box-Cox transforms methods

krigeTg(formula, locations, newdata, model = NULL, ..., nmax = Inf, nmin = 0, maxdist = Inf, block = numeric(0), nsim = 0, na.action = na.pass, debug.level = 1, lambda = 1.0)

Arguments

  • formula: formula that defines the dependent variable as a linear model of independent variables; suppose the dependent variable has name z, for ordinary and use a formula like z~1; the dependent variable should be NOT transformed.
  • locations: object of class Spatial, with observations
  • newdata: Spatial object with prediction/simulation locations; the coordinates should have names as defined in locations
  • model: variogram model of the TRANSFORMED dependent variable, see vgm , or fit.variogram
  • nmax: for local kriging: the number of nearest observations that should be used for a kriging prediction or simulation, where nearest is defined in terms of the space of the spatial locations. By default, all observations are used
  • nmin: for local kriging: if the number of nearest observations within distance maxdist is less than nmin, a missing value will be generated; see maxdist
  • maxdist: for local kriging: only observations within a distance of maxdist from the prediction location are used for prediction or simulation; if combined with nmax, both criteria apply
  • block: does not function correctly, afaik
  • nsim: does not function correctly, afaik
  • na.action: function determining what should be done with missing values in 'newdata'. The default is to predict 'NA'. Missing values in coordinates and predictors are both dealt with.
  • lambda: value for the Box-Cox transform
  • debug.level: debug level, passed to predict ; use -1 to see progress in percentage, and 0 to suppress all printed information
  • ...: other arguments that will be passed to gstat

Details

Function krigeTg uses transGaussian kriging as explained in https://www.math.umd.edu/~bnk/bak/Splus/kriging.html.

As it uses the R/gstat krige function to derive everything, it needs in addition to ordinary kriging on the transformed scale a simple kriging step to find m from the difference between the OK and SK prediction variance, and a kriging/BLUE estimation step to obtain the estimate of mumu.

For further details, see krige and predict .

Returns

an SpatialPointsDataFrame object containing the fields: m for the m (Lagrange) parameter for each location; var1SK.pred the c0Cinvc0 Cinv correction obtained by muhat for the mean estimate at each location; var1SK.var the simple kriging variance; var1.pred the OK prediction on the transformed scale; var1.var the OK kriging variance on the transformed scale; var1TG.pred the transGaussian kriging predictor; var1TG.var the transGaussian kriging variance, obtained by phi(muhat,lambda)2var1.varphi'(muhat, lambda)^2 * var1.var

References

N.A.C. Cressie, 1993, Statistics for Spatial Data, Wiley.

http://www.gstat.org/

Author(s)

Edzer Pebesma

See Also

gstat , predict

Examples

library(sp) data(meuse) coordinates(meuse) = ~x+y data(meuse.grid) gridded(meuse.grid) = ~x+y v = vgm(1, "Exp", 300) x1 = krigeTg(zinc~1,meuse,meuse.grid,v, lambda=1) # no transform x2 = krige(zinc~1,meuse,meuse.grid,v) summary(x2$var1.var-x1$var1TG.var) summary(x2$var1.pred-x1$var1TG.pred) lambda = -0.25 m = fit.variogram(variogram((zinc^lambda-1)/lambda ~ 1,meuse), vgm(1, "Exp", 300)) x = krigeTg(zinc~1,meuse,meuse.grid,m,lambda=-.25) spplot(x["var1TG.pred"], col.regions=bpy.colors()) summary(meuse$zinc) summary(x$var1TG.pred)