Title: | Convolution-Based Nonstationary Spatial Modeling |
---|---|
Description: | Fits convolution-based nonstationary Gaussian process models to point-referenced spatial data. The nonstationary covariance function allows the user to specify the underlying correlation structure and which spatial dependence parameters should be allowed to vary over space: the anisotropy, nugget variance, and process variance. The parameters are estimated via maximum likelihood, using a local likelihood approach. Also provided are functions to fit stationary spatial models for comparison, calculate the Kriging predictor and standard errors, and create various plots to visualize nonstationarity. |
Authors: | Mark D. Risser [aut, cre] |
Maintainer: | Mark D. Risser <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.2.5 |
Built: | 2025-02-12 04:41:05 UTC |
Source: | https://github.com/markdrisser/convospat |
Aniso_fit
estimates the parameters of the stationary spatial model.
Required inputs are the observed data and locations (a geoR object
with $coords and $data). Optional inputs include the covariance model
(exponential is the default).
Aniso_fit(geodata = NULL, sp.SPDF = NULL, coords = geodata$coords, data = geodata$data, cov.model = "exponential", mean.model = data ~ 1, fixed.nugg2.var = NULL, method = "reml", fix.tausq = FALSE, tausq = 0, fix.kappa = FALSE, kappa = 0.5, local.pars.LB = NULL, local.pars.UB = NULL, local.ini.pars = NULL)
Aniso_fit(geodata = NULL, sp.SPDF = NULL, coords = geodata$coords, data = geodata$data, cov.model = "exponential", mean.model = data ~ 1, fixed.nugg2.var = NULL, method = "reml", fix.tausq = FALSE, tausq = 0, fix.kappa = FALSE, kappa = 0.5, local.pars.LB = NULL, local.pars.UB = NULL, local.ini.pars = NULL)
geodata |
A list containing elements |
sp.SPDF |
A " |
coords |
An N x 2 matrix where each row has the two-dimensional
coordinates of the N data locations. By default, it takes the |
data |
A vector or matrix with N rows, containing the data values. Inputting a vector corresponds to a single replicate of data, while inputting a matrix corresponds to replicates. In the case of replicates, the model assumes the replicates are independent and identically distributed. |
cov.model |
A string specifying the model for the correlation
function; following |
mean.model |
An object of class |
fixed.nugg2.var |
Optional; describes the variance/covariance for a fixed (second) nugget term (represents a known error term). Either a vector of length N containing a station-specific variances (implying independent error) or an NxN covariance matrix (implying dependent error). Defaults to zero. |
method |
Indicates the estimation method, either maximum likelihood
( |
fix.tausq |
Logical; indicates whether the default nugget term
(tau^2) should be fixed ( |
tausq |
Scalar; fixed value for the nugget variance (when
|
fix.kappa |
Logical; indicates if the kappa parameter should be
fixed ( |
kappa |
Scalar; value of the kappa parameter. Only used if
|
local.pars.LB , local.pars.UB
|
Optional vectors of lower and upper
bounds, respectively, used by the |
local.ini.pars |
Optional vector of initial values used by the
|
A list with the following components:
MLEs.save |
Table of local maximum likelihood estimates for each mixture component location. |
data |
Observed data values. |
beta.GLS |
Vector of generalized least squares estimates of beta, the mean coefficients. |
beta.cov |
Covariance matrix of the generalized least squares estimate of beta. |
Mean.coefs |
"Regression table" for the mean coefficient estimates, listing the estimate, standard error, and t-value. |
Cov.mat |
Estimated covariance matrix ( |
Cov.mat.chol |
Cholesky of |
aniso.pars |
Vector of MLEs for the anisotropy parameters lam1, lam2, eta. |
aniso.mat |
2 x 2 anisotropy matrix, calculated from
|
tausq.est |
Scalar maximum likelihood estimate of tausq (nugget variance). |
sigmasq.est |
Scalar maximum likelihood estimate of sigmasq (process variance). |
kappa.MLE |
Scalar maximum likelihood estimate for kappa (when applicable). |
fixed.nugg2.var |
N x N matrix with the fixed variance/covariance for the second (measurement error) nugget term (defaults to zero). |
cov.model |
String; the correlation model used for estimation. |
coords |
N x 2 matrix of observation locations. |
global.loglik |
Scalar value of the maximized likelihood from the global optimization (if available). |
Xmat |
Design matrix, obtained from using |
fix.kappa |
Logical, indicating if kappa was fixed ( |
kappa |
Scalar; fixed value of kappa. |
# Using iid standard Gaussian data aniso.fit <- Aniso_fit( coords = cbind(runif(100), runif(100)), data = rnorm(100) )
# Using iid standard Gaussian data aniso.fit <- Aniso_fit( coords = cbind(runif(100), runif(100)), data = rnorm(100) )
Calculate three evaluation criteria – continuous rank probability score (CRPS), prediction mean square deviation ratio (pMSDR), and mean squared prediction error (MSPE) – comparing hold-out data and predictions.
evaluate_CV(holdout.data, pred.mean, pred.SDs)
evaluate_CV(holdout.data, pred.mean, pred.SDs)
holdout.data |
Observed/true data that has been held out for model comparison. |
pred.mean |
Predicted mean values corresponding to the hold-out locations. |
pred.SDs |
Predicted standard errors corresponding to the hold-out locations. |
A list with the following components:
CRPS |
The CRPS averaged over all hold-out locations. |
MSPE |
The mean squared prediction error. |
pMSDR |
The prediction mean square deviation ratio. |
## Not run: evaluate_CV( holdout.data = simdata$sim.data[holdout.index], pred.mean = pred.NS$pred.means, pred.SDs = pred.NS$pred.SDs ) ## End(Not run)
## Not run: evaluate_CV( holdout.data = simdata$sim.data[holdout.index], pred.mean = pred.NS$pred.means, pred.SDs = pred.NS$pred.SDs ) ## End(Not run)
f_mc_kernels
calculates spatially-varying mixture component kernels using
generalized linear models for each of the eigenvalues (lam1 and lam2) and
the angle of rotation (eta).
f_mc_kernels(y.min = 0, y.max = 5, x.min = 0, x.max = 5, N.mc = 3^2, lam1.coef = c(-1.3, 0.5, -0.6), lam2.coef = c(-1.4, -0.1, 0.2), logit.eta.coef = c(0, -0.15, 0.15))
f_mc_kernels(y.min = 0, y.max = 5, x.min = 0, x.max = 5, N.mc = 3^2, lam1.coef = c(-1.3, 0.5, -0.6), lam2.coef = c(-1.4, -0.1, 0.2), logit.eta.coef = c(0, -0.15, 0.15))
y.min |
Lower bound for the y-coordinate axis. |
y.max |
Upper bound for the y-coordinate axis. |
x.min |
Lower bound for the y-coordinate axis. |
x.max |
Upper bound for the y-coordinate axis. |
N.mc |
Number of mixture component locations. |
lam1.coef |
Log-linear regression coefficients for lam1; the coefficients correspond to the intercept, longitude, and latitude. |
lam2.coef |
Log-linear regression coefficients for lam2; the coefficients correspond to the intercept, longitude, and latitude. |
logit.eta.coef |
Scaled logit regression coefficients for eta; the coefficients correspond to the intercept, longitude, and latitude. |
A list with the following components:
mc.locations |
A |
mc.kernels |
A |
f_mc_kernels( y.min = 0, y.max = 5, x.min = 0, x.max = 5, N.mc = 3^2, lam1.coef = c(-1.3, 0.5, -0.6), lam2.coef = c(-1.4, -0.1, 0.2), logit.eta.coef = c(0, -0.15, 0.15) )
f_mc_kernels( y.min = 0, y.max = 5, x.min = 0, x.max = 5, N.mc = 3^2, lam1.coef = c(-1.3, 0.5, -0.6), lam2.coef = c(-1.4, -0.1, 0.2), logit.eta.coef = c(0, -0.15, 0.15) )
kernel_cov
calculates a 2 x 2 matrix based on the eigendecomposition
components (two eigenvalues and angle of rotation).
kernel_cov(params)
kernel_cov(params)
params |
A vector of three parameters, corresponding to (lam1, lam2, eta). The eigenvalues (lam1 and lam2) must be positive. |
A 2 x 2 kernel covariance matrix.
kernel_cov(c(1, 2, pi/3))
kernel_cov(c(1, 2, pi/3))
This function generates another function to be used within optim
to
obtain maximum likelihood estimates of
global variance parameters tausq, sigmasq with a fixed correlation
matrix (smoothness is fixed).
make_global_loglik1(data, Xmat, Corr, nugg2.var)
make_global_loglik1(data, Xmat, Corr, nugg2.var)
data |
A vector or matrix of data to use in the likelihood calculation. |
Xmat |
The design matrix for the mean model. |
Corr |
The correlation matrix. |
nugg2.var |
Fixed values for the covariance of the second nugget term. |
This function returns another function for use in optim
.
## Not run: make_global_loglik1( data, Xmat, Corr, nugg2.var ) ## End(Not run)
## Not run: make_global_loglik1( data, Xmat, Corr, nugg2.var ) ## End(Not run)
This function generates another function to be used within optim
to
obtain maximum likelihood estimates of
global variance parameters tausq, sigmasq, and nu.
make_global_loglik1_kappa(data, Xmat, cov.model, Scalemat, Distmat, nugg2.var)
make_global_loglik1_kappa(data, Xmat, cov.model, Scalemat, Distmat, nugg2.var)
data |
A vector or matrix of data to use in the likelihood calculation. |
Xmat |
The design matrix for the mean model. |
cov.model |
String; the covariance model. |
Scalemat |
Matrix; contains the scaling quantities from the covariance function. |
Distmat |
Matrix; contains the scaled distances. |
nugg2.var |
Fixed values for the covariance of the second nugget term. |
This function returns another function for use in optim
.
## Not run: make_global_loglik1_kappa( data, Xmat, cov.model, Scalemat, Distmat, nugg2.var ) ## End(Not run)
## Not run: make_global_loglik1_kappa( data, Xmat, cov.model, Scalemat, Distmat, nugg2.var ) ## End(Not run)
This function generates another function to be used within optim
to
obtain maximum likelihood estimates of
global variance parameter sigmasq with a fixed correlation
matrix (smoothness is fixed). The nugget variance is taken
to be spatially-varing.
make_global_loglik2(data, Xmat, Corr, obs.nuggets, nugg2.var)
make_global_loglik2(data, Xmat, Corr, obs.nuggets, nugg2.var)
data |
A vector or matrix of data to use in the likelihood calculation. |
Xmat |
The design matrix for the mean model. |
Corr |
The correlation matrix. |
obs.nuggets |
A vector containing the spatially-varying nuggets corresponding to each data location. |
nugg2.var |
Fixed values for the covariance of the second nugget term. |
This function returns another function for use in optim
.
## Not run: make_global_loglik2( data, Xmat, Corr, obs.nuggets, nugg2.var ) ## End(Not run)
## Not run: make_global_loglik2( data, Xmat, Corr, obs.nuggets, nugg2.var ) ## End(Not run)
This function generates another function to be used within optim
to
obtain maximum likelihood estimates of
global variance parameters sigmasq and nu. The nugget variance is
taken to be spatially-varying.
make_global_loglik2_kappa(data, Xmat, cov.model, Scalemat, Distmat, obs.nuggets, nugg2.var)
make_global_loglik2_kappa(data, Xmat, cov.model, Scalemat, Distmat, obs.nuggets, nugg2.var)
data |
A vector or matrix of data to use in the likelihood calculation. |
Xmat |
The design matrix for the mean model. |
cov.model |
String; the covariance model. |
Scalemat |
Matrix; contains the scaling quantities from the covariance function. |
Distmat |
Matrix; contains the scaled distances. |
obs.nuggets |
A vector containing the spatially-varying nuggets corresponding to each data location. |
nugg2.var |
Fixed values for the covariance of the second nugget term. |
This function returns another function for use in optim
.
## Not run: make_global_loglik2_kappa( data, Xmat, cov.model, Scalemat, Distmat, obs.nuggets, nugg2.var ) ## End(Not run)
## Not run: make_global_loglik2_kappa( data, Xmat, cov.model, Scalemat, Distmat, obs.nuggets, nugg2.var ) ## End(Not run)
This function generates another function to be used within optim
to
obtain maximum likelihood estimates of
global variance parameter tausq with a fixed correlation
matrix (smoothness is fixed). The process variance is taken
to be spatially-varing.
make_global_loglik3(data, Xmat, Corr, obs.variance, nugg2.var)
make_global_loglik3(data, Xmat, Corr, obs.variance, nugg2.var)
data |
A vector or matrix of data to use in the likelihood calculation. |
Xmat |
The design matrix for the mean model. |
Corr |
The correlation matrix matrix. |
obs.variance |
A vector containing the spatially-varying variance corresponding to each data location. |
nugg2.var |
Fixed values for the covariance of the second nugget term. |
This function returns another function for use in optim
.
## Not run: make_global_loglik3( data, Xmat, Corr, obs.variance, nugg2.var ) ## End(Not run)
## Not run: make_global_loglik3( data, Xmat, Corr, obs.variance, nugg2.var ) ## End(Not run)
This function generates another function to be used within optim
to
obtain maximum likelihood estimates of
global variance parameters tausq and nu. The process variance is
taken to be spatially-varying.
make_global_loglik3_kappa(data, Xmat, cov.model, Scalemat, Distmat, obs.variance, nugg2.var)
make_global_loglik3_kappa(data, Xmat, cov.model, Scalemat, Distmat, obs.variance, nugg2.var)
data |
A vector or matrix of data to use in the likelihood calculation. |
Xmat |
The design matrix for the mean model. |
cov.model |
String; the covariance model. |
Scalemat |
Matrix; contains the scaling quantities from the covariance function. |
Distmat |
Matrix; contains the scaled distances. |
obs.variance |
A vector containing the spatially-varying variance corresponding to each data location. |
nugg2.var |
Fixed values for the covariance of the second nugget term. |
This function returns another function for use in optim
.
## Not run: make_global_loglik3_kappa( data, Xmat, cov.model, Scalemat, Distmat, obs.variance, nugg2.var ) ## End(Not run)
## Not run: make_global_loglik3_kappa( data, Xmat, cov.model, Scalemat, Distmat, obs.variance, nugg2.var ) ## End(Not run)
This function generates another function to be used within optim
to
obtain maximum likelihood estimates of
global variance parameters nu. The process variance
and nugget variance are taken to be spatially-varying.
make_global_loglik4_kappa(data, Xmat, cov.model, Scalemat, Distmat, obs.variance, obs.nuggets, nugg2.var)
make_global_loglik4_kappa(data, Xmat, cov.model, Scalemat, Distmat, obs.variance, obs.nuggets, nugg2.var)
data |
A vector or matrix of data to use in the likelihood calculation. |
Xmat |
The design matrix for the mean model. |
cov.model |
String; the covariance model. |
Scalemat |
Matrix; contains the scaling quantities from the covariance function. |
Distmat |
Matrix; contains the scaled distances. |
obs.variance |
A vector containing the spatially-varying variance corresponding to each data location. |
obs.nuggets |
A vector containing the spatially-varying nuggets corresponding to each data location. |
nugg2.var |
Fixed values for the covariance of the second nugget term. |
This function returns another function for use in optim
.
## Not run: make_global_loglik4_kappa( data, Xmat, cov.model, Scalemat, Distmat, obs.variance, obs.nuggets, nugg2.var ) ## End(Not run)
## Not run: make_global_loglik4_kappa( data, Xmat, cov.model, Scalemat, Distmat, obs.variance, obs.nuggets, nugg2.var ) ## End(Not run)
This function generates another function to be used within optim
to
obtain maximum likelihood estimates of covariance (and possibly mean) parameters.
The function includes options for
(1) maximum likelihood ("ml"
) vs. restricted maximum likelihood
("reml"
),
(2) smoothness (kappa
): models without smoothness vs. estimating the
smoothness vs. using fixed smoothness,
(3) locally isotropic vs. locally anisotropic, and
(4) fixed nugget variance (tausq
): fixed vs. estimated.
make_local_lik(locations, cov.model, data, Xmat, nugg2.var = matrix(0, nrow(locations), nrow(locations)), tausq = 0, kappa = 0.5, fixed = rep(FALSE, 6), method = "reml", local.aniso = TRUE, fix.tausq = FALSE, fix.kappa = FALSE)
make_local_lik(locations, cov.model, data, Xmat, nugg2.var = matrix(0, nrow(locations), nrow(locations)), tausq = 0, kappa = 0.5, fixed = rep(FALSE, 6), method = "reml", local.aniso = TRUE, fix.tausq = FALSE, fix.kappa = FALSE)
locations |
A matrix of locations. |
cov.model |
String; the covariance model. |
data |
A vector or matrix of data to use in the likelihood calculation. |
Xmat |
The design matrix for the mean model. |
nugg2.var |
Fixed values for the variance/covariance of the second nugget term; defaults to a matrix of zeros. |
tausq |
Scalar; fixed value for the nugget variance (when
|
kappa |
Scalar; fixed value for the smoothness (when |
fixed |
Logical vector of |
method |
Indicates the estimation method, either maximum likelihood ( |
local.aniso |
Logical; indicates if the local covariance should be
anisotropic ( |
fix.tausq |
Logical; indicates whether the default nugget term
(tau^2) should be fixed ( |
fix.kappa |
Logical; indicates if the kappa parameter should be
fixed ( |
This function returns another function for use in optim
.
## Not run: make_local_lik( locations, cov.model, data, Xmat ) ## End(Not run)
## Not run: make_local_lik( locations, cov.model, data, Xmat ) ## End(Not run)
mc_N
calculates the number of observations (sample size) that
fall within a certain fit radius for each mixture component location.
mc_N(coords, mc.locations, fit.radius)
mc_N(coords, mc.locations, fit.radius)
coords |
A matrix of observation locations. |
mc.locations |
A matrix of the mixture component locations to use in the model fitting. |
fit.radius |
Scalar; defines the fitting radius for local likelihood estimation. |
A vector mc.N.fit
, which summarizes the number of
observation locations in coords
that fall within the fit radius
for each mixture component location.
## Not run: mc_N( coords = simdata$sim.locations, mc.locations = simdata$mc.locations, fit.radius = 1 ) ## End(Not run)
## Not run: mc_N( coords = simdata$sim.locations, mc.locations = simdata$mc.locations, fit.radius = 1 ) ## End(Not run)
NSconvo_fit
estimates the parameters of the nonstationary
convolution-based spatial model. Required inputs are the observed data and
locations (a geoR object with $coords and $data).
Optional inputs include mixture component locations (if not provided,
the number of mixture component locations are required), the fit radius,
the covariance model (exponential is the default), and whether or not the
nugget and process variance will be spatially-varying.
NSconvo_fit(geodata = NULL, sp.SPDF = NULL, coords = geodata$coords, data = geodata$data, cov.model = "exponential", mean.model = data ~ 1, mc.locations = NULL, N.mc = NULL, lambda.w = NULL, fixed.nugg2.var = NULL, mean.model.df = NULL, mc.kernels = NULL, fit.radius = NULL, ns.nugget = FALSE, ns.variance = FALSE, ns.mean = FALSE, local.aniso = TRUE, fix.tausq = FALSE, tausq = 0, fix.kappa = FALSE, kappa = 0.5, method = "reml", print.progress = TRUE, local.pars.LB = NULL, local.pars.UB = NULL, global.pars.LB = NULL, global.pars.UB = NULL, local.ini.pars = NULL, global.ini.pars = NULL)
NSconvo_fit(geodata = NULL, sp.SPDF = NULL, coords = geodata$coords, data = geodata$data, cov.model = "exponential", mean.model = data ~ 1, mc.locations = NULL, N.mc = NULL, lambda.w = NULL, fixed.nugg2.var = NULL, mean.model.df = NULL, mc.kernels = NULL, fit.radius = NULL, ns.nugget = FALSE, ns.variance = FALSE, ns.mean = FALSE, local.aniso = TRUE, fix.tausq = FALSE, tausq = 0, fix.kappa = FALSE, kappa = 0.5, method = "reml", print.progress = TRUE, local.pars.LB = NULL, local.pars.UB = NULL, global.pars.LB = NULL, global.pars.UB = NULL, local.ini.pars = NULL, global.ini.pars = NULL)
geodata |
A list containing elements |
sp.SPDF |
A " |
coords |
An N x 2 matrix where each row has the two-dimensional
coordinates of the N data locations. By default, it takes the |
data |
A vector or matrix with N rows, containing the data values. Inputting a vector corresponds to a single replicate of data, while inputting a matrix corresponds to replicates. In the case of replicates, the model assumes the replicates are independent and identically distributed. |
cov.model |
A string specifying the model for the correlation
function; following |
mean.model |
An object of class |
mc.locations |
Optional; matrix of mixture component locations. |
N.mc |
Optional; if |
lambda.w |
Scalar; tuning parameter for the weight function. Defaults to be the square of one-half of the minimum distance between mixture component locations. |
fixed.nugg2.var |
Optional; describes the variance/covariance for a fixed (second) nugget term (represents a known error term). Either a vector of length N containing a station-specific variances (implying independent error) or an NxN covariance matrix (implying dependent error). Defaults to zero. |
mean.model.df |
Optional data frame; refers to the variables used
in |
mc.kernels |
Optional specification of mixture component kernel matrices (based on expert opinion, etc.). |
fit.radius |
Scalar; specifies the fit radius or neighborhood size for the local likelihood estimation. |
ns.nugget |
Logical; indicates if the nugget variance (tausq) should
be spatially-varying ( |
ns.variance |
Logical; indicates if the process variance (sigmasq)
should be spatially-varying ( |
ns.mean |
Logical; indicates if the mean coefficeints (beta)
should be spatially-varying ( |
local.aniso |
Logical; indicates if the local covariance should be
anisotropic ( |
fix.tausq |
Logical; indicates whether the default nugget term
(tau^2) should be fixed ( |
tausq |
Scalar; fixed value for the nugget variance (when
|
fix.kappa |
Logical; indicates if the kappa parameter should be
fixed ( |
kappa |
Scalar; value of the kappa parameter. Only used if
|
method |
Indicates the estimation method, either maximum likelihood
( |
print.progress |
Logical; if |
local.pars.LB , local.pars.UB
|
Optional vectors of lower and upper
bounds, respectively, used by the |
global.pars.LB , global.pars.UB
|
Optional vectors of lower and upper
bounds, respectively, used by the |
local.ini.pars |
Optional vector of initial values used by the
|
global.ini.pars |
Optional vector of initial values used by the
|
A "NSconvo" object, with the following components:
mc.locations |
Mixture component locations used for the simulated data. |
mc.kernels |
Mixture component kernel matrices used for the simulated data. |
MLEs.save |
Table of local maximum likelihood estimates for each mixture component location. |
kernel.ellipses |
|
data |
Observed data values. |
beta.GLS |
Generalized least squares estimates of beta,
the mean coefficients. For |
beta.cov |
Covariance matrix of the generalized least squares
estimate of beta. For |
Mean.coefs |
"Regression table" for the mean coefficient estimates,
listing the estimate, standard error, and t-value (for |
tausq.est |
Estimate of tausq (nugget variance), either scalar (when
|
sigmasq.est |
Estimate of sigmasq (process variance), either scalar
(when |
beta.est |
Estimate of beta (mean coefficients), either a vector
(when |
kappa.MLE |
Scalar maximum likelihood estimate for kappa (when applicable). |
Cov.mat |
Estimated covariance matrix ( |
Cov.mat.chol |
Cholesky of |
cov.model |
String; the correlation model used for estimation. |
ns.nugget |
Logical, indicating if the nugget variance was estimated
as spatially-varing ( |
ns.variance |
Logical, indicating if the process variance was
estimated as spatially-varying ( |
fixed.nugg2.var |
N x N matrix with the fixed variance/covariance for the second (measurement error) nugget term (defaults to zero). |
coords |
N x 2 matrix of observation locations. |
global.loglik |
Scalar value of the maximized likelihood from the global optimization (if available). |
Xmat |
Design matrix, obtained from using |
lambda.w |
Tuning parameter for the weight function. |
fix.kappa |
Logical, indicating if kappa was fixed ( |
kappa |
Scalar; fixed value of kappa. |
# Using white noise data fit.model <- NSconvo_fit( coords = cbind( runif(100), runif(100)), data = rnorm(100), fit.radius = 0.4, N.mc = 4 )
# Using white noise data fit.model <- NSconvo_fit( coords = cbind( runif(100), runif(100)), data = rnorm(100), fit.radius = 0.4, N.mc = 4 )
NSconvo_sim
simulates data from the nonstationary model, given
mixture component kernel matrices. The function requires either a mixture
component kernel object, from the function f.mc.kernels(), or a direct
specification of the mixture component locations and mixture component
kernels.
NSconvo_sim(grid = TRUE, y.min = 0, y.max = 5, x.min = 0, x.max = 5, N.obs = 20^2, sim.locations = NULL, mc.kernels.obj = NULL, mc.kernels = NULL, mc.locations = NULL, lambda.w = NULL, tausq = 0.1, sigmasq = 1, beta.coefs = 4, kappa = NULL, covariates = rep(1, N.obs), cov.model = "exponential")
NSconvo_sim(grid = TRUE, y.min = 0, y.max = 5, x.min = 0, x.max = 5, N.obs = 20^2, sim.locations = NULL, mc.kernels.obj = NULL, mc.kernels = NULL, mc.locations = NULL, lambda.w = NULL, tausq = 0.1, sigmasq = 1, beta.coefs = 4, kappa = NULL, covariates = rep(1, N.obs), cov.model = "exponential")
grid |
Logical; indicates of the simulated data should fall on a
grid ( |
y.min |
Lower bound for the y-coordinate axis. |
y.max |
Upper bound for the y-coordinate axis. |
x.min |
Lower bound for the y-coordinate axis. |
x.max |
Upper bound for the y-coordinate axis. |
N.obs |
Number of simulated data values. |
sim.locations |
Optional |
mc.kernels.obj |
Object from the |
mc.kernels |
Optional specification of mixture component kernel matrices. |
mc.locations |
Optional specification of mixture component locations. |
lambda.w |
Scalar; tuning parameter for the weight function. |
tausq |
Scalar; true nugget variance. |
sigmasq |
Scalar; true process variance. |
beta.coefs |
Vector of true regression coefficients. Length must
match the number of columns in |
kappa |
Scalar; true smoothness. |
covariates |
Matrix with |
cov.model |
A string specifying the correlation function. Options
available in this package are: " |
A list with the following components:
sim.locations |
Matrix of locations for the simulated values. |
mc.locations |
Mixture component locations used for the simulated data. |
mc.kernels |
Mixture component kernel matrices used for the simulated data. |
kernel.ellipses |
|
Cov.mat |
True covariance matrix ( |
sim.data |
Simulated data values. |
lambda.w |
Tuning parameter for the weight function. |
## Not run: NSconvo_sim( grid = TRUE, y.min = 0, y.max = 5, x.min = 0, x.max = 5, N.obs = 20^2, sim.locations = NULL, mc.kernels.obj = NULL, mc.kernels = NULL, mc.locations = NULL, lambda.w = NULL, tausq = 0.1, sigmasq = 1, beta.coefs = 4, kappa = NULL, covariates = rep(1,N.obs), cov.model = "exponential" ) ## End(Not run)
## Not run: NSconvo_sim( grid = TRUE, y.min = 0, y.max = 5, x.min = 0, x.max = 5, N.obs = 20^2, sim.locations = NULL, mc.kernels.obj = NULL, mc.kernels = NULL, mc.locations = NULL, lambda.w = NULL, tausq = 0.1, sigmasq = 1, beta.coefs = 4, kappa = NULL, covariates = rep(1,N.obs), cov.model = "exponential" ) ## End(Not run)
This function plots the estimated correlation between a reference point and all other prediction locations.
## S3 method for class 'Aniso' plot(x, ref.loc = NULL, all.pred.locs = NULL, grid = TRUE, ...)
## S3 method for class 'Aniso' plot(x, ref.loc = NULL, all.pred.locs = NULL, grid = TRUE, ...)
x |
An "Aniso" object, from Aniso_fit(). |
ref.loc |
Vector of length 2; the reference location. |
all.pred.locs |
A matrix of all prediction locations. |
grid |
Logical; indicates if the |
... |
Arguments passed to |
A plot of either the estimated ellipses or estimated correlation is printed.
## Not run: plot.Aniso( Aniso.object ) ## End(Not run)
## Not run: plot.Aniso( Aniso.object ) ## End(Not run)
This function plots either the estimated anisotropy ellipses for each of the mixture component locations or the estimated correlation between a reference point and all other prediction locations.
## S3 method for class 'NSconvo' plot(x, plot.ellipses = TRUE, fit.radius = NULL, aniso.mat = NULL, true.mc = NULL, ref.loc = NULL, all.pred.locs = NULL, grid = TRUE, true.col = 1, aniso.col = 4, ns.col = 2, plot.mc.locs = TRUE, ...)
## S3 method for class 'NSconvo' plot(x, plot.ellipses = TRUE, fit.radius = NULL, aniso.mat = NULL, true.mc = NULL, ref.loc = NULL, all.pred.locs = NULL, grid = TRUE, true.col = 1, aniso.col = 4, ns.col = 2, plot.mc.locs = TRUE, ...)
x |
A "NSconvo" object, from NSconvo_fit(). |
plot.ellipses |
Logical; indicates whether the estimated
ellipses should be plotted ( |
fit.radius |
Scalar; defines the fit radius used for the local likelihood estimation. |
aniso.mat |
2 x 2 matrix; contains the estimated anisotropy ellipse from the stationary model (for comparison). |
true.mc |
The true mixture component ellipses, if known. |
ref.loc |
Vector of length 2; the reference location. |
all.pred.locs |
A matrix of all prediction locations. |
grid |
Logical; indicates if the |
true.col |
Color value for the true mixture component ellipses (if plotted). |
aniso.col |
Color value for the anisotropy ellipse (if plotted). |
ns.col |
Color value for the mixture component ellipses. |
plot.mc.locs |
Logical; indicates whether the mixture
component locations should be plotted ( |
... |
Other options passed to |
A plot of either the estimated ellipses or estimated correlation is printed.
## Not run: plot.NSconvo( NSconvo.object ) ## End(Not run)
## Not run: plot.NSconvo( NSconvo.object ) ## End(Not run)
predict.Aniso
calculates the kriging predictor and corresponding
standard errors at unmonitored sites.
## S3 method for class 'Aniso' predict(object, pred.coords, pred.covariates = NULL, pred.fixed.nugg2.var = NULL, ...)
## S3 method for class 'Aniso' predict(object, pred.coords, pred.covariates = NULL, pred.fixed.nugg2.var = NULL, ...)
object |
An "Aniso" object, from |
pred.coords |
Matrix of locations where predictions are required. |
pred.covariates |
Matrix of covariates for the prediction locations,
NOT including an intercept. The number of columns for this matrix must
match the design matrix from |
pred.fixed.nugg2.var |
An optional vector or matrix describing the
the variance/covariance a fixed second nugget term (corresponds to
|
... |
additional arguments affecting the predictions produced. |
A list with the following components:
pred.means |
Vector of the kriging predictor, for each location in
|
pred.SDs |
Vector of the kriging standard errors, for each location
in |
## Not run: pred.S <- predict( Aniso.obj, pred.coords = cbind(runif(300),runif(300)) ) ## End(Not run)
## Not run: pred.S <- predict( Aniso.obj, pred.coords = cbind(runif(300),runif(300)) ) ## End(Not run)
predict.NSconvo
calculates the kriging predictor and corresponding
standard errors at unmonitored sites.
## S3 method for class 'NSconvo' predict(object, pred.coords, pred.covariates = NULL, pred.fixed.nugg2.var = NULL, ...)
## S3 method for class 'NSconvo' predict(object, pred.coords, pred.covariates = NULL, pred.fixed.nugg2.var = NULL, ...)
object |
A "NSconvo" object, from |
pred.coords |
Matrix of locations where predictions are required. |
pred.covariates |
Matrix of covariates for the prediction locations,
NOT including an intercept. The number of columns for this matrix must
match the design matrix from |
pred.fixed.nugg2.var |
An optional vector or matrix describing the
the variance/covariance a fixed second nugget term (corresponds to
|
... |
additional arguments affecting the predictions produced. |
A list with the following components:
pred.means |
Vector of the kriging predictor, for each location in
|
pred.SDs |
Vector of the kriging standard errors, for each location
in |
## Not run: pred.NS <- predict( NSconvo.obj, pred.coords = matrix(c(1,1), ncol=2), pred.covariates = matrix(c(1,1), ncol=2) ) ## End(Not run)
## Not run: pred.NS <- predict( NSconvo.obj, pred.coords = matrix(c(1,1), ncol=2), pred.covariates = matrix(c(1,1), ncol=2) ) ## End(Not run)
A data set containing the necessary components to fit the nonstationary spatial model, simulated from the true model.
simdata
simdata
A list with the following objects:
A matrix of longitude/latitude coordinates of the simulated locations.
A matrix of longitude/latitude coordinates of the mixture component locations.
A three-dimensional array, containing the true 2 x 2 kernel covariance matrices for each mixture component location.
A three-dimensional array, containing the true 2 x 2 kernel covariance matrices for each simulated location.
A matrix of the simulated data; each of the ten columns correspond to an independent and identically distribured replicate.
Scalar; the value of the tuning parameter used in the weight function.
Vector; indicates which of the simulated locations should be used in the hold-out sample.
summary.Aniso
prints relevant output from the model fitting
procedure.
## S3 method for class 'Aniso' summary(object, ...)
## S3 method for class 'Aniso' summary(object, ...)
object |
An "Aniso" object, from |
... |
additional arguments affecting the summary produced. |
Text containing the model fitting results.
## Not run: summary.Aniso( Aniso.object ) ## End(Not run)
## Not run: summary.Aniso( Aniso.object ) ## End(Not run)
summary.NSconvo
prints relevant output from the model fitting
procedure.
## S3 method for class 'NSconvo' summary(object, ...)
## S3 method for class 'NSconvo' summary(object, ...)
object |
A "NSconvo" object, from |
... |
additional arguments affecting the summary produced. |
Text containing the model fitting results.
## Not run: summary.NSconvo( NSconvo.object ) ## End(Not run)
## Not run: summary.NSconvo( NSconvo.object ) ## End(Not run)
A list of two mixture component grids for fitting the nonstationary model to the western United States precipitation data.
US.mc.grids
US.mc.grids
A list with two elements:
Coarse mixture component grid.
Fine mixture component grid.
A matrix with two columns containing a fine grid of locations for which to make a filled-in prediction map for the western United States.
US.prediction.locs
US.prediction.locs
A matrix with two columns:
Longitude of the prediction grid.
Latitude of the prediction grid.
A data set containing the annual precipitation for 1270 locations in the western United States.
USprecip97
USprecip97
A data frame with the following variables:
Longitude of the monitoring site.
Latitude of the monitoring site.
Annual precipitation for the monitoring site, in millimeters.
Annual precipitation for the monitoring site, in log millimeters.
http://www.image.ucar.edu/GSP/Data/US.monthly.met/