learn_bipartite_k_component_graph {spectralGraphTopology} | R Documentation |
Learns a bipartite k-component graph Jointly learns the Laplacian and Adjacency matrices of a graph on the basis of an observed data matrix
Description
Learns a bipartite k-component graph
Jointly learns the Laplacian and Adjacency matrices of a graph on the basis of an observed data matrix
Usage
learn_bipartite_k_component_graph(
S,
is_data_matrix = FALSE,
z = 0,
k = 1,
w0 = "naive",
m = 7,
alpha = 0,
beta = 10000,
rho = 0.01,
fix_beta = TRUE,
beta_max = 1e+06,
nu = 10000,
lb = 0,
ub = 10000,
maxiter = 10000,
abstol = 1e-06,
reltol = 1e-04,
eigtol = 1e-09,
record_weights = FALSE,
record_objective = FALSE,
verbose = TRUE
)
Arguments
S |
either a pxp sample covariance/correlation matrix, or a pxn data matrix, where p is the number of nodes and n is the number of features (or data points per node) |
is_data_matrix |
whether the matrix S should be treated as data matrix or sample covariance matrix |
z |
the number of zero eigenvalues for the Adjancecy matrix |
k |
the number of components of the graph |
w0 |
initial estimate for the weight vector the graph or a string selecting an appropriate method. Available methods are: "qp": finds w0 that minimizes ||ginv(S) - L(w0)||_F, w0 >= 0; "naive": takes w0 as the negative of the off-diagonal elements of the pseudo inverse, setting to 0 any elements s.t. w0 < 0 |
m |
in case is_data_matrix = TRUE, then we build an affinity matrix based on Nie et. al. 2017, where m is the maximum number of possible connections for a given node |
alpha |
L1 regularization hyperparameter |
beta |
regularization hyperparameter for the term ||L(w) - U Lambda U'||^2_F |
rho |
how much to increase (decrease) beta in case fix_beta = FALSE |
fix_beta |
whether or not to fix the value of beta. In case this parameter is set to false, then beta will increase (decrease) depending whether the number of zero eigenvalues is lesser (greater) than k |
beta_max |
maximum allowed value for beta |
nu |
regularization hyperparameter for the term ||A(w) - V Psi V'||^2_F |
lb |
lower bound for the eigenvalues of the Laplacian matrix |
ub |
upper bound for the eigenvalues of the Laplacian matrix |
maxiter |
the maximum number of iterations |
abstol |
absolute tolerance on the weight vector w |
reltol |
relative tolerance on the weight vector w |
eigtol |
value below which eigenvalues are considered to be zero |
record_weights |
whether to record the edge values at each iteration |
record_objective |
whether to record the objective function values at each iteration |
verbose |
whether to output a progress bar showing the evolution of the iterations |
Value
A list containing possibly the following elements:
laplacian |
the estimated Laplacian Matrix |
adjacency |
the estimated Adjacency Matrix |
w |
the estimated weight vector |
psi |
optimization variable accounting for the eigenvalues of the Adjacency matrix |
lambda |
optimization variable accounting for the eigenvalues of the Laplacian matrix |
V |
eigenvectors of the estimated Adjacency matrix |
U |
eigenvectors of the estimated Laplacian matrix |
elapsed_time |
elapsed time recorded at every iteration |
beta_seq |
sequence of values taken by beta in case fix_beta = FALSE |
convergence |
boolean flag to indicate whether or not the optimization converged |
obj_fun |
values of the objective function at every iteration in case record_objective = TRUE |
negloglike |
values of the negative loglikelihood at every iteration in case record_objective = TRUE |
w_seq |
sequence of weight vectors at every iteration in case record_weights = TRUE |
Author(s)
Ze Vinicius and Daniel Palomar
References
S. Kumar, J. Ying, J. V. M. Cardoso, D. P. Palomar. A unified framework for structured graph learning via spectral constraints. Journal of Machine Learning Research, 2020. http://jmlr.org/papers/v21/19-276.html
Examples
library(spectralGraphTopology)
library(igraph)
library(viridis)
library(corrplot)
set.seed(42)
w <- c(1, 0, 0, 1, 0, 1) * runif(6)
Laplacian <- block_diag(L(w), L(w))
Atrue <- diag(diag(Laplacian)) - Laplacian
bipartite <- graph_from_adjacency_matrix(Atrue, mode = "undirected", weighted = TRUE)
n <- ncol(Laplacian)
Y <- MASS::mvrnorm(40 * n, rep(0, n), MASS::ginv(Laplacian))
graph <- learn_bipartite_k_component_graph(cov(Y), k = 2, beta = 1e2, nu = 1e2, verbose = FALSE)
graph$adjacency[graph$adjacency < 1e-2] <- 0
# Plot Adjacency matrices: true, noisy, and estimated
corrplot(Atrue / max(Atrue), is.corr = FALSE, method = "square", addgrid.col = NA, tl.pos = "n",
cl.cex = 1.25)
corrplot(graph$adjacency / max(graph$adjacency), is.corr = FALSE, method = "square",
addgrid.col = NA, tl.pos = "n", cl.cex = 1.25)
# Plot networks
estimated_bipartite <- graph_from_adjacency_matrix(graph$adjacency, mode = "undirected",
weighted = TRUE)
V(bipartite)$type <- rep(c(TRUE, FALSE), 4)
V(estimated_bipartite)$type <- rep(c(TRUE, FALSE), 4)
la = layout_as_bipartite(estimated_bipartite)
colors <- viridis(20, begin = 0, end = 1, direction = -1)
c_scale <- colorRamp(colors)
E(estimated_bipartite)$color = apply(
c_scale(E(estimated_bipartite)$weight / max(E(estimated_bipartite)$weight)), 1,
function(x) rgb(x[1]/255, x[2]/255, x[3]/255))
E(bipartite)$color = apply(c_scale(E(bipartite)$weight / max(E(bipartite)$weight)), 1,
function(x) rgb(x[1]/255, x[2]/255, x[3]/255))
la = la[, c(2, 1)]
# Plot networks: true and estimated
plot(bipartite, layout = la,
vertex.color = c("red","black")[V(bipartite)$type + 1],
vertex.shape = c("square", "circle")[V(bipartite)$type + 1],
vertex.label = NA, vertex.size = 5)
plot(estimated_bipartite, layout = la,
vertex.color = c("red","black")[V(estimated_bipartite)$type + 1],
vertex.shape = c("square", "circle")[V(estimated_bipartite)$type + 1],
vertex.label = NA, vertex.size = 5)