Title: | Spatial Misalignment: Interpolation, Linkage, and Estimation |
Version: | 1.0.5 |
Date: | 2024-06-13 |
Description: | Provides functions to estimate, predict and interpolate areal data. For estimation and prediction we assume areal data is an average of an underlying continuous spatial process as in Moraga et al. (2017) <doi:10.1016/j.spasta.2017.04.006>, Johnson et al. (2020) <doi:10.1186/s12942-020-00200-w>, and Wilson and Wakefield (2020) <doi:10.1093/biostatistics/kxy041>. The interpolation methodology is (mostly) based on Goodchild and Lam (1980, ISSN:01652273). |
License: | GPL-3 |
Encoding: | UTF-8 |
LazyData: | true |
RoxygenNote: | 7.3.1 |
SystemRequirements: | GDAL (>= 2.0.1), GEOS (>= 3.4.0), PROJ (>= 4.8.0) |
LinkingTo: | Rcpp, RcppArmadillo |
Imports: | numDeriv, Rcpp, sf, mvtnorm, stats, parallel, Matrix |
Depends: | R (≥ 4.0) |
URL: | https://lcgodoy.me/smile/, https://github.com/lcgodoy/smile/ |
BugReports: | https://github.com/lcgodoy/smile/issues/ |
Suggests: | knitr, rmarkdown, ggplot2, graphics, spelling |
VignetteBuilder: | knitr |
Language: | en-US |
NeedsCompilation: | yes |
Packaged: | 2024-06-13 22:55:44 UTC; lcgodoy |
Author: | Lucas da Cunha Godoy
|
Maintainer: | Lucas da Cunha Godoy <lcgodoy@duck.com> |
Repository: | CRAN |
Date/Publication: | 2024-06-14 08:00:02 UTC |
Areal Interpolation
Description
This function estimates variables observed at a "source" region into a "target" region. "Source" and "target" regions represent two different ways to divide a city, for example. For more details, see https://lcgodoy.me/smile/articles/sai.html.
Usage
ai(source, target, vars)
ai_var(source, target, vars, vars_var, sc_vars = FALSE, var_method = "CS")
Arguments
source |
a |
target |
a |
vars |
a |
vars_var |
a scalar of type |
sc_vars |
boolean indicating whether |
var_method |
a |
Value
the target (of type sf
) with estimates of the variables
observed at the source data.
Internal use only
Description
aggregate.sf taken from sf
.
Usage
aggregate_aux(
x,
by,
FUN,
...,
do_union = TRUE,
simplify = TRUE,
join = st_intersects
)
Arguments
x |
internal use |
by |
internal use |
FUN |
internal use |
... |
internal use |
do_union |
internal use |
simplify |
internal use |
join |
internal use |
Mean of a (Cubic spline) covariance function (Internal use)
Description
This is an auxiliary function for internal use. It helps to numerically integrate a covariance function evaluated at a grid of points within a polygon and speed-up the computations.
Usage
aux_cs(dist, sigsq, phi)
Arguments
dist |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
Value
The mean of spher_cov(dist, sigsq, phi)
.
Mean of a (Gaussian) covariance function (Internal use)
Description
This is an auxiliary function for internal use. It helps to numerically integrate a covariance function evaluated at a grid of points within a polygon and speed-up the computations.
Usage
aux_gauss(dist, sigsq, phi)
Arguments
dist |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
Value
The mean of gauss_cov(dist, sigsq, phi)
.
Mean of a Generalized Wendland covariance function (Internal use)
Description
This is an auxiliary function for internal use. It helps to numerically integrate a covariance function evaluated at a grid of points within a polygon and speed-up the computations.
Usage
aux_gw(dist, sigsq, phi, kappa, mu)
Arguments
dist |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
kappa |
|
mu |
|
Value
The mean of mat_cov(dist, sigsq, phi, kappa)
.
Internal use only
Description
Internal use only
Usage
single_dists(mat_list)
single_dists_tr(mat_list, tr_vec, tr_inp)
mult_dists(mat_list1, mat_list2, return_single)
pred_cdist(mat_list, pred_mat)
get_grid_list(x_to_list, by)
dist_from_grids(y_grid, by)
dist_from_grids_tr(y_grid, by, tr_vec, tr_inp)
mult_dist_from_grids(y_grid, x_grid, by)
Arguments
mat_list |
internal use |
tr_vec |
index for distance truncation |
tr_inp |
truncation input |
mat_list1 |
internal use |
mat_list2 |
internal use |
return_single |
internal use |
pred_mat |
internal use |
x_to_list |
internal use |
by |
internal use |
y_grid |
internal use |
x_grid |
internal use |
Mean of a (Matern) covariance function (Internal use)
Description
This is an auxiliary function for internal use. It helps to numerically integrate a covariance function evaluated at a grid of points within a polygon and speed-up the computations.
Usage
aux_matern(dist, sigsq, phi, nu)
Arguments
dist |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
nu |
the |
Value
The mean of mat_cov(dist, sigsq, phi, nu)
.
Mean of a (Powered Exponential) covariance function (Internal use)
Description
This is an auxiliary function for internal use. It helps to numerically integrate a covariance function evaluated at a grid of points within a polygon and speed-up the computations.
Usage
aux_pexp(dist, sigsq, phi, nu)
Arguments
dist |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
nu |
the |
Value
The mean of pexp_cov(dist, sigsq, phi, nu)
.
Mean of a (Spherical) covariance function (Internal use)
Description
This is an auxiliary function for internal use. It helps to numerically integrate a covariance function evaluated at a grid of points within a polygon and speed-up the computations.
Usage
aux_spher(dist, sigsq, phi)
Arguments
dist |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
Value
The mean of spher_cov(dist, sigsq, phi)
.
Mean of a (Matern - Wendland-1 tapper) covariance function (Internal use)
Description
This is an auxiliary function for internal use. It helps to numerically integrate a covariance function evaluated at a grid of points within a polygon and speed-up the computations.
Usage
aux_tapmat(dist, sigsq, phi, nu, theta)
Arguments
dist |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
nu |
smoothness parameter |
theta |
the |
Value
The mean of mat_cov(dist, sigsq, phi, nu)
.
Cubic spline covariance function for a polygons.
Description
Computing the Spherical covariance function between polygons.
Usage
comp_cs_cov(cross_dists, n, n2, sigsq, phi)
Arguments
cross_dists |
a |
n |
an integer representing number of polygons (note that, this is
different than the size of the list |
n2 |
usually, equal to |
sigsq |
the |
phi |
the |
Value
The spherical covariance matrix associated with a set of polygons.
See Also
single_exp
, single_matern
,
mat_cov
Gaussian covariance function for a polygons.
Description
Computing the Gaussian covariance function between polygons.
Usage
comp_gauss_cov(cross_dists, n, n2, sigsq, phi)
Arguments
cross_dists |
a |
n |
an integer representing number of polygons (note that, this is
different than the size of the list |
n2 |
usually, equal to |
sigsq |
the |
phi |
the |
Value
The gaussian covariance matrix associated with a set of polygons.
See Also
single_exp
, single_matern
,
mat_cov
Generalized Wendland covariance function for a polygons.
Description
Computing the Matern covariance function between polygons.
Usage
comp_gw_cov(cross_dists, n, n2, sigsq, phi, kappa, mu)
Arguments
cross_dists |
a |
n |
an integer representing number of polygons (note that, this is
different than the size of the list |
n2 |
usually, equal to |
sigsq |
the |
phi |
the |
kappa |
|
mu |
|
Value
The wendland-1 covariance matrix associated with a set of polygons.
Matern covariance function for a polygons.
Description
Computing the Matern covariance function between polygons.
Usage
comp_mat_cov(cross_dists, n, n2, sigsq, phi, nu)
Arguments
cross_dists |
a |
n |
an integer representing number of polygons (note that, this is
different than the size of the list |
n2 |
usually, equal to |
sigsq |
the |
phi |
the |
nu |
the |
Value
The matern covariance matrix associated with a set of polygons.
See Also
single_exp
, single_matern
,
mat_cov
Powered Exponential covariance function for a polygons.
Description
Computing the Powered Exponential covariance function between polygons.
Usage
comp_pexp_cov(cross_dists, n, n2, sigsq, phi, nu)
Arguments
cross_dists |
a |
n |
an integer representing number of polygons (note that, this is
different than the size of the list |
n2 |
usually, equal to |
sigsq |
the |
phi |
the |
nu |
the |
Value
The powered exponential covariance matrix associated with a set of polygons.
See Also
single_exp
, single_matern
,
mat_cov
Spherical covariance function for a polygons.
Description
Computing the Spherical covariance function between polygons.
Usage
comp_spher_cov(cross_dists, n, n2, sigsq, phi)
Arguments
cross_dists |
a |
n |
an integer representing number of polygons (note that, this is
different than the size of the list |
n2 |
usually, equal to |
sigsq |
the |
phi |
the |
Value
The spherical covariance matrix associated with a set of polygons.
See Also
single_exp
, single_matern
,
mat_cov
Wendland-1 covariance function for a polygons.
Description
Computing the Matern covariance function between polygons.
Usage
comp_tapmat_cov(cross_dists, n, n2, sigsq, phi, nu, theta)
Arguments
cross_dists |
a |
n |
an integer representing number of polygons (note that, this is
different than the size of the list |
n2 |
usually, equal to |
sigsq |
the |
phi |
the |
nu |
the smoothness parameter |
theta |
the taper distance. |
Value
The wendland-1 covariance matrix associated with a set of polygons.
Pairwise distances between matrices
Description
Internal use.
Usage
crossdist(m1, m2)
Arguments
m1 |
a matrix representing a grid of points within a polygon. |
m2 |
a matrix representing a grid of points within a polygon. |
Computing the Cubic spline covariance function for a single distance measure.
Description
Computing the Cubic spline covariance function for a single distance measure.
Usage
cs_cov(dists, sigsq, phi)
Arguments
dists |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
Value
The Spherical covariance function (for a stationary and
isotropic process) associated with the provided distances (dists
)
and the given set of parameters.
Creating a distance matrix
Description
Internal use. For now it only supports euclidean distance. (May import parallelDist in the future).
Usage
distmat(my_mat)
Arguments
my_mat |
a matrix representing a grid of points. |
MLEs for fixed V.
Description
MLEs for fixed V.
Usage
est_mle(y, Vinv)
Arguments
y |
response variable. |
Vinv |
inverse of |
Details
internal use.
Find phi parameter for the Exponential spatial auto-correlation function
Description
Function designed to find the phi parameter such that the
correlation between points within a given distance d
is at most a
given value.
Usage
find_phi(
d,
nu,
kappa,
mu2,
family = "matern",
range = c(1e-04, 1000),
cut = 0.05
)
Arguments
d |
maximum distance for spatial dependence equal to |
nu |
smoothness parameter associated with the Matern cov. function. |
kappa |
one of the smoothness parameters associated with the Generalized Wendland covariance function |
mu2 |
one of the smoothness parameters associated with the Generalized Wendland covariance function |
family |
covariance function family, the options are |
range |
Minimum and maximum distance to be considered. The default is
|
cut |
desired spatial correlation at a distance |
Value
a numeric
value indicating the range parameter such that the
spatial correlation between two points at distance d
is
cut
.
Fitting an underlying continuous process to areal data
Description
Fitting an underlying continuous process to areal data
Usage
fit_spm(x, ...)
## S3 method for class 'spm'
fit_spm(
x,
model,
theta_st,
nu = NULL,
tr = NULL,
kappa = 1,
mu2 = 1.5,
apply_exp = FALSE,
opt_method = "Nelder-Mead",
control_opt = list(),
comp_hess = TRUE,
...
)
fit_spm2(
x,
model,
nu,
tr,
kappa = 1,
mu2 = 1.5,
comp_hess = TRUE,
phi_min,
phi_max,
nphi = 10,
cores = getOption("mc.cores", 1L)
)
Arguments
x |
an object of type |
... |
additional parameters, either passed to stats::optim. |
model |
a |
theta_st |
a |
nu |
a |
tr |
tapper range |
kappa |
|
mu2 |
the smoothness parameter |
apply_exp |
a |
opt_method |
a |
control_opt |
a named |
comp_hess |
a |
phi_min |
a |
phi_max |
a |
nphi |
a |
cores |
a |
Details
This function uses the stats::optim function optimization
algorithms to find the Maximum Likelihood estimators, and their standard
errors, from a model adapted from. The function allows the user to input
the control parameters from the stats::optim function through the argument
control_opt
, which is a named list. Additionally, the one can
input lower and upper boundaries for the optimization problem, as well
as the preferred optimization algorithm (as long as it is available for
stats::optim). The preferred algorithm is selected by the argument
opt_method
. In addition to the control of the optimization, the
user can select a covariance function among the following: Matern,
Exponential, Powered Exponential, Gaussian, and Spherical. The parameter
apply_exp
is a logical
scalar such that, if set to
TRUE
, the \exp
function is applied to the nonnegative
parameters, allowing the optimization algorithm to search for all the
parameters over the real numbers.
The model assumes \deqn{Y(\mathbf{s}) = \mu + S(\mathbf{s})} at the point level. Where \eqn{S ~ GP(0, \sigma^2 C(\lVert \mathbf{s} - \mathbf{s}_2 \rVert; \theta))}. Further, the observed data is supposed to be \eqn{Y(B) = \lvert B \rvert^{-1} \int_{B} Y(\mathbf{s}) \, \textrm{d} \mathbf{s}}.
Value
a spm_fit
object containing the information about the
estimation of the model parameters.
Examples
data(liv_lsoa) ## loading the LSOA data
msoa_spm <- sf_to_spm(sf_obj = liv_msoa, n_pts = 500,
type = "regular", by_polygon = FALSE,
poly_ids = "msoa11cd",
var_ids = "leb_est")
## fitting model
theta_st_msoa <- c("phi" = 1) # initial value for the range parameter
fit_msoa <-
fit_spm(x = msoa_spm,
theta_st = theta_st_msoa,
model = "matern",
nu = .5,
apply_exp = TRUE,
opt_method = "L-BFGS-B",
control = list(maxit = 500))
AIC(fit_msoa)
summary_spm_fit(fit_msoa, sig = .05)
Computing the Gaussian covariance function for a single distance measure.
Description
Computing the Gaussian covariance function for a single distance measure.
Usage
gauss_cov(dists, sigsq, phi)
Arguments
dists |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
Value
The Gaussian covariance function (for a stationary and
isotropic process) associated with the provided distances (dists
)
and the given set of parameters.
See Also
Akaike's (and Bayesian) An Information Criterion for spm_fit
objects.
Description
Akaike's (and Bayesian) An Information Criterion for spm_fit
objects.
Usage
## S3 method for class 'spm_fit'
AIC(object, ..., k = 2)
## S3 method for class 'spm_fit'
BIC(object, ...)
Arguments
object |
a |
... |
optionally more fitted model objects. |
k |
|
Value
a numeric
scalar corresponding to the goodness of fit
measure.
Generalized Wendland covariance function for a given distance matrix.
Description
Computing the Matern covariance function for a matrix of distances.
Usage
gw_cov(dists, sigsq, phi, kappa, mu)
Arguments
dists |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
|
kappa |
|
mu |
|
Value
The GW (isotropic) covariance function associated with the provided
distances (dists
) and the given set of parameters.
Liverpool Lower Super Output Area.
Description
A dataset containing containing the LSOA's for Liverpool along with estimates for Index of Multiple Deprivation. Data taken from Johnson et al. 2020
Usage
liv_lsoa
Format
A sf
data frame with 298 rows and 6 variables:
- lsoa11cd
LSOA code
- lsoa11cd
LSOA name
- male
Male population
- female
Female population
- imdscore
Index of Multiple Deprivation
- area
LMSOA area, in
km^2
Details
The data was projected to EPSG 27700 and units changed to km
Source
https://ij-healthgeographics.biomedcentral.com/articles/10.1186/s12942-020-00200-w
Liverpool Middle Super Output Area.
Description
A dataset containing containing the MSOA's for Liverpool along with estimates for Life Expectancy at Birth. Data taken from Johnson et al. 2020
Usage
liv_msoa
Format
A sf
data frame with 61 rows and 4 variables:
- msoa11cd
MSOA code
- msoa11cd
MSOA name
- lev_est
Estimated life expectancy at birth, in years
- area
MSOA area, in
km^2
Details
The data was projected to EPSG 27700 and units changed to km
Source
https://ij-healthgeographics.biomedcentral.com/articles/10.1186/s12942-020-00200-w
Matern covariance function for a given distance matrix.
Description
Computing the Matern covariance function for a matrix of distances.
Usage
mat_cov(dists, sigsq, phi, nu)
Arguments
dists |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
nu |
the |
Value
The matern covariance function (for a stationary and isotropic
process) associated with the provided distances (dists
) and the
given set of parameters.
See Also
Calculates the (global) Moran's I
Description
Calculates the (global) Moran's I
Usage
morans_i(sf_dt, variable)
Arguments
sf_dt |
a |
variable |
a |
Value
a numeric
scalar.
Nova Lima census tracts
Description
A dataset containing containing the census tracts for the city of Nova Lima in Minas Gerais - Brazil.
Usage
nl_ct
Format
A sf
data frame with 113 rows and 14 variables:
- cd_setor
unique identifier
- hh_density
average household density
- var_hhd
variance of the household density
- avg_income
average income per household
- var_income
variance of the income per household
- pop
population in the census tract
- avg_age
average age of the inhabitants in the census tract
- var_age
variance of the variable age in the census tract
- prop_women
proportion of women
- prop_elder
proportion of people with 55 years of age or older
- illit_rate
illiteracy rate
- prop_white
proportion of self-declared white people
- prop_black
proportion of self-declared black people
- prop_native
proportion of self-declared native people
Details
The data is project using the SIRGAS 2000.
Powered Exponential covariance function for a given distance matrix.
Description
Computing the Matern covariance function for a matrix of distances.
Usage
pexp_cov(dists, sigsq, phi, nu)
Arguments
dists |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
nu |
the |
Value
The powered exponential covariance function (for a stationary and
isotropic process) associated with the provided distances (dists
)
and the given set of parameters.
See Also
Prediction over the same or a different set of regions (or points).
Description
Realizes predictions that can be useful when researchers are interested in predict a variable observed in one political division of a city (or state) on another division of the same region.
Usage
predict_spm(x, ...)
## S3 method for class 'spm_fit'
predict_spm(x, .aggregate = TRUE, ...)
## S3 method for class 'sf'
predict_spm(x, spm_obj, n_pts, type, outer_poly = NULL, id_var, ...)
Arguments
x |
a |
... |
additional parameters |
.aggregate |
|
spm_obj |
an object of either class |
n_pts |
a |
type |
|
outer_poly |
(object) |
id_var |
if |
Value
a list
of size 4 belonging to the class spm_pred
. This
list contains the predicted values and the mean and covariance matrix
associated with the conditional distribution used to compute the
predictions.
Examples
data(liv_lsoa) ## loading the LSOA data
data(liv_msoa) ## loading the MSOA data
msoa_spm <- sf_to_spm(sf_obj = liv_msoa, n_pts = 500,
type = "regular", by_polygon = FALSE,
poly_ids = "msoa11cd",
var_ids = "leb_est")
## fitting model
theta_st_msoa <- c("phi" = 1) # initial value for the range parameter
fit_msoa <-
fit_spm(x = msoa_spm,
theta_st = theta_st_msoa,
model = "matern",
nu = .5,
apply_exp = TRUE,
opt_method = "L-BFGS-B",
control = list(maxit = 500))
pred_lsoa <- predict_spm(x = liv_lsoa, spm_obj = fit_msoa, id_var = "lsoa11cd")
single sf
to spm
Description
Transforming a sf
into a spm
object (Internal use)
Usage
single_sf_to_spm(
sf_obj,
n_pts,
type = "regular",
by_polygon = FALSE,
poly_ids = NULL,
var_ids = NULL,
trunc_d = NULL
)
sf_to_spm(
sf_obj,
n_pts,
type = "regular",
by_polygon = FALSE,
poly_ids = NULL,
var_ids = NULL,
trunc_d = NULL
)
Arguments
sf_obj |
a |
n_pts |
a |
type |
a |
by_polygon |
a |
poly_ids |
a |
var_ids |
a scalar or vector of type |
trunc_d |
truncation distance for grid points. Consider using half of the maximum distance between polygons |
Value
a named list
of size 6 belonging to the class
spm
. This list stores all the objects necessary to fit models
using the fit_spm
.
Examples
data(liv_lsoa) # loading the LSOA data
msoa_spm <- sf_to_spm(sf_obj = liv_msoa, n_pts = 1000,
type = "regular", by_polygon = FALSE,
poly_ids = "msoa11cd",
var_ids = "leb_est")
Evaluate log-lik
Description
Evaluate the log-likelihood for a given set of parameters
Usage
singl_ll_nn_hess(
theta,
.dt,
dists,
npix,
model,
nu = NULL,
tr = NULL,
kappa = 1,
mu2 = 1.5,
apply_exp = FALSE
)
Arguments
theta |
a |
.dt |
a |
dists |
a |
npix |
a |
model |
a |
nu |
|
tr |
taper range |
kappa |
|
mu2 |
the smoothness parameter |
apply_exp |
a |
Details
Internal use.
Value
a scalar representing -log.lik
.
Evaluate log-lik
Description
Evaluate the log-likelihood for a given set of parameters
Usage
singl_log_lik(
theta,
.dt,
dists,
npix,
model,
nu = NULL,
tr = NULL,
kappa = 1,
mu2 = 1.5,
apply_exp = FALSE
)
Arguments
theta |
a |
.dt |
a |
dists |
a |
npix |
a |
model |
a |
nu |
|
tr |
|
kappa |
|
mu2 |
the smoothness parameter |
apply_exp |
a |
Details
Internal use.
Value
a scalar representing -log.lik
.
Evaluate log-lik
Description
Evaluate the log-likelihood for a given set of parameters - No nugget + profile likelihood
Usage
singl_log_lik_nn(
theta,
.dt,
dists,
npix,
model,
nu = NULL,
tr = NULL,
kappa = 1,
mu2 = 1.5,
apply_exp = FALSE
)
Arguments
theta |
a scalar for the |
.dt |
a |
dists |
a |
npix |
a |
model |
a |
nu |
|
kappa |
|
mu2 |
the smoothness parameter |
apply_exp |
a |
Details
Internal use.
Value
a scalar representing -log.lik
.
Evaluate log-lik
Description
Evaluate the log-likelihood for a given set of parameters - New parametrization + profile likelihood
Usage
singl_log_plik(
theta,
.dt,
dists,
npix,
model,
nu = NULL,
tr = NULL,
kappa = 1,
mu2 = 1.5,
apply_exp = FALSE
)
Arguments
theta |
a |
.dt |
a |
dists |
a |
npix |
a |
model |
a |
nu |
|
tr |
|
kappa |
|
mu2 |
the smoothness parameter |
apply_exp |
a |
Details
Internal use.
Value
a scalar representing -log.lik
.
Cubic spline covariance function (scalar)
Description
Computing the Spherical covariance function for a single distance measure.
Usage
single_cs(d, sigsq, phi)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
Value
a scalar representing the (gaussian) covariance between two
observations d
apart of each other.
See Also
single_exp
, single_matern
,
single_matern3
, single_matern5
,
mat_cov
Exponential covariance function (scalar)
Description
Computing the Exponential covariance function for a single distance measure.
Usage
single_exp(d, sigsq, phi)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
Value
a scalar representing the (exponential) covariance between two
observations d
apart of each other.
See Also
single_pexp
, single_matern
,
single_matern3
, single_matern5
,
mat_cov
Gaussian covariance function (scalar)
Description
Computing the Gaussian covariance function for a single distance measure.
Usage
single_gauss(d, sigsq, phi)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
Value
a scalar representing the (gaussian) covariance between two
observations d
apart of each other.
See Also
single_exp
, single_matern
,
single_matern3
, single_matern5
,
mat_cov
Matern Generalized Wendland (GW) covariance function with kappa = 0 (scalar - generic)
Description
adapted from Bevilacqua et al. 2019.
Usage
single_gw(d, sigsq, phi, kappa, mu)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
kappa |
|
mu |
a parameter that controls the smoothness of the covariance
function. Note that, |
Value
a scalar representing the GW covariance between two
observations d
apart of each other.
Matern Generalized Wendland (GW) covariance function with kappa = 0 (scalar - generic)
Description
adapted from Bevilacqua et al. 2019.
Usage
single_gw0(d, sigsq, phi, mu)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
mu |
a parameter that controls the smoothness of the covariance
function. Note that, |
Value
a scalar representing the GW covariance between two
observations d
apart of each other.
Matern Generalized Wendland (GW) covariance function with kappa = 1 (scalar - generic)
Description
adapted from Bevilacqua et al. 2019.
Usage
single_gw1(d, sigsq, phi, mu)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
mu |
a parameter that controls the smoothness of the covariance
function. Note that, |
Value
a scalar representing the GW covariance between two
observations d
apart of each other.
Matern Generalized Wendland (GW) covariance function with kappa = 2 (scalar - generic)
Description
adapted from Bevilacqua et al. 2019.
Usage
single_gw2(d, sigsq, phi, mu)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
mu |
a parameter that controls the smoothness of the covariance
function. Note that, |
Value
a scalar representing the GW covariance between two
observations d
apart of each other.
Matern Generalized Wendland (GW) covariance function with kappa = 3 (scalar - generic)
Description
adapted from Bevilacqua et al. 2019.
Usage
single_gw3(d, sigsq, phi, mu)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
mu |
a parameter that controls the smoothness of the covariance
function. Note that, |
Value
a scalar representing the GW covariance between two
observations d
apart of each other.
Matern covariance function (scalar - generic)
Description
Computing the Matern covariance function for a single distance
measure, adapted from geoR
.
Usage
single_matern(d, sigsq, phi, nu)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
nu |
the |
Value
a scalar representing the (matern) covariance between two
observations d
apart of each other.
See Also
single_matern
, single_matern5
single_exp
, mat_cov
Matern covariance function (scalar - nu = 3/2)
Description
Computing the Matern covariance function for a single distance
measure, with \nu = 3/2
.
Usage
single_matern3(d, sigsq, phi)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
Value
a scalar representing the (matern) covariance between two
observations d
apart of each other.
See Also
single_matern
, single_matern5
single_exp
, mat_cov
Matern covariance function (scalar - nu = 5/2)
Description
Computing the Matern covariance function for a single distance
measure, with \nu = 5/2
.
Usage
single_matern5(d, sigsq, phi)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
Value
a scalar representing the (matern) covariance between two
observations d
apart of each other.
See Also
single_matern
, single_matern3
single_exp
, mat_cov
Powered Exponential covariance function (scalar)
Description
Computing the Powered Exponential covariance function for a single distance measure.
Usage
single_pexp(d, sigsq, phi, nu)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
nu |
the |
Value
a scalar representing the (exponential) covariance between two
observations d
apart of each other.
See Also
single_exp
, single_matern
,
single_matern3
, single_matern5
,
mat_cov
Spherical covariance function (scalar)
Description
Computing the Spherical covariance function for a single distance measure.
Usage
single_spher(d, sigsq, phi)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
Value
a scalar representing the (gaussian) covariance between two
observations d
apart of each other.
See Also
single_exp
, single_matern
,
single_matern3
, single_matern5
,
mat_cov
Matern (tapered) covariance function (scalar - generic)
Description
Computing the Matern covariance function for a single distance
measure, adapted from geoR
using Wendland-1 as a tapper.
Usage
single_tapmat(d, sigsq, phi, nu, theta)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
nu |
the smoothness parameter |
theta |
the |
Value
a scalar representing the (tapered matern) covariance between two
observations d
apart of each other.
See Also
single_matern
, single_matern5
single_exp
, mat_cov
tapered Matern covariance function (scalar - nu = 1/2)
Description
Computing the Matern covariance function for a single distance
measure, with \nu = 1/2
.
Usage
single_tapmat1(d, sigsq, phi, theta)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
theta |
|
Value
a scalar representing the (tapered matern) covariance between two
observations d
apart of each other.
See Also
single_matern
, single_matern5
single_exp
, mat_cov
Tapered Matern covariance function (scalar - nu = 3/2)
Description
Computing the Matern covariance function for a single distance
measure, with \nu = 3/2
.
Usage
single_tapmat3(d, sigsq, phi, theta)
Arguments
d |
a scalar representing the distance on which it is desired to evaluate the covariance function. |
sigsq |
the |
phi |
the |
theta |
|
Value
a scalar representing the (tapered matern) covariance between two
observations d
apart of each other.
See Also
single_matern
, single_matern5
single_exp
, mat_cov
smile: Spatial MIsaLignment Estimation
Description
smile: Spatial MIsaLignment Estimation
Computing the Spherical covariance function for a single distance measure.
Description
Computing the Spherical covariance function for a single distance measure.
Usage
spher_cov(dists, sigsq, phi)
Arguments
dists |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
Value
The Spherical covariance function (for a stationary and
isotropic process) associated with the provided distances (dists
)
and the given set of parameters.
See Also
Remove holes from a sfc
POLYGON
Description
internal use. Taken from https://github.com/michaeldorman/nngeo/
Usage
st_remove_holes(x)
Arguments
x |
a |
Value
a sf
or sfc
polygon.
Summarizing spm_fit
Description
Provides a data.frame
with point estimates and
confidence intervals for the parameters of the model fitted using the
spm_fit
function.
Usage
summary_spm_fit(x, sig = 0.05)
Arguments
x |
a |
sig |
a real number between 0 and 1 indicating significance level to be used to compute the confidence intervals for the parameter estimates. |
Value
a data.frame
summarizing the parameters estimated by the
fit_spm
function.
Tapered Matern covariance function for a given distance matrix.
Description
Computing the tapered Matern covariance function for a matrix of distances.
Usage
tapmat_cov(dists, sigsq, phi, nu, theta)
Arguments
dists |
a numeric matrix representing the distance between spatial entities. |
sigsq |
the |
phi |
the |
nu |
smoothness parameter |
theta |
|
Value
The tapered matern covariance function (for a stationary and isotropic
process) associated with the provided distances (dists
) and the
given set of parameters.
Transform method for sf
objects
Description
Internal usage.
Usage
## S3 method for class 'sf'
transform(`_data`, ...)
Arguments
_data |
a |
... |
additional options |
Value
an sf
object.
Voronoi Data Linkage
Description
Reminder, have to create an example. This will be exported after we submit the paper for publication.
Usage
vdl(coords_sf, areal_sf, vars, buff)
Arguments
coords_sf |
|
areal_sf |
|
vars |
a |
buff |
scalar |
Value
a sf
object for the coords_sf
spatial data set.
Voronoi Data Linkage - Single variable and variance
Description
Reminder, have to create an example. This will be exported after we submit the paper for publication.
Usage
vdl_var(coords_sf, areal_sf, res_var, variance, var_method = "CS", buff)
Arguments
coords_sf |
|
areal_sf |
|
res_var |
a |
variance |
a |
var_method |
a |
buff |
scalar |
Value
a sf
object, containing the id_coords
variable and the
list_vars
for the coords_sf
spatial data set.
Voronoi Tessellation inside a polygon
Description
voronoi tessellation of a given a set of points inside a polygon. This is an internal use function.
Usage
vor_build(points_sf, poly_sf)
Arguments
points_sf |
|
poly_sf |
a |
Value
a sf
object containing the polygons associated with the
voronoi tessellation of points_sf
the polygon poly_sf
Building weight matrix W for Areal Interpolation
Description
internal use. W_{ij} = | A_i \, cap \, B_j |
.
Usage
w_col(source_unit, target)
build_w(source, target)
est_w(W, source_dt, target)
var_w(W, var_vec, target, method = "CS", rho_mi)
Arguments
source_unit |
a single |
target |
a |
source |
a |
W |
the weight matrix. |
source_dt |
a |
var_vec |
a |
method |
a |
rho_mi |
|
Value
A n \times m
numeric
matrix. Where n
is the
number of observations in the target and m
is the sample size in
the source dataset.