Spatial (bivariate) and spatio temporal optimal linear local prediction for Gaussian and non Gaussian RFs.
Spatial (bivariate) and spatio temporal optimal linear local prediction for Gaussian and non Gaussian RFs.
For a given set of spatial location sites (and temporal instants), the function computes optmal local linear prediction and the associated mean squared error for the Gaussian and non Gaussian case using a spatial (temporal) neighborhood computed using the function GeoNeighborhood
UTF-8
estobj: An object of class Geofit that includes information about data, model and estimates.
data: A d-dimensional vector (a single spatial realisation) or a (dxd)-matrix (a single spatial realisation on regular grid) or a (txd)-matrix (a single spatial-temporal realisation) or an (dxdxt)-array (a single spatial-temporal realisation on regular grid) giving the data used for prediction.
coordx: A numeric (dx2)-matrix or (dx3)-matrix Coordinates on a sphere for a fixed radius radius
are passed in lon/lat format expressed in decimal degrees.
coordy: A numeric vector giving 1-dimension of spatial coordinates; Optional argument, the default is NULL.
coordz: A numeric vector giving 1-dimension of spatial coordinates; Optional argument, the default is NULL.
coordt: A numeric vector giving 1-dimension of temporal coordinates used for prediction; the default is NULL
then a spatial random field is expected.
coordx_dyn: A list of m numeric (dx2)-matrices containing dynamical (in time) spatial coordinates. Optional argument, the default is NULL
corrmodel: String; the name of a correlation model, for the description see the Section Details .
distance: String; the name of the spatial distance. The default is Eucl, the euclidean distance. See the Section Details of GeoFit.
grid: Logical; if FALSE (the default) the data used for prediction are interpreted as spatial or spatial-temporal realisations on a set of non-equispaced spatial sites (irregular grid).
loc: A numeric (nx2)-matrix (where n is the number of spatial sites) giving 2-dimensions of spatial coordinates to be predicted.
neighb: Numeric; an optional positive integer indicating the order of the neighborhood.
maxdist: Numeric; an optional positive value indicating the distance in the spatial neighborhood.
maxtime: Numeric; an optional positive integer value indicating the order of the temporal neighborhood.
method: String; the type of matrix decomposition used in the simulation. Default is cholesky. The other possible choices is svd.
n: Numeric; the number of trials in a binomial random fields. Default is 1.
nloc: Numeric; the number of trials of the locations sites to be predicted in the binomial random field. If missing then a rounded mean of n is considered.
mse: Logical; if TRUE (the default) MSE of the kriging predictor is computed
model: String; the type of RF and therefore the densities associated to the likelihood objects. Gaussian is the default, see the Section Details .
param: A list of parameter values required for the correlation model.See the Section Details .
anisopars: A list of two elements: "angle" and "ratio" i.e. the anisotropy angle and the anisotropy ratio, respectively.
radius: Numeric: the radius of the sphere if coordinates are passed in lon/lat format;
sparse: Logical; if TRUE kriging is computed with sparse matrices algorithms using spam package. Default is FALSE. It should be used with compactly supported covariances.
time: A numeric (mx1) vector (where m is the number of temporal instants) giving the temporal instants to be predicted; the default is NULL
then only spatial prediction is performed.
type: String; if Standard then standard kriging is performed;if Tapering
then kriging with covariance tapering is performed;if Pairwise then pairwise kriging is performed
type_mse: String; if Theoretical then theoretical MSE pairwise kriging is computed. If SubSamp then an estimation based on subsampling is computed.
type_krig: String; the type of kriging. If Simple (the default) then simple kriging is performed. If Optim then optimal kriging is performed for some non-Gaussian RFs
weigthed: Logical; if TRUE then decreasing weigths coming from a compactly supported correlation function with compact support maxdist (maxtime)are used in the pairwise kriging.
which: Numeric; In the case of bivariate (tapered) cokriging it indicates which variable to predict. It can be 1 or 2
copula: String; the type of copula. It can be "Clayton" or "Gaussian"
X: Numeric; Matrix of spatio(temporal)covariates in the linear mean specification.
Xloc: Numeric; Matrix of spatio(temporal)covariates in the linear mean specification associated to predicted locations.
Mloc: Numeric; Vector of spatio(temporal) estimated means associated to predicted locations.
spobj: An object of class sp or spacetime
spdata: Character:The name of data in the sp or spacetime object
parallel: Logical; if TRUE then the prediction computation is parallelized
ncores: Numeric; number of cores involved in parallelization.
Details
This function use the GeoKrig with a spatial or spatio-temporal neighborhood computed using the function GeoNeighborhood. The neighborhood is specified with the option maxdist and maxtime.
Returns
Returns an object of class Kg. An object of class Kg is a list containing at most the following components: - bivariate: TRUE if spatial bivariate cokriging is performed, otherwise FALSE;
coordx: A d-dimensional vector of spatial coordinates used for prediction;
coordy: A d-dimensional vector of spatial coordinates used for prediction;
coordz: A d-dimensional vector of spatial coordinates used for prediction;
coordt: A t-dimensional vector of temporal coordinates used for prediction;
corrmodel: String: the correlation model;
covmatrix: The covariance matrix if type is Standard. An object of class spam if type is Tapering
data: The vector or matrix or array of data used for prediction
distance: String: the type of spatial distance;
grid: TRUE if the spatial data used for prediction are observed in a regular grid, otherwise FALSE;
loc: A (nx2)-matrix of spatial locations to be predicted.
n: The number of trial for Binomial RFs
nozero: In the case of tapered simple kriging the percentage of non zero values in the covariance matrix. Otherwise is NULL.
numcoord: Numeric:he number d of spatial coordinates used for prediction;
numloc: Numeric: the number n of spatial coordinates to be predicted;
numtime: Numeric: the number d of the temporal instants used for prediction;
numt: Numeric: the number m of the temporal instants to be predicted;
model: The type of RF, see GeoFit.
param: Numeric: The covariance parameters;
pred: A (nxm)-matrix of spatio or spatio temporal kriging prediction;
radius: Numeric: the radius of the sphere if coordinates are pssed in lon/lat format;
spacetime: TRUE if spatio-temporal kriging and FALSE if spatial kriging;
tapmod: String: the taper model if type is Tapering. Otherwise is NULL.
time: A m-dimensional vector of temporal coordinates to be predicted;
type: String: the type of kriging (Standard or Tapering).
type_krig: String: the type of kriging.
mse: A (nxm)-matrix of spatio or spatio temporal mean square error kriging prediction;
References
Gaetan, C. and Guyon, X. (2010) Spatial Statistics and Modelling. Spring Verlang, New York. Furrer R., Genton, M.G. and Nychka D. (2006). Covariance Tapering for Interpolation of Large Spatial Datasets. Journal of Computational and Graphical Statistics, 15-3 , 502--523.
############################################################################### Examples of Spatial local kriging #############################################################################require(GeoModels)####model="Gaussian"# Define the spatial-coordinates of the points:set.seed(759)x = runif(1000,0,1)y = runif(1000,0,1)coords=cbind(x,y)# Set the exponential cov parameters:corrmodel ="GenWend"mean=0; sill=1nugget=0; scale=0.2param=list(mean=mean,sill=sill,nugget=nugget,smooth=0,scale=scale,power2=4)# Simulation of the spatial Gaussian random field:data = GeoSim(coordx=coords, corrmodel=corrmodel, param=param)$data
# Maximum pairwise likelihood fitting of the space time random field:start=list(scale=scale,sill=sill,mean=mean)fixed=list(power2=4,smooth=0,nugget=0)fit = GeoFit(data, coordx=coords, corrmodel=corrmodel, start=start,fixed=fixed, likelihood='Conditional', type='Pairwise', neighb=3)# locations to predictloc_to_pred=matrix(runif(8),4,2)###################################################################### Example 1. Comparing spatial kriging with local kriging for### a Gaussian random field with GenWend correlation.### ###############################################################param=append(fit$param,fit$fixed)pr=GeoKrig(fit,loc=loc_to_pred,mse=TRUE)pr_loc=GeoKrigloc(fit,loc=loc_to_pred,neighb=100,mse=TRUE)pr$pred;pr_loc$pred
################################################################ Example: spatio temporal Gaussian local kriging ##################################################################require(GeoModels)require(fields)set.seed(78)coords=cbind(runif(100),runif(100))coordt=seq(0,5,0.25)corrmodel="Matern_Matern"param=list(nugget=0,mean=0,scale_s=0.2/3,scale_t=0.25/3,sill=2, smooth_s=0.5,smooth_t=0.5)data = GeoSim(coordx=coords, coordt=coordt, corrmodel=corrmodel, param=param)$data
# Maximum pairwise likelihood fitting of the space time random field:start = list(scale_s=0.2/3,scale_t=0.25,sill=2,mean=0)fixed = list(smooth_s=0.5,smooth_t=0.5,nugget=0)I=Inflower=list(scale_s=0,scale_t=0,sill=0,mean=-I)upper=list(scale_s=I,scale_t=I,sill=I,mean=I)fit = GeoFit(data, coordx=coords, coordt=coordt, model=model, corrmodel=corrmodel, likelihood='Conditional', type='Pairwise',start=start,fixed=fixed, optimizer="nlminb",lower=lower,upper=upper, neighb=3,maxtime=1)## four location to predictloc_to_pred=matrix(runif(8),4,2)## three temporal instants to predicttime=c(0.5,1.5,3.5)pr=GeoKrig(fit,loc=loc_to_pred,time=time,mse=TRUE)pr_loc=GeoKrigloc(fit,loc=loc_to_pred,time=time, neigh=25,maxtime=1, mse=TRUE)## full and local prediction pr$pred
pr_loc$pred
################################################################ Example: spatio bivariate Gaussian local cokriging ###################################################################set.seed(6)#NN=1500 # number of spatial locations #x = runif(NN, 0, 1); #y = runif(NN, 0, 1) #coords=cbind(x,y) ## setting parameters#mean_1 = 2; mean_2= -1#nugget_1 =0;nugget_2=0#sill_1 =0.5; sill_2 =1; ### correlation parameters#CorrParam("Bi_Matern")#scale_1=0.2/3; scale_2=0.15/3; scale_12=0.5*(scale_2+scale_1) #smooth_1=smooth_2=smooth_12=0.5#pcol = -0.4#param= list(nugget_1=nugget_1,nugget_2=nugget_2,# sill_1=sill_1,sill_2=sill_2,# mean_1=mean_1,mean_2=mean_2,# smooth_1=smooth_1, smooth_2=smooth_2,smooth_12=smooth_12,# scale_1=scale_1, scale_2=scale_2,scale_12=scale_12,# pcol=pcol)## simulation#data = GeoSim(coordx=coords, corrmodel="Bi_Matern",model=model,param=param)$data#fixed=list(mean_1=mean_1,mean_2=mean_2, nugget_1=nugget_1,nugget_2=nugget_2, # smooth_1=smooth_1, smooth_2=smooth_2,smooth_12=smooth_12)#start=list( sill_1=sill_1,sill_2=sill_2,# scale_1=scale_1,scale_2=scale_2,scale_12=scale_12, pcol=pcol)## estimation with maximum likelihood #fit = GeoFit(data=data,coordx=coords, corrmodel="Bi_Matern",# likelihood="Marginal",type="Pairwise",optimizer="BFGS",neighb=5,#start=start,fixed=fixed)###### co-kriging for the fist component ###############xx=seq(0,1,0.022)#loc_to_pred=as.matrix(expand.grid(xx,xx))#pr1 = GeoKrigloc(fit,which=1,mse=TRUE,loc=loc_to_pred,neighb=100)#opar=par(no.readonly = TRUE)#par(mfrow=c(1,2))#zlim=c(-2.5,2.5)#colour = rainbow(100)#quilt.plot(coords,data[1,] ,col=colour,main = paste(" Fist component")) #quilt.plot(loc_to_pred,pr1$pred,col=colour,# main = paste(" Kriging first component"),ylab="")#par(opar)