Title: | Bayesian Shrinkage Estimators for Precision Matrices in Gaussian Graphical Models |
---|---|
Description: | This R package offers block Gibbs samplers for the Bayesian (adaptive) graphical lasso, ridge, and naive elastic net priors. These samplers facilitate the simulation of the posterior distribution of precision matrices for Gaussian distributed data and were originally proposed by: Wang (2012) <doi:10.1214/12-BA729>; Smith et al. (2022) <doi:10.48550/arXiv.2210.16290> and Smith et al. (2023) <doi:10.48550/arXiv.2306.14199>, respectively. |
Authors: | Jarod Smith [aut, cre] , Mohammad Arashi [aut] , Andriette Bekker [aut] |
Maintainer: | Jarod Smith <[email protected]> |
License: | GPL (>= 3) |
Version: | 0.3.0 |
Built: | 2024-11-05 06:28:59 UTC |
Source: | CRAN |
Implements the Type I naive Bayesian adaptive graphical elastic net block Gibbs sampler to simulate the posterior distribution of the precision matrix for Gaussian graphical models.
blockBAGENI( X, burnin, iterations, verbose = TRUE, r = 0.001, s = 0.01, a = 0.001, b = 0.1 )
blockBAGENI( X, burnin, iterations, verbose = TRUE, r = 0.001, s = 0.01, a = 0.001, b = 0.1 )
X |
A numeric matrix, assumed to be generated from a multivariate Gaussian distribution. |
burnin |
An integer specifying the number of burn-in iterations. |
iterations |
An integer specifying the length of the Markov chain after the burn-in iterations. |
verbose |
A logical determining whether the progress of the MCMC sampler should be displayed. |
r |
A double specifying the value of the shape parameter for the gamma prior associated with the Bayesian graphical lasso penalty term. |
s |
A double specifying the value of the scale parameter for the gamma prior associated with the Bayesian graphical lasso penalty term. |
a |
A double specifying the value of the shape parameter for the inverse gamma prior associated with the Bayesian graphical ridge penalty term. |
b |
A double specifying the value of the scale parameter for the inverse gamma prior associated with the Bayesian graphical ridge penalty term. |
A list containing precision 'Omega' and covariance 'Sigma' matrices from the Markov chains.
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBAGENI(X, iterations = 1000, burnin = 500) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBAGENI(X, iterations = 1000, burnin = 500) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
Implements the Type II naive Bayesian adaptive graphical elastic net block Gibbs sampler to simulate the posterior distribution of the precision matrix for Gaussian graphical models.
blockBAGENII(X, burnin, iterations, verbose = TRUE, s = 0.1, b = 0.001)
blockBAGENII(X, burnin, iterations, verbose = TRUE, s = 0.1, b = 0.001)
X |
A numeric matrix, assumed to be generated from a multivariate Gaussian distribution. |
burnin |
An integer specifying the number of burn-in iterations. |
iterations |
An integer specifying the length of the Markov chain after the burn-in iterations. |
verbose |
A logical determining whether the progress of the MCMC sampler should be displayed. |
s |
A double specifying the value of the rate parameter for the exponential prior associated with the Bayesian graphical lasso penalty term. |
b |
A double specifying the value of the rate parameter for the exponential prior associated with the Bayesian graphical ridge penalty term. |
A list containing precision 'Omega' and covariance 'Sigma' matrices from the Markov chains.
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBAGENII(X, iterations = 1000, burnin = 500) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBAGENII(X, iterations = 1000, burnin = 500) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
Implements a Bayesian adaptive graphical lasso block Gibbs sampler to simulate the posterior distribution of the precision matrix for Gaussian graphical models.
blockBAGL(X, burnin, iterations, verbose = TRUE, r = 0.01, s = 1e-06)
blockBAGL(X, burnin, iterations, verbose = TRUE, r = 0.01, s = 1e-06)
X |
A numeric matrix, assumed to be generated from a multivariate Gaussian distribution. |
burnin |
An integer specifying the number of burn-in iterations. |
iterations |
An integer specifying the length of the Markov chain after the burn-in iterations. |
verbose |
A logical determining whether the progress of the MCMC sampler should be displayed. |
r |
A double specifying the value of the shape parameter for the gamma prior. |
s |
A double specifying the value of the scale parameter for the gamma prior. |
A list containing precision 'Omega' and covariance 'Sigma' matrices from the Markov chains.
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBAGL(X, iterations = 1000, burnin = 500) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBAGL(X, iterations = 1000, burnin = 500) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
Implements a Bayesian adaptive graphical ridge block Gibbs sampler to simulate the posterior distribution of the precision matrix for Gaussian graphical models.
blockBAGR(X, burnin, iterations, verbose = TRUE, a = 1, b = 0.01)
blockBAGR(X, burnin, iterations, verbose = TRUE, a = 1, b = 0.01)
X |
A numeric matrix, assumed to be generated from a multivariate Gaussian distribution. |
burnin |
An integer specifying the number of burn-in iterations. |
iterations |
An integer specifying the length of the Markov chain after the burn-in iterations. |
verbose |
A logical determining whether the progress of the MCMC sampler should be displayed. |
a |
A double specifying the value of the shape parameter for the inverse gamma prior. |
b |
A double specifying the value of the scale parameter for the inverse gamma prior. |
A list containing precision 'Omega' and covariance 'Sigma' matrices from the Markov chains.
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBAGR(X, iterations = 1000, burnin = 500) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBAGR(X, iterations = 1000, burnin = 500) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
Implements the Bayesian graphical elastic net block Gibbs sampler to simulate the posterior distribution of the precision matrix for Gaussian graphical models.
blockBGEN(X, burnin, iterations, lambda = 1, sig = 1, verbose = TRUE)
blockBGEN(X, burnin, iterations, lambda = 1, sig = 1, verbose = TRUE)
X |
A numeric matrix, assumed to be generated from a multivariate Gaussian distribution. |
burnin |
An integer specifying the number of burn-in iterations. |
iterations |
An integer specifying the length of the Markov chain after the burn-in iterations. |
lambda |
A numeric value representing the rate parameter for the double exponential and exponential prior associated with the Bayesian graphical lasso penalty term. |
sig |
A numeric value representing the standard deviation parameter for the double Gaussian and truncated Gaussian prior associated with the Bayesian graphical ridge penalty term. |
verbose |
A logical determining whether the progress of the MCMC sampler should be displayed. |
A list containing precision 'Omega' and covariance 'Sigma' matrices from the Markov chains.
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBGEN(X, iterations = 1000, burnin = 500, lambda = 1, sig = 1) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBGEN(X, iterations = 1000, burnin = 500, lambda = 1, sig = 1) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
Implements a Bayesian graphical lasso block Gibbs sampler to simulate the posterior distribution of the precision matrix for Gaussian graphical models.
blockBGL(X, burnin, iterations, lambda = 1, verbose = TRUE)
blockBGL(X, burnin, iterations, lambda = 1, verbose = TRUE)
X |
A numeric matrix, assumed to be generated from a multivariate Gaussian distribution. |
burnin |
An integer representing the number of burn-in iterations. |
iterations |
An integer representing the length of the Markov chain post burn-in. |
lambda |
A numeric value representing the rate parameter for the double exponential and exponential prior. |
verbose |
A logical indicating if the MCMC sampler progress should be printed. |
A list containing precision 'Omega' and covariance 'Sigma' matrices from the Markov chains.
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBGL(X, iterations = 1000, burnin = 500, lambda = 0.5) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBGL(X, iterations = 1000, burnin = 500, lambda = 0.5) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
Implements a Bayesian graphical ridge block Gibbs sampler to simulate the posterior distribution of the precision matrix for Gaussian graphical models.
blockBGR(X, burnin, iterations, sig = 1, verbose = TRUE)
blockBGR(X, burnin, iterations, sig = 1, verbose = TRUE)
X |
A numeric matrix, assumed to be generated from a multivariate Gaussian distribution. |
burnin |
An integer representing the number of burn-in iterations. |
iterations |
An integer representing the length of the Markov chain post burn-in. |
sig |
A numeric value representing the standard deviation parameter for the double Gaussian and truncated Gaussian prior. |
verbose |
A logical indicating if the MCMC sampler progress should be printed. |
A list containing precision 'Omega' and covariance 'Sigma' matrices from the Markov chains.
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBGR(X, iterations = 1000, burnin = 500, sig = 0.5) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)
# Generate true precision matrix: p <- 10 n <- 500 OmegaTrue <- pracma::Toeplitz(c(0.7^rep(1:p-1))) SigTrue <- pracma::inv(OmegaTrue) # Generate expected value vector: mu <- rep(0,p) # Generate multivariate normal distribution: set.seed(123) X <- MASS::mvrnorm(n, mu = mu, Sigma = SigTrue) # Generate posterior distribution: posterior <- blockBGR(X, iterations = 1000, burnin = 500, sig = 0.5) # Estimated precision matrix using the mean of the posterior: OmegaEst <- apply(simplify2array(posterior$Omega), 1:2, mean)