Title: | Covariate-Based Covariance Functions for Nonstationary Spatial Modeling |
---|---|
Description: | Estimation and prediction of nonstationary Gaussian process with modular covariate-based covariance functions. An induced compact-supported nonstationary covariance function is provided to speed up computations when handling densly sampled domains. |
Authors: | Federico Blasi [aut, cre] , Reinhard Furrer [ctb] |
Maintainer: | Federico Blasi <[email protected]> |
License: | GPL (>= 3) |
Version: | 0.1.3 |
Built: | 2024-10-10 09:25:40 UTC |
Source: | https://github.com/blasif/cocons |
Provides routines and methods for estimating and predicting nonstationary Gaussian process models with modular covariate-based covariance functions. Several sources of nonstationarity can be modeled based on spatial information, including a trend, marginal standard deviation, local geometric anisotropy, local nugget, and spatially varying smoothness. Each of these components is modeled separately. An induced compact-supported nonstationary covariance function is provided to speed up computations when handling densly sampled domains. Models are estimated via maximum likelihood (and flavours of it, such as penalized and profile maximum likelihood). A variety of functions are also included to compute prediction metrics and to visualize, simulate, and summarize these types of models. Details of the models can be found in the vignette and in help(coco).
This package is provided "as is" without warranty of any kind, either express or implied. Backwards compatibility will not be offered until later versions.
Federico Blasi [aut, cre], [email protected]
## Not run: vignette("cocons", package = "cocons") methods(class = "coco") ## End(Not run)
## Not run: vignette("cocons", package = "cocons") methods(class = "coco") ## End(Not run)
Creates an S4 object of class coco
, which is the centerpiece of the cocons package. The function provides a set of consistency checks for ensuring the suitability of the different objects involved.
coco(type, data, locs, z, model.list, info, output = list())
coco(type, data, locs, z, model.list, info, output = list())
type |
( |
data |
( |
locs |
( |
z |
( |
model.list |
( |
info |
( |
output |
( |
Two types of coco
objects are available, each assuming a different type of covariance matrix for the Gaussian process.
Type "dense"
builds dense covariance matrices (non zero elements), while type "sparse"
builds sparse covariance
matrices by tapering the dense covariance matrix with a compact isotropic compact-supported correlation matrix [1].
Type "sparse"
allows a set of efficient algorithms, thus making it more suitable for large sample sizes.
An important component of the coco
S4 class is the model.list
specification, involving individual formulas provided as a list, where each of them specifies a covariate-based parametric model for a specific source of nonstationarity.
It involves "trend"
for the spatial trend, the "std.dev"
for the marginal standard deviation,
"scale"
, "aniso"
and "tilt"
, each of them shaping specific aspects of the local spatial geometrically anisotropy structure,
"smooth"
handling local smoothness, and "nugget"
handling the local nugget effect. The models are defined as:
Source | Related to | Description | Model |
mean | |
Mean function | |
std.dev | |
Marginal standard deviation | |
scale | |
Local scale | |
aniso | |
Local geometric anisotropy | |
tilt | |
(Restricted) local tilt | |
smooth | |
Local smoothness | |
nugget | |
Local micro-scale variability | |
where ,
,
,
,
,
, and
are the parameter vectors of each model,
, and
are the lower and upper bounds limiting the range of variation of the spatially-varying smoothness, and where
relates to a specific design matrix defined by the specific models for each of the source of nonstationarity.
Lastly, arguments for the "info"
list argument involve:
"lambda"
: (numeric
) a positive scalar specifying the regularization parameter.
"smooth.limits"
: (numeric vector
) specifying the allowed range of variation for the spatially varying smoothness.
"taper"
: (numeric
) specifying the desired taper function from the spam package (only for "sparse" coco objects).
"delta"
: (numeric
) specifying the taper range/scale (only for "sparse" coco objects).
"skip.scale"
: (integer vector
) By default, all covariates are scaled. skip.scale
allows to specify the index of those variables in data
that should not be scaled during the optimization.
(S4
) An S4 object of class coco
.
Federico Blasi
[1] Furrer, Reinhard, Marc G. Genton, and Douglas Nychka. "Covariance tapering for interpolation of large spatial datasets." Journal of Computational and Graphical Statistics 15.3 (2006): 502-523.
## Not run: locs <- expand.grid(seq(0,1,length.out = 10), seq(0,1,length.out = 10)) toydata <- data.frame('x' = locs[,1]) set.seed(1) z <- rnorm(100) model.list <- list('mean' = 0, 'std.dev' = formula( ~ 1), 'scale' = formula( ~ 1 + x), 'aniso' = 0, 'tilt' = 0, 'smooth' = 3/2, 'nugget' = -Inf) coco_object <- coco(type = 'dense', data = toydata, locs = as.matrix(locs), z = z, model.list = model.list) coco_object ## End(Not run)
## Not run: locs <- expand.grid(seq(0,1,length.out = 10), seq(0,1,length.out = 10)) toydata <- data.frame('x' = locs[,1]) set.seed(1) z <- rnorm(100) model.list <- list('mean' = 0, 'std.dev' = formula( ~ 1), 'scale' = formula( ~ 1 + x), 'aniso' = 0, 'tilt' = 0, 'smooth' = 3/2, 'nugget' = -Inf) coco_object <- coco(type = 'dense', data = toydata, locs = as.matrix(locs), z = z, model.list = model.list) coco_object ## End(Not run)
An S4 class to store information
type
(character
) One of two available types "dense" or "sparse". See description.
data
(data.frame
) A data.frame with covariates information, where colnames(data) matches model.list specification
locs
(numeric matrix
) a matrix with locs matching data
z
(numeric matrix
) A matrix of dimension n x p with response values
model.list
(list
) A list specifying a model for each aspect of the spatial structure.
info
(list
) a list with information about the coco object
output
(list
) if building an already fitted coco
object (not the standard approach), then requires an output from Optimparallel output, including as well boundaries, etc.
Federico Blasi
This function estimates the spatial model parameters using the L-BFGS-B optimizer [1].
cocoOptim(coco.object, boundaries = list(), ncores = "auto", optim.type, safe, optim.control)
cocoOptim(coco.object, boundaries = list(), ncores = "auto", optim.type, safe, optim.control)
coco.object |
( |
boundaries |
( |
ncores |
( |
optim.type |
(
|
safe |
( |
optim.control |
( |
(S4
) An optimized S4 object of class coco
.
Federico Blasi
[1] Byrd, Richard H., et al. "A limited memory algorithm for bound constrained optimization." SIAM Journal on scientific computing 16.5 (1995): 1190-1208.
[2] Gerber, Florian, and Reinhard Furrer. "optimParallel: An R package providing a parallel version of the L-BFGS-B optimization method." R Journal 11.1 (2019): 352-358.
## Not run: model.list <- list('mean' = 0, 'std.dev' = formula( ~ 1 + cov_x + cov_y), 'scale' = formula( ~ 1 + cov_x + cov_y), 'aniso' = 0, 'tilt' = 0, 'smooth' = 3/2, 'nugget' = -Inf) coco_object <- coco(type = 'dense', data = holes[[1]][1:100,], locs = as.matrix(holes[[1]][1:100,1:2]), z = holes[[1]][1:100,]$z, model.list = model.list) optim_coco <- cocoOptim(coco_object, boundaries = getBoundaries(coco_object, lower.value = -3, 3)) plotOptimInfo(optim_coco) plot(optim_coco) plot(optim_coco, type = 'ellipse') plot(optim_coco, type = 'correlations', index = c(2,3,5)) summary(optim_coco) getEstims(optim_coco) ## End(Not run)
## Not run: model.list <- list('mean' = 0, 'std.dev' = formula( ~ 1 + cov_x + cov_y), 'scale' = formula( ~ 1 + cov_x + cov_y), 'aniso' = 0, 'tilt' = 0, 'smooth' = 3/2, 'nugget' = -Inf) coco_object <- coco(type = 'dense', data = holes[[1]][1:100,], locs = as.matrix(holes[[1]][1:100,1:2]), z = holes[[1]][1:100,]$z, model.list = model.list) optim_coco <- cocoOptim(coco_object, boundaries = getBoundaries(coco_object, lower.value = -3, 3)) plotOptimInfo(optim_coco) plot(optim_coco) plot(optim_coco, type = 'ellipse') plot(optim_coco, type = 'correlations', index = c(2,3,5)) summary(optim_coco) getEstims(optim_coco) ## End(Not run)
Computes point predictions and standard errors based on conditional Gaussian distributions for nonstationary spatial models.
cocoPredict(coco.object, newdataset, newlocs, type = 'mean', ...)
cocoPredict(coco.object, newdataset, newlocs, type = 'mean', ...)
coco.object |
( |
newdataset |
( |
newlocs |
( |
type |
( |
... |
Additional arguments. If |
A list containing:
trend
: The systematic large-scale variability.
mean
: The stochastic mean.
sd.pred
: The standard errors, when type = 'pred'
is specified.
Federico Blasi
## Not run: # Stationary model model.list_stat <- list('mean' = 0, 'std.dev' = formula( ~ 1), 'scale' = formula( ~ 1), 'aniso' = 0, 'tilt' = 0, 'smooth' = 3/2, 'nugget' = -Inf) model.list_ns <- list('mean' = 0, 'std.dev' = formula( ~ 1 + cov_x + cov_y), 'scale' = formula( ~ 1 + cov_x + cov_y), 'aniso' = 0, 'tilt' = 0, 'smooth' = 3/2, 'nugget' = -Inf) coco_object <- coco(type = 'dense', data = holes[[1]][1:100, ], locs = as.matrix(holes[[1]][1:100, 1:2]), z = holes[[1]][1:100, ]$z, model.list = model.list_stat) optim_coco_stat <- cocoOptim(coco_object, boundaries = getBoundaries(coco_object, lower.value = -3, 3)) coco_preds_stat <- cocoPredict(optim_coco_stat, newdataset = holes[[2]], newlocs = as.matrix(holes[[2]][, 1:2]), type = "pred") # Update model [email protected] <- model.list_ns optim_coco_ns <- cocoOptim(coco_object, boundaries = getBoundaries(coco_object, lower.value = -3, 3)) coco_preds_ns <- cocoPredict(optim_coco_ns, newdataset = holes[[2]], newlocs = as.matrix(holes[[2]][, 1:2]), type = "pred") par(mfrow = c(1, 3)) fields::quilt.plot(main = "full data", holes[[1]][, 1:2], holes[[1]]$z, xlim = c(-1, 1), ylim = c(-1, 1)) fields::quilt.plot(main = "stationary se", holes[[2]][, 1:2], coco_preds_stat$sd.pred, xlim = c(-1, 1), ylim = c(-1, 1)) fields::quilt.plot(main = "nonstationary se", holes[[2]][, 1:2], coco_preds_ns$sd.pred, xlim = c(-1, 1), ylim = c(-1, 1)) ## End(Not run)
## Not run: # Stationary model model.list_stat <- list('mean' = 0, 'std.dev' = formula( ~ 1), 'scale' = formula( ~ 1), 'aniso' = 0, 'tilt' = 0, 'smooth' = 3/2, 'nugget' = -Inf) model.list_ns <- list('mean' = 0, 'std.dev' = formula( ~ 1 + cov_x + cov_y), 'scale' = formula( ~ 1 + cov_x + cov_y), 'aniso' = 0, 'tilt' = 0, 'smooth' = 3/2, 'nugget' = -Inf) coco_object <- coco(type = 'dense', data = holes[[1]][1:100, ], locs = as.matrix(holes[[1]][1:100, 1:2]), z = holes[[1]][1:100, ]$z, model.list = model.list_stat) optim_coco_stat <- cocoOptim(coco_object, boundaries = getBoundaries(coco_object, lower.value = -3, 3)) coco_preds_stat <- cocoPredict(optim_coco_stat, newdataset = holes[[2]], newlocs = as.matrix(holes[[2]][, 1:2]), type = "pred") # Update model coco_object@model.list <- model.list_ns optim_coco_ns <- cocoOptim(coco_object, boundaries = getBoundaries(coco_object, lower.value = -3, 3)) coco_preds_ns <- cocoPredict(optim_coco_ns, newdataset = holes[[2]], newlocs = as.matrix(holes[[2]][, 1:2]), type = "pred") par(mfrow = c(1, 3)) fields::quilt.plot(main = "full data", holes[[1]][, 1:2], holes[[1]]$z, xlim = c(-1, 1), ylim = c(-1, 1)) fields::quilt.plot(main = "stationary se", holes[[2]][, 1:2], coco_preds_stat$sd.pred, xlim = c(-1, 1), ylim = c(-1, 1)) fields::quilt.plot(main = "nonstationary se", holes[[2]][, 1:2], coco_preds_ns$sd.pred, xlim = c(-1, 1), ylim = c(-1, 1)) ## End(Not run)
draw realizations of stationary and nonstationary Gaussian processes with covariate-based covariance functions.
cocoSim(coco.object, pars, n, seed, standardize, type = 'classic', sim.type = NULL, cond.info = NULL)
cocoSim(coco.object, pars, n, seed, standardize, type = 'classic', sim.type = NULL, cond.info = NULL)
coco.object |
( |
pars |
( |
n |
( |
seed |
( |
standardize |
( |
type |
( |
sim.type |
( |
cond.info |
( |
#' The argument sim.type = 'cond'
specifies a conditional simulation, requiring cond.info
to be provided.
cond.info
is a list including newdataset
, a data.frame containing covariates present in model.list
at the simulation locations, and newlocs
,
a matrix specifying the locations corresponding to the simulation, with indexing that matches newdataset
.
The argument type = 'classic'
assumes a simplified parameterization for the covariance function, with log-parameterizations applied to the parameters std.dev
,
scale
, and smooth
.
(matrix
) a matrix dim(data)[1] x n.
Federico Blasi
## Not run: model.list <- list('mean' = 0, 'std.dev' = formula( ~ 1 + cov_x + cov_y), 'scale' = formula( ~ 1 + cov_x + cov_y), 'aniso' = 0, 'tilt' = 0, 'smooth' = 0.5, 'nugget' = -Inf) coco_object <- coco(type = 'dense', data = holes[[1]][1:1000,], locs = as.matrix(holes[[1]][1:1000,1:2]), z = holes[[1]][1:1000,]$z, model.list = model.list) coco_sim <- cocoSim(coco.object = coco_object, pars = c(0,0.25,0.25, # pars related to std.dev log(0.25),1,-1), # pars related to scale n = 1, standardize = TRUE) fields::quilt.plot(coco_object@locs,coco_sim) ## End(Not run)
## Not run: model.list <- list('mean' = 0, 'std.dev' = formula( ~ 1 + cov_x + cov_y), 'scale' = formula( ~ 1 + cov_x + cov_y), 'aniso' = 0, 'tilt' = 0, 'smooth' = 0.5, 'nugget' = -Inf) coco_object <- coco(type = 'dense', data = holes[[1]][1:1000,], locs = as.matrix(holes[[1]][1:1000,1:2]), z = holes[[1]][1:1000,]$z, model.list = model.list) coco_sim <- cocoSim(coco.object = coco_object, pars = c(0,0.25,0.25, # pars related to std.dev log(0.25),1,-1), # pars related to scale n = 1, standardize = TRUE) fields::quilt.plot(coco_object@locs,coco_sim) ## End(Not run)
Dense covariance function (difference parameterization)
cov_rns(theta, locs, x_covariates, smooth_limits)
cov_rns(theta, locs, x_covariates, smooth_limits)
theta |
vector of parameters |
locs |
a matrix with locations |
x_covariates |
design data.frame |
smooth_limits |
smooth limits |
dense covariance matrix
Dense covariance function (classic parameterization)
cov_rns_classic(theta, locs, x_covariates)
cov_rns_classic(theta, locs, x_covariates)
theta |
vector of parameters |
locs |
a matrix with locations |
x_covariates |
design data.frame |
dense covariance matrix with classic parameterization
Dense covariance function
cov_rns_pred( theta, locs, locs_pred, x_covariates, x_covariates_pred, smooth_limits )
cov_rns_pred( theta, locs, locs_pred, x_covariates, x_covariates_pred, smooth_limits )
theta |
vector of parameters |
locs |
a matrix with locations |
locs_pred |
a matrix with prediction locations |
x_covariates |
design data.frame |
x_covariates_pred |
design data.frame at prediction locations |
smooth_limits |
smooth limits |
dense covariance matrix
Sparse covariance function
cov_rns_taper( theta, locs, x_covariates, colindices, rowpointers, smooth_limits )
cov_rns_taper( theta, locs, x_covariates, colindices, rowpointers, smooth_limits )
theta |
vector of parameters |
locs |
a matrix with locations |
x_covariates |
design data.frame |
colindices |
from spam object |
rowpointers |
from spam object |
smooth_limits |
smooth limits |
sparse covariance matrix between locs and pred_locs
Sparse covariance function
cov_rns_taper_pred( theta, locs, locs_pred, x_covariates, x_covariates_pred, colindices, rowpointers, smooth_limits )
cov_rns_taper_pred( theta, locs, locs_pred, x_covariates, x_covariates_pred, colindices, rowpointers, smooth_limits )
theta |
vector of parameters |
locs |
a matrix with locations |
locs_pred |
a matrix with prediction locations |
x_covariates |
design data.frame |
x_covariates_pred |
design data.frame at prediction locations |
colindices |
from spam object |
rowpointers |
from spam object |
smooth_limits |
smooth limits |
sparse covariance matrix at locs
Retrieve the Akaike information criterion from a fitted coco object.
getAIC(coco.object)
getAIC(coco.object)
coco.object |
|
(numeric
) the associated AIC value
Federico Blasi
Retrieve BIC from a fitted coco object.
getBIC(coco.object)
getBIC(coco.object)
coco.object |
|
(numeric
) the associated BIC value
Federico Blasi
provides a generic set of upper and lower bounds for the L-BFGS-B routine
getBoundaries(x, lower.value, upper.value)
getBoundaries(x, lower.value, upper.value)
x |
|
lower.value |
|
upper.value |
|
(list
) a list with boundaries and simple init values for the optim L-BFGS-B routine
Federico Blasi
provides a generic set of upper and lower bounds for the L-BFGS-B routine
getBoundariesV2(coco.object, mean.limits, std.dev.limits, scale.limits, aniso.limits, tilt.limits, smooth.limits, nugget.limits)
getBoundariesV2(coco.object, mean.limits, std.dev.limits, scale.limits, aniso.limits, tilt.limits, smooth.limits, nugget.limits)
coco.object |
|
mean.limits |
|
std.dev.limits |
|
scale.limits |
|
aniso.limits |
|
tilt.limits |
|
smooth.limits |
|
nugget.limits |
|
(list
) a list with boundaries for the optim L-BFGS-B routine
Federico Blasi
provides a generic set of upper and lower bounds for the L-BFGS-B routine
getBoundariesV3(coco.object, mean.limits, global.lower, std.dev.max.effects, scale.max.effects, aniso.max.effects, tilt.max.effects, smooth.max.effects, nugget.max.effects)
getBoundariesV3(coco.object, mean.limits, global.lower, std.dev.max.effects, scale.max.effects, aniso.max.effects, tilt.max.effects, smooth.max.effects, nugget.max.effects)
coco.object |
|
mean.limits |
|
global.lower |
|
std.dev.max.effects |
|
scale.max.effects |
|
aniso.max.effects |
|
tilt.max.effects |
|
smooth.max.effects |
|
nugget.max.effects |
|
(list
) a list with boundaries for the optim L-BFGS-B routine
Federico Blasi
Compute approximate confidence intervals for a (fitted) coco object.
getCIs(coco.object, inv.hess, alpha = 0.05)
getCIs(coco.object, inv.hess, alpha = 0.05)
coco.object |
|
inv.hess |
|
alpha |
|
(numeric matrix
) a matrix with approximate confidence intervals for each parameter in the model.
Federico Blasi
Compute the covariance matrix of coco.object
.
getCovMatrix(coco.object, type = 'global', index = NULL)
getCovMatrix(coco.object, type = 'global', index = NULL)
coco.object |
|
type |
|
index |
|
(matrix
or S4
) a n x n covariance matrix (for 'dense' coco objects) or a S4 spam object (for 'sparse' coco objects).
Federico Blasi
## Not run: model.list <- list('mean' = 0, 'std.dev' = formula( ~ 1 + cov_x + cov_y), 'scale' = formula( ~ 1 + cov_x + cov_y), 'aniso' = 0, 'tilt' = 0, 'smooth' = 3/2, 'nugget' = -Inf) coco_object <- coco(type = 'dense', data = holes[[1]][1:100,], locs = as.matrix(holes[[1]][1:100,1:2]), z = holes[[1]][1:100,]$z, model.list = model.list) optim_coco <- cocoOptim(coco_object, boundaries = getBoundaries(coco_object, lower.value = -3, 3)) getCovMatrix(optim_coco) ## End(Not run)
## Not run: model.list <- list('mean' = 0, 'std.dev' = formula( ~ 1 + cov_x + cov_y), 'scale' = formula( ~ 1 + cov_x + cov_y), 'aniso' = 0, 'tilt' = 0, 'smooth' = 3/2, 'nugget' = -Inf) coco_object <- coco(type = 'dense', data = holes[[1]][1:100,], locs = as.matrix(holes[[1]][1:100,1:2]), z = holes[[1]][1:100,]$z, model.list = model.list) optim_coco <- cocoOptim(coco_object, boundaries = getBoundaries(coco_object, lower.value = -3, 3)) getCovMatrix(optim_coco) ## End(Not run)
Retrieves the Continuous Ranked Probability Score (CRPS) [1].
getCRPS(z.pred, mean.pred, sd.pred)
getCRPS(z.pred, mean.pred, sd.pred)
z.pred |
|
mean.pred |
|
sd.pred |
|
(numeric vector
) retrieves CRPS.
Federico Blasi
[1] Gneiting, Tilmann, and Adrian E. Raftery. "Strictly proper scoring rules, prediction, and estimation." Journal of the American statistical Association 102.477 (2007): 359-378.
Creates a unique design matrix based on model specification for each of the different potentially spatially varying aspects.
getDesignMatrix(model.list, data)
getDesignMatrix(model.list, data)
model.list |
|
data |
|
(list
) a list with two elements: a design matrix of dimension
(n x p), and a par.pos object, indexing columns of the design matrix to each of the spatially-varying functions.
Federico Blasi
Retrieve estimates from a fitted coco object.
getEstims(coco.object)
getEstims(coco.object)
coco.object |
|
(list
) a list with the estimates parameters for the different aspects
Federico Blasi
numerically approximate the Hessian. Hessians of parameters based on "pmle" are based on full likelihoods.
getHessian(coco.object, ncores = parallel::detectCores() - 1, eps = .Machine$double.eps^(1/4))
getHessian(coco.object, ncores = parallel::detectCores() - 1, eps = .Machine$double.eps^(1/4))
coco.object |
|
ncores |
|
eps |
|
(numeric matrix
) a symmetric matrix pxp of the approximated (observed) Hessian
Federico Blasi
Retrieve the loglikelihood value from a fitted coco object.
getLoglik(coco.object)
getLoglik(coco.object)
coco.object |
|
(numeric
) wrap for value from a OptimParallel object
Federico Blasi
Computes the Log-Score [1].
getLogScore(z.pred, mean.pred, sd.pred)
getLogScore(z.pred, mean.pred, sd.pred)
z.pred |
|
mean.pred |
|
sd.pred |
|
(numeric vector
) retrieves Log-Score.
Federico Blasi
[1] Gneiting, Tilmann, and Adrian E. Raftery. "Strictly proper scoring rules, prediction, and estimation." Journal of the American statistical Association 102.477 (2007): 359-378.
Returns a list of parameter vectors for each of the aspects.
getModelLists(theta, par.pos, type = 'diff')
getModelLists(theta, par.pos, type = 'diff')
theta |
|
par.pos |
|
type |
|
(list
) a list of different spatial aspects and mean required for the cov.rns functions
Federico Blasi
Based on the inverse of the Hessian (based on the difference parameterization for the std.dev and scale parameters), retrieves the modified inverse of the hessian (i.e. std.dev and scale).
getModHess(coco.object, inv.hess)
getModHess(coco.object, inv.hess)
coco.object |
|
inv.hess |
|
(numeric matrix
) the modified inverse of the hessian matrix
Federico Blasi
compute the negative 2 log likelihood based on theta
GetNeg2loglikelihood(theta, par.pos, locs, x_covariates, smooth.limits, z, n, lambda, safe = TRUE)
GetNeg2loglikelihood(theta, par.pos, locs, x_covariates, smooth.limits, z, n, lambda, safe = TRUE)
theta |
|
par.pos |
|
locs |
|
x_covariates |
|
smooth.limits |
|
z |
|
n |
|
lambda |
|
safe |
|
value
Federico Blasi
compute the negative 2 log likelihood based on theta
GetNeg2loglikelihoodProfile(theta, par.pos, locs, x_covariates, smooth.limits, z, n, x_betas,lambda, safe = TRUE)
GetNeg2loglikelihoodProfile(theta, par.pos, locs, x_covariates, smooth.limits, z, n, x_betas,lambda, safe = TRUE)
theta |
|
par.pos |
|
locs |
|
x_covariates |
|
smooth.limits |
|
z |
|
n |
|
x_betas |
|
lambda |
|
safe |
|
value
Federico Blasi
compute the negative 2 log likelihood based on theta
GetNeg2loglikelihoodTaper(theta, par.pos, ref_taper, locs, x_covariates, smooth.limits, cholS, z, n, lambda, safe = TRUE)
GetNeg2loglikelihoodTaper(theta, par.pos, ref_taper, locs, x_covariates, smooth.limits, cholS, z, n, lambda, safe = TRUE)
theta |
|
par.pos |
|
ref_taper |
|
locs |
|
x_covariates |
|
smooth.limits |
|
cholS |
|
z |
|
n |
|
lambda |
|
safe |
|
value
Federico Blasi
compute the negative 2 log likelihood based on theta
GetNeg2loglikelihoodTaperProfile(theta, par.pos, ref_taper, locs, x_covariates, smooth.limits, cholS, z, n, lambda, safe = TRUE)
GetNeg2loglikelihoodTaperProfile(theta, par.pos, ref_taper, locs, x_covariates, smooth.limits, cholS, z, n, lambda, safe = TRUE)
theta |
|
par.pos |
|
ref_taper |
|
locs |
|
x_covariates |
|
smooth.limits |
|
cholS |
|
z |
|
n |
|
lambda |
|
safe |
|
(numeric)
Federico Blasi
Centers and scale the design matrix.
getScale(x, mean.vector = NULL, sd.vector = NULL)
getScale(x, mean.vector = NULL, sd.vector = NULL)
x |
|
mean.vector |
|
sd.vector |
|
(list
) a list with a scaled design matrix of dimension n x (p+1), and a set of mean and sd vectors
employed to scale the matrix
Federico Blasi
Evaluates the spatially-varying functions of the nonstationary spatial structure.
getSpatEffects(coco.object)
getSpatEffects(coco.object)
coco.object |
|
(list
) a list with the different estimated surfaces.
Federico Blasi
Compute the trend of the (fitted) coco object.
getTrend(coco.object)
getTrend(coco.object)
coco.object |
|
(numeric vector
) a vector with the adjusted trend.
Federico Blasi
The synthetic "holes" provides a set of training and test data.frame of a Gaussian process realization with a (inherently dense) nonstationary covariance function. Four holes are present in the training dataset, and the task is to predict them.
holes
holes
A list with training and test data.frame with rows and variables:
first spatial coordinate
second spatial coordinate
first spatial characteristic
second spatial characteristic
response variable
Source of the data
data(holes)
data(holes)
The synthetic "holes_bm" provides a set of training and test data.frame of a Gaussian process realization with a (inherently dense) nonstationary covariance function. Four holes are present in the training dataset, and the task is to predict them. This version provides ten independent realizations of the process, as well as considers a spatial mean effect.
holes_bm
holes_bm
A list with training, training.z, test, and test.z data.frames with rows and variables:
first spatial coordinate
second spatial coordinate
first spatial characteristic
second spatial characteristic
third spatial characteristic
i-th response variable
Source of the data
data(holes_bm)
data(holes_bm)
check whether an R object is a formula
is.formula(x)
is.formula(x)
x |
(ANY) an R object. |
TRUE/FALSE
Federico Blasi
This method plots objects of class coco
.
## S4 method for signature 'coco,missing' plot(x, y, type = NULL, index = NULL, factr = 0.1, ...)
## S4 method for signature 'coco,missing' plot(x, y, type = NULL, index = NULL, factr = 0.1, ...)
x |
( |
y |
Not used. |
type |
( |
index |
( |
factr |
( |
... |
Additional arguments passed to quilt.plot. |
Several plots are created.
Federico Blasi
plot output of optim
plotOptimInfo(coco.object, ...)
plotOptimInfo(coco.object, ...)
coco.object |
an optimized coco.object |
... |
arguments for par() |
Outputs a sequence of plots detailing parameters during the optimization routine
Federico Blasi
This method show objects of class 'coco'.
## S4 method for signature 'coco' show(object)
## S4 method for signature 'coco' show(object)
object |
( |
A plot is created.
Federico Blasi
The synthetic "stripes" provides a set of training and test data.frame of a Gaussian process realization with a (inherently sparse) nonstationary covariance function. Several stripes are present in the training dataset, and the task is to predict them.
stripes
stripes
A list with training and test data.frame with rows and variables:
first spatial coordinate
second spatial coordinate
first spatial characteristic
second spatial characteristic
third spatial characteristic
response variable
Source of the data
data(stripes)
data(stripes)
method summary for objects of class 'coco'.
## S4 method for signature 'coco' summary(object, inv.hess = NULL)
## S4 method for signature 'coco' summary(object, inv.hess = NULL)
object |
( |
inv.hess |
( |
summary the coco object
Federico Blasi