Title: | Covariate-varying Networks |
---|---|
Description: | A package for inferring high-dimensional Gaussian graphical networks that change with multiple discrete covariates |
Authors: | Louis Dijkstra [aut, cre] |
Maintainer: | Louis Dijkstra <[email protected]> |
License: | GPL (>= 3) |
Version: | 1.1 |
Built: | 2025-01-04 06:12:55 UTC |
Source: | https://github.com/bips-hb/CVN |
Checks whether the input for the function CVN
is valid.
This function does not return anything. The execution of the function
halts when an issue has been detected.
check_correctness_input(raw_data, W, lambda1, lambda2, gamma1, gamma2, rho)
check_correctness_input(raw_data, W, lambda1, lambda2, gamma1, gamma2, rho)
raw_data |
A list with matrices. The number of columns should be the same for each matrix |
W |
The |
lambda1 |
A vector of |
lambda2 |
A vector of |
gamma1 |
A vector of |
gamma2 |
A vector of |
rho |
The |
data.frame
for the Edges for visNetwork
In order to visualize a graph, we need to create a
data.frame
that can be used by the visNetwork
package.
This function returns the needed data.frame
given a
adjacency matrix.
create_edges_visnetwork(adj_matrix)
create_edges_visnetwork(adj_matrix)
adj_matrix |
A symmetric adjacency matrix |
Data frame that be used as input for visNetwork
to be used for the Generalized LASSOGenerates a matrix to be used for the generalized LASSO.
We solve a generalized LASSO problem for each edge
for each update step for
.
create_matrix_D(W, lambda1, lambda2, rho = 1, remove_zero_row = TRUE)
create_matrix_D(W, lambda1, lambda2, rho = 1, remove_zero_row = TRUE)
W |
The |
lambda1 |
The |
lambda2 |
The |
rho |
The |
remove_zero_row |
If |
A -dimensional matrix
Tibshirani, R. J., & Taylor, J. (2011). The solution path of the generalized lasso. Annals of Statistics, 39(3), 1335–1371. https://doi.org/10.1214/11-AOS878
m <- 4 # number of graphs W <- matrix(1, nrow = m, ncol = m) # penalty terms: lambda1 <- .2 lambda2 <- .4 rho <- 1 CVN::create_matrix_D(W, lambda1, lambda2, rho) # [,1] [,2] [,3] [,4] # [1,] 0.2 0.0 0.0 0.0 # [2,] 0.0 0.2 0.0 0.0 # [3,] 0.0 0.0 0.2 0.0 # [4,] 0.0 0.0 0.0 0.2 # [5,] 0.4 -0.4 0.0 0.0 # [6,] 0.4 0.0 -0.4 0.0 # [7,] 0.4 0.0 0.0 -0.4 # [8,] 0.0 0.4 -0.4 0.0 # [9,] 0.0 0.4 0.0 -0.4 # [10,] 0.0 0.0 0.4 -0.4
m <- 4 # number of graphs W <- matrix(1, nrow = m, ncol = m) # penalty terms: lambda1 <- .2 lambda2 <- .4 rho <- 1 CVN::create_matrix_D(W, lambda1, lambda2, rho) # [,1] [,2] [,3] [,4] # [1,] 0.2 0.0 0.0 0.0 # [2,] 0.0 0.2 0.0 0.0 # [3,] 0.0 0.0 0.2 0.0 # [4,] 0.0 0.0 0.0 0.2 # [5,] 0.4 -0.4 0.0 0.0 # [6,] 0.4 0.0 -0.4 0.0 # [7,] 0.4 0.0 0.0 -0.4 # [8,] 0.0 0.4 -0.4 0.0 # [9,] 0.0 0.4 0.0 -0.4 # [10,] 0.0 0.0 0.4 -0.4
visNetwork
packageCreates a data frame that can be used for the
visNetwork
package.
create_nodes_visnetwork(n_nodes, labels = 1:n_nodes)
create_nodes_visnetwork(n_nodes, labels = 1:n_nodes)
n_nodes |
Number of nodes in the graph |
labels |
The labels for the individual nodes
(Default: |
Data frame with two columns: id
and title
This function generates different weight matrices, .
There are several types:
full
All graphs are fully connected with weight 1
glasso
All graphs are disconnected with weight 0. This
mimicks the GLASSO, where each graph is estimated independently
grid
A weight matrix for a grid
uniform-random
Fully-connected, but the entries are drawn
from a uniform distribution
create_weight_matrix( type = c("full", "glasso", "grid", "uniform-random"), m = 9 )
create_weight_matrix( type = c("full", "glasso", "grid", "uniform-random"), m = 9 )
type |
The type of weight matrix |
m |
Number of graphs |
Weight matrix
Estimates a covariate-varying network model (CVN), i.e.,
Gaussian graphical models that change with (multiple) external covariate(s).
The smoothing between the graphs is specified by the
-dimensional
weight matrix
. The function returns the estimated precision matrices
for each graph.
CVN( data, W, lambda1 = 1:2, lambda2 = 1:2, gamma1 = NULL, gamma2 = NULL, rho = 1, eps = 1e-04, maxiter = 100, truncate = 1e-05, rho_genlasso = 1, eps_genlasso = 1e-10, maxiter_genlasso = 100, truncate_genlasso = 1e-04, n_cores = min(length(lambda1) * length(lambda2), parallel::detectCores() - 1), normalized = FALSE, warmstart = TRUE, minimal = FALSE, verbose = TRUE )
CVN( data, W, lambda1 = 1:2, lambda2 = 1:2, gamma1 = NULL, gamma2 = NULL, rho = 1, eps = 1e-04, maxiter = 100, truncate = 1e-05, rho_genlasso = 1, eps_genlasso = 1e-10, maxiter_genlasso = 100, truncate_genlasso = 1e-04, n_cores = min(length(lambda1) * length(lambda2), parallel::detectCores() - 1), normalized = FALSE, warmstart = TRUE, minimal = FALSE, verbose = TRUE )
data |
A list with matrices, each entry associated with a single graph. The number of columns should be the same for each matrix. Number of observations can differ |
W |
The |
lambda1 |
Vector with different |
lambda2 |
Vector with different |
gamma1 |
A vector of |
gamma2 |
A vector of |
rho |
The |
eps |
If the relative difference between two update steps is
smaller than |
maxiter |
Maximum number of iterations (Default: |
truncate |
All values of the final |
rho_genlasso |
The |
eps_genlasso |
If the relative difference between two update steps is
smaller than |
maxiter_genlasso |
Maximum number of iterations for solving
the generalized LASSO problem (Default: |
truncate_genlasso |
All values of the final |
n_cores |
Number of cores used (Default: max. number of cores - 1, or the total number penalty term pairs if that is less) |
normalized |
Data is normalized if |
warmstart |
If |
minimal |
If |
verbose |
Verbose (Default: |
A CVN
object containing the estimates for all the graphs
for each different value of . General results for
the different values of
can be found in the data frame
results
. It consists of multiple columns, namely:
lambda1 |
|
lambda2 |
|
converged |
whether algorithm converged or not |
value |
value of the negative log-likelihood function |
n_iterations |
number of iterations of the ADMM |
aic |
Aikake information criterion |
gamma1 |
|
gamma2 |
|
id |
The id. This corresponds to the indices of the lists |
bic |
Bayesian information criterion |
The estimates of the precision matrices and the corresponding adjacency matrices
for the different values of can be found
Theta |
A list with the estimated precision matrices |
adj_matrices |
A list with the estimated adjacency matrices corresponding to the
estimated precision matrices in |
In addition, the input given to the CVN function is stored in the object as well:
Sigma |
Empirical covariance matrices |
m |
Number of graphs |
p |
Number of variables |
n_obs |
Vector of length |
data |
The |
W |
The |
maxiter |
Maximum number of iterations for the ADMM |
rho |
The |
eps |
The stopping criterion |
truncate |
Truncation value for |
maxiter_genlasso |
Maximum number of iterations for the generarlzed LASSO |
rho_genlasso |
The |
eps_genlasso |
The stopping criterion |
truncate_genlasso |
Truncation value for |
n_lambda_values |
Total number of |
normalized |
If |
warmstart |
If |
minimal |
If |
hits_border_aic |
If |
hits_border_bic |
If |
When estimating the graph for different values of
and
, we use the graph estimated (if available)
for other
and
values closest to them.
data(grid) m <- 9 # must be 9 for this example #' Choice of the weight matrix W. #' (uniform random) W <- matrix(runif(m*m), ncol = m) W <- W %*% t(W) W <- W / max(W) diag(W) <- 0 # lambdas: lambda1 = 1:2 lambda2 = 1:2 (cvn <- CVN::CVN(grid, W, lambda1 = lambda1, lambda2 = lambda2, eps = 1e-3, maxiter = 1000, verbose = TRUE))
data(grid) m <- 9 # must be 9 for this example #' Choice of the weight matrix W. #' (uniform random) W <- matrix(runif(m*m), ncol = m) W <- W %*% t(W) W <- W / max(W) diag(W) <- 0 # lambdas: lambda1 = 1:2 lambda2 = 1:2 (cvn <- CVN::CVN(grid, W, lambda1 = lambda1, lambda2 = lambda2, eps = 1e-3, maxiter = 1000, verbose = TRUE))
cvn
objectDetermines a given information criteria for a cvn
object,
see CVN
.
determine_information_criterion_cvn(cvn, type = c("AIC", "BIC"))
determine_information_criterion_cvn(cvn, type = c("AIC", "BIC"))
A function for estimating a CVN for a single -value.
See for more details
CVN
estimate( m, p, W, Theta0, Z0, Y0, a, eta1, eta2, Sigma, n_obs, rho, rho_genlasso, eps, eps_genlasso, maxiter, maxiter_genlasso, truncate, truncate_genlasso, verbose = FALSE )
estimate( m, p, W, Theta0, Z0, Y0, a, eta1, eta2, Sigma, n_obs, rho, rho_genlasso, eps, eps_genlasso, maxiter, maxiter_genlasso, truncate, truncate_genlasso, verbose = FALSE )
Finds the 'core graph', i.e., those edges that are present in all of the estimated graphs
find_core_graph(cvn)
find_core_graph(cvn)
cvn |
A |
A list of adjacency matrix, one for each value of (lambda1, lambda2)
genlassoRcpp
See for details genlassoRcpp
genlasso_wrapper(y, W, m, c, eta1, eta2, a, rho, max_iter, eps, truncate)
genlasso_wrapper(y, W, m, c, eta1, eta2, a, rho, max_iter, eps, truncate)
Solves efficiently the generalized LASSO problem of the form
where and
are
-dimensional vectors and
is a
-matrix where
.
We solve this optimization problem using an adaption of the ADMM
algorithm presented in Zhu (2017).
genlassoRcpp(Y, W, m, eta1, eta2, a, rho, max_iter, eps, truncate)
genlassoRcpp(Y, W, m, eta1, eta2, a, rho, max_iter, eps, truncate)
W |
The weight matrix |
m |
The number of graphs |
eta1 |
Equals |
eta2 |
Equals |
a |
Value added to the diagonal of |
rho |
The ADMM's parameter |
max_iter |
Maximum number of iterations |
eps |
Stopping criterion. If differences
are smaller than |
truncate |
Values below |
y |
The |
The estimated vector
Zhu, Y. (2017). An Augmented ADMM Algorithm With Application to the Generalized Lasso Problem. Journal of Computational and Graphical Statistics, 26(1), 195–204. https://doi.org/10.1080/10618600.2015.1114491
A wrapper for the GLASSO in the context of CVNs. Each graph
is estimated individually. There is NO smoothing between the graphs.
This function relies completely on the glasso
package.
The output is, therefore, slightly different than for the
CVN
function.
glasso( data, lambda1 = 1:2, eps = 1e-04, maxiter = 10000, n_cores = min(length(lambda1), parallel::detectCores() - 1), normalized = FALSE, verbose = TRUE )
glasso( data, lambda1 = 1:2, eps = 1e-04, maxiter = 10000, n_cores = min(length(lambda1), parallel::detectCores() - 1), normalized = FALSE, verbose = TRUE )
data |
A list with matrices, each entry associated with a single graph. The number of columns should be the same for each matrix. Number of observations can differ |
lambda1 |
Vector with different |
eps |
Threshold for convergence (Default: |
maxiter |
Maximum number of iterations (Default: 10,000) |
n_cores |
Number of cores used (Default: max. number of cores - 1, or the total number penalty term pairs if that is less) |
normalized |
Data is normalized if |
verbose |
Verbose (Default: |
A CVN
object containing the estimates for all the graphs
for different value of . General results for
the different value of
can be found in the data frame
results
. It consists of multiple columns, namely:
lambda1 |
|
value |
value of the negative log-likelihood function |
aic |
Aikake information criteration |
id |
The id. This corresponds to the indices of the lists |
The estimates of the precision matrices and the corresponding adjacency matrices
for the different values of can be found
Theta |
A list with the estimated precision matrices |
adj_matrices |
A list with the estimated adjacency matrices corresponding to the
estimated precision matrices in |
In addition, the input given to this function is stored in the object as well:
Sigma |
Empirical covariance matrices |
m |
Number of graphs |
p |
Number of variables |
n_obs |
Vector of length |
data |
The |
maxiter |
Maximum number of iterations for the ADMM |
eps |
The stopping criterion |
n_lambda_values |
Total number of |
normalized |
If |
data(grid) m <- 9 # must be 9 for this example #' Choice of the weight matrix W. #' (uniform random) W <- matrix(runif(m*m), ncol = m) W <- W %*% t(W) W <- W / max(W) diag(W) <- 0 # lambdas: lambda1 = 1:4 (glasso_est <- CVN::glasso(grid, lambda1 = lambda1))
data(grid) m <- 9 # must be 9 for this example #' Choice of the weight matrix W. #' (uniform random) W <- matrix(runif(m*m), ncol = m) W <- W %*% t(W) W <- W / max(W) diag(W) <- 0 # lambdas: lambda1 = 1:4 (glasso_est <- CVN::glasso(grid, lambda1 = lambda1))
Data generated for 9 graphs in total, organized in a grid of
(3x3). See the package CVNSim
for more information
on how the grid is constructed: https://github.com/bips-hb/CVNSim
data(grid)
data(grid)
List
https://github.com/bips-hb/CVNSim
cvn
ObjectReturns the structural Hamming distances
hamming_distance(cvn, verbose = TRUE)
hamming_distance(cvn, verbose = TRUE)
cvn |
A |
verbose |
If |
A list of symmetric matrices. Each matrix contains the structural
Hamming distances between the different graphs. Each item in the
list corresponds to one pair
Returns the structural Hamming distance between multiple graphs
hamming_distance_adj_matrices(adj_matrices)
hamming_distance_adj_matrices(adj_matrices)
adj_matrices |
A list of adjacency matrices |
Matrix of Hamming distances
intervalOne often selected the optimal model for the -values
based on the AIC and BIC.
This function checks and warns when the optimal value lies on the border of the
values
takes.
hits_end_lambda_intervals(results)
hits_end_lambda_intervals(results)
results |
Results of the |
List with two values:
hits_border_aic |
If |
hits_border_bic |
If |
Estimates a graph for which there are no observation based on a previously fitted CVN model
interpolate(cvn, weights, truncate = NULL)
interpolate(cvn, weights, truncate = NULL)
cvn |
A CVN fit with |
weights |
A vector of length |
truncate |
Truncation value. When a value in the precision matrix is
considered 0. If |
A list with
adj_matrices |
A list of adjacency matrix. One for each pair
of |
m |
Number of graphs |
p |
Number of variables |
weights |
The weights used for interpolation |
truncate |
Truncation value |
n_lambda_values |
Total number of |
results
. It consists of two columns:
lambda1 |
|
lambda2 |
|
for inner-ADMM for the
-update stepThe -update step, see
updateZ
, requires
us to solve a special Generalized LASSO problem of the form
where and
are
-dimensional vectors and
is a
-matrix where
.
We solve this optimization problem using an adaption of the ADMM
algorithm presented in Zhu (2017).
This algorithm requires the choice of a matrix
such that
is positive semidefinite. In order to optimize the ADMM,
we choose the matrix
to be diagonal with a fixed value
.
This function determines the smallest value of
such that
is indeed positive semidefinite. We do this be determining
the largest eigenvalue
matrix_A_inner_ADMM(W, eta1, eta2)
matrix_A_inner_ADMM(W, eta1, eta2)
W |
Weight matrix |
eta1 , eta2
|
The values |
Value of
Zhu, Y. (2017). An Augmented ADMM Algorithm With Application to the Generalized Lasso Problem. Journal of Computational and Graphical Statistics, 26(1), 195–204. https://doi.org/10.1080/10618600.2015.1114491
Returns a heat map of the distance matrix for a particular CVN
plot_hamming_distances( distance_matrix, absolute = TRUE, limits = c(NA, NA), title = "", legend_label = "Hamming Distance", add_counts_to_cells = TRUE, add_ticks_labels = TRUE, t = -6, r = -8 )
plot_hamming_distances( distance_matrix, absolute = TRUE, limits = c(NA, NA), title = "", legend_label = "Hamming Distance", add_counts_to_cells = TRUE, add_ticks_labels = TRUE, t = -6, r = -8 )
distance_matrix |
Symmetric matrix with distances |
absolute |
If |
limits |
The limits for the values of the Hamming distance |
title |
Title plot (Default is none) |
legend_label |
Title of the legend (Default: "Hamming Distance") |
add_counts_to_cells |
If |
add_ticks_labels |
If |
t |
Distance between tick labels and x-axis (Default: -6) |
r |
Distance between tick labels and y-axis (Default: -8) |
A heatmap plot
Creates all the heatmaps for a CVN, a heatmap for each
pair of
plot_hamming_distances_cvn( cvn, absolute = TRUE, same_range = TRUE, titles = rep("", cvn$n_lambda_values), legend_label = "Hamming Distance", add_counts_to_cells = TRUE, add_ticks_labels = TRUE, t = -6, r = -8, verbose = TRUE )
plot_hamming_distances_cvn( cvn, absolute = TRUE, same_range = TRUE, titles = rep("", cvn$n_lambda_values), legend_label = "Hamming Distance", add_counts_to_cells = TRUE, add_ticks_labels = TRUE, t = -6, r = -8, verbose = TRUE )
cvn |
A |
absolute |
If |
same_range |
If |
titles |
Title of the plots (Default is none) |
legend_label |
Title of the legend (Default: "Hamming Distance") |
add_counts_to_cells |
If |
add_ticks_labels |
If |
t |
Distance between tick labels and x-axis (Default: -6) |
r |
Distance between tick labels and y-axis (Default: -8) |
verbose |
If |
List of plots
Returns a heat map of the AIC or BIC for a fitted CVN
plot_information_criterion( cvn, criterion = c("bic", "aic"), use_gammas = TRUE, show_minimum = TRUE, title = "", xlabel = NULL, ylabel = NULL, legend_label = NULL, limits = c(NA, NA) )
plot_information_criterion( cvn, criterion = c("bic", "aic"), use_gammas = TRUE, show_minimum = TRUE, title = "", xlabel = NULL, ylabel = NULL, legend_label = NULL, limits = c(NA, NA) )
cvn |
Fitted CVN, see |
criterion |
The information criterion, must be either |
use_gammas |
If |
show_minimum |
If |
title |
Title plot (Default is none) |
xlabel |
Label for the |
ylabel |
Label for the |
legend_label |
Title for the legend. Default depends on |
limits |
The limits for the values of the Hamming distance |
A heatmap plot
Returns a heat map of a weight matrix
plot_weight_matrix( W, title = "", legend_label = "weight", add_counts_to_cells = TRUE, add_ticks_labels = TRUE, t = -6, r = -8 )
plot_weight_matrix( W, title = "", legend_label = "weight", add_counts_to_cells = TRUE, add_ticks_labels = TRUE, t = -6, r = -8 )
W |
Symmetric weight matrix |
title |
Title plot (Default is none) |
legend_label |
Title of the legend (Default: "weight") |
add_counts_to_cells |
If |
add_ticks_labels |
If |
t |
Distance between tick labels and x-axis (Default: -6) |
r |
Distance between tick labels and y-axis (Default: -8) |
A heatmap plot
Plot Function for CVN Object Class
## S3 method for class 'cvn' plot(x, ...)
## S3 method for class 'cvn' plot(x, ...)
Print Function for the CVN Object Class
## S3 method for class 'cvn' print(x, ...)
## S3 method for class 'cvn' print(x, ...)
Returns the relative difference between precision matrix
(parameter
Theta_new
) and (parameter
Theta_old
).
relative_difference_precision_matrices(Theta_new, Theta_old)
relative_difference_precision_matrices(Theta_new, Theta_old)
Theta_new |
A list with matrices with the updated values of |
Theta_old |
A list with matrices with the old values of |
This is used for checking whether the stopping condition has been met.
The relative difference between and
visNetwork
A subset of edges can be assign a different thickness or color.
set_attributes_to_edges_visnetwork( edges, subset_edges, width = c(NA, NA), color = c(NULL, NULL) )
set_attributes_to_edges_visnetwork( edges, subset_edges, width = c(NA, NA), color = c(NULL, NULL) )
edges |
A data.frame create by |
subset_edges |
A list with the elements |
width |
Vector with two values. The first is assigned to the
edges in the subset given by |
color |
Vector with two values. The first is assigned to the
edges in the subset given by |
A data frame that can be used by the visNetwork
package
Function that removes most of the items to make the CVN object more memory sufficient. This is especially important when the graphs are rather larger
strip_cvn(cvn)
strip_cvn(cvn)
cvn |
|
Reduced cvn where Theta
, data
and Sigma
are removed
-update step for the ADMMupdateTheta
returns the updated value of for the
ADMM given the previously updated values of
and
updateTheta(m, Z, Y, Sigma, n_obs, rho = 1)
updateTheta(m, Z, Y, Sigma, n_obs, rho = 1)
m |
Number of graphs |
Z |
A list with matrices with the current values of |
Y |
A list with matrices with the current values of |
Sigma |
A list with empirical covariance matrices |
n_obs |
A |
rho |
The |
A list with matrices with the new values of
-update step for the ADMMReturns the updated value of for the ADMM given the previously
updated values of
and
updateY(Theta, Z, Y)
updateY(Theta, Z, Y)
Theta |
A list with matrices with the current values of |
Z |
A list with matrices with the current values of |
Y |
A list with matrices with the current values of |
A list with matrices with the new values of
-update Step for the ADMMA wrapper for the C
function that returns the
updated value of for the ADMM given the previously
updated values of
and
updateZ_wrapper( m, p, nrow_D, Theta, Y, W, eta1, eta2, a, rho_genlasso, maxiter_genlasso, eps_genlasso, truncate_genlasso )
updateZ_wrapper( m, p, nrow_D, Theta, Y, W, eta1, eta2, a, rho_genlasso, maxiter_genlasso, eps_genlasso, truncate_genlasso )
m |
Number of graphs |
p |
Number of variables |
nrow_D |
Number of rows of the |
Theta |
A list with matrices with the current values of |
Y |
A list with matrices with the current values of |
W |
Weight matrix |
eta1 |
TODO |
rho_genlasso |
The |
maxiter_genlasso |
Maximum number of iterations for solving the generalized LASSO problem |
eps_genlasso |
If the relative difference between two update steps is
smaller than |
truncate_genlasso |
All values of the final |
A list with matrices with the new values of
-update StepA C
implementation of the -update step. We
solve a generalized LASSO problem repeatedly for each of the
individual edges
updateZRcpp(m, p, Theta, Y, W, eta1, eta2, a, rho, max_iter, eps, truncate)
updateZRcpp(m, p, Theta, Y, W, eta1, eta2, a, rho, max_iter, eps, truncate)
m |
The number of graphs |
p |
The number of variables |
Theta |
A list of matrices with the |
Y |
A list of matrices with the |
eta1 |
Equals |
eta2 |
Equals |
a |
Value added to the diagonal of |
rho |
The ADMM's parameter |
max_iter |
Maximum number of iterations |
eps |
Stopping criterion. If differences
are smaller than |
truncate |
Values below |
The estimated vector
visNetwork
plotCreates a visNetwork
plot given a list of
nodes and edges. The nodes data frame can be
created with create_nodes_visnetwork
;
the edges with create_edges_visnetwork
.
In order to highlight edges, you can use
set_attributes_to_edges_visnetwork
.
visnetwork( nodes, edges, node_titles = nodes$id, title = "", igraph_layout = "layout_in_circle" )
visnetwork( nodes, edges, node_titles = nodes$id, title = "", igraph_layout = "layout_in_circle" )
A visNetwork
plot
nodes <- CVN::create_nodes_visnetwork(n_nodes = 5, labels = LETTERS[1:5]) adj_matrix <- matrix(c(0, 1, 0, 1, 0, 1, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1, 0), ncol = 5) edges <- CVN::create_edges_visnetwork(adj_matrix) edges <- set_attributes_to_edges_visnetwork(edges, subset_edges = list(from = c(1, 2), to = c(4, 3)), width = c(3, .5), color = c("red", "blue")) CVN::visnetwork(nodes, edges)
nodes <- CVN::create_nodes_visnetwork(n_nodes = 5, labels = LETTERS[1:5]) adj_matrix <- matrix(c(0, 1, 0, 1, 0, 1, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1, 0), ncol = 5) edges <- CVN::create_edges_visnetwork(adj_matrix) edges <- set_attributes_to_edges_visnetwork(edges, subset_edges = list(from = c(1, 2), to = c(4, 3)), width = c(3, .5), color = c("red", "blue")) CVN::visnetwork(nodes, edges)
visNetwork
plots for a CVN objectCreates all visNetwork
plots, see CVN::visnetwork
,
for all graphs in a cvn
object
visnetwork_cvn( cvn, node_titles = 1:cvn$p, titles = lapply(1:cvn$n_lambda_values, function(i) sapply(1:cvn$m, function(j) "")), show_core_graph = TRUE, width = c(3, 1), color = c("red", "blue"), igraph_layout = "layout_in_circle", verbose = TRUE )
visnetwork_cvn( cvn, node_titles = 1:cvn$p, titles = lapply(1:cvn$n_lambda_values, function(i) sapply(1:cvn$m, function(j) "")), show_core_graph = TRUE, width = c(3, 1), color = c("red", "blue"), igraph_layout = "layout_in_circle", verbose = TRUE )
cvn |
A |
node_titles |
Vector with title of the nodes (Default: |
titles |
A list with |
show_core_graph , width , color
|
Show the core graph using the width and colors |
List