Title: | Multiscale Graphical Lasso |
---|---|
Description: | Inference of Multiscale graphical models with neighborhood selection approach. The method is based on solving a convex optimization problem combining a Lasso and fused-group Lasso penalties. This allows to infer simultaneously a conditional independence graph and a clustering partition. The optimization is based on the Continuation with Nesterov smoothing in a Shrinkage-Thresholding Algorithm solver (Hadj-Selem et al. 2018) <doi:10.1109/TMI.2018.2829802> implemented in python. |
Authors: | Edmond Sanou [aut, cre], Tung Le [ctb], Christophe Ambroise [ths], Geneviève Robin [ths] |
Maintainer: | Edmond Sanou <[email protected]> |
License: | MIT + file LICENSE |
Version: | 0.1.3 |
Built: | 2025-03-10 05:51:35 UTC |
Source: | https://github.com/desanou/mglasso |
Adjacency matrix
adj_mat(mat, sym_rule = "and")
adj_mat(mat, sym_rule = "and")
mat |
matrix of regression coefficients |
sym_rule |
symmetrization rule, either AND or OR |
adjacency matrix
Init Beta 1 matrix
Init Beta 1 matrix
beta_idty(p) beta_idty(p)
beta_idty(p) beta_idty(p)
Init Beta via OLS
Init Beta via OLS
Initialize regression matrix
beta_ols(X) beta_ols(X) beta_ols(X)
beta_ols(X) beta_ols(X) beta_ols(X)
X |
data |
A zero-diagonal matrix of regression vectors.
vectorize beta matrix
vectorize beta matrix
Transform a matrix of regression coefficients to vector removing the diagonal
beta_to_vector(beta_mat) beta_to_vector(beta_mat) beta_to_vector(beta_mat)
beta_to_vector(beta_mat) beta_to_vector(beta_mat) beta_to_vector(beta_mat)
beta_mat |
matrix of regressions vectors |
A numeric vector of all regression coefficients.
return precision matrix
bloc_diag(n_vars, connectivity_mat, prop_clusters, rho)
bloc_diag(n_vars, connectivity_mat, prop_clusters, rho)
cah_glasso
cah_glasso(num_clusters, data, lam1, hclust_obj)
cah_glasso(num_clusters, data, lam1, hclust_obj)
Solve the MGLasso optimization problem using CONESTA algorithm. Interface to the pylearn.parsimony python library.
conesta( X, lam1, lam2, beta_warm = c(0), type_ = "initial", W_ = NULL, mean_ = FALSE, max_iter_ = 10000, prec_ = 0.01 )
conesta( X, lam1, lam2, beta_warm = c(0), type_ = "initial", W_ = NULL, mean_ = FALSE, max_iter_ = 10000, prec_ = 0.01 )
X |
Data matrix nxp. |
lam1 |
Sparsity penalty. |
lam2 |
Total variation penalty. |
beta_warm |
Warm initialization vector. |
type_ |
Character scalar. By default set to initial version which doesn't use weights |
W_ |
Weights matrix for total variation penalties. |
mean_ |
Logical scalar. If TRUE weights the optimization function by the inverse of sample size. |
max_iter_ |
Numeric scalar. Maximum number of iterations. |
prec_ |
Numeric scalar. Tolerance for the stopping criterion (duality gap). |
COntinuation with NEsterov smoothing in a Shrinkage-Thresholding Algorithm (CONESTA, Hadj-Selem et al. 2018) doi:10.1109/TMI.2018.2829802 is an algorithm design for solving optimization problems including group-wise penalties. This function is an interface with the python solver. The MGLasso problem is first reformulated in a problem of the form
where vector is
the vectorized form of matrix
.
Numeric matrix of size pxp. Line k
of the matrix represents
the coefficients obtained from the L1-L2 penalized regression of variable
k
on the others.
mglasso()
for the MGLasso model estimation.
## Not run: # because of installation of external packages during checks mglasso::install_pylearn_parsimony(envname = "rmglasso", method = "conda") reticulate::use_condaenv("rmglasso", required = TRUE) reticulate::py_config() n = 30 K = 2 p = 4 rho = 0.85 blocs <- list() for (j in 1:K) { bloc <- matrix(rho, nrow = p/K, ncol = p/K) for(i in 1:(p/K)) { bloc[i,i] <- 1 } blocs[[j]] <- bloc } mat.covariance <- Matrix::bdiag(blocs) mat.covariance set.seed(11) X <- mvtnorm::rmvnorm(n, mean = rep(0,p), sigma = as.matrix(mat.covariance)) X <- scale(X) res <- conesta(X, 0.1, 0.1) ## End(Not run)
## Not run: # because of installation of external packages during checks mglasso::install_pylearn_parsimony(envname = "rmglasso", method = "conda") reticulate::use_condaenv("rmglasso", required = TRUE) reticulate::py_config() n = 30 K = 2 p = 4 rho = 0.85 blocs <- list() for (j in 1:K) { bloc <- matrix(rho, nrow = p/K, ncol = p/K) for(i in 1:(p/K)) { bloc[i,i] <- 1 } blocs[[j]] <- bloc } mat.covariance <- Matrix::bdiag(blocs) mat.covariance set.seed(11) X <- mvtnorm::rmvnorm(n, mean = rep(0,p), sigma = as.matrix(mat.covariance)) X <- scale(X) res <- conesta(X, 0.1, 0.1) ## End(Not run)
CONESTA solver for numerical experiments.
conesta_rwrapper( X, lam1, lam2, beta_warm = c(0), type_ = "initial", W_ = NULL, mean_ = FALSE, max_iter_ = 10000, prec_ = 0.01 )
conesta_rwrapper( X, lam1, lam2, beta_warm = c(0), type_ = "initial", W_ = NULL, mean_ = FALSE, max_iter_ = 10000, prec_ = 0.01 )
## Not run: # because of installation of external packages during checks mglasso::install_pylearn_parsimony(envname = "rmglasso", method = "conda") reticulate::use_condaenv("rmglasso", required = TRUE) reticulate::py_config() n = 30 K = 2 p = 4 rho = 0.85 blocs <- list() for (j in 1:K) { bloc <- matrix(rho, nrow = p/K, ncol = p/K) for(i in 1:(p/K)) { bloc[i,i] <- 1 } blocs[[j]] <- bloc } mat.covariance <- Matrix::bdiag(blocs) mat.covariance set.seed(11) X <- mvtnorm::rmvnorm(n, mean = rep(0,p), sigma = as.matrix(mat.covariance)) X <- scale(X) res <- conesta_rwrapper(X, 0.1, 0.1) ## End(Not run)
## Not run: # because of installation of external packages during checks mglasso::install_pylearn_parsimony(envname = "rmglasso", method = "conda") reticulate::use_condaenv("rmglasso", required = TRUE) reticulate::py_config() n = 30 K = 2 p = 4 rho = 0.85 blocs <- list() for (j in 1:K) { bloc <- matrix(rho, nrow = p/K, ncol = p/K) for(i in 1:(p/K)) { bloc[i,i] <- 1 } blocs[[j]] <- bloc } mat.covariance <- Matrix::bdiag(blocs) mat.covariance set.seed(11) X <- mvtnorm::rmvnorm(n, mean = rep(0,p), sigma = as.matrix(mat.covariance)) X <- scale(X) res <- conesta_rwrapper(X, 0.1, 0.1) ## End(Not run)
cost
computes the cost function of Mglasso
method.
cost(beta, x, lambda1 = 0, lambda2 = 0) cost(beta, x, lambda1 = 0, lambda2 = 0) cost(beta, x, lambda1 = 0, lambda2 = 0)
cost(beta, x, lambda1 = 0, lambda2 = 0) cost(beta, x, lambda1 = 0, lambda2 = 0) cost(beta, x, lambda1 = 0, lambda2 = 0)
beta |
p by p numeric matrix. In rows, regression vectors coefficients after node-wise regression. |
x |
n by p numeric matrix. Data with variables in columns. |
lambda1 |
numeric scalar. Lasso penalization parameter. |
lambda2 |
numeric scalar. Fused-group Lasso penalization parameter. |
numeric scalar. The cost.
distances Beta
Compute distance matrix between regression vectors
dist_beta(beta, distance = "euclidean") dist_beta(beta, distance = "euclidean")
dist_beta(beta, distance = "euclidean") dist_beta(beta, distance = "euclidean")
beta |
matrix of regression vectors |
distance |
euclidean or relative distance |
A numeric matrix of distances.
Mean error from classical regression
error(Theta, X)
error(Theta, X)
Formula from Huge paper
error_huge(Theta, X)
error_huge(Theta, X)
TO DO: Fill upper triangular matrix then sum up with the transpose to have full matrix
expand_beta(beta_level, clusters)
expand_beta(beta_level, clusters)
doesn't work when dealing with matrix where diagonal of zero should be adjusted
expand_beta_deprecated(beta_level, clusters)
expand_beta_deprecated(beta_level, clusters)
extracts meta-variables indices
extract_meta(full_graph = NULL, clusters)
extract_meta(full_graph = NULL, clusters)
fun_lines
applies function fun
to regression vectors while reordering the coefficients,
such that the j
-th coefficient in beta[j, ]
is permuted with the i
-th coefficient.
fun_lines(i, j, beta, fun = `-`, ni = 1, nj = 1)
fun_lines(i, j, beta, fun = `-`, ni = 1, nj = 1)
i |
integer scalar. Index of the first vector. |
j |
integer scalar. Index of the second vector. |
beta |
p by p numeric matrix. In rows, regression vectors coefficients after node-wise regression. |
fun |
function. Applied on lines. |
ni |
integer scalar. Weight for vector |
nj |
integer scalar. Weight for vector |
numeric vector
beta <- matrix(round(rnorm(9),2), ncol = 3) diag(beta) <- 0 beta fun_lines(1, 2, beta) fun_lines(2, 1, beta)
beta <- matrix(round(rnorm(9),2), ncol = 3) diag(beta) <- 0 beta fun_lines(1, 2, beta) fun_lines(2, 1, beta)
Plot ROC curve and calculate AUC
get_auc(omega_hat_list, omega, to = to_)
get_auc(omega_hat_list, omega, to = to_)
type |
Classical ROC curve tpr = f(FPR) or TPR = f(precision) adjusted version. Compute AUC and partial AUC |
compute range of number of clusters from ROC outputs take in parameter an object from reorder_mglasso_roc_calculations
get_range_nclusters(out, thresh_fuse = 1e-06, p = 40)
get_range_nclusters(out, thresh_fuse = 1e-06, p = 40)
Title
ggplot_roc( omega_hat_list, omega, type = c("classical", "precision_recall"), main = NULL )
ggplot_roc( omega_hat_list, omega, type = c("classical", "precision_recall"), main = NULL )
main |
Neighborhood selection estimate
graph_estimate(rho, X)
graph_estimate(rho, X)
Plot the image of a matrix
image_sparse(matrix, main_ = "", sub_ = "", col_names = FALSE)
image_sparse(matrix, main_ = "", sub_ = "", col_names = FALSE)
matrix |
matrix of regression coefficients |
main_ |
title |
sub_ |
subtitle |
col_names |
columns names |
No return value.
pylearn-parsimony contains the solver CONESTA used for the mglasso problem and is available on github at https://github.com/neurospin/pylearn-parsimony It is advised to use a python version ">=3.7,<3.10". Indeed, the latest version of scipy under which mglasso was developped is scipy 1.7.1 which is based on python ">=3.7,<3.10". In turn, this version of scipy can only be associated with a version of numpy ">=1.16.5,<1.23.0"
install_pylearn_parsimony( method = c("auto", "virtualenv", "conda"), conda = "auto", extra_pack = c("scipy == 1.7.1", "scikit-learn", "numpy == 1.22.4", "six", "matplotlib"), python_version = "3.8", restart_session = TRUE, envname = NULL, ... )
install_pylearn_parsimony( method = c("auto", "virtualenv", "conda"), conda = "auto", extra_pack = c("scipy == 1.7.1", "scikit-learn", "numpy == 1.22.4", "six", "matplotlib"), python_version = "3.8", restart_session = TRUE, envname = NULL, ... )
method |
Installation method. By default, "auto" automatically finds a method that will work in the local environment. Change the default to force a specific installation method. Note that the "virtualenv" method is not available on Windows. |
conda |
The path to a |
extra_pack |
Character vector. Extra-packages to be installed. |
python_version |
The requested Python version. Ignored when attempting to install with a Python virtual environment. |
restart_session |
Restart R session after installing (note this will only occur within RStudio) |
envname |
The name, or full path, of the environment in which Python
packages are to be installed. When |
... |
additionnal arguments passed to |
No return value.
Beta and X must have the same number of variables
Beta and X must have the same number of variables
lagrangian(Beta, X, lambda1 = 0, lambda2 = 0) lagrangian(Beta, X, lambda1 = 0, lambda2 = 0)
lagrangian(Beta, X, lambda1 = 0, lambda2 = 0) lagrangian(Beta, X, lambda1 = 0, lambda2 = 0)
Beta |
numeric matrix. In rows, regression vectors coefficients following of node-wise regression. diag(Beta) = 0 |
X |
numeric matrix. Data with variables in columns. |
lambda1 |
numeric scalar. Lasso penalization parameter. |
lambda2 |
numeric scalar. Fused-group Lasso penalization parameter. |
numeric scalar. The lagrangian
numeric scalar. The lagrangian
## Generation of K block partitions n = 50 K = 3 p = 6 rho = 0.85 blocs <- list() for (j in 1:K) { bloc <- matrix(rho, nrow = p/K, ncol = p/K) for(i in 1:(p/K)) { bloc[i,i] <- 1 } blocs[[j]] <- bloc } mat.covariance <- bdiag(blocs) set.seed(11) X <- rmvnorm(n, mean = rep(0,p), sigma = as.matrix(mat.covariance)) X <- scale(X) ## Initialization for Beta Beta1 <- matrix(0, nrow = p, ncol = p) for(i in 1:p){ Beta1[i,-i] <- solve(t(X[,-i])%*%X[,-i]) %*% t(X[,-i]) %*% X[,i] } lagrangian(Beta, X, 0, 0) ## Generation of K block partitions n = 50 K = 3 p = 6 rho = 0.85 blocs <- list() for (j in 1:K) { bloc <- matrix(rho, nrow = p/K, ncol = p/K) for(i in 1:(p/K)) { bloc[i,i] <- 1 } blocs[[j]] <- bloc } mat.covariance <- bdiag(blocs) set.seed(11) X <- rmvnorm(n, mean = rep(0,p), sigma = as.matrix(mat.covariance)) X <- scale(X) ## Initialization for Beta Beta1 <- matrix(0, nrow = p, ncol = p) for(i in 1:p){ Beta1[i,-i] <- solve(t(X[,-i])%*%X[,-i]) %*% t(X[,-i]) %*% X[,i] } lagrangian(Beta, X, 0, 0)
## Generation of K block partitions n = 50 K = 3 p = 6 rho = 0.85 blocs <- list() for (j in 1:K) { bloc <- matrix(rho, nrow = p/K, ncol = p/K) for(i in 1:(p/K)) { bloc[i,i] <- 1 } blocs[[j]] <- bloc } mat.covariance <- bdiag(blocs) set.seed(11) X <- rmvnorm(n, mean = rep(0,p), sigma = as.matrix(mat.covariance)) X <- scale(X) ## Initialization for Beta Beta1 <- matrix(0, nrow = p, ncol = p) for(i in 1:p){ Beta1[i,-i] <- solve(t(X[,-i])%*%X[,-i]) %*% t(X[,-i]) %*% X[,i] } lagrangian(Beta, X, 0, 0) ## Generation of K block partitions n = 50 K = 3 p = 6 rho = 0.85 blocs <- list() for (j in 1:K) { bloc <- matrix(rho, nrow = p/K, ncol = p/K) for(i in 1:(p/K)) { bloc[i,i] <- 1 } blocs[[j]] <- bloc } mat.covariance <- bdiag(blocs) set.seed(11) X <- rmvnorm(n, mean = rep(0,p), sigma = as.matrix(mat.covariance)) X <- scale(X) ## Initialization for Beta Beta1 <- matrix(0, nrow = p, ncol = p) for(i in 1:p){ Beta1[i,-i] <- solve(t(X[,-i])%*%X[,-i]) %*% t(X[,-i]) %*% X[,i] } lagrangian(Beta, X, 0, 0)
Check first estimate coeffs with glm
lasso_estimate(response_variable_number, penalty_value)
lasso_estimate(response_variable_number, penalty_value)
mean of randomly simulated precision matrices in the same configuration
mean_prec_mat(nrep = 10, config = config_)
mean_prec_mat(nrep = 10, config = config_)
Merge Beta Different types of merging and their effect
merge_beta(Beta, pair_to_merge, clusters)
merge_beta(Beta, pair_to_merge, clusters)
compute clusters partition from pairs of variables to merge
merge_clusters(pairs_to_merge, clusters)
merge_clusters(pairs_to_merge, clusters)
pairs_to_merge |
table of the indices of variables to be merge |
clusters |
numeric vector. By default 1:p where p is the number of variables |
A numeric vector.
Merge labels
merge_labels(merged_pair, labels, level)
merge_labels(merged_pair, labels, level)
merge clusters from table
merge_proc( pairs_to_merge, clusters, X, Beta, level, gain_level, gains, labels, merge )
merge_proc( pairs_to_merge, clusters, X, Beta, level, gain_level, gains, labels, merge )
weighted mean
mergeX(X, pair_to_merge, clusters)
mergeX(X, pair_to_merge, clusters)
Cluster variables using L2 fusion penalty and simultaneously estimates a gaussian graphical model structure with the addition of L1 sparsity penalty.
mglasso( x, lambda1 = 0, fuse_thresh = 0.001, maxit = NULL, distance = c("euclidean", "relative"), lambda2_start = 1e-04, lambda2_factor = 1.5, precision = 0.01, weights_ = NULL, type = c("initial"), compact = TRUE, verbose = FALSE )
mglasso( x, lambda1 = 0, fuse_thresh = 0.001, maxit = NULL, distance = c("euclidean", "relative"), lambda2_start = 1e-04, lambda2_factor = 1.5, precision = 0.01, weights_ = NULL, type = c("initial"), compact = TRUE, verbose = FALSE )
x |
Numeric matrix ( |
lambda1 |
Positive numeric scalar. Lasso penalty. |
fuse_thresh |
Positive numeric scalar. Threshold for clusters fusion. |
maxit |
Integer scalar. Maximum number of iterations. |
distance |
Character. Distance between regression vectors with permutation on symmetric coefficients. |
lambda2_start |
Numeric scalar. Starting value for fused-group Lasso penalty (clustering penalty). |
lambda2_factor |
Numeric scalar. Step used to update fused-group Lasso penalty in a multiplicative way.. |
precision |
Tolerance for the stopping criterion (duality gap). |
weights_ |
Matrix of weights. |
type |
If "initial" use classical version of MGLasso without weights. |
compact |
Logical scalar. If TRUE, only save results when previous clusters are different from current. |
verbose |
Logical scalar. Print trace. Default value is FALSE. |
Estimates a gaussian graphical model structure while hierarchically grouping
variables by optimizing a pseudo-likelihood function combining Lasso and
fuse-group Lasso penalties. The problem is solved via the COntinuation
with NEsterov smoothing in a Shrinkage-Thresholding Algorithm (Hadj-Selem et
al. 2018). Varying the fusion penalty in a multiplicative
fashion allow to uncover a seemingly hierarchical clustering structure. For
, the approach is theoretically equivalent to the
Meinshausen-Bühlmann (2006) neighborhood selection as resuming to the
optimization of pseudo-likelihood function with
penalty
(Rocha et al., 2008). The algorithm stops when all the variables have merged
into one cluster. The criterion used to merge clusters is the
-norm distance between regression vectors.
For each iteration of the algorithm, the following function is optimized :
where is the vector of coefficients obtained after regression
on the others and
is a permutation exchanging
with
.
A list-like object of class mglasso
is returned.
out |
List of lists. Each element of the list corresponds to a
clustering level. An element named |
l1 |
the sparsity penalty |
conesta()
for the problem solver,
plot_mglasso()
for plotting the output of mglasso
.
## Not run: mglasso::install_pylearn_parsimony(envname = "rmglasso", method = "conda") reticulate::use_condaenv("rmglasso", required = TRUE) reticulate::py_config() n = 50 K = 3 p = 9 rho = 0.85 blocs <- list() for (j in 1:K) { bloc <- matrix(rho, nrow = p/K, ncol = p/K) for(i in 1:(p/K)) { bloc[i,i] <- 1 } blocs[[j]] <- bloc } mat.covariance <- Matrix::bdiag(blocs) mat.covariance set.seed(11) X <- mvtnorm::rmvnorm(n, mean = rep(0,p), sigma = as.matrix(mat.covariance)) X <- scale(X) res <- mglasso(X, 0.1, lambda2_start = 0.1) res$out[[1]]$clusters res$out[[1]]$beta ## End(Not run)
## Not run: mglasso::install_pylearn_parsimony(envname = "rmglasso", method = "conda") reticulate::use_condaenv("rmglasso", required = TRUE) reticulate::py_config() n = 50 K = 3 p = 9 rho = 0.85 blocs <- list() for (j in 1:K) { bloc <- matrix(rho, nrow = p/K, ncol = p/K) for(i in 1:(p/K)) { bloc[i,i] <- 1 } blocs[[j]] <- bloc } mat.covariance <- Matrix::bdiag(blocs) mat.covariance set.seed(11) X <- mvtnorm::rmvnorm(n, mean = rep(0,p), sigma = as.matrix(mat.covariance)) X <- scale(X) res <- mglasso(X, 0.1, lambda2_start = 0.1) res$out[[1]]$clusters res$out[[1]]$beta ## End(Not run)
neighbor_select
neighbor_select( data = data$X, config, lambda_min_ratio = 0.01, nlambda = 10, nresamples = 20, lambdas = NULL, model = NULL, verbose = FALSE, estim_var = NULL )
neighbor_select( data = data$X, config, lambda_min_ratio = 0.01, nlambda = 10, nresamples = 20, lambdas = NULL, model = NULL, verbose = FALSE, estim_var = NULL )
One simulation configuration
one_config(n, p, pi, alpha, rho)
one_config(n, p, pi, alpha, rho)
compute TPR, FPR, SHD given estimated and true precision matrices
perf_one(omega_hat, omega)
perf_one(omega_hat, omega)
SHD: structural hamming distance
get performances from list of estimations
perf_vec(omega_hat_list, omega)
perf_vec(omega_hat_list, omega)
Plot MGLasso Clusterpath
plot_clusterpath(X, mglasso_res, colnames_ = NULL)
plot_clusterpath(X, mglasso_res, colnames_ = NULL)
X |
numeric matrix |
mglasso_res |
object of class |
colnames_ |
columns labels |
This function plot the clustering path of mglasso method on the 2 principal components axis of X. As the centroids matrices are not of the same dimension as X, we choose to plot the predicted X matrix path.
no return value.
Plot the object returned by the mglasso
function.
plot_mglasso(mglasso_, levels_ = NULL) plot_mglasso(mglasso_, levels_ = NULL)
plot_mglasso(mglasso_, levels_ = NULL) plot_mglasso(mglasso_, levels_ = NULL)
mglasso_ |
Object of class |
levels_ |
Character vector. Selected levels for which estimated matrices will be plot. If NULL plot all levels. |
No return value.
Compute precision matrix from regression vectors
precision_to_regression(K)
precision_to_regression(K)
K |
precision matrix |
A numeric matrix.
pareil pour les clusters ACP dimensions ?
repart(cor_)
repart(cor_)
EBIC
select_ebic_weighted(Thetas, ploglik, n_edges, n, p, gam = 0.5, pen_params)
select_ebic_weighted(Thetas, ploglik, n_edges, n, p, gam = 0.5, pen_params)
K-fold cross validation neghborhood lasso selection
select_kfold( X, Thetas, lambdas = NULL, n_lambda, K_fold = 10, criterion = "ploglik", verbose = TRUE )
select_kfold( X, Thetas, lambdas = NULL, n_lambda, K_fold = 10, criterion = "ploglik", verbose = TRUE )
criterion |
c("ploglik", "rmse") |
if criterion = "ploglik" use pseudo-log likelihood formula from Huge paper with matrix approach if "rmse" use mean squarred prediction error
K-fold cross validation mglasso
select_kfold_mglasso( X, lambda1s = NULL, lambda2s = NULL, K_fold = 5, nl1 = 1, nl2 = 1, lam1_min_ratio, verbose = TRUE )
select_kfold_mglasso( X, lambda1s = NULL, lambda2s = NULL, K_fold = 5, nl1 = 1, nl2 = 1, lam1_min_ratio, verbose = TRUE )
Finds the optimal number of clusters using slope heuristic
select_partition(gains)
select_partition(gains)
integer scalar. The indice of the selected model.
stability selection mglasso
select_stab_mglasso(X, l1_, l2_, subsample_ratio, nrep, stab_thresh)
select_stab_mglasso(X, l1_, l2_, subsample_ratio, nrep, stab_thresh)
stability selection mglasso II stars way
select_stars_mglasso( X, lambda1s = NULL, lambda2s = NULL, subsample_ratio = NULL, nrep = 1, stars_thresh = 0.1, nl1 = 1, nl2 = 1 )
select_stars_mglasso( X, lambda1s = NULL, lambda2s = NULL, subsample_ratio = NULL, nrep = 1, stars_thresh = 0.1, nl1 = 1, nl2 = 1 )
def sequences for lambda1s and lambda2s not sure if max of lambda1 s still the same as in the lasso case. But if find better equivalence will update this part
seq_l1l2( X, nlam1 = 2, nlam2 = 2, logscale = TRUE, mean = FALSE, lambda1_min_ratio = 0.01, require_non_list = FALSE, l2_max = NULL )
seq_l1l2( X, nlam1 = 2, nlam2 = 2, logscale = TRUE, mean = FALSE, lambda1_min_ratio = 0.01, require_non_list = FALSE, l2_max = NULL )
mean |
in conesta_rwrapper is the mean criterion used ie averaged by np |
simulate data with given graph structure
sim_data( p = 20, np_ratio = 2, structure = c("block_diagonal", "hub", "scale_free", "erdos"), alpha, prob_mat, rho, g, inter_cluster_edge_prob = 0.01, p_erdos = 0.1, verbose = FALSE )
sim_data( p = 20, np_ratio = 2, structure = c("block_diagonal", "hub", "scale_free", "erdos"), alpha, prob_mat, rho, g, inter_cluster_edge_prob = 0.01, p_erdos = 0.1, verbose = FALSE )
verbose |
A list: graph : precision
symmetrize matrix of regression vectors pxp
Apply symmetrization on estimated graph
symmetrize(mat, rule = "and") symmetrize(mat, rule = "and")
symmetrize(mat, rule = "and") symmetrize(mat, rule = "and")
mat |
graph or precision matrix |
rule |
"and" or "or" rule |
A numeric matrix.