| GeoFit2 | R Documentation |
Maximum weighted composite-likelihood fitting for Gaussian and some Non-Gaussian univariate spatial, spatio-temporal and bivariate spatial RFs. A first preliminary estimation is performed using independence composite-likelihood for the marginal parameters of the model. The estimates are then used as starting values in the second final estimation step. The function allows fixing any of the parameters and setting upper/lower bounds in the optimization.
GeoFit2(data, coordx, coordy=NULL, coordz=NULL, coordt=NULL,
coordx_dyn=NULL, copula=NULL, corrmodel, distance="Eucl",
fixed=NULL, anisopars=NULL, est.aniso=c(FALSE,FALSE),
grid=FALSE, likelihood="Marginal", lower=NULL, maxdist=Inf,
neighb=NULL, p_neighb=1, maxtime=Inf, memdist=TRUE,
method="cholesky", model="Gaussian", n=1, onlyvar=FALSE,
optimizer="Nelder-Mead", radius=1, score=FALSE,
sensitivity=FALSE, sparse=FALSE, start=NULL,
thin_method="iid", type="Pairwise", upper=NULL,
varest=FALSE, weighted=FALSE, X=NULL,
spobj=NULL, spdata=NULL)
data |
A |
coordx |
A numeric ( |
coordy |
A numeric vector giving 1-dimension of spatial coordinates; optional argument, default is |
coordz |
A numeric vector giving 1-dimension of spatial coordinates; optional argument, default is |
coordt |
A numeric vector assigning 1-dimension of temporal coordinates. Optional argument, default is |
coordx_dyn |
A list of |
copula |
String; the type of copula. It can be |
corrmodel |
String; the name of a correlation model, for the description see the Section Details. |
distance |
String; the name of the spatial distance. Default is |
fixed |
An optional named list giving the values of the parameters that will be considered as known values. The listed parameters for a given correlation function will not be estimated. |
anisopars |
A list of two elements: |
est.aniso |
A bivariate logical vector providing which anisotropic parameters must be estimated. |
grid |
Logical; if |
likelihood |
String; the configuration of the composite likelihood. |
lower |
An optional named list giving lower bounds for parameters when the optimizer is |
maxdist |
Numeric; an optional positive value indicating the maximum spatial distance considered in the composite likelihood computation. See the Section Details for more information. |
neighb |
Numeric; an optional positive integer indicating the order of neighborhood in the composite likelihood computation. See the Section Details for more information. |
p_neighb |
Numeric; a value in |
maxtime |
Numeric; an optional positive integer indicating the order of temporal neighborhood in the composite likelihood computation. |
memdist |
Logical; if |
method |
String; the type of matrix decomposition used in the likelihood computation. Default is |
model |
String; the type of RF and therefore the densities associated to the likelihood objects.
|
n |
Numeric; number of trials in a binomial RF; number of successes in a negative binomial RF. |
onlyvar |
Logical; if |
optimizer |
String; the optimization algorithm (see |
radius |
Numeric; the radius of the sphere in the case of lon-lat coordinates. Default value is |
score |
Logical; should score function be computed? Default is |
sensitivity |
Logical; if |
sparse |
Logical; if |
start |
An optional named list with initial values for parameters used by the numerical routines in the maximization procedure.
Default is |
thin_method |
String; thinning scheme used when |
type |
String; the type of the likelihood objects. If |
upper |
An optional named list giving upper bounds for parameters when the optimizer is |
varest |
Logical; if |
weighted |
Logical; if |
X |
Numeric; matrix of spatio(temporal) covariates in the linear mean specification. |
spobj |
An object of class |
spdata |
Character; the name of data in the |
The function GeoFit2 is similar to GeoFit.
However, GeoFit2 performs a preliminary estimation using maximum independence composite likelihood
for the marginal parameters of the model and then uses the obtained estimates as starting values in the final
weighted composite likelihood estimation (that includes both marginal and dependence parameters).
This provides robust starting values for the marginal parameters in the optimization algorithm.
Stochastic thinning of nearest-neighbor pairs can be enabled via p_neighb<1. The argument thin_method
controls the thinning scheme (default "iid").
Returns an object of class GeoFit.
An object of class GeoFit is a list containing at most the following components:
bivariate |
Logical: |
clic |
The composite information criterion; if the full likelihood is considered then it coincides with AIC. |
coordx |
A |
coordy |
A |
coordt |
A |
coordx_dyn |
A list of dynamical (in time) spatial coordinates. |
conf.int |
Confidence intervals for standard maximum likelihood estimation. |
convergence |
A string that denotes if convergence is reached. |
copula |
The type of copula. |
corrmodel |
The correlation model. |
data |
The vector or matrix or array (or list) of data. |
distance |
The type of spatial distance. |
fixed |
A list of fixed parameters. |
iterations |
The number of iterations used by the numerical routine. |
likelihood |
The configuration of the composite likelihood. |
logCompLik |
The value of the log composite-likelihood at the maximum. |
maxdist |
The maximum spatial distance used in the weighted composite likelihood (or |
maxtime |
The order of temporal neighborhood in the composite likelihood computation. |
message |
Extra message passed from the numerical routines. |
model |
The density associated to the likelihood objects. |
missp |
|
n |
The number of trials in a binomial RF; the number of successes in a negative binomial RF. |
neighb |
The order of spatial neighborhood in the composite likelihood computation. |
ns |
The number of (different) location sites in the bivariate case. |
numcoord |
The number of spatial coordinates. |
numtime |
The number of temporal realisations of the RF. |
param |
A list of parameter estimates. |
radius |
The radius of the sphere in the case of great-circle distance. |
stderr |
Standard errors for standard maximum likelihood estimation. |
sensmat |
The sensitivity matrix. |
varcov |
The variance-covariance matrix of the estimates. |
type |
The type of the likelihood objects. |
X |
The matrix of covariates. |
Moreno Bevilacqua, moreno.bevilacqua89@gmail.com, https://sites.google.com/view/moreno-bevilacqua/home, Víctor Morales Oñate, victor.morales@uv.cl, https://sites.google.com/site/moralesonatevictor/, Christian Caamaño-Carrillo, chcaaman@ubiobio.cl, https://www.researchgate.net/profile/Christian-Caamano
library(GeoModels)
###############################################################
############ Examples of spatial Gaussian RFs ################
###############################################################
################################################################
###
### Example 1 : Maximum pairwise conditional likelihood fitting
### of a Gaussian RF with Matern correlation
###
###############################################################
model="Gaussian"
# Define the spatial-coordinates of the points:
set.seed(3)
N=400 # number of location sites
x <- runif(N, 0, 1)
set.seed(6)
y <- runif(N, 0, 1)
coords <- cbind(x,y)
# Define spatial matrix covariates
X=cbind(rep(1,N),runif(N))
# Set the covariance model's parameters:
corrmodel <- "Matern"
mean <- 0.2
mean1 <- -0.5
sill <- 1
nugget <- 0
scale <- 0.2/3
smooth=0.5
param<-list(mean=mean,mean1=mean1,sill=sill,nugget=nugget,scale=scale,smooth=smooth)
# Simulation of the spatial Gaussian RF:
data <- GeoSim(coordx=coords,model=model,corrmodel=corrmodel, param=param,X=X)$data
fixed<-list(nugget=nugget,smooth=smooth)
start<-list(mean=mean,mean1=mean1,scale=scale,sill=sill)
################################################################
###
### Maximum pairwise likelihood fitting of
### Gaussian RFs with exponential correlation.
###
###############################################################
fit1 <- GeoFit2(data=data,coordx=coords,corrmodel=corrmodel,
optimizer="BFGS",neighb=3,likelihood="Conditional",
type="Pairwise", start=start,fixed=fixed,X=X)
print(fit1)
###############################################################
############ Examples of spatial non-Gaussian RFs #############
###############################################################
################################################################
###
### Example 2. Maximum pairwise likelihood fitting of
### a LogGaussian RF with Generalized Wendland correlation
###
###############################################################
set.seed(524)
# Define the spatial-coordinates of the points:
N=500
x <- runif(N, 0, 1)
y <- runif(N, 0, 1)
coords <- cbind(x,y)
X=cbind(rep(1,N),runif(N))
mean=1; mean1=2 # regression parameters
nugget=0
sill=0.5
scale=0.2
smooth=0
model="LogGaussian"
corrmodel="GenWend"
param=list(mean=mean,mean1=mean1,sill=sill,scale=scale,
nugget=nugget,power2=4,smooth=smooth)
# Simulation of a non stationary LogGaussian RF:
data <- GeoSim(coordx=coords, corrmodel=corrmodel,model=model,X=X,
param=param)$data
fixed<-list(nugget=nugget,power2=4,smooth=smooth)
start<-list(mean=mean,mean1=mean1,scale=scale,sill=sill)
I=Inf
lower<-list(mean=-I,mean1=-I,scale=0,sill=0)
upper<-list(mean= I,mean1= I,scale=I,sill=I)
# Maximum pairwise composite-likelihood fitting of the RF:
fit <- GeoFit2(data=data,coordx=coords,corrmodel=corrmodel, model=model,
neighb=3,likelihood="Conditional",type="Pairwise",X=X,
optimizer="nlminb",lower=lower,upper=upper,
start=start,fixed=fixed)
print(unlist(fit$param))
################################################################
###
### Example 3. Maximum pairwise likelihood fitting of
### SinhAsinh RFs with Wendland0 correlation
###
###############################################################
set.seed(261)
model="SinhAsinh"
# Define the spatial-coordinates of the points:
x <- runif(500, 0, 1)
y <- runif(500, 0, 1)
coords <- cbind(x,y)
corrmodel="Wend0"
mean=0;nugget=0
sill=1
skew=-0.5
tail=1.5
power2=4
c_supp=0.2
# model parameters
param=list(power2=power2,skew=skew,tail=tail,
mean=mean,sill=sill,scale=c_supp,nugget=nugget)
data <- GeoSim(coordx=coords, corrmodel=corrmodel,model=model, param=param)$data
plot(density(data))
fixed=list(power2=power2,nugget=nugget)
start=list(scale=c_supp,skew=skew,tail=tail,mean=mean,sill=sill)
# Maximum pairwise likelihood:
fit1 <- GeoFit2(data=data,coordx=coords,corrmodel=corrmodel, model=model,
neighb=3,likelihood="Marginal",type="Pairwise",
start=start,fixed=fixed)
print(unlist(fit1$param))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.