Title: | Draw Samples of Truncated Multivariate Normal Distributions |
---|---|
Description: | Draw samples from truncated multivariate normal distribution using the sequential nearest neighbor (SNN) method introduced in "Scalable Sampling of Truncated Multivariate Normals Using Sequential Nearest-Neighbor Approximation" <doi:10.48550/arXiv.2406.17307>. |
Authors: | Jian Cao [aut, cre], Matthias Katzfuss [aut] |
Maintainer: | Jian Cao <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.2.0 |
Built: | 2025-02-08 01:15:48 UTC |
Source: | https://github.com/jcatwood/nntmvn |
Find ordered nearest neighbors based on correlation, assuming the absolute value of the correlation is monotonically decreasing with distance. Returns an n X (m + 1) matrix, each row indicating the m + 1 nearest neighbors including itself.
corr_nn(covmat, m)
corr_nn(covmat, m)
covmat |
the covariance matrix |
m |
the number of nearest neighbors |
an n X (m + 1) matrix
library(RANN) library(nntmvn) set.seed(123) d <- 3 n <- 100 locs <- matrix(runif(d * n), n, d) covparms <- c(2, 0.01, 0) covmat <- GpGp::matern15_isotropic(covparms, locs) m <- 10 NNarray_test <- RANN::nn2(locs, k = m + 1)[[1]] NNarray <- nntmvn::corr_nn(covmat, m) cat("Number of mismatch is", sum(NNarray != NNarray_test, na.rm = TRUE))
library(RANN) library(nntmvn) set.seed(123) d <- 3 n <- 100 locs <- matrix(runif(d * n), n, d) covparms <- c(2, 0.01, 0) covmat <- GpGp::matern15_isotropic(covparms, locs) m <- 10 NNarray_test <- RANN::nn2(locs, k = m + 1)[[1]] NNarray <- nntmvn::corr_nn(covmat, m) cat("Number of mismatch is", sum(NNarray != NNarray_test, na.rm = TRUE))
Draw one sample of the underlying GP responses for a partially censored Gaussian process using sequential nearest neighbor (SNN) method
rptmvn( y, cens_lb, cens_ub, mask_cens, m = 30, covmat = NULL, locs = NULL, cov_name = NULL, cov_parm = NULL, NN = NULL, ordering = 0, seed = NULL )
rptmvn( y, cens_lb, cens_ub, mask_cens, m = 30, covmat = NULL, locs = NULL, cov_name = NULL, cov_parm = NULL, NN = NULL, ordering = 0, seed = NULL )
y |
uncensored responses of length n, where n is the number of all responses |
cens_lb |
lower bound vector for TMVN of length n |
cens_ub |
upper bound vector for TMVN of length n |
mask_cens |
mask for censored responses (also locations) of length n |
m |
positive integer for the number of nearest neighbors used |
covmat |
n-by-n dense covariance matrix, either |
locs |
location matrix n X d |
cov_name |
covariance function name from the |
cov_parm |
parameters for the covariance function from the |
NN |
n X m matrix for nearest neighbors. i-th row is the nearest neighbor indices of y_i. |
ordering |
|
seed |
set seed for reproducibility |
a vector of length n representing the underlying GP responses
library(GpGp) library(RANN) library(nntmvn) set.seed(123) x <- matrix(seq(from = 0, to = 1, length.out = 51), ncol = 1) cov_name <- "matern15_isotropic" cov_parm <- c(1.0, 0.1, 0.001) #' variance, range, nugget cov_func <- getFromNamespace(cov_name, "GpGp") covmat <- cov_func(cov_parm, x) y <- t(chol(covmat)) %*% rnorm(length(x)) mask <- y < 0.3 y_cens <- y y_cens[mask] <- NA lb <- rep(-Inf, 100) ub <- rep(0.3, 100) m <- 10 y_samp_mtd1 <- rptmvn(y_cens, lb, ub, mask, m = m, locs = x, cov_name = cov_name, cov_parm = cov_parm, seed = 123 ) y_samp_mtd2 <- rptmvn(y_cens, lb, ub, mask, m = m, covmat = covmat, seed = 123 ) plot(x, y_cens, ylim = range(y)) points(x[mask, ], y[mask], col = "blue") plot(x, y_cens, ylim = range(y)) points(x[mask, ], y_samp_mtd1[mask], col = "red") plot(x, y_cens, ylim = range(y)) points(x[mask, ], y_samp_mtd2[mask], col = "brown")
library(GpGp) library(RANN) library(nntmvn) set.seed(123) x <- matrix(seq(from = 0, to = 1, length.out = 51), ncol = 1) cov_name <- "matern15_isotropic" cov_parm <- c(1.0, 0.1, 0.001) #' variance, range, nugget cov_func <- getFromNamespace(cov_name, "GpGp") covmat <- cov_func(cov_parm, x) y <- t(chol(covmat)) %*% rnorm(length(x)) mask <- y < 0.3 y_cens <- y y_cens[mask] <- NA lb <- rep(-Inf, 100) ub <- rep(0.3, 100) m <- 10 y_samp_mtd1 <- rptmvn(y_cens, lb, ub, mask, m = m, locs = x, cov_name = cov_name, cov_parm = cov_parm, seed = 123 ) y_samp_mtd2 <- rptmvn(y_cens, lb, ub, mask, m = m, covmat = covmat, seed = 123 ) plot(x, y_cens, ylim = range(y)) points(x[mask, ], y[mask], col = "blue") plot(x, y_cens, ylim = range(y)) points(x[mask, ], y_samp_mtd1[mask], col = "red") plot(x, y_cens, ylim = range(y)) points(x[mask, ], y_samp_mtd2[mask], col = "brown")
Draw one sample from a truncated multivariate normal (TMVN) distribution using sequential nearest neighbor (SNN) method
rtmvn( cens_lb, cens_ub, m = 30, covmat = NULL, locs = NULL, cov_name = NULL, cov_parm = NULL, NN = NULL, ordering = 0, seed = NULL )
rtmvn( cens_lb, cens_ub, m = 30, covmat = NULL, locs = NULL, cov_name = NULL, cov_parm = NULL, NN = NULL, ordering = 0, seed = NULL )
cens_lb |
lower bound vector for TMVN of length n |
cens_ub |
upper bound vector for TMVN of length n |
m |
positive integer for the number of nearest neighbors used |
covmat |
n-by-n dense covariance matrix, either |
locs |
location matrix n X d |
cov_name |
covariance function name from the |
cov_parm |
parameters for the covariance function from the |
NN |
n X m matrix for nearest neighbors. i-th row is the nearest neighbor indices of y_i. |
ordering |
|
seed |
set seed for reproducibility |
a vector of length n representing the underlying GP responses
library(nntmvn) library(TruncatedNormal) set.seed(123) x <- matrix(seq(from = 0, to = 1, length.out = 51), ncol = 1) cov_name <- "matern15_isotropic" cov_parm <- c(1.0, 0.1, 0.001) #'' variance, range, nugget cov_func <- getFromNamespace(cov_name, "GpGp") covmat <- cov_func(cov_parm, x) lb <- rep(-Inf, nrow(x)) ub <- rep(-1, nrow(x)) m <- 30 samp_SNN <- matrix(NA, 3, nrow(x)) for (i in 1:3) { samp_SNN[i, ] <- nntmvn::rtmvn(lb, ub, m = m, covmat = covmat, locs = x, ordering = 0) } samp_TN <- TruncatedNormal::rtmvnorm(3, rep(0, nrow(x)), covmat, lb, ub) qqplot(samp_SNN, samp_TN, xlim = range(samp_SNN, samp_TN), ylim = range(samp_SNN, samp_TN)) abline(a = 0, b = 1, lty = "dashed", col = "red")
library(nntmvn) library(TruncatedNormal) set.seed(123) x <- matrix(seq(from = 0, to = 1, length.out = 51), ncol = 1) cov_name <- "matern15_isotropic" cov_parm <- c(1.0, 0.1, 0.001) #'' variance, range, nugget cov_func <- getFromNamespace(cov_name, "GpGp") covmat <- cov_func(cov_parm, x) lb <- rep(-Inf, nrow(x)) ub <- rep(-1, nrow(x)) m <- 30 samp_SNN <- matrix(NA, 3, nrow(x)) for (i in 1:3) { samp_SNN[i, ] <- nntmvn::rtmvn(lb, ub, m = m, covmat = covmat, locs = x, ordering = 0) } samp_TN <- TruncatedNormal::rtmvnorm(3, rep(0, nrow(x)), covmat, lb, ub) qqplot(samp_SNN, samp_TN, xlim = range(samp_SNN, samp_TN), ylim = range(samp_SNN, samp_TN)) abline(a = 0, b = 1, lty = "dashed", col = "red")