The SCG functions provide a framework, called Spectral Coarse Graining (SCG), for reducing large graphs while preserving their spectral-related features, that is features closely related with the eigenvalues and eigenvectors of a graph matrix (which for now can be the adjacency, the stochastic, or the Laplacian matrix).
Common examples of such features comprise the first-passage-time of random walkers on Markovian graphs, thermodynamic properties of lattice models in statistical physics (e.g. Ising model), and the epidemic threshold of epidemic network models (SIR and SIS models).
SCG differs from traditional clustering schemes by producing a coarse-grained graph (not just a partition of the vertices), representative of the original one. As shown in [1], Principal Component Analysis can be viewed as a particular SCG, called exact SCG, where the matrix to be coarse-grained is the covariance matrix of some data set.
SCG should be of interest to practitioners of various fields dealing with problems where matrix eigenpairs play an important role, as for instance is the case of dynamical processes on networks.
The main idea of SCG is to operate on a matrix a shrinkage operation specifically designed to preserve some of the matrix eigenpairs while not altering other important matrix features (such as its structure). Mathematically, this idea was expressed as follows. Consider a (complex) n x n matrix M and form the product
M'=LMR*,
where n' < n and L, R are from C[n'xn]} and are such that LR*=I[n'] (R* denotes the conjugate transpose of R). Under these assumptions, it can be shown that P=R*L is an n'-rank projector and that, if (lambda, v) is a (right) eigenpair of M (i.e. Mv=lambda v} and P is orthogonal, there exists an eigenvalue lambda' of M' such that
|lambda-lambda'| <= const ||e[P](v)|| [1+O(||e[P](v)||2)],
where ||e[P](v)||=||v-Pv||. Hence, if P (or equivalently L, R) is chosen so as to make ||e[P](v)|| as small as possible, one can preserve to any desired level the original eigenvalue lambda in the coarse-grained matrix M'; under extra assumptions on M, this result can be generalized to eigenvectors [1]. This leads to the following generic definition of a SCG problem.
Given M (C[nxn]) and (lambda, v), a (right) eigenpair of M to be preserved by the coarse graining, the problem is to find a projector P' solving
min(||e[P](v)||, p in Omega),
where Omega is a set of projectors in C[nxn] described by some ad hoc constraints c[1], ..., c[r] (e.g. c[1]: P in R[nxn], c[2]: P=t(P), c[3]: P[i,j] >= 0}, etc).
Choosing pertinent constraints to solve the SCG problem is of great importance in applications. For instance, in the absence of constraints the SCG problem is solved trivially by P'=vv* (v is assumed normalized). We have designed a particular constraint, called homogeneous mixing, which ensures that vertices belonging to the same group are merged consistently from a physical point of view (see [1] for details). Under this constraint the SCG problem reduces to finding the partition of 1, ..., n (labeling the original vertices) minimizing
||e[P](v)||2 = sum([v(i)-(Pv)(i)]2; alpha=1,...,n', i in alpha),
where alpha denotes a group (i.e. a block) in a partition of {1, ..., n}, and |alpha| is the number of elements in alpha.
If M is symmetric or stochastic, for instance, then it may be desirable (or mandatory) to choose L, R so that M' is symmetric or stochastic as well. This structural constraint has led to the construction of particular semi-projectors for symmetric [1], stochastic [3] and Laplacian [2] matrices, that are made available.
In short, the coarse graining of matrices and graphs involves:
Retrieving a matrix or a graph matrix M from the problem.
Computing the eigenpairs of M to be preserved in the coarse-grained graph or matrix.
Setting some problem-specific constraints (e.g. dimension of the coarse-grained object).
Solving the constrained SCG problem, that is finding P'.
Computing from P' two semi-projectors L' and R' (e.g. following the method proposed in [1]).
Working out the product M'=L'MR'* and, if needed, defining from M' a coarse-grained graph.
The main functions are igraph_scg_adjacency()
, igraph_scg_laplacian()
and igraph_scg_stochastic()
.
These functions handle all the steps involved in the
Spectral Coarse Graining (SCG) of some particular matrices and graphs
as described above and in reference [1]. In more details,
they compute some prescribed eigenpairs of a matrix or a
graph matrix, (for now adjacency, Laplacian and stochastic matrices are
available), work out an optimal partition to preserve the eigenpairs,
and finally output a coarse-grained matrix or graph along with other
useful information.
These steps can also be carried out independently: (1) Use
igraph_get_adjacency()
, igraph_get_sparsemat()
,
igraph_laplacian()
, igraph_get_stochastic()
or igraph_get_stochastic_sparsemat()
to compute a matrix M.
(2) Work out some prescribed eigenpairs of M e.g. by
means of igraph_arpack_rssolve()
or igraph_arpack_rnsolve()
. (3) Invoke one the four
algorithms of the function igraph_scg_grouping()
to get a
partition that will preserve the eigenpairs in the coarse-grained
matrix. (4) Compute the semi-projectors L and R using
igraph_scg_semiprojectors()
and from there the coarse-grained
matrix M'=LMR*. If necessary, construct a coarse-grained graph from
M' (e.g. as in [1]).
[1] D. Morton de Lachapelle, D. Gfeller, and P. De Los Rios, Shrinking Matrices while Preserving their Eigenpairs with Application to the Spectral Coarse Graining of Graphs. Submitted to SIAM Journal on Matrix Analysis and Applications, 2008. http://people.epfl.ch/david.morton
[2] D. Gfeller, and P. De Los Rios, Spectral Coarse Graining and Synchronization in Oscillator Networks. Physical Review Letters, 100(17), 2008. http://arxiv.org/abs/0708.2055
[3] D. Gfeller, and P. De Los Rios, Spectral Coarse Graining of Complex Networks, Physical Review Letters, 99(3), 2007. http://arxiv.org/abs/0706.0812
igraph_scg_adjacency
— Spectral coarse graining, symmetric case.igraph_scg_stochastic
— Spectral coarse graining, stochastic case.igraph_scg_laplacian
— Spectral coarse graining, Laplacian case.igraph_scg_grouping
— SCG problem solver.igraph_scg_semiprojectors
— Compute SCG semi-projectors for a given partition.igraph_scg_norm_eps
— Calculate SCG residuals.
int igraph_scg_adjacency(const igraph_t *graph, const igraph_matrix_t *matrix, const igraph_sparsemat_t *sparsemat, const igraph_vector_t *ev, igraph_integer_t nt, const igraph_vector_t *nt_vec, igraph_scg_algorithm_t algo, igraph_vector_t *values, igraph_matrix_t *vectors, igraph_vector_t *groups, igraph_bool_t use_arpack, igraph_integer_t maxiter, igraph_t *scg_graph, igraph_matrix_t *scg_matrix, igraph_sparsemat_t *scg_sparsemat, igraph_matrix_t *L, igraph_matrix_t *R, igraph_sparsemat_t *Lsparse, igraph_sparsemat_t *Rsparse);
This function handles all the steps involved in the Spectral Coarse Graining (SCG) of some matrices and graphs as described in the reference below.
Arguments:
|
The input graph. Exactly one of |
|
The input matrix. Exactly one of |
|
The input sparse matrix. Exactly one of |
|
A vector of positive integers giving the indexes of the eigenpairs to be preserved. 1 designates the eigenvalue with largest algebraic value, 2 the one with second largest algebraic value, etc. |
|
Positive integer. When |
|
A numeric vector of length one or the length must
match the number of eigenvectors given in |
|
The algorithm to solve the SCG problem. Possible
values: |
|
If this is not |
|
If this is not |
|
If this is not |
|
Whether to use ARPACK for solving the eigenproblem. Currently ARPACK is not implemented. |
|
A positive integer giving the number of iterations
of the k-means algorithm when |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
Returns:
Error code. |
Time complexity: TODO.
See also:
|
Example 27.1. File examples/simple/scg.c
/* -*- mode: C -*- */ /* IGraph library. Copyright (C) 2011-2012 Gabor Csardi <csardi.gabor@gmail.com> 334 Harvard st, Cambridge MA, 02139 USA This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA */ #include <igraph.h> int main() { igraph_t g; igraph_vector_t ev; igraph_t scg_graph; igraph_matrix_t scg_matrix; igraph_sparsemat_t scg_sparsemat; igraph_matrix_t L, R; igraph_sparsemat_t Lsparse, Rsparse; igraph_matrix_t input_matrix; igraph_vector_t groups; igraph_vector_t eval; igraph_matrix_t evec; igraph_tree(&g, 10, /* children= */ 3, IGRAPH_TREE_UNDIRECTED); igraph_vector_init(&ev, 1); igraph_matrix_init(&L, 0, 0); igraph_matrix_init(&R, 0, 0); igraph_matrix_init(&scg_matrix, 0, 0); igraph_vector_init(&groups, 0); igraph_vector_init(&eval, 0); igraph_matrix_init(&evec, 0, 0); #define CALLSYM(algo) do { \ igraph_vector_clear(&eval); \ igraph_matrix_resize(&evec, 0, 0); \ igraph_scg_adjacency(&g, /*matrix=*/ 0, /*sparsemat=*/ 0, &ev, \ /* intervals= */ 3, /* intervals_vector= */ 0, \ /* algorithm= */ algo, &eval, &evec, \ /* groups= */ &groups, /* use_arpack= */ 0, \ /* maxiter= */ 0, &scg_graph, &scg_matrix, \ &scg_sparsemat, &L, &R, \ &Lsparse, &Rsparse); } while(0) #define PRINTRES() \ do { \ printf("------------------------------------\n"); \ igraph_write_graph_edgelist(&scg_graph, stdout); \ printf("---\n"); \ igraph_vector_print(&groups); \ printf("---\n"); \ igraph_vector_print(&eval); \ igraph_matrix_print(&evec); \ printf("---\n"); \ igraph_sparsemat_print(&scg_sparsemat, stdout); \ printf("---\n"); \ igraph_sparsemat_print(&Lsparse, stdout); \ printf("---\n"); \ igraph_sparsemat_print(&Rsparse, stdout); \ printf("---\n"); \ } while (0) VECTOR(ev)[0] = 1; CALLSYM(IGRAPH_SCG_EXACT); PRINTRES(); igraph_destroy(&scg_graph); igraph_sparsemat_destroy(&scg_sparsemat); igraph_sparsemat_destroy(&Lsparse); igraph_sparsemat_destroy(&Rsparse); VECTOR(ev)[0] = 3; CALLSYM(IGRAPH_SCG_EXACT); PRINTRES(); igraph_destroy(&scg_graph); igraph_sparsemat_destroy(&scg_sparsemat); igraph_sparsemat_destroy(&Lsparse); igraph_sparsemat_destroy(&Rsparse); igraph_vector_resize(&ev, 2); VECTOR(ev)[0] = 1; VECTOR(ev)[1] = 3; CALLSYM(IGRAPH_SCG_EXACT); PRINTRES(); igraph_destroy(&scg_graph); igraph_sparsemat_destroy(&scg_sparsemat); igraph_sparsemat_destroy(&Lsparse); igraph_sparsemat_destroy(&Rsparse); #define CALLSYM2(algo) do { \ igraph_vector_clear(&eval); \ igraph_matrix_resize(&evec, 0, 0); \ igraph_scg_adjacency(/* graph=*/ 0, &input_matrix, /*sparsemat=*/ 0, \ &ev, /* intervals= */ 3, \ /* intervals_vector= */ 0, \ /* algorithm= */ algo, &eval, &evec, \ /* groups= */ &groups, /* use_arpack= */ 0, \ /* maxiter= */ 0, &scg_graph, &scg_matrix, \ &scg_sparsemat, &L, &R, \ &Lsparse, &Rsparse); } while (0) igraph_matrix_init(&input_matrix, 0, 0); igraph_get_adjacency(&g, &input_matrix, IGRAPH_GET_ADJACENCY_BOTH, /* eids= */ 0); igraph_vector_resize(&ev, 1); VECTOR(ev)[0] = 1; CALLSYM2(IGRAPH_SCG_EXACT); PRINTRES(); igraph_destroy(&scg_graph); igraph_sparsemat_destroy(&scg_sparsemat); igraph_sparsemat_destroy(&Lsparse); igraph_sparsemat_destroy(&Rsparse); VECTOR(ev)[0] = 3; CALLSYM2(IGRAPH_SCG_EXACT); PRINTRES(); igraph_destroy(&scg_graph); igraph_sparsemat_destroy(&scg_sparsemat); igraph_sparsemat_destroy(&Lsparse); igraph_sparsemat_destroy(&Rsparse); igraph_vector_resize(&ev, 2); VECTOR(ev)[0] = 1; VECTOR(ev)[1] = 3; CALLSYM2(IGRAPH_SCG_EXACT); PRINTRES(); igraph_destroy(&scg_graph); igraph_sparsemat_destroy(&scg_sparsemat); igraph_sparsemat_destroy(&Lsparse); igraph_sparsemat_destroy(&Rsparse); igraph_matrix_destroy(&evec); igraph_vector_destroy(&eval); igraph_vector_destroy(&groups); igraph_matrix_destroy(&input_matrix); igraph_matrix_destroy(&scg_matrix); igraph_matrix_destroy(&L); igraph_matrix_destroy(&R); igraph_vector_destroy(&ev); igraph_destroy(&g); /* -------------------------------------------------------------------- */ return 0; }
int igraph_scg_stochastic(const igraph_t *graph, const igraph_matrix_t *matrix, const igraph_sparsemat_t *sparsemat, const igraph_vector_t *ev, igraph_integer_t nt, const igraph_vector_t *nt_vec, igraph_scg_algorithm_t algo, igraph_scg_norm_t norm, igraph_vector_complex_t *values, igraph_matrix_complex_t *vectors, igraph_vector_t *groups, igraph_vector_t *p, igraph_bool_t use_arpack, igraph_integer_t maxiter, igraph_t *scg_graph, igraph_matrix_t *scg_matrix, igraph_sparsemat_t *scg_sparsemat, igraph_matrix_t *L, igraph_matrix_t *R, igraph_sparsemat_t *Lsparse, igraph_sparsemat_t *Rsparse);
This function handles all the steps involved in the Spectral Coarse Graining (SCG) of some matrices and graphs as described in the reference below.
Arguments:
|
The input graph. Exactly one of |
|
The input matrix. Exactly one of |
|
The input sparse matrix. Exactly one of |
|
A vector of positive integers giving the indexes of the eigenpairs to be preserved. 1 designates the eigenvalue with largest magnitude, 2 the one with second largest magnitude, etc. |
|
Positive integer. When |
|
A numeric vector of length one or the length must
match the number of eigenvectors given in |
|
The algorithm to solve the SCG problem. Possible
values: |
|
Either |
|
If this is not |
|
If this is not |
|
If this is not |
|
If this is not |
|
Whether to use ARPACK for solving the eigenproblem. Currently ARPACK is not implemented. |
|
A positive integer giving the number of iterations
of the k-means algorithm when |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
Returns:
Error code. |
Time complexity: TODO.
See also:
|
int igraph_scg_laplacian(const igraph_t *graph, const igraph_matrix_t *matrix, const igraph_sparsemat_t *sparsemat, const igraph_vector_t *ev, igraph_integer_t nt, const igraph_vector_t *nt_vec, igraph_scg_algorithm_t algo, igraph_scg_norm_t norm, igraph_scg_direction_t direction, igraph_vector_complex_t *values, igraph_matrix_complex_t *vectors, igraph_vector_t *groups, igraph_bool_t use_arpack, igraph_integer_t maxiter, igraph_t *scg_graph, igraph_matrix_t *scg_matrix, igraph_sparsemat_t *scg_sparsemat, igraph_matrix_t *L, igraph_matrix_t *R, igraph_sparsemat_t *Lsparse, igraph_sparsemat_t *Rsparse);
This function handles all the steps involved in the Spectral Coarse Graining (SCG) of some matrices and graphs as described in the reference below.
Arguments:
|
The input graph. Exactly one of |
|
The input matrix. Exactly one of |
|
The input sparse matrix. Exactly one of |
|
A vector of positive integers giving the indexes of the eigenpairs to be preserved. 1 designates the eigenvalue with largest magnitude, 2 the one with second largest magnitude, etc. |
|
Positive integer. When |
|
A numeric vector of length one or the length must
match the number of eigenvectors given in |
|
The algorithm to solve the SCG problem. Possible
values: |
|
Either |
|
Whether to work with left or right eigenvectors.
Possible values: |
|
If this is not |
|
If this is not |
|
If this is not |
|
Whether to use ARPACK for solving the eigenproblem. Currently ARPACK is not implemented. |
|
A positive integer giving the number of iterations
of the k-means algorithm when |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
Returns:
Error code. |
Time complexity: TODO.
See also:
|
int igraph_scg_grouping(const igraph_matrix_t *V, igraph_vector_t *groups, igraph_integer_t nt, const igraph_vector_t *nt_vec, igraph_scg_matrix_t mtype, igraph_scg_algorithm_t algo, const igraph_vector_t *p, igraph_integer_t maxiter);
This function solves the Spectral Coarse Graining (SCG) problem; either exactly, or approximately but faster.
The algorithm IGRAPH_SCG_OPTIMUM
solves the SCG problem exactly
for each eigenvector in V
. The running time of this algorithm is
O(max(nt) m^2) for the symmetric and Laplacian matrix problems.
It is O(m^3) for the stochastic problem. Here m is the number
of rows in V
. In all three cases, the memory usage is O(m^2).
The algorithms IGRAPH_SCG_INTERV
and IGRAPH_SCG_INTERV_KM
solve
the SCG problem approximately by performing a (for now) constant
binning of the components of the eigenvectors, that is nt_vec[i]
constant-size bins are used to partition the i
th eigenvector in V
.
When algo
is IGRAPH_SCG_INTERV_KM
, the (Lloyd) k-means algorithm is
run on each partition obtained by IGRAPH_SCG_INTERV
to improve
accuracy.
Once a minimizing partition (either exact or approximate) has been
found for each eigenvector, the final grouping is worked out as
follows: two vertices are grouped together in the final partition if
they are grouped together in each minimizing partition. In general, the
size of the final partition is not known in advance when the number
of columns in V
is larger than one.
Finally, the algorithm IGRAPH_SCG_EXACT
groups the vertices with
equal components in each eigenvector. The last three algorithms
essentially have linear running time and memory load.
Arguments:
|
The matrix of eigenvectors to be preserved by coarse graining, each column is an eigenvector. |
|
Pointer to an initialized vector; the result of the SCG is stored here. |
|
Positive integer. When |
|
May be (1) a numeric vector of length one, or
(2) a vector of the same length as the number of eigenvectors given in |
|
The type of semi-projectors used in the SCG. Possible
values are |
|
The algorithm to solve the SCG problem. Possible
values: |
|
A probability vector, or |
|
A positive integer giving the number of iterations
of the k-means algorithm when |
Returns:
Error code. |
Time complexity: see description above.
See also:
Example 27.2. File examples/simple/igraph_scg_grouping.c
/* -*- mode: C -*- */ /* IGraph library. Copyright (C) 2011-2012 Gabor Csardi <csardi.gabor@gmail.com> 334 Harvard st, Cambridge, MA, 02138 USA This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA */ #include <igraph.h> #define SIZE (1000) int main() { igraph_matrix_t M, M2; igraph_vector_t lambda; igraph_matrix_t V; igraph_vector_t groups; igraph_vector_t ivec; int i, j; int n; igraph_rng_seed(igraph_rng_default(), 42); /* Symmetric matrix, exponentially distributed elements */ igraph_matrix_init(&M, SIZE, SIZE); n = igraph_matrix_nrow(&M); for (i = 0; i < n; i++) { for (j = 0; j < n; j++) { MATRIX(M, i, j) = igraph_rng_get_exp(igraph_rng_default(), 1); } } igraph_matrix_init(&M2, n, n); igraph_matrix_update(&M2, &M); igraph_matrix_transpose(&M2); igraph_matrix_add(&M, &M2); igraph_matrix_scale(&M, 0.5); igraph_matrix_destroy(&M2); /* Get first (most positive) two eigenvectors */ igraph_vector_init(&lambda, 0); igraph_matrix_init(&V, 0, 0); igraph_lapack_dsyevr(&M, IGRAPH_LAPACK_DSYEV_SELECT, /*vl=*/ 0, /*vu=*/ 0, /*vestimate=*/ 0, /*il=*/ n - 1, /*iu=*/ n, /*abstol=*/ 0.0, /*values=*/ &lambda, /*vectors=*/ &V, /*support=*/ 0); /* Grouping */ igraph_vector_init(&groups, 0); igraph_vector_init(&ivec, 2); VECTOR(ivec)[0] = 2; VECTOR(ivec)[1] = 3; igraph_scg_grouping(&V, &groups, /*invervals=*/ 0, /*intervals_vector=*/ &ivec, IGRAPH_SCG_SYMMETRIC, IGRAPH_SCG_OPTIMUM, /*p=*/ 0, /*maxiter=*/ 100); igraph_vector_print(&groups); igraph_vector_destroy(&ivec); igraph_vector_destroy(&groups); igraph_vector_destroy(&lambda); igraph_matrix_destroy(&V); igraph_matrix_destroy(&M); return 0; }
Example 27.3. File examples/simple/igraph_scg_grouping2.c
/* -*- mode: C -*- */ /* IGraph library. Copyright (C) 2011-2012 Gabor Csardi <csardi.gabor@gmail.com> 334 Harvard st, Cambridge MA, 02139 USA This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA */ #include <igraph.h> int main() { igraph_t g; igraph_matrix_t adj, V; igraph_vector_t groups; igraph_eigen_which_t which; igraph_matrix_init(&adj, 0, 0); igraph_matrix_init(&V, 0, 0); igraph_vector_init(&groups, 0); igraph_tree(&g, 10, /* children= */ 3, IGRAPH_TREE_UNDIRECTED); igraph_get_adjacency(&g, &adj, IGRAPH_GET_ADJACENCY_BOTH, /*eids=*/ 0); which.pos = IGRAPH_EIGEN_LM; which.howmany = 1; igraph_eigen_matrix_symmetric(&adj, /*sparsemat=*/ 0, /*fun=*/ 0, igraph_vcount(&g), /*extra=*/ 0, /*algorithm=*/ IGRAPH_EIGEN_LAPACK, &which, /*options=*/ 0, /*storage=*/ 0, /*values=*/ 0, &V); igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_SYMMETRIC, IGRAPH_SCG_OPTIMUM, /*p=*/ 0, /*maxiter=*/ 10000); igraph_vector_print(&groups); igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_SYMMETRIC, IGRAPH_SCG_INTERV_KM, /*p=*/ 0, /*maxiter=*/ 10000); igraph_vector_print(&groups); igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_SYMMETRIC, IGRAPH_SCG_INTERV, /*p=*/ 0, /*maxiter=*/ 10000); igraph_vector_print(&groups); igraph_scg_grouping(&V, &groups, /*(ignored) intervals=*/ 0, /*intervals_vector=*/ 0, IGRAPH_SCG_SYMMETRIC, IGRAPH_SCG_EXACT, /*p=*/ 0, /*maxiter=*/ 10000); igraph_vector_print(&groups); igraph_vector_destroy(&groups); igraph_matrix_destroy(&V); igraph_matrix_destroy(&adj); igraph_destroy(&g); return 0; }
Example 27.4. File examples/simple/igraph_scg_grouping3.c
/* -*- mode: C -*- */ /* IGraph library. Copyright (C) 2011-2012 Gabor Csardi <csardi.gabor@gmail.com> 334 Harvard st, Cambridge MA, 02139 USA This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA */ #include <igraph.h> int main() { const int nodes = 10; igraph_t g; igraph_matrix_t V, V3; igraph_matrix_complex_t V2; igraph_sparsemat_t stochastic; igraph_vector_t groups; igraph_eigen_which_t which; igraph_vector_t p, selcol; /* This is a 10-node tree with no non-trivial automorphisms. */ igraph_small(&g, nodes, IGRAPH_UNDIRECTED, 3, 5, 4, 5, 4, 9, 8, 9, 0, 9, 0, 6, 1, 6, 1, 2, 7, 8, -1); igraph_matrix_complex_init(&V2, 0, 0); igraph_matrix_init(&V, 0, 0); igraph_matrix_init(&V3, 0, 0); igraph_vector_init(&groups, 0); igraph_vector_init(&p, 0); igraph_vector_init(&selcol, 1); igraph_get_stochastic_sparsemat(&g, &stochastic, /*column-wise=*/ 0); /* p is always the eigenvector corresponding to the 1-eigenvalue. * Since the graph is undirected, p is proportional to the degree vector. */ igraph_degree(&g, &p, igraph_vss_all(), IGRAPH_ALL, IGRAPH_LOOPS); which.pos = IGRAPH_EIGEN_LR; which.howmany = 3; igraph_eigen_matrix(/*matrix=*/ 0, &stochastic, /*fun=*/ 0, nodes, /*extra=*/ 0, /*algorithm=*/ IGRAPH_EIGEN_LAPACK, &which, /*options=*/ 0, /*storage=*/ 0, /*values=*/ 0, &V2); igraph_matrix_complex_real(&V2, &V3); VECTOR(selcol)[0] = 2; igraph_matrix_select_cols(&V3, &V, &selcol); /* ------------ */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_STOCHASTIC, IGRAPH_SCG_OPTIMUM, &p, /*maxiter=*/ 10000); igraph_vector_print(&groups); /* ------------ */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_STOCHASTIC, IGRAPH_SCG_INTERV_KM, &p, /*maxiter=*/ 10000); igraph_vector_print(&groups); /* ------------ */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_STOCHASTIC, IGRAPH_SCG_INTERV, &p, /*maxiter=*/ 10000); igraph_vector_print(&groups); /* ------------ */ igraph_scg_grouping(&V, &groups, /*(ignored) intervals=*/ 0, /*intervals_vector=*/ 0, IGRAPH_SCG_STOCHASTIC, IGRAPH_SCG_EXACT, &p, /*maxiter=*/ 10000); igraph_vector_print(&groups); /* ------------ */ igraph_vector_destroy(&p); igraph_vector_destroy(&selcol); igraph_vector_destroy(&groups); igraph_matrix_destroy(&V); igraph_matrix_destroy(&V3); igraph_matrix_complex_destroy(&V2); igraph_sparsemat_destroy(&stochastic); igraph_destroy(&g); return 0; }
Example 27.5. File examples/simple/igraph_scg_grouping4.c
/* -*- mode: C -*- */ /* IGraph library. Copyright (C) 2011-2012 Gabor Csardi <csardi.gabor@gmail.com> 334 Harvard st, Cambridge MA, 02139 USA This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA */ #include <igraph.h> int main() { const int nodes = 10; igraph_t g; igraph_matrix_t V; igraph_matrix_complex_t V2; igraph_sparsemat_t laplacian; igraph_vector_t groups; igraph_eigen_which_t which; igraph_tree(&g, nodes, /* children= */ 3, IGRAPH_TREE_UNDIRECTED); igraph_matrix_complex_init(&V2, 0, 0); igraph_matrix_init(&V, 0, 0); igraph_vector_init(&groups, 0); igraph_sparsemat_init(&laplacian, 0, 0, 0); igraph_laplacian(&g, /*res=*/ 0, /*sparseres=*/ &laplacian, /*normalized=*/ 0, /*weights=*/ 0); which.pos = IGRAPH_EIGEN_LR; which.howmany = 1; igraph_eigen_matrix(/*matrix=*/ 0, &laplacian, /*fun=*/ 0, nodes, /*extra=*/ 0, /*algorithm=*/ IGRAPH_EIGEN_LAPACK, &which, /*options=*/ 0, /*storage=*/ 0, /*values=*/ 0, &V2); igraph_matrix_complex_real(&V2, &V); /* ------------ */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_LAPLACIAN, IGRAPH_SCG_OPTIMUM, /*p=*/ 0, /*maxiter=*/ 10000); igraph_vector_print(&groups); /* ------------ */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_LAPLACIAN, IGRAPH_SCG_INTERV_KM, /*p=*/ 0, /*maxiter=*/ 10000); igraph_vector_print(&groups); /* ------------ */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_LAPLACIAN, IGRAPH_SCG_INTERV, /*p=*/ 0, /*maxiter=*/ 10000); igraph_vector_print(&groups); /* ------------ */ igraph_scg_grouping(&V, &groups, /*(ignored) intervals=*/ 0, /*intervals_vector=*/ 0, IGRAPH_SCG_LAPLACIAN, IGRAPH_SCG_EXACT, /*p=*/ 0, /*maxiter=*/ 10000); igraph_vector_print(&groups); /* ------------ */ igraph_vector_destroy(&groups); igraph_matrix_destroy(&V); igraph_matrix_complex_destroy(&V2); igraph_sparsemat_destroy(&laplacian); igraph_destroy(&g); return 0; }
int igraph_scg_semiprojectors(const igraph_vector_t *groups, igraph_scg_matrix_t mtype, igraph_matrix_t *L, igraph_matrix_t *R, igraph_sparsemat_t *Lsparse, igraph_sparsemat_t *Rsparse, const igraph_vector_t *p, igraph_scg_norm_t norm);
The three types of semi-projectors are defined as follows. Let gamma(j) label the group of vertex j in a partition of all the vertices.
The symmetric semi-projectors are defined as
L[alpha,j] = R[alpha,j] = 1/sqrt(|alpha|) delta[alpha,gamma(j)],
the (row) Laplacian semi-projectors as
L[alpha,j] = 1/|alpha| delta[alpha,gamma(j)]
and
R[alpha,j] = delta[alpha,gamma(j)],
and the (row) stochastic semi-projectors as
L[alpha,j] = p[1][j] / sum(p[1][k]; k in gamma(j)) delta[alpha,gamma(j)]
and
R[alpha,j] = delta[alpha,gamma(j)],
where p[1] is the (left) eigenvector associated with the
one-eigenvalue of the stochastic matrix. L and R are
defined in a symmetric way when norm
is IGRAPH_SCG_NORM_COL
. All these semi-projectors verify various
properties described in the reference.
Arguments:
|
A vector of integers, giving the group label of every vertex in the partition. Group labels should start at zero and should be sequential. |
|
The type of semi-projectors. For now |
|
If not a |
|
If not a |
|
If not a |
|
If not a |
|
|
|
Either |
Returns:
Error code. |
Time complexity: TODO.
See also:
Example 27.6. File examples/simple/igraph_scg_semiprojectors.c
/* -*- mode: C -*- */ /* IGraph library. Copyright (C) 2011-2012 Gabor Csardi <csardi.gabor@gmail.com> 334 Harvard st, Cambridge MA, 02139 USA This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA */ #include <igraph.h> int main() { igraph_t g; igraph_matrix_t L, R; igraph_sparsemat_t Lsparse, Rsparse; igraph_matrix_t adj, V; igraph_vector_t groups; igraph_eigen_which_t which; igraph_matrix_init(&L, 0, 0); igraph_matrix_init(&R, 0, 0); igraph_matrix_init(&adj, 0, 0); igraph_matrix_init(&V, 0, 0); igraph_vector_init(&groups, 0); igraph_tree(&g, 10, /* children= */ 3, IGRAPH_TREE_UNDIRECTED); igraph_get_adjacency(&g, &adj, IGRAPH_GET_ADJACENCY_BOTH, /*eids=*/ 0); which.pos = IGRAPH_EIGEN_LM; which.howmany = 1; igraph_eigen_matrix_symmetric(&adj, /*sparsemat=*/ 0, /*fun=*/ 0, igraph_vcount(&g), /*extra=*/ 0, /*algorithm=*/ IGRAPH_EIGEN_LAPACK, &which, /*options=*/ 0, /*storage=*/ 0, /*values=*/ 0, &V); #define SEMI() \ do { \ igraph_scg_semiprojectors(&groups, IGRAPH_SCG_SYMMETRIC, &L, &R, \ &Lsparse, &Rsparse, /*p=*/ 0, \ IGRAPH_SCG_NORM_ROW); \ } while(0) #define PRINTRES() \ do { \ printf("----------------------\n"); \ igraph_matrix_print(&L); \ printf("---\n"); \ igraph_matrix_print(&R); \ printf("---\n"); \ igraph_sparsemat_destroy(&Lsparse); \ igraph_sparsemat_destroy(&Rsparse); \ } while (0) /* -------------- */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_SYMMETRIC, IGRAPH_SCG_OPTIMUM, /*p=*/ 0, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 2, /*intervals_vector=*/ 0, IGRAPH_SCG_SYMMETRIC, IGRAPH_SCG_INTERV_KM, /*p=*/ 0, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 2, /*intervals_vector=*/ 0, IGRAPH_SCG_SYMMETRIC, IGRAPH_SCG_INTERV, /*p=*/ 0, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_scg_grouping(&V, &groups, /*(ignored) intervals=*/ 0, /*intervals_vector=*/ 0, IGRAPH_SCG_SYMMETRIC, IGRAPH_SCG_EXACT, /*p=*/ 0, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_vector_destroy(&groups); igraph_matrix_destroy(&L); igraph_matrix_destroy(&R); igraph_matrix_destroy(&V); igraph_matrix_destroy(&adj); igraph_destroy(&g); return 0; }
Example 27.7. File examples/simple/igraph_scg_semiprojectors2.c
/* -*- mode: C -*- */ /* IGraph library. Copyright (C) 2011-2012 Gabor Csardi <csardi.gabor@gmail.com> 334 Harvard st, Cambridge MA, 02139 USA This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA */ #include <igraph.h> int main() { igraph_t g; igraph_matrix_t L, R; igraph_sparsemat_t Lsparse, Rsparse; igraph_matrix_t V, V3; igraph_matrix_complex_t V2; igraph_sparsemat_t stochastic; igraph_vector_t groups; igraph_eigen_which_t which; igraph_vector_t p, selcol; igraph_matrix_init(&L, 0, 0); igraph_matrix_init(&R, 0, 0); igraph_matrix_init(&V, 0, 0); igraph_matrix_init(&V3, 0, 0); igraph_vector_init(&groups, 0); igraph_vector_init(&selcol, 1); /* This is a 10-node tree with no non-trivial automorphisms. */ igraph_small(&g, 10, IGRAPH_UNDIRECTED, 3, 5, 4, 5, 4, 9, 8, 9, 0, 9, 0, 6, 1, 6, 1, 2, 7, 8, -1); igraph_matrix_complex_init(&V2, 0, 0); igraph_vector_init(&p, 0); igraph_get_stochastic_sparsemat(&g, &stochastic, /*column-wise=*/ 0); /* p is always the eigenvector corresponding to the 1-eigenvalue. * Since the graph is undirected, p is proportional to the degree vector. */ igraph_degree(&g, &p, igraph_vss_all(), IGRAPH_ALL, IGRAPH_LOOPS); which.pos = IGRAPH_EIGEN_LR; which.howmany = 3; igraph_eigen_matrix(/*matrix=*/ 0, &stochastic, /*fun=*/ 0, 10, /*extra=*/ 0, /*algorithm=*/ IGRAPH_EIGEN_LAPACK, &which, /*options=*/ 0, /*storage=*/ 0, /*values=*/ 0, &V2); igraph_matrix_complex_real(&V2, &V3); VECTOR(selcol)[0] = 2; igraph_matrix_select_cols(&V3, &V, &selcol); #define SEMI() \ do { \ igraph_scg_semiprojectors(&groups, IGRAPH_SCG_STOCHASTIC, &L, &R, \ &Lsparse, &Rsparse, &p, \ IGRAPH_SCG_NORM_ROW); \ } while(0) #define PRINTRES() \ do { \ printf("----------------------\n"); \ igraph_matrix_print(&L); \ printf("---\n"); \ igraph_matrix_print(&R); \ printf("---\n"); \ igraph_sparsemat_destroy(&Lsparse); \ igraph_sparsemat_destroy(&Rsparse); \ } while (0) /* -------------- */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_STOCHASTIC, IGRAPH_SCG_OPTIMUM, &p, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_STOCHASTIC, IGRAPH_SCG_INTERV_KM, &p, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_STOCHASTIC, IGRAPH_SCG_INTERV, &p, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_scg_grouping(&V, &groups, /*(ignored) intervals=*/ 0, /*intervals_vector=*/ 0, IGRAPH_SCG_STOCHASTIC, IGRAPH_SCG_EXACT, &p, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_vector_destroy(&p); igraph_vector_destroy(&selcol); igraph_vector_destroy(&groups); igraph_matrix_destroy(&L); igraph_matrix_destroy(&R); igraph_matrix_destroy(&V); igraph_matrix_destroy(&V3); igraph_matrix_complex_destroy(&V2); igraph_sparsemat_destroy(&stochastic); igraph_destroy(&g); return 0; }
Example 27.8. File examples/simple/igraph_scg_semiprojectors3.c
/* -*- mode: C -*- */ /* IGraph library. Copyright (C) 2011-2012 Gabor Csardi <csardi.gabor@gmail.com> 334 Harvard st, Cambridge MA, 02139 USA This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA */ #include <igraph.h> int main() { int nodes = 10; igraph_t g; igraph_matrix_t L, R; igraph_sparsemat_t Lsparse, Rsparse; igraph_matrix_t V; igraph_matrix_complex_t V2; igraph_sparsemat_t laplacian; igraph_vector_t groups; igraph_eigen_which_t which; igraph_matrix_init(&L, 0, 0); igraph_matrix_init(&R, 0, 0); igraph_matrix_init(&V, 0, 0); igraph_matrix_complex_init(&V2, 0, 0); igraph_vector_init(&groups, 0); igraph_tree(&g, 10, /* children= */ 3, IGRAPH_TREE_UNDIRECTED); igraph_sparsemat_init(&laplacian, nodes, nodes, igraph_ecount(&g) * 2); igraph_laplacian(&g, /*res=*/ 0, /*sparseres=*/ &laplacian, /*normalized=*/ 0, /*weights=*/ 0); which.pos = IGRAPH_EIGEN_LM; which.howmany = 1; igraph_eigen_matrix(/*matrix=*/ 0, &laplacian, /*fun=*/ 0, 10, /*extra=*/ 0, /*algorithm=*/ IGRAPH_EIGEN_LAPACK, &which, /*options=*/ 0, /*storage=*/ 0, /*values=*/ 0, &V2); igraph_matrix_complex_real(&V2, &V); #define SEMI() \ do { \ igraph_scg_semiprojectors(&groups, IGRAPH_SCG_LAPLACIAN, &L, &R, \ &Lsparse, &Rsparse, /*p=*/ 0, \ IGRAPH_SCG_NORM_ROW); \ } while(0) #define PRINTRES() \ do { \ printf("----------------------\n"); \ igraph_matrix_print(&L); \ printf("---\n"); \ igraph_matrix_print(&R); \ printf("---\n"); \ igraph_sparsemat_destroy(&Lsparse); \ igraph_sparsemat_destroy(&Rsparse); \ } while (0) /* -------------- */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 3, /*intervals_vector=*/ 0, IGRAPH_SCG_LAPLACIAN, IGRAPH_SCG_OPTIMUM, /*p=*/ 0, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 2, /*intervals_vector=*/ 0, IGRAPH_SCG_LAPLACIAN, IGRAPH_SCG_INTERV_KM, /*p=*/ 0, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_scg_grouping(&V, &groups, /*intervals=*/ 2, /*intervals_vector=*/ 0, IGRAPH_SCG_LAPLACIAN, IGRAPH_SCG_INTERV, /*p=*/ 0, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_scg_grouping(&V, &groups, /*(ignored) intervals=*/ 0, /*intervals_vector=*/ 0, IGRAPH_SCG_LAPLACIAN, IGRAPH_SCG_EXACT, /*p=*/ 0, /*maxiter=*/ 10000); SEMI(); PRINTRES(); /* -------------- */ igraph_matrix_destroy(&L); igraph_matrix_destroy(&R); igraph_matrix_destroy(&V); igraph_matrix_complex_destroy(&V2); igraph_vector_destroy(&groups); igraph_sparsemat_destroy(&laplacian); igraph_destroy(&g); return 0; }
int igraph_scg_norm_eps(const igraph_matrix_t *V, const igraph_vector_t *groups, igraph_vector_t *eps, igraph_scg_matrix_t mtype, const igraph_vector_t *p, igraph_scg_norm_t norm);
Computes |v[i]-Pv[i]|, where v[i] is the i-th eigenvector in V
and P is the projector corresponding to the mtype
argument.
Arguments:
|
The matrix of eigenvectors to be preserved by coarse graining, each column is an eigenvector. |
|
A vector of integers, giving the group label of every vertex in the partition. Group labels should start at zero and should be sequential. |
|
Pointer to a real value, the result is stored here. |
|
The type of semi-projectors. For now |
|
|
|
Either |
Returns:
Error code. |
Time complexity: TODO.
See also:
|
← Chapter 26. Hierarchical random graphs | Chapter 28. Embedding of graphs → |