Title: | Fast Inference for Large Spatial Datasets using BRISC |
---|---|
Description: | Fits bootstrap with univariate spatial regression models using Bootstrap for Rapid Inference on Spatial Covariances (BRISC) for large datasets using nearest neighbor Gaussian processes detailed in Saha and Datta (2018) <doi:10.1002/sta4.184>. |
Authors: | Arkajyoti Saha [aut, cre], Abhirup Datta [aut], Jorge Nocedal [ctb], Naoaki Okazaki [ctb], Lukas M. Weber [ctb] |
Maintainer: | Arkajyoti Saha <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.0.6 |
Built: | 2024-11-02 06:38:49 UTC |
Source: | CRAN |
The function BRISC_bootstrap
performs bootstrap to provide confidence intervals for parameters of univariate spatial
regression models using outputs of BRISC_estimation
. The details of the bootstrap method can be found in BRISC
(Saha & Datta, 2018). The optimization is performed with C library of limited-memory BFGS libLBFGS: a library of
Limited-memory Broyden-Fletcher-Goldfarb-Shanno (L-BFGS), http://www.chokkan.org/software/liblbfgs/ (Naoaki Okazaki).
For user convenience the soure codes of the package libLBFGS are provided in the package. Some code blocks are borrowed
from the R package: spNNGP: Spatial Regression Models for Large Datasets using Nearest Neighbor Gaussian Processes
https://CRAN.R-project.org/package=spNNGP .
BRISC_bootstrap(BRISC_Out, n_boot = 100, h = 1, n_omp = 1, init = "Initial", verbose = TRUE, nugget_status = 1)
BRISC_bootstrap(BRISC_Out, n_boot = 100, h = 1, n_omp = 1, init = "Initial", verbose = TRUE, nugget_status = 1)
BRISC_Out |
an object of class |
n_boot |
number of bootstrap samples. Default value is 100. |
h |
number of core to be used in parallel computing setup for bootstrap samples. If |
n_omp |
number of threads to be used, value can be more than 1 if source code is compiled with OpenMP support. Default is 1. |
init |
keyword that specifies initialization scheme to be used. Supported keywords are: |
verbose |
if |
nugget_status |
if |
A list comprising of the following:
boot.Theta |
estimates of spatial covariance parameters corresponding to bootstrap samples. |
boot.Beta |
estimates of beta corresponding to bootstrap samples. |
confidence.interval |
confidence intervals corresponding to the parameters. |
boot.time |
time (in seconds) required to perform the bootstrapping after preprocessing data in |
Arkajyoti Saha [email protected],
Abhirup Datta [email protected]
Saha, A., & Datta, A. (2018). BRISC: bootstrap for rapid inference on spatial covariances. Stat, e184, DOI: 10.1002/sta4.184.
Okazaki N. libLBFGS: a library of Limited-memory Broyden-Fletcher-Goldfarb-Shanno
(L-BFGS),
http://www.chokkan.org/software/liblbfgs/ .
Andrew Finley, Abhirup Datta and Sudipto Banerjee (2017). spNNGP: Spatial Regression Models for Large Datasets using Nearest Neighbor Gaussian Processes. R package version 0.1.1. https://CRAN.R-project.org/package=spNNGP
rmvn <- function(n, mu = 0, V = matrix(1)){ p <- length(mu) if(any(is.na(match(dim(V),p)))) stop("Dimension not right!") D <- chol(V) t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p))) } set.seed(1) n <- 300 coords <- cbind(runif(n,0,1), runif(n,0,1)) beta <- c(1,5) x <- cbind(rnorm(n), rnorm(n)) sigma.sq = 1 phi = 5 tau.sq = 0.1 B <- as.matrix(beta) D <- as.matrix(dist(coords)) R <- exp(-phi*D) w <- rmvn(1, rep(0,n), sigma.sq*R) y <- rnorm(n, x%*%B + w, sqrt(tau.sq)) estimation_result <- BRISC_estimation(coords, y, x) bootstrap_result <- BRISC_bootstrap(estimation_result, n_boot = 10)
rmvn <- function(n, mu = 0, V = matrix(1)){ p <- length(mu) if(any(is.na(match(dim(V),p)))) stop("Dimension not right!") D <- chol(V) t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p))) } set.seed(1) n <- 300 coords <- cbind(runif(n,0,1), runif(n,0,1)) beta <- c(1,5) x <- cbind(rnorm(n), rnorm(n)) sigma.sq = 1 phi = 5 tau.sq = 0.1 B <- as.matrix(beta) D <- as.matrix(dist(coords)) R <- exp(-phi*D) w <- rmvn(1, rep(0,n), sigma.sq*R) y <- rnorm(n, x%*%B + w, sqrt(tau.sq)) estimation_result <- BRISC_estimation(coords, y, x) bootstrap_result <- BRISC_bootstrap(estimation_result, n_boot = 10)
The function BRISC_correlation
creates correlated data (known structure) using Nearest Neighbor
Gaussian Processes (NNGP). BRISC_correlation
uses the sparse Cholesky representation of Vecchia’s
likelihood developed in Datta et al., 2016. Some code blocks are borrowed from the R package: spNNGP:
Spatial Regression Models for Large Datasets using Nearest Neighbor Gaussian Processes
https://CRAN.R-project.org/package=spNNGP .
BRISC_correlation(coords, sim, sigma.sq = 1, tau.sq = 0, phi = 1, nu = 1.5, n.neighbors = NULL, n_omp = 1, cov.model = "exponential", search.type = "tree", stabilization = NULL, verbose = TRUE, tol = 12)
BRISC_correlation(coords, sim, sigma.sq = 1, tau.sq = 0, phi = 1, nu = 1.5, n.neighbors = NULL, n_omp = 1, cov.model = "exponential", search.type = "tree", stabilization = NULL, verbose = TRUE, tol = 12)
coords |
an |
sim |
an |
sigma.sq |
value of sigma square. Default value is 1. |
tau.sq |
value of tau square. Default value is 0.1. |
phi |
value of phi. Default value is 1. |
nu |
value of nu, only required for matern covariance model. Default value is 1.5. |
n.neighbors |
number of neighbors used in the NNGP. Default value is |
n_omp |
number of threads to be used, value can be more than 1 if source code is compiled with OpenMP support. Default is 1. |
cov.model |
keyword that specifies the covariance function to be used in modelling the spatial dependence structure
among the observations. Supported keywords are: |
search.type |
keyword that specifies type of nearest neighbor search algorithm to be used. Supported keywords are:
|
stabilization |
when we use a very smooth covarince model (lower values of phi for spherical and Gaussian
covariance and low phi and high nu for Matern covarinace) in absence of a non-negligble nugget, the correlation process may fail
due to computational instability. If |
verbose |
if |
tol |
the input observation coordinates are rounded to this many places after the decimal. The default value is 12. |
Denote be the input
sim
. Let be the precision matrix associated with the covariance model determined by the
and model parameters. Then
BRISC_correlation
calculates , where
is given as follows:
where, is a sparse approximation of the cholesky factor
of the precision matrix
, obtained from NNGP.
A list comprising of the following:
coords |
the matrix |
n.neighbors |
the used value of |
cov.model |
the used covariance model. |
Theta |
parameters of covarinace model; accounts for |
input.data |
the matrix |
output.data |
the output matrix |
time |
time (in seconds) required after preprocessing data in R, |
Arkajyoti Saha [email protected],
Abhirup Datta [email protected]
Datta, A., S. Banerjee, A.O. Finley, and A.E. Gelfand. (2016) Hierarchical Nearest-Neighbor Gaussian process models for large geostatistical datasets. Journal of the American Statistical Association, 111:800-812.
Andrew Finley, Abhirup Datta and Sudipto Banerjee (2017). spNNGP: Spatial Regression Models for Large Datasets using Nearest Neighbor Gaussian Processes. R package version 0.1.1. https://CRAN.R-project.org/package=spNNGP
set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) sigma.sq = 1 phi = 1 set.seed(1) sim <- matrix(rnorm(3*n),n, 3) correlation_result <- BRISC_correlation(coords, sigma.sq = sigma.sq, phi = phi, sim = sim)
set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) sigma.sq = 1 phi = 1 set.seed(1) sim <- matrix(rnorm(3*n),n, 3) correlation_result <- BRISC_correlation(coords, sigma.sq = sigma.sq, phi = phi, sim = sim)
The function BRISC_decorrelation
is used to decorrelate data (known structure) using Nearest Neighbor
Gaussian Processes (NNGP). BRISC_decorrelation
uses the sparse Cholesky representation of Vecchia’s
likelihood developed in Datta et al., 2016. Some code blocks are borrowed from the R package: spNNGP:
Spatial Regression Models for Large Datasets using Nearest Neighbor Gaussian Processes
https://CRAN.R-project.org/package=spNNGP .
BRISC_decorrelation(coords, sim, sigma.sq = 1, tau.sq = 0, phi = 1, nu = 1.5, n.neighbors = NULL, n_omp = 1, cov.model = "exponential", search.type = "tree", stabilization = NULL, verbose = TRUE, tol = 12)
BRISC_decorrelation(coords, sim, sigma.sq = 1, tau.sq = 0, phi = 1, nu = 1.5, n.neighbors = NULL, n_omp = 1, cov.model = "exponential", search.type = "tree", stabilization = NULL, verbose = TRUE, tol = 12)
coords |
an |
sim |
an |
sigma.sq |
value of sigma square. Default value is 1. |
tau.sq |
value of tau square. Default value is 0.1. |
phi |
value of phi. Default value is 1. |
nu |
value of nu, only required for Matern covariance model. Default value is 1.5. |
n.neighbors |
number of neighbors used in the NNGP. Default value is |
n_omp |
number of threads to be used, value can be more than 1 if source code is compiled with OpenMP support. Default is 1. |
cov.model |
keyword that specifies the covariance function to be used in modelling the spatial dependence structure
among the observations. Supported keywords are: |
search.type |
keyword that specifies type of nearest neighbor search algorithm to be used. Supported keywords are:
|
stabilization |
when the correlated data are generated from a very smooth covarince model (lower values of phi for spherical and Gaussian
covariance and low phi and high nu for Matern covarinace), the decorrelation process may fail due to computational
instability. If |
verbose |
if |
tol |
the input observation coordinates are rounded to this many places after the decimal. The default value is 12. |
Denote be the input
sim
. Let be the covariance matrix associated with the covariance model determined by the
and model parameters. Then
BRISC_decorrelation
calculates , where
is given as follows:
where, is a sparse approximation of the cholesky factor
of the precision matrix
, obtained from NNGP.
A list comprising of the following:
coords |
the matrix |
n.neighbors |
the used value of |
cov.model |
the used covariance model. |
Theta |
parameters of covarinace model; accounts for |
input.data |
if |
output.data |
the output matrix |
time |
time (in seconds) required after preprocessing data in |
Arkajyoti Saha [email protected],
Abhirup Datta [email protected]
Datta, A., S. Banerjee, A.O. Finley, and A.E. Gelfand. (2016) Hierarchical Nearest-Neighbor Gaussian process models for large geostatistical datasets. Journal of the American Statistical Association, 111:800-812.
Andrew Finley, Abhirup Datta and Sudipto Banerjee (2017). spNNGP: Spatial Regression Models for Large Datasets using Nearest Neighbor Gaussian Processes. R package version 0.1.1. https://CRAN.R-project.org/package=spNNGP
rmvn <- function(n, mu = 0, V = matrix(1)){ p <- length(mu) if(any(is.na(match(dim(V),p)))) stop("Dimension not right!") D <- chol(V) t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p))) } set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) sigma.sq = 1 phi = 1 set.seed(1) D <- as.matrix(dist(coords)) R <- exp(-phi*D) sim <- rmvn(3, rep(0,n), sigma.sq*R) decorrelation_result <- BRISC_decorrelation(coords, sim = sim)
rmvn <- function(n, mu = 0, V = matrix(1)){ p <- length(mu) if(any(is.na(match(dim(V),p)))) stop("Dimension not right!") D <- chol(V) t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p))) } set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) sigma.sq = 1 phi = 1 set.seed(1) D <- as.matrix(dist(coords)) R <- exp(-phi*D) sim <- rmvn(3, rep(0,n), sigma.sq*R) decorrelation_result <- BRISC_decorrelation(coords, sim = sim)
The function BRISC_estimation
fits univariate spatial regression models for large
spatial data using Vecchia's approximate likelihood (Vecchia, 1988). BRISC_estimation
uses the sparse Cholesky representation of Vecchia’s likelihood developed in Datta et al., 2016.
The Maximum Likelihood Estimates (MLE) of the parameters are used later for calculating the
confidence interval via the BRISC_bootstrap
(BRISC, Saha & Datta, 2018).
We recommend
using BRISC_estimation
followed by BRISC_bootstrap
to obtain the confidence
intervals for the model parameters.
The optimization is performed with C library of limited-memory BFGS
libLBFGS: a library of Limited-memory Broyden-Fletcher-Goldfarb-Shanno (L-BFGS),
http://www.chokkan.org/software/liblbfgs/ (Naoaki Okazaki). For user convenience
the source codes of the package libLBFGS are provided in the package. The code for
the coordinate ordering method, approximate Maximum Minimum Distance (Guinness, 2018)
is available in
https://github.com/joeguinness/gp_reorder/tree/master/R and is adopted
with minor modification. Some code blocks are borrowed from the R package: spNNGP:
Spatial Regression Models for Large Datasets using Nearest Neighbor Gaussian Processes
https://CRAN.R-project.org/package=spNNGP .
BRISC_estimation(coords, y, x = NULL, sigma.sq = 1, tau.sq = 0.1, phi = 1, nu = 1.5, n.neighbors = 15, n_omp = 1, order = "Sum_coords", cov.model = "exponential", search.type = "tree", stabilization = NULL, pred.stabilization = 1e-8, verbose = TRUE, eps = 2e-05, nugget_status = 1, ordering = NULL, neighbor = NULL, tol = 12)
BRISC_estimation(coords, y, x = NULL, sigma.sq = 1, tau.sq = 0.1, phi = 1, nu = 1.5, n.neighbors = 15, n_omp = 1, order = "Sum_coords", cov.model = "exponential", search.type = "tree", stabilization = NULL, pred.stabilization = 1e-8, verbose = TRUE, eps = 2e-05, nugget_status = 1, ordering = NULL, neighbor = NULL, tol = 12)
coords |
an |
y |
an |
x |
an |
sigma.sq |
starting value of sigma square. Default value is 1. |
tau.sq |
starting value of tau square. Default value is 0.1. |
phi |
starting value of phi. Default value is 1. |
nu |
starting value of nu, only required for matern covariance model. Default value is 1.5. |
n.neighbors |
number of neighbors used in the NNGP. Default value is 15. |
n_omp |
number of threads to be used, value can be more than 1 if source code is compiled with OpenMP support. Default is 1. |
order |
keyword that specifies the ordering scheme to be used in ordering the observations. Supported keywords are:
|
cov.model |
keyword that specifies the covariance function to be used in modelling the spatial dependence structure
among the observations. Supported keywords are: |
search.type |
keyword that specifies type of nearest neighbor search algorithm to be used. Supported keywords are:
|
stabilization |
when the spatial errors are generated from a very smooth covarince model (lower values of phi for spherical and Gaussian
covariance and low phi and high nu for Matern covarinace), the estimation process may fail due to computational
instability. If |
pred.stabilization |
if not |
verbose |
if |
eps |
tolerance to be used in centred finite difference approximation of derivatives. Default value is 2e-05. |
nugget_status |
if |
ordering |
if |
neighbor |
if |
tol |
the input observation coordinates, response and the covariates are rounded to this many places after the decimal. The default value is 12. |
An object of class BRISC_Out
, which is a list comprising:
ord |
the vector of indices used to order data necessary for fitting the NNGP model. |
coords |
the matrix |
y |
If |
X |
the matrix |
n.neighbors |
the used value of |
cov.model |
the used covariance model. |
eps |
value of used |
init |
initial values of the parameters of the covariance model; |
Beta |
estimate of beta. |
Theta |
estimate of parameters of covarinace model. |
log_likelihood |
value of Vecchia’s approximate log likelihood with parameter estimates. |
estimation.time |
time (in seconds) required to perform the model fitting after ordering and preprocessing data in |
BRISC_Object |
object required for bootstrap and prediction. |
Arkajyoti Saha [email protected],
Abhirup Datta [email protected]
Saha, A., & Datta, A. (2018). BRISC: bootstrap for rapid inference on spatial covariances. Stat, e184, DOI: 10.1002/sta4.184.
Datta, A., S. Banerjee, A.O. Finley, and A.E. Gelfand. (2016) Hierarchical Nearest-Neighbor Gaussian process models for large geostatistical datasets. Journal of the American Statistical Association, 111:800-812.
Guinness, J. (2018) Permutation and Grouping Methods for
Sharpening Gaussian Process Approximations, Technometrics,
DOI: 10.1080/00401706.2018.1437476,
https://github.com/joeguinness/gp_reorder/tree/master/R .
Vecchia, A. V. (1988) Estimation and model identification for continuous spatial processes. Journal of the Royal Statistical Society. Series B (Methodological), 297-312.
Okazaki N. libLBFGS: a library of Limited-memory Broyden-Fletcher-Goldfarb-Shanno
(L-BFGS),
http://www.chokkan.org/software/liblbfgs/ .
Andrew Finley, Abhirup Datta and Sudipto Banerjee (2020). spNNGP: Spatial Regression Models for Large Datasets using Nearest Neighbor Gaussian Processes. R package version 0.1.4. https://CRAN.R-project.org/package=spNNGP
Ra, S. W., & Kim, J. K. (1993). A fast mean-distance-ordered partial codebook search algorithm for image vector quantization. IEEE Transactions on Circuits and Systems II: Analog and Digital Signal Processing, 40(9), 576-579.
rmvn <- function(n, mu = 0, V = matrix(1)){ p <- length(mu) if(any(is.na(match(dim(V),p)))) stop("Dimension not right!") D <- chol(V) t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p))) } set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) beta <- c(1,5) x <- cbind(rnorm(n), rnorm(n)) sigma.sq = 1 phi = 1 tau.sq = 0.1 B <- as.matrix(beta) D <- as.matrix(dist(coords)) R <- exp(-phi*D) w <- rmvn(1, rep(0,n), sigma.sq*R) y <- rnorm(n, x%*%B + w, sqrt(tau.sq)) estimation_result <- BRISC_estimation(coords, y, x) estimation_result$Theta ##Estimates of covariance model parameters. estimation_result$Beta ##Estimates of Beta
rmvn <- function(n, mu = 0, V = matrix(1)){ p <- length(mu) if(any(is.na(match(dim(V),p)))) stop("Dimension not right!") D <- chol(V) t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p))) } set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) beta <- c(1,5) x <- cbind(rnorm(n), rnorm(n)) sigma.sq = 1 phi = 1 tau.sq = 0.1 B <- as.matrix(beta) D <- as.matrix(dist(coords)) R <- exp(-phi*D) w <- rmvn(1, rep(0,n), sigma.sq*R) y <- rnorm(n, x%*%B + w, sqrt(tau.sq)) estimation_result <- BRISC_estimation(coords, y, x) estimation_result$Theta ##Estimates of covariance model parameters. estimation_result$Beta ##Estimates of Beta
The function BRISC_neighbor
creates the set of nearest neighbors for a given set of coordinates, which can be used as an input
for "neighbor"
argument in BRISC_estimation
. This is especially useful for avoiding often computationally
intensive nearest neighbor finding scheme in case of multiple application of BRISC_estimation
on a fixed set of
coordinates.
BRISC_neighbor(coords, n.neighbors = 15, n_omp = 1, order = "Sum_coords", search.type = "tree", verbose = TRUE, ordering = NULL, tol = 12 )
BRISC_neighbor(coords, n.neighbors = 15, n_omp = 1, order = "Sum_coords", search.type = "tree", verbose = TRUE, ordering = NULL, tol = 12 )
coords |
an |
n.neighbors |
number of neighbors used in the NNGP. Default value is 15. |
n_omp |
number of threads to be used, value can be more than 1 if source code is compiled with OpenMP support. Default is 1. |
order |
keyword that specifies the ordering scheme to be used in ordering the observations. Supported keywords are:
|
search.type |
keyword that specifies type of nearest neighbor search algorithm to be used. Supported keywords are:
|
verbose |
if |
ordering |
if not |
tol |
the input observation coordinates, response and the covariates are rounded to this many places after the decimal. The default value is 12. |
A list containing information regarding nearest neighbors which can be used as an input
for "neighbor"
argument in BRISC_estimation
.
Arkajyoti Saha [email protected],
Abhirup Datta [email protected]
Saha, A., & Datta, A. (2018). BRISC: bootstrap for rapid inference on spatial covariances. Stat, e184, DOI: 10.1002/sta4.184.
Datta, A., S. Banerjee, A.O. Finley, and A.E. Gelfand. (2016) Hierarchical Nearest-Neighbor Gaussian process models for large geostatistical datasets. Journal of the American Statistical Association, 111:800-812.
Guinness, J. (2018) Permutation and Grouping Methods for
Sharpening Gaussian Process Approximations, Technometrics,
DOI: 10.1080/00401706.2018.1437476,
https://github.com/joeguinness/gp_reorder/tree/master/R .
set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) neighbor_result <- BRISC_neighbor(coords)
set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) neighbor_result <- BRISC_neighbor(coords)
The function BRISC_order
outputs the ordering for a set of coordinates, which can be used as an input for
"ordering"
argument in BRISC_estimation
. This is especially useful for avoiding often computationally
intensive location ordering scheme in case of multiple application of BRISC_estimation
on a fixed set of
coordinates.
BRISC_order(coords, order = "Sum_coords", verbose = TRUE)
BRISC_order(coords, order = "Sum_coords", verbose = TRUE)
coords |
an |
order |
keyword that specifies the ordering scheme to be used in ordering the observations. Supported keywords are:
|
verbose |
if |
An integer vector of ordering of the input coordinates which can be used as an input for
"ordering"
argument in BRISC_estimation
.
Arkajyoti Saha [email protected],
Abhirup Datta [email protected]
Saha, A., & Datta, A. (2018). BRISC: bootstrap for rapid inference on spatial covariances. Stat, e184, DOI: 10.1002/sta4.184.
Datta, A., S. Banerjee, A.O. Finley, and A.E. Gelfand. (2016) Hierarchical Nearest-Neighbor Gaussian process models for large geostatistical datasets. Journal of the American Statistical Association, 111:800-812.
Guinness, J. (2018) Permutation and Grouping Methods for
Sharpening Gaussian Process Approximations, Technometrics,
DOI: 10.1080/00401706.2018.1437476,
https://github.com/joeguinness/gp_reorder/tree/master/R .
set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) ordering_result <- BRISC_order(coords)
set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) ordering_result <- BRISC_order(coords)
The function BRISC_prediction
performs fast prediction on a set of new locations with univariate
spatial regression models using Nearest Neighbor Gaussian Processes (NNGP) (Datta et al., 2016).
BRISC_prediction uses the parameter estimates from BRISC_estimation
for the prediction. Some
code blocks are borrowed from the R package: spNNGP: Spatial Regression Models for Large Datasets
using Nearest Neighbor Gaussian Processes
https://CRAN.R-project.org/package=spNNGP .
BRISC_prediction(BRISC_Out, coords.0, X.0 = NULL, n_omp = 1, verbose = TRUE, tol = 12)
BRISC_prediction(BRISC_Out, coords.0, X.0 = NULL, n_omp = 1, verbose = TRUE, tol = 12)
BRISC_Out |
an object of class |
coords.0 |
the spatial coordinates corresponding to prediction locations. Its structure should be same as that of coords
in |
X.0 |
the covariates for prediction locations. Its Structure should be identical (including intercept) with that of
covariates provided for estimation purpose in |
n_omp |
number of threads to be used, value can be more than 1 if source code is compiled with OpenMP support. Default is 1. |
verbose |
if |
tol |
the coordinates and the covariates corresponding to the prediction locations are rounded to this many places after the decimal. The default value is 12. |
A list comprising of the following:
prediction |
predicted response corresponding to X.0 and coords.0. |
prediction.ci |
confidence intervals corresponding to the predictions. |
prediction.time |
time (in seconds) required to perform the prediction after preprocessing data in |
Arkajyoti Saha [email protected],
Abhirup Datta [email protected]
Datta, A., S. Banerjee, A.O. Finley, and A.E. Gelfand. (2016) Hierarchical Nearest-Neighbor Gaussian process models for large geostatistical datasets. Journal of the American Statistical Association, 111:800-812.
Andrew Finley, Abhirup Datta and Sudipto Banerjee (2017). spNNGP: Spatial Regression Models for Large Datasets using Nearest Neighbor Gaussian Processes. R package version 0.1.1. https://CRAN.R-project.org/package=spNNGP
rmvn <- function(n, mu = 0, V = matrix(1)){ p <- length(mu) if(any(is.na(match(dim(V),p)))) stop("Dimension not right!") D <- chol(V) t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p))) } set.seed(1) n <- 500 coords <- cbind(runif(n,0,1), runif(n,0,1)) beta <- c(1,5) x <- cbind(rnorm(n), rnorm(n)) sigma.sq = 1 phi = 1 tau.sq = 0.1 B <- as.matrix(beta) D <- as.matrix(dist(coords)) R <- exp(-phi*D) w <- rmvn(1, rep(0,n), sigma.sq*R) y <- rnorm(n, x%*%B + w, sqrt(tau.sq)) estimation_result <- BRISC_estimation(coords[1:400,], y[1:400], x[1:400,]) prediction_result <- BRISC_prediction(estimation_result, coords[401:500,], x[401:500,])
rmvn <- function(n, mu = 0, V = matrix(1)){ p <- length(mu) if(any(is.na(match(dim(V),p)))) stop("Dimension not right!") D <- chol(V) t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p))) } set.seed(1) n <- 500 coords <- cbind(runif(n,0,1), runif(n,0,1)) beta <- c(1,5) x <- cbind(rnorm(n), rnorm(n)) sigma.sq = 1 phi = 1 tau.sq = 0.1 B <- as.matrix(beta) D <- as.matrix(dist(coords)) R <- exp(-phi*D) w <- rmvn(1, rep(0,n), sigma.sq*R) y <- rnorm(n, x%*%B + w, sqrt(tau.sq)) estimation_result <- BRISC_estimation(coords[1:400,], y[1:400], x[1:400,]) prediction_result <- BRISC_prediction(estimation_result, coords[401:500,], x[401:500,])
The function BRISC_simulation
simulates correlated data (known structure) using Nearest Neighbor
Gaussian Processes (NNGP). BRISC_simulation
uses the sparse Cholesky representation of Vecchia’s
likelihood developed in Datta et al., 2016. BRISC_simulation
uses BRISC_correlation
for this purpose.
BRISC_simulation(coords, sim_number = 1, seeds = NULL, sigma.sq = 1, tau.sq = 0, phi = 1, nu = 1.5, n.neighbors = NULL, n_omp = 1, cov.model = "exponential", search.type = "tree", stabilization = NULL, verbose = TRUE, tol = 12)
BRISC_simulation(coords, sim_number = 1, seeds = NULL, sigma.sq = 1, tau.sq = 0, phi = 1, nu = 1.5, n.neighbors = NULL, n_omp = 1, cov.model = "exponential", search.type = "tree", stabilization = NULL, verbose = TRUE, tol = 12)
coords |
an |
sim_number |
number of simulations. Default value is 1. |
seeds |
seeds which are used in generation of the initial independent data. Default value is |
sigma.sq |
value of sigma square. Default value is 1. |
tau.sq |
value of tau square. Default value is 0.1. |
phi |
value of phi. Default value is 1. |
nu |
starting value of nu, only required for matern covariance model. Default value is 1.5. |
n.neighbors |
number of neighbors used in the NNGP. Default value is 15. |
n_omp |
number of threads to be used, value can be more than 1 if source code is compiled with OpenMP support. Default is 1. |
cov.model |
keyword that specifies the covariance function to be used in modelling the spatial dependence structure
among the observations. Supported keywords are: |
search.type |
keyword that specifies type of nearest neighbor search algorithm to be used. Supported keywords are:
|
stabilization |
when we use a very smooth covarince model (lower values of phi for spherical and Gaussian
covariance and low phi and high nu for Matern covarinace) in absence of a non-negligble nugget, the correlation process may fail
due to computational instability. If |
verbose |
if |
tol |
the input observation coordinates are rounded to this many places after the decimal. The default value is 12. |
A list comprising of the following:
coords |
the matrix |
n.neighbors |
the used value of |
cov.model |
the used covariance model. |
Theta |
parameters of covarinace model; accounts for |
input.data |
the |
output.data |
the |
time |
time (in seconds) required after preprocessing data in |
Arkajyoti Saha [email protected],
Abhirup Datta [email protected]
set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) sigma.sq = 1 phi = 1 simulation_result <- BRISC_simulation(coords, sim_number = 3)
set.seed(1) n <- 1000 coords <- cbind(runif(n,0,1), runif(n,0,1)) sigma.sq = 1 phi = 1 simulation_result <- BRISC_simulation(coords, sim_number = 3)
The function BRISC_variogram.ci
plots estimated Variogram and associated confience region.
BRISC_variogram.ci
uses the parameter estimates from BRISC_estimation
and associated confidence interval from BRISC_bootstrap
.
BRISC_variogram.ci(BRISC_Out, confidence_est, plot.variogram = FALSE)
BRISC_variogram.ci(BRISC_Out, confidence_est, plot.variogram = FALSE)
BRISC_Out |
an object of class |
confidence_est |
bootstrp sample of the Theta parameters, obtained from |
plot.variogram |
if |
A list comprising of the following:
variogram |
variogram and associated confidence region corresponding to lag ranging from 0 to 20, evaluated at 0.01 frequency. |
Plot |
plots the Variogram and associated confidence region with legends. |
Arkajyoti Saha [email protected],
Abhirup Datta [email protected]
rmvn <- function(n, mu = 0, V = matrix(1)){ p <- length(mu) if(any(is.na(match(dim(V),p)))) stop("Dimension not right!") D <- chol(V) t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p))) } set.seed(1) n <- 300 coords <- cbind(runif(n,0,1), runif(n,0,1)) beta <- c(1,5) x <- cbind(rnorm(n), rnorm(n)) sigma.sq = 1 phi = 5 tau.sq = 0.1 B <- as.matrix(beta) D <- as.matrix(dist(coords)) R <- exp(-phi*D) w <- rmvn(1, rep(0,n), sigma.sq*R) y <- rnorm(n, x%*%B + w, sqrt(tau.sq)) estimation_result <- BRISC_estimation(coords, y, x) bootstrap_result <- BRISC_bootstrap(estimation_result, n_boot = 10) varg <- BRISC_variogram.ci(estimation_result, bootstrap_result$boot.Theta, plot.variogram = TRUE)
rmvn <- function(n, mu = 0, V = matrix(1)){ p <- length(mu) if(any(is.na(match(dim(V),p)))) stop("Dimension not right!") D <- chol(V) t(matrix(rnorm(n*p), ncol=p)%*%D + rep(mu,rep(n,p))) } set.seed(1) n <- 300 coords <- cbind(runif(n,0,1), runif(n,0,1)) beta <- c(1,5) x <- cbind(rnorm(n), rnorm(n)) sigma.sq = 1 phi = 5 tau.sq = 0.1 B <- as.matrix(beta) D <- as.matrix(dist(coords)) R <- exp(-phi*D) w <- rmvn(1, rep(0,n), sigma.sq*R) y <- rnorm(n, x%*%B + w, sqrt(tau.sq)) estimation_result <- BRISC_estimation(coords, y, x) bootstrap_result <- BRISC_bootstrap(estimation_result, n_boot = 10) varg <- BRISC_variogram.ci(estimation_result, bootstrap_result$boot.Theta, plot.variogram = TRUE)