igraph Reference Manual

For using the igraph C library

Search the manual:

Chapter 24. Detecting community structure

1. Common functions related to community structure

1.1. igraph_modularity — Calculate the modularity of a graph with respect to some clusters or vertex types.

int igraph_modularity(const igraph_t *graph,
                      const igraph_vector_t *membership,
                      const igraph_vector_t *weights,
                      const igraph_real_t resolution,
                      const igraph_bool_t directed,
                      igraph_real_t *modularity);

The modularity of a graph with respect to some clustering of the vertices (or assignment of vertex types) measures how strongly separated the different clusters are from each other compared to a random null model. It is defined as

Q = 1/(2m) sum_ij (A_ij - gamma * k_i * k_j / (2m)) * d(c_i,c_j),

where m is the number of edges, A_ij is the adjacency matrix, k_i is the degree of vertex i, c_i is the cluster that vertex i belongs to (or its vertex type), d(i,j)=1 if i=j and 0 otherwise, and the sum goes over all i, j pairs of vertices.

The resolution parameter gamma allows weighting the random null model, which might be useful when finding partitions with a high modularity. Maximizing modularity with higher values of the resolution parameter typically results in more, smaller clusters when finding partitions with a high modularity. Lower values typically results in fewer, larger clusters. The original definition of modularity is retrieved when setting gamma=1.

Modularity can also be calculated on directed graphs. This only requires a relatively modest change

Q = 1/(m) sum_ij (A_ij - gamma * k^out_i * k^in_j / m) * d(c_i,c_j),

where k^out_i is the out-degree of node i and k^in_j is the in-degree of node j.

Modularity on weighted graphs is also meaningful. When taking edge weights into account, A_ij equals the weight of the corresponding edge (or 0 if there is no edge), k_i is the strength (i.e. the weighted degree) of vertex i, with similar counterparts for a directed graph, and m is the total weight of all edges.

Note that the modularity is not well-defined for graphs with no edges. igraph returns NaN for graphs with no edges; see https://github.com/igraph/igraph/issues/1539 for a detailed discussion.

For the original definition of modularity, see Newman, M. E. J., and Girvan, M. (2004). Finding and evaluating community structure in networks. Physical Review E 69, 026113. https://doi.org/10.1103/PhysRevE.69.026113

For the directed definition of modularity, see Leicht, E. A., and Newman, M. E. J. (2008). Community Structure in Directed Networks. Physical Review Letters 100, 118703. https://doi.org/10.1103/PhysRevLett.100.118703

For the introduction of the resolution parameter, see Reichardt, J., and Bornholdt, S. (2006). Statistical mechanics of community detection. Physical Review E 74, 016110. https://doi.org/10.1103/PhysRevE.74.016110

Arguments: 

graph:

The input graph.

membership:

Numeric vector of integer values which gives the type of each vertex, i.e. the cluster to which it belongs. It does not have to be consecutive, i.e. empty communities are allowed.

weights:

Weight vector or NULL if no weights are specified.

resolution:

Resolution parameter. Must be greater than or equal to 0. Set it to 1 to use the classical definition of modularity.

directed:

Whether to use the directed or undirected version of modularity. Ignored for undirected graphs.

modularity:

Pointer to a real number, the result will be stored here.

Returns: 

Error code.

See also: 

Time complexity: O(|V|+|E|), the number of vertices plus the number of edges.

1.2. igraph_modularity_matrix — Calculate the modularity matrix

int igraph_modularity_matrix(const igraph_t *graph,
                             const igraph_vector_t *weights,
                             const igraph_real_t resolution,
                             igraph_matrix_t *modmat,
                             igraph_bool_t directed);

This function returns the modularity matrix defined as

B_ij = A_ij - gamma * k_i * k_j / (2m)

for undirected graphs, where A_ij is the adjacency matrix, gamma is the resolution parameter, k_i is the degree of vertex i, and m is the number of edges in the graph. When there are no edges, or the weights add up to zero, the result is undefined.

For directed graphs the modularity matrix is changed to

B_ij = A_ij - gamma * k^out_i * k^in_j / m where k^out_i is the out-degree of node i and k^in_j is the in-degree of node j.

Note that self-loops in undirected graphs are multiplied by 2 in this implementation. If weights are specified, the weighted counterparts are used.

Arguments: 

graph:

The input graph.

weights:

Edge weights, pointer to a vector. If this is a null pointer then every edge is assumed to have a weight of 1.

resolution:

Resolution parameter. Must be greater than or equal to 0. Default is 1. Lower values favor fewer, larger communities; higher values favor more, smaller communities.

modmat:

Pointer to an initialized matrix in which the modularity matrix is stored.

directed:

For directed graphs: if the edges should be treated as undirected. For undirected graphs this is ignored.

See also: 

1.3. igraph_community_optimal_modularity — Calculate the community structure with the highest modularity value

int igraph_community_optimal_modularity(const igraph_t *graph,
                                        igraph_real_t *modularity,
                                        igraph_vector_t *membership,
                                        const igraph_vector_t *weights);

This function calculates the optimal community structure for a graph, in terms of maximal modularity score.

The calculation is done by transforming the modularity maximization into an integer programming problem, and then calling the GLPK library to solve that. Please see Ulrik Brandes et al.: On Modularity Clustering, IEEE Transactions on Knowledge and Data Engineering 20(2):172-188, 2008.

Note that modularity optimization is an NP-complete problem, and all known algorithms for it have exponential time complexity. This means that you probably don't want to run this function on larger graphs. Graphs with up to fifty vertices should be fine, graphs with a couple of hundred vertices might be possible.

Arguments: 

graph:

The input graph. It is always treated as undirected.

modularity:

Pointer to a real number, or a null pointer. If it is not a null pointer, then a optimal modularity value is returned here.

membership:

Pointer to a vector, or a null pointer. If not a null pointer, then the membership vector of the optimal community structure is stored here.

weights:

Vector giving the weights of the edges. If it is NULL then each edge is supposed to have the same weight.

Returns: 

Error code.

See also: 

igraph_modularity(), igraph_community_fastgreedy() for an algorithm that finds a local optimum in a greedy way.

Time complexity: exponential in the number of vertices.

Example 24.1.  File examples/simple/igraph_community_optimal_modularity.c

/* -*- mode: C -*-  */
/*
   IGraph library.
   Copyright (C) 2010-2012  Gabor Csardi <csardi.gabor@gmail.com>
   334 Harvard street, 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>

void prepare_weights_vector(igraph_vector_t* weights, const igraph_t* graph) {
    int i, n = igraph_ecount(graph);
    igraph_vector_resize(weights, n);
    for (i = 0; i < n; i++) {
        VECTOR(*weights)[i] = i % 5;
    }
}

int main() {
    igraph_t graph;

    igraph_vector_t v;
    igraph_real_t edges[] = { 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0, 6, 0, 7, 0, 8,
                              0, 10, 0, 11, 0, 12, 0, 13, 0, 17, 0, 19, 0, 21, 0, 31,
                              1, 2, 1, 3, 1, 7, 1, 13, 1, 17, 1, 19, 1, 21, 1, 30,
                              2, 3, 2, 7, 2, 27, 2, 28, 2, 32, 2, 9, 2, 8, 2, 13,
                              3, 7, 3, 12, 3, 13, 4, 6, 4, 10, 5, 6, 5, 10, 5, 16,
                              6, 16, 8, 30, 8, 32, 8, 33, 9, 33, 13, 33, 14, 32, 14, 33,
                              15, 32, 15, 33, 18, 32, 18, 33, 19, 33, 20, 32, 20, 33,
                              22, 32, 22, 33, 23, 25, 23, 27, 23, 32, 23, 33, 23, 29,
                              24, 25, 24, 27, 24, 31, 25, 31, 26, 29, 26, 33, 27, 33,
                              28, 31, 28, 33, 29, 32, 29, 33, 30, 32, 30, 33, 31, 32, 31, 33,
                              32, 33
                            };

    igraph_vector_t membership;
    igraph_vector_t weights;
    igraph_real_t modularity;
    igraph_bool_t simple;
    int retval;

    igraph_vector_view(&v, edges, sizeof(edges) / sizeof(double));
    igraph_create(&graph, &v, 0, IGRAPH_UNDIRECTED);

    igraph_vector_init(&weights, 0);

    igraph_is_simple(&graph, &simple);
    if (!simple) {
        return 1;
    }

    igraph_vector_init(&membership, 0);

    igraph_set_error_handler(&igraph_error_handler_printignore);

    /* Zachary karate club, unweighted */
    retval = igraph_community_optimal_modularity(&graph, &modularity,
             &membership, 0);
    if (retval == IGRAPH_UNIMPLEMENTED) {
        return 77;
    }
    if (fabs(modularity - 0.4197896) > 0.0000001) {
        return 2;
    }
    /* Zachary karate club, weighted */
    prepare_weights_vector(&weights, &graph);
    igraph_community_optimal_modularity(&graph, &modularity,
                                        &membership, &weights);
    if (fabs(modularity - 0.5115767) > 0.0000001) {
        return 4;
    }
    igraph_destroy(&graph);

    /* simple graph with loop edges, unweighted */
    igraph_small(&graph, 6, IGRAPH_UNDIRECTED,
                 0, 1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 0, 0, 0, 2, 2, -1);
    igraph_community_optimal_modularity(&graph, &modularity,
                                        &membership, 0);
    if (fabs(modularity - 0.28125) > 0.00001) {
        return 3;
    }
    /* simple graph with loop edges, weighted */
    prepare_weights_vector(&weights, &graph);
    igraph_community_optimal_modularity(&graph, &modularity,
                                        &membership, &weights);
    if (fabs(modularity - 0.36686) > 0.00001) {
        return 5;
    }
    igraph_destroy(&graph);

    igraph_vector_destroy(&membership);
    igraph_vector_destroy(&weights);

    return 0;
}


1.4. igraph_community_to_membership — Create membership vector from community structure dendrogram

int igraph_community_to_membership(const igraph_matrix_t *merges,
                                   igraph_integer_t nodes,
                                   igraph_integer_t steps,
                                   igraph_vector_t *membership,
                                   igraph_vector_t *csize);

This function creates a membership vector from a community structure dendrogram. A membership vector contains for each vertex the id of its graph component, the graph components are numbered from zero, see the same argument of igraph_clusters() for an example of a membership vector.

Many community detection algorithms return with a merges matrix, igraph_community_walktrap() and igraph_community_edge_betweenness() are two examples. The matrix contains the merge operations performed while mapping the hierarchical structure of a network. If the matrix has n-1 rows, where n is the number of vertices in the graph, then it contains the hierarchical structure of the whole network and it is called a dendrogram.

This function performs steps merge operations as prescribed by the merges matrix and returns the current state of the network.

If merges is not a complete dendrogram, it is possible to take steps steps if steps is not bigger than the number lines in merges.

Arguments: 

merges:

The two-column matrix containing the merge operations. See igraph_community_walktrap() for the detailed syntax.

nodes:

The number of leaf nodes in the dendrogram.

steps:

Integer constant, the number of steps to take.

membership:

Pointer to an initialized vector, the membership results will be stored here, if not NULL. The vector will be resized as needed.

csize:

Pointer to an initialized vector, or NULL. If not NULL then the sizes of the components will be stored here, the vector will be resized as needed.

See also: 

Time complexity: O(|V|), the number of vertices in the graph.

1.5. igraph_reindex_membership — Makes the IDs in a membership vector continuous

int igraph_reindex_membership(igraph_vector_t *membership,
                              igraph_vector_t *new_to_old,
                              igraph_integer_t *nb_clusters);

This function reindexes component IDs in a membership vector in a way that the new IDs start from zero and go up to C-1, where C is the number of unique component IDs in the original vector. The supplied membership is expected to fall in the range 0, ..., n - 1.

Arguments: 

membership:

Numeric vector which gives the type of each vertex, i.e. the component to which it belongs. The vector will be altered in-place.

new_to_old:

Pointer to a vector which will contain the old component ID for each new one, or NULL, in which case it is not returned. The vector will be resized as needed.

nb_clusters:

Pointer to an integer for the number of distinct clusters. If not NULL, this will be updated to reflect the number of distinct clusters found in membership.

Time complexity: should be O(n) for n elements.

1.6. igraph_compare_communities — Compares community structures using various metrics

int igraph_compare_communities(const igraph_vector_t *comm1,
                               const igraph_vector_t *comm2, igraph_real_t* result,
                               igraph_community_comparison_t method);

This function assesses the distance between two community structures using the variation of information (VI) metric of Meila (2003), the normalized mutual information (NMI) of Danon et al (2005), the split-join distance of van Dongen (2000), the Rand index of Rand (1971) or the adjusted Rand index of Hubert and Arabie (1985).

References:

Meila M: Comparing clusterings by the variation of information. In: Schölkopf B, Warmuth MK (eds.). Learning Theory and Kernel Machines: 16th Annual Conference on Computational Learning Theory and 7th Kernel Workshop, COLT/Kernel 2003, Washington, DC, USA. Lecture Notes in Computer Science, vol. 2777, Springer, 2003. ISBN: 978-3-540-40720-1.

Danon L, Diaz-Guilera A, Duch J, Arenas A: Comparing community structure identification. J Stat Mech P09008, 2005.

van Dongen S: Performance criteria for graph clustering and Markov cluster experiments. Technical Report INS-R0012, National Research Institute for Mathematics and Computer Science in the Netherlands, Amsterdam, May 2000.

Rand WM: Objective criteria for the evaluation of clustering methods. J Am Stat Assoc 66(336):846-850, 1971.

Hubert L and Arabie P: Comparing partitions. Journal of Classification 2:193-218, 1985.

Arguments: 

comm1:

the membership vector of the first community structure

comm2:

the membership vector of the second community structure

result:

the result is stored here.

method:

the comparison method to use. IGRAPH_COMMCMP_VI selects the variation of information (VI) metric of Meila (2003), IGRAPH_COMMCMP_NMI selects the normalized mutual information measure proposed by Danon et al (2005), IGRAPH_COMMCMP_SPLIT_JOIN selects the split-join distance of van Dongen (2000), IGRAPH_COMMCMP_RAND selects the unadjusted Rand index (1971) and IGRAPH_COMMCMP_ADJUSTED_RAND selects the adjusted Rand index.

Returns: 

Error code.

Time complexity: O(n log(n)).

1.7. igraph_split_join_distance — Calculates the split-join distance of two community structures

int igraph_split_join_distance(const igraph_vector_t *comm1,
                               const igraph_vector_t *comm2, igraph_integer_t *distance12,
                               igraph_integer_t *distance21);

The split-join distance between partitions A and B is the sum of the projection distance of A from B and the projection distance of B from A. The projection distance is an asymmetric measure and it is defined as follows:

First, each set in partition A is evaluated against all sets in partition B. For each set in partition A, the best matching set in partition B is found and the overlap size is calculated. (Matching is quantified by the size of the overlap between the two sets). Then, the maximal overlap sizes for each set in A are summed together and subtracted from the number of elements in A.

The split-join distance will be returned in two arguments, distance12 will contain the projection distance of the first partition from the second, while distance21 will be the projection distance of the second partition from the first. This makes it easier to detect whether a partition is a subpartition of the other, since in this case, the corresponding distance will be zero.

Reference:

van Dongen S: Performance criteria for graph clustering and Markov cluster experiments. Technical Report INS-R0012, National Research Institute for Mathematics and Computer Science in the Netherlands, Amsterdam, May 2000.

Arguments: 

comm1:

the membership vector of the first community structure

comm2:

the membership vector of the second community structure

distance12:

pointer to an igraph_integer_t, the projection distance of the first community structure from the second one will be returned here.

distance21:

pointer to an igraph_integer_t, the projection distance of the second community structure from the first one will be returned here.

Returns: 

Error code.

\see igraph_compare_communities() with the IGRAPH_COMMCMP_SPLIT_JOIN method if you are not interested in the individual distances but only the sum of them. Time complexity: O(n log(n)).

2. Community structure based on statistical mechanics

2.1. igraph_community_spinglass — Community detection based on statistical mechanics

int igraph_community_spinglass(const igraph_t *graph,
                               const igraph_vector_t *weights,
                               igraph_real_t *modularity,
                               igraph_real_t *temperature,
                               igraph_vector_t *membership,
                               igraph_vector_t *csize,
                               igraph_integer_t spins,
                               igraph_bool_t parupdate,
                               igraph_real_t starttemp,
                               igraph_real_t stoptemp,
                               igraph_real_t coolfact,
                               igraph_spincomm_update_t update_rule,
                               igraph_real_t gamma,
                               /* the rest is for the NegSpin implementation */
                               igraph_spinglass_implementation_t implementation,
                               /*                 igraph_matrix_t *adhesion, */
                               /*                 igraph_matrix_t *normalised_adhesion, */
                               /*                 igraph_real_t *polarization, */
                               igraph_real_t gamma_minus);

This function implements the community structure detection algorithm proposed by Joerg Reichardt and Stefan Bornholdt. The algorithm is described in their paper: Statistical Mechanics of Community Detection, http://arxiv.org/abs/cond-mat/0603718 .

From version 0.6, igraph also supports an extension to the algorithm that allows negative edge weights. This is described in V. A. Traag and Jeroen Bruggeman: Community detection in networks with positive and negative links, http://arxiv.org/abs/0811.2329 .

Arguments: 

graph:

The input graph, it may be directed but the direction of the edges is not used in the algorithm.

weights:

The vector giving the edge weights, it may be NULL, in which case all edges are weighted equally. The edge weights must be positive unless using the IGRAPH_SPINCOMM_IMP_NEG implementation. This condition is not verified by the function.

modularity:

Pointer to a real number, if not NULL then the modularity score of the solution will be stored here. This is the gereralized modularity that simplifies to the one defined in M. E. J. Newman and M. Girvan, Phys. Rev. E 69, 026113 (2004), if the gamma parameter is one.

temperature:

Pointer to a real number, if not NULL then the temperature at the end of the algorithm will be stored here.

membership:

Pointer to an initialized vector or NULL. If not NULL then the result of the clustering will be stored here. For each vertex, the number of its cluster is given, with the first cluster numbered zero. The vector will be resized as needed.

csize:

Pointer to an initialized vector or NULL. If not NULL then the sizes of the clusters will stored here in cluster number order. The vector will be resized as needed.

spins:

Integer giving the number of spins, i.e. the maximum number of clusters. Usually it is not a program to give a high number here, the default was 25 in the original code. Even if the number of spins is high the number of clusters in the result might be small.

parupdate:

A logical constant, whether to update all spins in parallel. The default for this argument was FALSE (i.e. 0) in the original code. It is not implemented in the IGRAPH_SPINCOMM_INP_NEG implementation.

starttemp:

Real number, the temperature at the start. The value of this argument was 1.0 in the original code.

stoptemp:

Real number, the algorithm stops at this temperature. The default was 0.01 in the original code.

coolfact:

Real number, the cooling factor for the simulated annealing. The default was 0.99 in the original code.

update_rule:

The type of the update rule. Possible values: IGRAPH_SPINCOMM_UPDATE_SIMPLE and IGRAPH_SPINCOMM_UPDATE_CONFIG. Basically this parameter defines the null model based on which the actual clustering is done. If this is IGRAPH_SPINCOMM_UPDATE_SIMPLE then the random graph (i.e. G(n,p)), if it is IGRAPH_SPINCOMM_UPDATE then the configuration model is used. The configuration means that the baseline for the clustering is a random graph with the same degree distribution as the input graph.

gamma:

Real number. The gamma parameter of the algorithm. This defines the weight of the missing and existing links in the quality function for the clustering. The default value in the original code was 1.0, which is equal weight to missing and existing edges. Smaller values make the existing links contibute more to the energy function which is minimized in the algorithm. Bigger values make the missing links more important. (If my understanding is correct.)

implementation:

Constant, chooses between the two implementations of the spin-glass algorithm that are included in igraph. IGRAPH_SPINCOMM_IMP_ORIG selects the original implementation, this is faster, IGRAPH_SPINCOMM_INP_NEG selects a new implementation by Vincent Traag that allows negative edge weights.

gamma_minus:

Real number. Parameter for the IGRAPH_SPINCOMM_IMP_NEG implementation. This specifies the balance between the importance of present and non-present negative weighted edges in a community. Smaller values of gamma_minus lead to communities with lesser negative intra-connectivity. If this argument is set to zero, the algorithm reduces to a graph coloring algorithm, using the number of spins as the number of colors.

Returns: 

Error code.

See also: 

igraph_community_spinglass_single() for calculating the community of a single vertex.

Time complexity: TODO.

2.2. igraph_community_spinglass_single — Community of a single node based on statistical mechanics

int igraph_community_spinglass_single(const igraph_t *graph,
                                      const igraph_vector_t *weights,
                                      igraph_integer_t vertex,
                                      igraph_vector_t *community,
                                      igraph_real_t *cohesion,
                                      igraph_real_t *adhesion,
                                      igraph_integer_t *inner_links,
                                      igraph_integer_t *outer_links,
                                      igraph_integer_t spins,
                                      igraph_spincomm_update_t update_rule,
                                      igraph_real_t gamma);

This function implements the community structure detection algorithm proposed by Joerg Reichardt and Stefan Bornholdt. It is described in their paper: Statistical Mechanics of Community Detection, http://arxiv.org/abs/cond-mat/0603718 .

This function calculates the community of a single vertex without calculating all the communities in the graph.

Arguments: 

graph:

The input graph, it may be directed but the direction of the edges is not used in the algorithm.

weights:

Pointer to a vector with the weights of the edges. Alternatively NULL can be supplied to have the same weight for every edge.

vertex:

The vertex id of the vertex of which ths community is calculated.

community:

Pointer to an initialized vector, the result, the ids of the vertices in the community of the input vertex will be stored here. The vector will be resized as needed.

cohesion:

Pointer to a real variable, if not NULL the cohesion index of the community will be stored here.

adhesion:

Pointer to a real variable, if not NULL the adhesion index of the community will be stored here.

inner_links:

Pointer to an integer, if not NULL the number of edges within the community is stored here.

outer_links:

Pointer to an integer, if not NULL the number of edges between the community and the rest of the graph will be stored here.

spins:

The number of spins to use, this can be higher than the actual number of clusters in the network, in which case some clusters will contain zero vertices.

update_rule:

The type of the update rule. Possible values: IGRAPH_SPINCOMM_UPDATE_SIMPLE and IGRAPH_SPINCOMM_UPDATE_CONFIG. Basically this parameter defined the null model based on which the actual clustering is done. If this is IGRAPH_SPINCOMM_UPDATE_SIMPLE then the random graph (ie. G(n,p)), if it is IGRAPH_SPINCOMM_UPDATE then the configuration model is used. The configuration means that the baseline for the clustering is a random graph with the same degree distribution as the input graph.

gamma:

Real number. The gamma parameter of the algorithm. This defined the weight of the missing and existing links in the quality function for the clustering. The default value in the original code was 1.0, which is equal weight to missing and existing edges. Smaller values make the existing links contibute more to the energy function which is minimized in the algorithm. Bigger values make the missing links more important. (If my understanding is correct.)

Returns: 

Error code.

See also: 

igraph_community_spinglass() for the traditional version of the algorithm.

Time complexity: TODO.

3. Community structure based on eigenvectors of matrices

The function documented in these section implements the leading eigenvector method developed by Mark Newman and published in MEJ Newman: Finding community structure using the eigenvectors of matrices, Phys Rev E 74:036104 (2006).

The heart of the method is the definition of the modularity matrix, B, which is B=A-P, A being the adjacency matrix of the (undirected) network, and P contains the probability that certain edges are present according to the configuration model In other words, a Pij element of P is the probability that there is an edge between vertices i and j in a random network in which the degrees of all vertices are the same as in the input graph.

The leading eigenvector method works by calculating the eigenvector of the modularity matrix for the largest positive eigenvalue and then separating vertices into two community based on the sign of the corresponding element in the eigenvector. If all elements in the eigenvector are of the same sign that means that the network has no underlying community structure. Check Newman's paper to understand why this is a good method for detecting community structure.

The leading eigenvector community structure detection method is implemented in igraph_community_leading_eigenvector(). After the initial split, the following splits are done in a way to optimize modularity regarding to the original network. Note that any further refinement, for example using Kernighan-Lin, as proposed in Section V.A of Newman (2006), is not implemented here.

Example 24.2.  File examples/simple/igraph_community_leading_eigenvector.c

/* -*- mode: C -*-  */
/*
   IGraph library.
   Copyright (C) 2007-2012  Gabor Csardi <csardi.gabor@gmail.com>
   334 Harvard street, 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 print_vector(const igraph_vector_t *v) {
    long int i, n = igraph_vector_size(v);
    for (i = 0; i < n; i++) {
        printf("%.2g", (double)VECTOR(*v)[i]);
        if (i != n - 1) {
            printf(" ");
        }
    }
    printf("\n");
    return 0;
}

int print_matrix(const igraph_matrix_t *m) {
    long int i, j, nrow = igraph_matrix_nrow(m), ncol = igraph_matrix_ncol(m);
    for (i = 0; i < nrow; i++) {
        for (j = 0; j < ncol; j++) {
            printf("%.2g", (double)MATRIX(*m, i, j));
            if (j != ncol - 1) {
                printf(" ");
            }
        }
        printf("\n");
    }
    return 0;
}

int main() {

    igraph_t g;
    igraph_matrix_t merges;
    igraph_vector_t membership;
    igraph_vector_t x;
    igraph_arpack_options_t options;

    /* Zachary Karate club */
    igraph_small(&g, 0, IGRAPH_UNDIRECTED,
                 0,  1,  0,  2,  0,  3,  0,  4,  0,  5,
                 0,  6,  0,  7,  0,  8,  0, 10,  0, 11,
                 0, 12,  0, 13,  0, 17,  0, 19,  0, 21,
                 0, 31,  1,  2,  1,  3,  1,  7,  1, 13,
                 1, 17,  1, 19,  1, 21,  1, 30,  2,  3,
                 2,  7,  2,  8,  2,  9,  2, 13,  2, 27,
                 2, 28,  2, 32,  3,  7,  3, 12,  3, 13,
                 4,  6,  4, 10,  5,  6,  5, 10,  5, 16,
                 6, 16,  8, 30,  8, 32,  8, 33,  9, 33,
                 13, 33, 14, 32, 14, 33, 15, 32, 15, 33,
                 18, 32, 18, 33, 19, 33, 20, 32, 20, 33,
                 22, 32, 22, 33, 23, 25, 23, 27, 23, 29,
                 23, 32, 23, 33, 24, 25, 24, 27, 24, 31,
                 25, 31, 26, 29, 26, 33, 27, 33, 28, 31,
                 28, 33, 29, 32, 29, 33, 30, 32, 30, 33,
                 31, 32, 31, 33, 32, 33,
                 -1);

    igraph_matrix_init(&merges, 0, 0);
    igraph_vector_init(&membership, 0);
    igraph_vector_init(&x, 0);
    igraph_arpack_options_init(&options);

    igraph_community_leading_eigenvector(&g, /*weights=*/ 0, &merges,
                                         &membership, 1,
                                         &options, /*modularity=*/ 0,
                                         /*start=*/ 0, /*eigenvalues=*/ 0,
                                         /*eigenvectors=*/ 0, /*history=*/ 0,
                                         /*callback=*/ 0,
                                         /*callback_extra=*/ 0);

    print_matrix(&merges);
    print_vector(&membership);

    printf("\n");

    /* Make all the steps */
    igraph_community_leading_eigenvector(&g, /*weights=*/ 0, &merges,
                                         &membership, igraph_vcount(&g),
                                         &options, /*modularity=*/ 0,
                                         /*start=*/ 0, /*eigenvalues=*/ 0,
                                         /*eigenvectors=*/ 0, /*history=*/ 0,
                                         /*callback=*/ 0,
                                         /*callback_extra=*/ 0);

    print_matrix(&merges);
    print_vector(&membership);

    igraph_vector_destroy(&x);
    igraph_vector_destroy(&membership);
    igraph_matrix_destroy(&merges);
    igraph_destroy(&g);

    return 0;
}


3.1. igraph_community_leading_eigenvector — Leading eigenvector community finding (proper version).

int igraph_community_leading_eigenvector(const igraph_t *graph,
        const igraph_vector_t *weights,
        igraph_matrix_t *merges,
        igraph_vector_t *membership,
        igraph_integer_t steps,
        igraph_arpack_options_t *options,
        igraph_real_t *modularity,
        igraph_bool_t start,
        igraph_vector_t *eigenvalues,
        igraph_vector_ptr_t *eigenvectors,
        igraph_vector_t *history,
        igraph_community_leading_eigenvector_callback_t *callback,
        void *callback_extra);

Newman's leading eigenvector method for detecting community structure. This is the proper implementation of the recursive, divisive algorithm: each split is done by maximizing the modularity regarding the original network, see MEJ Newman: Finding community structure in networks using the eigenvectors of matrices, Phys Rev E 74:036104 (2006).

Arguments: 

graph:

The undirected input graph.

weights:

The weights of the edges, or a null pointer for unweighted graphs.

merges:

The result of the algorithm, a matrix containing the information about the splits performed. The matrix is built in the opposite way however, it is like the result of an agglomerative algorithm. If at the end of the algorithm (after steps steps was done) there are p communities, then these are numbered from zero to p-1. The first line of the matrix contains the first merge (which is in reality the last split) of two communities into community p, the merge in the second line forms community p+1, etc. The matrix should be initialized before calling and will be resized as needed. This argument is ignored of it is NULL.

membership:

The membership of the vertices after all the splits were performed will be stored here. The vector must be initialized before calling and will be resized as needed. This argument is ignored if it is NULL. This argument can also be used to supply a starting configuration for the community finding, in the format of a membership vector. In this case the start argument must be set to 1.

steps:

The maximum number of steps to perform. It might happen that some component (or the whole network) has no underlying community structure and no further steps can be done. If you want as many steps as possible then supply the number of vertices in the network here.

options:

The options for ARPACK. n is always overwritten. ncv is set to at least 4.

modularity:

If not a null pointer, then it must be a pointer to a real number and the modularity score of the final division is stored here.

start:

Boolean, whether to use the community structure given in the membership argument as a starting point.

eigenvalues:

Pointer to an initialized vector or a null pointer. If not a null pointer, then the eigenvalues calculated along the community structure detection are stored here. The non-positive eigenvalues, that do not result a split, are stored as well.

eigenvectors:

If not a null pointer, then the eigenvectors that are calculated in each step of the algorithm, are stored here, in a pointer vector. Each eigenvector is stored in an igraph_vector_t object. The user is responsible of deallocating the memory that belongs to the individual vectors, by calling first igraph_vector_destroy(), and then igraph_free() on them.

history:

Pointer to an initialized vector or a null pointer. If not a null pointer, then a trace of the algorithm is stored here, encoded numerically. The various operations:

IGRAPH_LEVC_HIST_START_FULL

Start the algorithm from an initial state where each connected component is a separate community.

IGRAPH_LEVC_HIST_START_GIVEN

Start the algorithm from a given community structure. The next value in the vector contains the initial number of communities.

IGRAPH_LEVC_HIST_SPLIT

Split a community into two communities. The id of the splitted community is given in the next element of the history vector. The id of the first new community is the same as the id of the splitted community. The id of the second community equals to the number of communities before the split.

IGRAPH_LEVC_HIST_FAILED

Tried to split a community, but it was not worth it, as it does not result in a bigger modularity value. The id of the community is given in the next element of the vector.

callback:

A null pointer or a function of type igraph_community_leading_eigenvector_callback_t. If given, this callback function is called after each eigenvector/eigenvalue calculation. If the callback returns a non-zero value, then the community finding algorithm stops. See the arguments passed to the callback at the documentation of igraph_community_leading_eigenvector_callback_t.

callback_extra:

Extra argument to pass to the callback function.

Returns: 

Error code.

See also: 

igraph_community_walktrap() and igraph_community_spinglass() for other community structure detection methods.

Time complexity: O(|E|+|V|^2*steps), |V| is the number of vertices, |E| the number of edges, steps the number of splits performed.

3.2. igraph_community_leading_eigenvector_callback_t — Callback for the leading eigenvector community finding method.

typedef int igraph_community_leading_eigenvector_callback_t(
    const igraph_vector_t *membership,
    long int comm,
    igraph_real_t eigenvalue,
    const igraph_vector_t *eigenvector,
    igraph_arpack_function_t *arpack_multiplier,
    void *arpack_extra,
    void *extra);

The leading eigenvector community finding implementation in igraph is able to call a callback function, after each eigenvalue calculation. This callback function must be of igraph_community_leading_eigenvector_callback_t type. The following arguments are passed to the callback:

Arguments: 

membership:

The actual membership vector, before recording the potential change implied by the newly found eigenvalue.

comm:

The id of the community that the algorithm tried to split in the last iteration. The community ids are indexed from zero here!

eigenvalue:

The eigenvalue the algorithm has just found.

eigenvector:

The eigenvector corresponding to the eigenvalue the algorithm just found.

arpack_multiplier:

A function that was passed to igraph_arpack_rssolve() to solve the last eigenproblem.

arpack_extra:

The extra argument that was passed to the ARPACK solver.

extra:

Extra argument that as passed to igraph_community_leading_eigenvector().

See also: 

3.3. igraph_le_community_to_membership — Vertex membership from the leading eigenvector community structure

int igraph_le_community_to_membership(const igraph_matrix_t *merges,
                                      igraph_integer_t steps,
                                      igraph_vector_t *membership,
                                      igraph_vector_t *csize);

This function creates a membership vector from the result of igraph_community_leading_eigenvector(), It takes membership and performs steps merges, according to the supplied merges matrix.

Arguments: 

merges:

The two-column matrix containing the merge operations. See igraph_community_walktrap() for the detailed syntax. This is usually from the output of the leading eigenvector community structure detection routines.

steps:

The number of steps to make according to merges.

membership:

Initially the starting membership vector, on output the resulting membership vector, after performing steps merges.

csize:

Optionally the sizes of the communities is stored here, if this is not a null pointer, but an initialized vector.

Returns: 

Error code.

Time complexity: O(|V|), the number of vertices.

4. Walktrap: Community structure based on random walks

4.1. igraph_community_walktrap — This function is the implementation of the Walktrap community

int igraph_community_walktrap(const igraph_t *graph,
                              const igraph_vector_t *weights,
                              int steps,
                              igraph_matrix_t *merges,
                              igraph_vector_t *modularity,
                              igraph_vector_t *membership);

finding algorithm, see Pascal Pons, Matthieu Latapy: Computing communities in large networks using random walks, https://arxiv.org/abs/physics/0512106

Currently the original C++ implementation is used in igraph, see https://www-complexnetworks.lip6.fr/~latapy/PP/walktrap.html We are grateful to Matthieu Latapy and Pascal Pons for providing this source code.

In contrast to the original implementation, isolated vertices are allowed in the graph and they are assumed to have a single incident loop edge with weight 1.

Arguments: 

graph:

The input graph, edge directions are ignored.

weights:

Numeric vector giving the weights of the edges. If it is a NULL pointer then all edges will have equal weights. The weights are expected to be positive.

steps:

Integer constant, the length of the random walks.

merges:

Pointer to a matrix, the merges performed by the algorithm will be stored here (if not NULL). Each merge is a row in a two-column matrix and contains the ids of the merged clusters. Clusters are numbered from zero and cluster numbers smaller than the number of nodes in the network belong to the individual vertices as singleton clusters. In each step a new cluster is created from two other clusters and its id will be one larger than the largest cluster id so far. This means that before the first merge we have n clusters (the number of vertices in the graph) numbered from zero to n-1. The first merge creates cluster n, the second cluster n+1, etc.

modularity:

Pointer to a vector. If not NULL then the modularity score of the current clustering is stored here after each merge operation.

membership:

Pointer to a vector. If not a NULL pointer, then the membership vector corresponding to the maximal modularity score is stored here. If it is not a NULL pointer, then neither modularity nor merges may be NULL.

Returns: 

Error code.

See also: 

Time complexity: O(|E||V|^2) in the worst case, O(|V|^2 log|V|) typically, |V| is the number of vertices, |E| is the number of edges.

Example 24.3.  File examples/simple/walktrap.c

/* -*- mode: C -*-  */
/*
   IGraph library.
   Copyright (C) 2007-2012  Gabor Csardi <csardi.gabor@gmail.com>
   334 Harvard street, 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 merges;
    igraph_vector_t modularity;
    long int no_of_nodes;
    long int i;

    igraph_rng_seed(igraph_rng_default(), 42);

    igraph_small(&g, 5, IGRAPH_UNDIRECTED,
                 0, 1, 0, 2, 0, 3, 0, 4, 1, 2, 1, 3, 1, 4, 2, 3, 2, 4, 3, 4,
                 5, 6, 5, 7, 5, 8, 5, 9, 6, 7, 6, 8, 6, 9, 7, 8, 7, 9, 8, 9, 0, 5, -1);
    igraph_vector_init(&modularity, 0);
    igraph_matrix_init(&merges, 0, 0);

    igraph_community_walktrap(&g, 0 /* no weights */,
                              4 /* steps */,
                              &merges, &modularity,
                              /* membership=*/ 0);

    no_of_nodes = igraph_vcount(&g);
    printf("Merges:\n");
    for (i = 0; i < igraph_matrix_nrow(&merges); i++) {
        printf("%2.1li + %2.li -> %2.li (modularity %4.2f)\n",
               (long int)MATRIX(merges, i, 0),
               (long int)MATRIX(merges, i, 1),
               no_of_nodes + i,
               VECTOR(modularity)[i]);
    }

    igraph_destroy(&g);

    /* isolated vertices */
    igraph_small(&g, 5, IGRAPH_UNDIRECTED, -1);
    if (igraph_community_walktrap(&g, 0 /* no weights */, 4 /* steps */, &merges,
                                  &modularity, /* membership = */ 0)) {
        return 1;
    }
    if (igraph_vector_min(&modularity) != 0 || igraph_vector_max(&modularity) != 0) {
        return 2;
    }
    igraph_destroy(&g);

    igraph_matrix_destroy(&merges);
    igraph_vector_destroy(&modularity);
    return 0;
}


5. Edge betweenness based community detection

5.1. igraph_community_edge_betweenness — Community finding based on edge betweenness.

int igraph_community_edge_betweenness(const igraph_t *graph,
                                      igraph_vector_t *result,
                                      igraph_vector_t *edge_betweenness,
                                      igraph_matrix_t *merges,
                                      igraph_vector_t *bridges,
                                      igraph_vector_t *modularity,
                                      igraph_vector_t *membership,
                                      igraph_bool_t directed,
                                      const igraph_vector_t *weights);

Community structure detection based on the betweenness of the edges in the network. The algorithm was invented by M. Girvan and M. Newman, see: M. Girvan and M. E. J. Newman: Community structure in social and biological networks, Proc. Nat. Acad. Sci. USA 99, 7821-7826 (2002).

The idea is that the betweenness of the edges connecting two communities is typically high, as many of the shortest paths between nodes in separate communities go through them. So we gradually remove the edge with highest betweenness from the network, and recalculate edge betweenness after every removal. This way sooner or later the network falls off to two components, then after a while one of these components falls off to two smaller components, etc. until all edges are removed. This is a divisive hierarchical approach, the result is a dendrogram.

Arguments: 

graph:

The input graph.

result:

Pointer to an initialized vector, the result will be stored here, the ids of the removed edges in the order of their removal. It will be resized as needed. It may be NULL if the edge IDs are not needed by the caller.

edge_betweenness:

Pointer to an initialized vector or NULL. In the former case the edge betweenness of the removed edge is stored here. The vector will be resized as needed.

merges:

Pointer to an initialized matrix or NULL. If not NULL then merges performed by the algorithm are stored here. Even if this is a divisive algorithm, we can replay it backwards and note which two clusters were merged. Clusters are numbered from zero, see the merges argument of igraph_community_walktrap() for details. The matrix will be resized as needed.

bridges:

Pointer to an initialized vector of NULL. If not NULL then all edge removals which separated the network into more components are marked here.

modularity:

If not a null pointer, then the modularity values of the different divisions are stored here, in the order corresponding to the merge matrix. The modularity values will take weights into account if weights is not null.

membership:

If not a null pointer, then the membership vector, corresponding to the highest modularity value, is stored here.

directed:

Logical constant, whether to calculate directed betweenness (i.e. directed paths) for directed graphs. It is ignored for undirected graphs.

weights:

An optional vector containing edge weights. If null, the unweighted edge betweenness scores will be calculated and used. If not null, the weighted edge betweenness scores will be calculated and used.

Returns: 

Error code.

See also: 

Time complexity: O(|V||E|^2), as the betweenness calculation requires O(|V||E|) and we do it |E|-1 times.

Example 24.4.  File examples/simple/igraph_community_edge_betweenness.c

/* -*- mode: C -*-  */
/*
   IGraph library.
   Copyright (C) 2007-2012  Gabor Csardi <csardi.gabor@gmail.com>
   334 Harvard street, 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 igraph_vector_between(const igraph_vector_t* v, const igraph_vector_t* lo,
                          const igraph_vector_t* hi) {
    return igraph_vector_all_le(lo, v) && igraph_vector_all_ge(hi, v);
}

void test_unweighted() {
    igraph_t g;
    igraph_vector_t edges, eb;
    long int i;
    long int no_of_edges;

    /* Zachary Karate club */
    igraph_small(&g, 0, IGRAPH_UNDIRECTED,
                 0,  1,  0,  2,  0,  3,  0,  4,  0,  5,
                 0,  6,  0,  7,  0,  8,  0, 10,  0, 11,
                 0, 12,  0, 13,  0, 17,  0, 19,  0, 21,
                 0, 31,  1,  2,  1,  3,  1,  7,  1, 13,
                 1, 17,  1, 19,  1, 21,  1, 30,  2,  3,
                 2,  7,  2,  8,  2,  9,  2, 13,  2, 27,
                 2, 28,  2, 32,  3,  7,  3, 12,  3, 13,
                 4,  6,  4, 10,  5,  6,  5, 10,  5, 16,
                 6, 16,  8, 30,  8, 32,  8, 33,  9, 33,
                 13, 33, 14, 32, 14, 33, 15, 32, 15, 33,
                 18, 32, 18, 33, 19, 33, 20, 32, 20, 33,
                 22, 32, 22, 33, 23, 25, 23, 27, 23, 29,
                 23, 32, 23, 33, 24, 25, 24, 27, 24, 31,
                 25, 31, 26, 29, 26, 33, 27, 33, 28, 31,
                 28, 33, 29, 32, 29, 33, 30, 32, 30, 33,
                 31, 32, 31, 33, 32, 33,
                 -1);

    igraph_vector_init(&edges, 0);
    igraph_vector_init(&eb, 0);
    igraph_community_edge_betweenness(&g, &edges, &eb, 0 /*merges */,
                                      0 /*bridges */, /*modularity=*/ 0,
                                      /*membership=*/ 0,
                                      IGRAPH_UNDIRECTED,
                                      /*weights=*/ 0);

    no_of_edges = igraph_ecount(&g);
    for (i = 0; i < no_of_edges; i++) {
        printf("%li ", (long int)VECTOR(edges)[i]);
    }
    printf("\n");

    for (i = 0; i < no_of_edges; i++) {
        printf("%.2f ", VECTOR(eb)[i]);
    }
    printf("\n");

    /* Try it once again without storage space for edges */
    igraph_community_edge_betweenness(&g, 0, &eb, 0 /*merges */,
                                      0 /*bridges */, /*modularity=*/ 0,
                                      /*membership=*/ 0,
                                      IGRAPH_UNDIRECTED,
                                      /*weights=*/ 0);
    for (i = 0; i < no_of_edges; i++) {
        printf("%.2f ", VECTOR(eb)[i]);
    }
    printf("\n");

    igraph_vector_destroy(&eb);
    igraph_vector_destroy(&edges);
    igraph_destroy(&g);
}

#define EPS 1e-4

void test_weighted() {
    igraph_t g;
    igraph_vector_t edges, eb, weights;
    igraph_real_t weights_array[] = { 4, 1, 3, 2, 5, 8, 6, 7 };

    igraph_real_t edges_array1[] = { 2, 3, 0, 1, 4, 7, 5, 6 };
    igraph_real_t edges_array2[] = { 2, 3, 6, 5, 0, 1, 4, 7 };
    igraph_real_t eb_array1_lo[] = { 4, 5, 3 + 1 / 3.0 - EPS, 4, 2.5, 4, 1, 1 };
    igraph_real_t eb_array1_hi[] = { 4, 5, 3 + 1 / 3.0 + EPS, 4, 2.5, 4, 1, 1 };
    igraph_real_t eb_array2_lo[] = { 4, 5, 3 + 1 / 3.0 - EPS, 6, 1.5, 2, 1, 1 };
    igraph_real_t eb_array2_hi[] = { 4, 5, 3 + 1 / 3.0 + EPS, 6, 1.5, 2, 1, 1 };

    igraph_vector_t edges_sol1, edges_sol2, eb_sol1_lo, eb_sol1_hi, eb_sol2_lo, eb_sol2_hi;

    igraph_vector_view(&edges_sol1, edges_array1,
                       sizeof(edges_array1) / sizeof(double));
    igraph_vector_view(&edges_sol2, edges_array2,
                       sizeof(edges_array2) / sizeof(double));
    igraph_vector_view(&eb_sol1_lo, eb_array1_lo, sizeof(eb_array1_lo) / sizeof(double));
    igraph_vector_view(&eb_sol2_lo, eb_array2_lo, sizeof(eb_array2_lo) / sizeof(double));
    igraph_vector_view(&eb_sol1_hi, eb_array1_hi, sizeof(eb_array1_hi) / sizeof(double));
    igraph_vector_view(&eb_sol2_hi, eb_array2_hi, sizeof(eb_array2_hi) / sizeof(double));

    /* Small graph as follows: A--B--C--A, A--D--E--A, B--D, C--E */
    igraph_small(&g, 0, IGRAPH_UNDIRECTED,
                 0, 1, 0, 2, 0, 3, 0, 4, 1, 2, 1, 3, 2, 4, 3, 4, -1);
    igraph_vector_view(&weights, weights_array, igraph_ecount(&g));

    igraph_vector_init(&edges, 0);
    igraph_vector_init(&eb, 0);
    igraph_community_edge_betweenness(&g, &edges, &eb, 0 /*merges */,
                                      0 /*bridges */, /*modularity=*/ 0,
                                      /*membership=*/ 0,
                                      IGRAPH_UNDIRECTED,
                                      &weights);

    if (!igraph_vector_all_e(&edges_sol1, &edges) &&
        !igraph_vector_all_e(&edges_sol2, &edges)) {
        printf("Error, edges vector was: \n");
        igraph_vector_print(&edges);
        exit(2);
    }
    if (!igraph_vector_between(&eb, &eb_sol1_lo, &eb_sol1_hi) &&
        !igraph_vector_between(&eb, &eb_sol2_lo, &eb_sol2_hi)) {
        printf("Error, eb vector was: \n");
        igraph_vector_print(&eb);
        exit(2);
    }

    /* Try it once again without storage space for edges */
    igraph_community_edge_betweenness(&g, 0, &eb, 0 /*merges */,
                                      0 /*bridges */, /*modularity=*/ 0,
                                      /*membership=*/ 0,
                                      IGRAPH_UNDIRECTED,
                                      &weights);

    if (!igraph_vector_between(&eb, &eb_sol1_lo, &eb_sol1_hi) &&
        !igraph_vector_between(&eb, &eb_sol2_lo, &eb_sol2_hi)) {
        printf("Error, eb vector was: \n");
        igraph_vector_print(&eb);
        exit(2);
    }

    igraph_vector_destroy(&eb);
    igraph_vector_destroy(&edges);
    igraph_destroy(&g);
}

void test_zero_edge_graph() {
    igraph_t g;
    igraph_vector_t eb;
    igraph_vector_t res;

    igraph_full(&g, 1, 0, 0);
    igraph_vector_init(&res, igraph_ecount(&g));
    igraph_vector_init(&eb, igraph_ecount(&g));

    igraph_community_edge_betweenness(&g,
        &res, // result
        &eb, // edge_betweenness result
        NULL, // merges result
        NULL, // bridges
        NULL, // modularity
        NULL, // membership
        IGRAPH_UNDIRECTED, // directed
        NULL // weights
        );

    igraph_vector_destroy(&eb);
    printf("No crash\n");
    igraph_vector_destroy(&res);
    igraph_destroy(&g);
}

int main() {
    test_unweighted();
    test_weighted();
    test_zero_edge_graph();
    return 0;
}


5.2. igraph_community_eb_get_merges — Calculating the merges, i.e. the dendrogram for an edge betweenness community structure.

int igraph_community_eb_get_merges(const igraph_t *graph,
                                   const igraph_bool_t directed,
                                   const igraph_vector_t *edges,
                                   const igraph_vector_t *weights,
                                   igraph_matrix_t *res,
                                   igraph_vector_t *bridges,
                                   igraph_vector_t *modularity,
                                   igraph_vector_t *membership);

This function is handy if you have a sequence of edges which are gradually removed from the network and you would like to know how the network falls apart into separate components. The edge sequence may come from the igraph_community_edge_betweenness() function, but this is not necessary. Note that igraph_community_edge_betweenness() can also calculate the dendrogram, via its merges argument.

Arguments: 

graph:

The input graph.

edges:

Vector containing the edges to be removed from the network, all edges are expected to appear exactly once in the vector.

directed:

Whether to use the directed or undirected version of modularity. Will be ignored for undirected graphs.

weights:

An optional vector containing edge weights. If null, the unweighted modularity scores will be calculated. If not null, the weighted modularity scores will be calculated. Ignored if both modularity and membership are NULL pointers.

res:

Pointer to an initialized matrix, if not NULL then the dendrogram will be stored here, in the same form as for the igraph_community_walktrap() function: the matrix has two columns and each line is a merge given by the ids of the merged components. The component ids are numbered from zero and component ids smaller than the number of vertices in the graph belong to individual vertices. The non-trivial components containing at least two vertices are numbered from n, where n is the number of vertices in the graph. So if the first line contains a and b that means that components a and b are merged into component n, the second line creates component n+1, etc. The matrix will be resized as needed.

bridges:

Pointer to an initialized vector or NULL. If not null then the index of the edge removals which split the network will be stored here. The vector will be resized as needed.

modularity:

If not a null pointer, then the modularity values for the different divisions, corresponding to the merges matrix, will be stored here.

membership:

If not a null pointer, then the membership vector for the best division (in terms of modularity) will be stored here.

Returns: 

Error code.

See also: 

Time complexity: O(|E|+|V|log|V|), |V| is the number of vertices, |E| is the number of edges.

6. Community structure based on the optimization of modularity

6.1. igraph_community_fastgreedy — Finding community structure by greedy optimization of modularity.

int igraph_community_fastgreedy(const igraph_t *graph,
                                const igraph_vector_t *weights,
                                igraph_matrix_t *merges,
                                igraph_vector_t *modularity,
                                igraph_vector_t *membership);

This function implements the fast greedy modularity optimization algorithm for finding community structure, see A Clauset, MEJ Newman, C Moore: Finding community structure in very large networks, http://www.arxiv.org/abs/cond-mat/0408187 for the details.

Some improvements proposed in K Wakita, T Tsurumi: Finding community structure in mega-scale social networks, http://www.arxiv.org/abs/cs.CY/0702048v1 have also been implemented.

Arguments: 

graph:

The input graph. It must be a graph without multiple edges. This is checked and an error message is given for graphs with multiple edges.

weights:

Potentially a numeric vector containing edge weights. Supply a null pointer here for unweighted graphs. The weights are expected to be non-negative.

merges:

Pointer to an initialized matrix or NULL, the result of the computation is stored here. The matrix has two columns and each merge corresponds to one merge, the ids of the two merged components are stored. The component ids are numbered from zero and the first n components are the individual vertices, n is the number of vertices in the graph. Component n is created in the first merge, component n+1 in the second merge, etc. The matrix will be resized as needed. If this argument is NULL then it is ignored completely.

modularity:

Pointer to an initialized vector or NULL pointer, in the former case the modularity scores along the stages of the computation are recorded here. The vector will be resized as needed.

membership:

Pointer to a vector. If not a null pointer, then the membership vector corresponding to the best split (in terms of modularity) is stored here.

Returns: 

Error code.

See also: 

igraph_community_walktrap(), igraph_community_edge_betweenness() for other community detection algorithms, igraph_community_to_membership() to convert the dendrogram to a membership vector.

Time complexity: O(|E||V|log|V|) in the worst case, O(|E|+|V|log^2|V|) typically, |V| is the number of vertices, |E| is the number of edges.

Example 24.5.  File examples/simple/igraph_community_fastgreedy.c

/* -*- mode: C -*-  */
/*
   IGraph library.
   Copyright (C) 2006-2012  Gabor Csardi <csardi.gabor@gmail.com>
   334 Harvard street, 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>

void show_results(igraph_t *g, igraph_vector_t *mod, igraph_matrix_t *merges,
                  igraph_vector_t *membership, FILE* f) {
    long int i = 0;
    igraph_vector_t our_membership;

    igraph_vector_init(&our_membership, 0);

    if (mod != 0) {
        i = igraph_vector_which_max(mod);
        fprintf(f, "Modularity:  %f\n", VECTOR(*mod)[i]);
    } else {
        fprintf(f, "Modularity:  ---\n");
    }

    if (membership != 0) {
        igraph_vector_update(&our_membership, membership);
    } else if (merges != 0) {
        igraph_community_to_membership(merges, igraph_vcount(g), i, &our_membership, 0);
    }

    printf("Membership: ");
    for (i = 0; i < igraph_vector_size(&our_membership); i++) {
        printf("%li ", (long int)VECTOR(our_membership)[i]);
    }
    printf("\n");

    igraph_vector_destroy(&our_membership);
}

int main() {
    igraph_t g;
    igraph_vector_t modularity, weights, membership;
    igraph_matrix_t merges;

    igraph_vector_init(&modularity, 0);
    igraph_matrix_init(&merges, 0, 0);
    igraph_vector_init(&weights, 0);
    igraph_vector_init(&membership, 0);

    /* Simple unweighted graph */
    igraph_small(&g, 10, IGRAPH_UNDIRECTED,
                 0, 1, 0, 2, 0, 3, 0, 4, 1, 2, 1, 3, 1, 4, 2, 3, 2, 4, 3, 4,
                 5, 6, 5, 7, 5, 8, 5, 9, 6, 7, 6, 8, 6, 9, 7, 8, 7, 9, 8, 9,
                 0, 5, -1);
    igraph_community_fastgreedy(&g, 0, &merges, &modularity, /*membership=*/ 0);
    show_results(&g, &modularity, &merges, 0, stdout);

    /* Same simple graph, with uniform edge weights */
    igraph_vector_resize(&weights, igraph_ecount(&g));
    igraph_vector_fill(&weights, 2);
    igraph_community_fastgreedy(&g, &weights, &merges, &modularity,
                                /*membership=*/ 0);
    show_results(&g, &modularity, &merges, 0, stdout);
    igraph_destroy(&g);

    /* Simple nonuniform weighted graph, with and without weights */
    igraph_small(&g, 6, IGRAPH_UNDIRECTED,
                 0, 1, 1, 2, 2, 3, 2, 4, 2, 5, 3, 4, 3, 5, 4, 5, -1);
    igraph_vector_resize(&weights, 8);
    igraph_vector_fill(&weights, 1);
    VECTOR(weights)[0] = 10;
    VECTOR(weights)[1] = 10;
    igraph_community_fastgreedy(&g, 0, &merges, &modularity, /*membership=*/ 0);
    show_results(&g, &modularity, &merges, 0, stdout);
    igraph_community_fastgreedy(&g, &weights, &merges, &modularity,
                                /*membership=*/ 0);
    show_results(&g, &modularity, &merges, 0, stdout);
    igraph_destroy(&g);

    /* Zachary Karate club */
    igraph_small(&g, 0, IGRAPH_UNDIRECTED,
                 0,  1,  0,  2,  0,  3,  0,  4,  0,  5,
                 0,  6,  0,  7,  0,  8,  0, 10,  0, 11,
                 0, 12,  0, 13,  0, 17,  0, 19,  0, 21,
                 0, 31,  1,  2,  1,  3,  1,  7,  1, 13,
                 1, 17,  1, 19,  1, 21,  1, 30,  2,  3,
                 2,  7,  2,  8,  2,  9,  2, 13,  2, 27,
                 2, 28,  2, 32,  3,  7,  3, 12,  3, 13,
                 4,  6,  4, 10,  5,  6,  5, 10,  5, 16,
                 6, 16,  8, 30,  8, 32,  8, 33,  9, 33,
                 13, 33, 14, 32, 14, 33, 15, 32, 15, 33,
                 18, 32, 18, 33, 19, 33, 20, 32, 20, 33,
                 22, 32, 22, 33, 23, 25, 23, 27, 23, 29,
                 23, 32, 23, 33, 24, 25, 24, 27, 24, 31,
                 25, 31, 26, 29, 26, 33, 27, 33, 28, 31,
                 28, 33, 29, 32, 29, 33, 30, 32, 30, 33,
                 31, 32, 31, 33, 32, 33,
                 -1);
    igraph_community_fastgreedy(&g, 0, &merges, &modularity,
                                /*membership=*/ 0);
    show_results(&g, &modularity, &merges, 0, stdout);
    igraph_destroy(&g);

    /* Simple disconnected graph with isolates */
    igraph_small(&g, 9, IGRAPH_UNDIRECTED,
                 0,  1,  0,  2,  0,  3,  1,  2,  1,  3,  2,  3,
                 4,  5,  4,  6,  4,  7,  5,  6,  5,  7,  6,  7,
                 -1);
    igraph_community_fastgreedy(&g, 0, &merges, &modularity, /*membership=*/ 0);
    show_results(&g, &modularity, &merges, 0, stdout);
    igraph_destroy(&g);

    /* Disjoint union of two rings */
    igraph_small(&g, 20, IGRAPH_UNDIRECTED,
                 0, 1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 6, 6, 7, 7, 8, 8, 9, 0, 9,
                 10, 11, 11, 12, 12, 13, 13, 14, 14, 15, 15, 16, 16, 17, 17, 18, 18, 19, 10, 19, -1);
    igraph_community_fastgreedy(&g, 0, &merges, &modularity, /*membership=*/ 0);
    show_results(&g, &modularity, &merges, 0, stdout);
    igraph_destroy(&g);

    /* Completely empty graph */
    igraph_small(&g, 10, IGRAPH_UNDIRECTED, -1);
    igraph_community_fastgreedy(&g, 0, &merges, &modularity, /*membership=*/ 0);
    show_results(&g, &modularity, &merges, 0, stdout);
    igraph_destroy(&g);

    /* Ring graph with loop edges */
    igraph_small(&g, 6, IGRAPH_UNDIRECTED,
                 0, 1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 0, 0, 0, 2, 2, -1);
    igraph_community_fastgreedy(&g, 0, &merges, &modularity, /*membership=*/ 0);
    show_results(&g, &modularity, &merges, 0, stdout);
    igraph_destroy(&g);

    /* Regression test -- graph with two vertices and two edges */
    igraph_small(&g, 2, IGRAPH_UNDIRECTED, 0, 0, 1, 1, -1);
    igraph_community_fastgreedy(&g, 0, &merges, &modularity, /*membership=*/ 0);
    show_results(&g, &modularity, &merges, 0, stdout);
    igraph_destroy(&g);

    /* Regression test -- asking for optimal membership vector but not
     * providing a modularity vector */
    igraph_small(&g, 10, IGRAPH_UNDIRECTED,
                 0, 1, 0, 2, 0, 3, 0, 4, 1, 2, 1, 3, 1, 4, 2, 3, 2, 4, 3, 4,
                 5, 6, 5, 7, 5, 8, 5, 9, 6, 7, 6, 8, 6, 9, 7, 8, 7, 9, 8, 9,
                 0, 5, -1);
    igraph_community_fastgreedy(&g, 0, &merges, 0, &membership);
    show_results(&g, 0, &merges, &membership, stdout);
    igraph_destroy(&g);

    /* Regression test -- asking for optimal membership vector but not
     * providing a merge matrix */
    igraph_small(&g, 10, IGRAPH_UNDIRECTED,
                 0, 1, 0, 2, 0, 3, 0, 4, 1, 2, 1, 3, 1, 4, 2, 3, 2, 4, 3, 4,
                 5, 6, 5, 7, 5, 8, 5, 9, 6, 7, 6, 8, 6, 9, 7, 8, 7, 9, 8, 9,
                 0, 5, -1);
    igraph_community_fastgreedy(&g, 0, 0, &modularity, &membership);
    show_results(&g, &modularity, 0, &membership, stdout);

    /* Regression test -- asking for optimal membership vector but not
     * providing a merge matrix or a modularity vector */
    igraph_community_fastgreedy(&g, 0, 0, 0, &membership);
    show_results(&g, 0, 0, &membership, stdout);
    igraph_destroy(&g);

    igraph_vector_destroy(&membership);
    igraph_vector_destroy(&modularity);
    igraph_vector_destroy(&weights);
    igraph_matrix_destroy(&merges);

    return 0;
}


6.2. igraph_community_multilevel — Finding community structure by multi-level optimization of modularity.

int igraph_community_multilevel(const igraph_t *graph,
                                const igraph_vector_t *weights,
                                const igraph_real_t resolution,
                                igraph_vector_t *membership,
                                igraph_matrix_t *memberships, igraph_vector_t *modularity);

This function implements the multi-level modularity optimization algorithm for finding community structure, see Blondel, V. D., Guillaume, J.-L., Lambiotte, R., & Lefebvre, E. (2008). Fast unfolding of communities in large networks. Journal of Statistical Mechanics: Theory and Experiment, 10008(10), 6. https://doi.org/10.1088/1742-5468/2008/10/P10008 for the details (preprint: http://arxiv.org/abs/0803.0476). The algorithm is sometimes known as the "Louvain" algorithm.

The algorithm is based on the modularity measure and a hierarchical approach. Initially, each vertex is assigned to a community on its own. In every step, vertices are re-assigned to communities in a local, greedy way: in a random order, each vertex is moved to the community with which it achieves the highest contribution to modularity. When no vertices can be reassigned, each community is considered a vertex on its own, and the process starts again with the merged communities. The process stops when there is only a single vertex left or when the modularity cannot be increased any more in a step.

The resolution parameter gamma allows finding communities at different resolutions. Higher values of the resolution parameter typically result in more, smaller communities. Lower values typically result in fewer, larger communities. The original definition of modularity is retrieved when setting gamma=1. Note that the returned modularity value is calculated using the indicated resolution parameter. See igraph_modularity() for more details. This function was contributed by Tom Gregorovic.

Arguments: 

graph:

The input graph. It must be an undirected graph.

weights:

Numeric vector containing edge weights. If NULL, every edge has equal weight. The weights are expected to be non-negative.

resolution:

Resolution parameter. Must be greater than or equal to 0. Lower values favor fewer, larger communities; higher values favor more, smaller communities. Set it to 1 to use the classical definition of modularity.

membership:

The membership vector, the result is returned here. For each vertex it gives the ID of its community. The vector must be initialized and it will be resized accordingly.

memberships:

Numeric matrix that will contain the membership vector after each level, if not NULL. It must be initialized and it will be resized accordingly.

modularity:

Numeric vector that will contain the modularity score after each level, if not NULL. It must be initialized and it will be resized accordingly.

Returns: 

Error code.

Time complexity: in average near linear on sparse graphs.

Example 24.6.  File examples/simple/igraph_community_multilevel.c

/* -*- mode: C -*-  */
/* vim:set ts=4 sts=4 sw=4 et: */
/*
   IGraph library.
   Copyright (C) 2006-2012  Gabor Csardi <csardi.gabor@gmail.com>
   334 Harvard street, 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>

void show_results(igraph_t *g, igraph_vector_t *membership, igraph_matrix_t *memberships, igraph_vector_t *modularity, FILE* f) {
    long int i, j, no_of_nodes = igraph_vcount(g);

    j = igraph_vector_which_max(modularity);
    for (i = 0; i < igraph_vector_size(membership); i++) {
        if (VECTOR(*membership)[i] != MATRIX(*memberships, j, i)) {
            fprintf(f, "WARNING: best membership vector element %li does not match the best one in the membership matrix\n", i);
        }
    }

    fprintf(f, "Modularities:\n");
    igraph_vector_print(modularity);

    for (i = 0; i < igraph_matrix_nrow(memberships); i++) {
        for (j = 0; j < no_of_nodes; j++) {
            fprintf(f, "%ld ", (long int)MATRIX(*memberships, i, j));
        }
        fprintf(f, "\n");
    }

    fprintf(f, "\n");
}

int main() {
    igraph_t g;
    igraph_vector_t modularity, membership, edges;
    igraph_matrix_t memberships;
    int i, j, k;

    igraph_vector_init(&modularity, 0);
    igraph_vector_init(&membership, 0);
    igraph_matrix_init(&memberships, 0, 0);

    igraph_rng_seed(igraph_rng_default(), 42);

    /* Unweighted test graph from the paper of Blondel et al */
    igraph_small(&g, 16, IGRAPH_UNDIRECTED,
                 0, 2, 0, 3, 0, 4, 0, 5,
                 1, 2, 1, 4, 1, 7,
                 2, 4, 2, 5, 2, 6,
                 3, 7,
                 4, 10,
                 5, 7, 5, 11,
                 6, 7, 6, 11,
                 8, 9, 8, 10, 8, 11, 8, 14, 8, 15,
                 9, 12, 9, 14,
                 10, 11, 10, 12, 10, 13, 10, 14,
                 11, 13,
                 -1);
    igraph_community_multilevel(&g, 0, 1, &membership, &memberships, &modularity);
    show_results(&g, &membership, &memberships, &modularity, stdout);

    /* Higher resolution */
    igraph_community_multilevel(&g, 0, 1.5, &membership, &memberships, &modularity);
    show_results(&g, &membership, &memberships, &modularity, stdout);

    igraph_destroy(&g);

    /* Ring of 30 cliques */
    igraph_vector_init(&edges, 0);
    for (i = 0; i < 30; i++) {
        for (j = 0; j < 5; j++) {
            for (k = j + 1; k < 5; k++) {
                igraph_vector_push_back(&edges, i * 5 + j);
                igraph_vector_push_back(&edges, i * 5 + k);
            }
        }
    }
    for (i = 0; i < 30; i++) {
        igraph_vector_push_back(&edges, i * 5 % 150);
        igraph_vector_push_back(&edges, (i * 5 + 6) % 150);
    }
    igraph_create(&g, &edges, 150, 0);
    igraph_community_multilevel(&g, 0, 1, &membership, &memberships, &modularity);
    show_results(&g, &membership, &memberships, &modularity, stdout);
    igraph_destroy(&g);

    igraph_vector_destroy(&modularity);
    igraph_vector_destroy(&membership);
    igraph_vector_destroy(&edges);
    igraph_matrix_destroy(&memberships);

    return 0;
}


6.3. igraph_community_leiden — Finding community structure using the Leiden algorithm.

int igraph_community_leiden(const igraph_t *graph,
                            const igraph_vector_t *edge_weights, const igraph_vector_t *node_weights,
                            const igraph_real_t resolution_parameter, const igraph_real_t beta, const igraph_bool_t start,
                            igraph_vector_t *membership, igraph_integer_t *nb_clusters, igraph_real_t *quality);

This function implements the Leiden algorithm for finding community structure, see Traag, V. A., Waltman, L., & van Eck, N. J. (2019). From Louvain to Leiden: guaranteeing well-connected communities. Scientific reports, 9(1), 5233. http://dx.doi.org/10.1038/s41598-019-41695-z

It is similar to the multilevel algorithm, often called the Louvain algorithm, but it is faster and yields higher quality solutions. It can optimize both modularity and the Constant Potts Model, which does not suffer from the resolution-limit (see preprint http://arxiv.org/abs/1104.3083).

The Leiden algorithm consists of three phases: (1) local moving of nodes, (2) refinement of the partition and (3) aggregation of the network based on the refined partition, using the non-refined partition to create an initial partition for the aggregate network. In the local move procedure in the Leiden algorithm, only nodes whose neighborhood has changed are visited. The refinement is done by restarting from a singleton partition within each cluster and gradually merging the subclusters. When aggregating, a single cluster may then be represented by several nodes (which are the subclusters identified in the refinement).

The Leiden algorithm provides several guarantees. The Leiden algorithm is typically iterated: the output of one iteration is used as the input for the next iteration. At each iteration all clusters are guaranteed to be connected and well-separated. After an iteration in which nothing has changed, all nodes and some parts are guaranteed to be locally optimally assigned. Finally, asymptotically, all subsets of all clusters are guaranteed to be locally optimally assigned. For more details, please see Traag, Waltman & van Eck (2019).

The objective function being optimized is

1 / 2m sum_ij (A_ij - gamma n_i n_j)d(s_i, s_j)

where m is the total edge weight, A_ij is the weight of edge (i, j), gamma is the so-called resolution parameter, n_i is the node weight of node i, s_i is the cluster of node i and d(x, y) = 1 if and only if x = y and 0 otherwise. By setting n_i = k_i, the degree of node i, and dividing gamma by 2m, you effectively obtain an expression for modularity. Hence, the standard modularity will be optimized when you supply the degrees as node_weights and by supplying as a resolution parameter 1.0/(2*m), with m the number of edges.

Arguments: 

graph:

The input graph. It must be an undirected graph.

edge_weights:

Numeric vector containing edge weights. If NULL, every edge has equal weight of 1. The weights need not be non-negative.

node_weights:

Numeric vector containing node weights.

resolution_parameter:

The resolution parameter used, which is represented by gamma in the objective function mentioned in the documentation.

beta:

The randomness used in the refinement step when merging. A small amount of randomness (beta = 0.01) typically works well.

start:

Start from membership vector. If this is true, the optimization will start from the provided membership vector. If this is false, the optimization will start from a singleton partition.

membership:

The membership vector. This is both used as the initial membership from which optimisation starts and is updated in place. It must hence be properly initialized. When finding clusters from scratch it is typically started using a singleton clustering. This can be achieved using igraph_vector_init_seq.

nb_clusters:

The number of clusters contained in membership. Must not be a NULL pointer.

quality:

The quality of the partition, in terms of the objective function as included in the documentation. If NULL the quality will not be calculated.

Returns: 

Error code.

Time complexity: near linear on sparse graphs.

Example 24.7.  File examples/simple/igraph_community_leiden.c

/* -*- mode: C -*-  */
/*
   IGraph library.
   Copyright (C) 2006-2012  Gabor Csardi <csardi.gabor@gmail.com>
   334 Harvard street, 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 graph;
    igraph_vector_t membership, degree;
    igraph_integer_t nb_clusters;
    igraph_real_t quality;

    /* Set default seed to get reproducible results */
    igraph_rng_seed(igraph_rng_default(), 0);

    /* Simple unweighted graph */
    igraph_small(&graph, 10, IGRAPH_UNDIRECTED,
                 0, 1, 0, 2, 0, 3, 0, 4, 1, 2, 1, 3, 1, 4, 2, 3, 2, 4, 3, 4,
                 5, 6, 5, 7, 5, 8, 5, 9, 6, 7, 6, 8, 6, 9, 7, 8, 7, 9, 8, 9,
                 0, 5, -1);

    /* Perform Leiden algorithm using CPM */
    igraph_vector_init(&membership, igraph_vcount(&graph));
    igraph_community_leiden(&graph, NULL, NULL, 0.05, 0.01, 0, &membership, &nb_clusters, &quality);

    printf("Leiden found %" IGRAPH_PRId " clusters using CPM (resolution parameter 0.05), quality is %.4f.\n", nb_clusters, quality);
    printf("Membership: ");
    igraph_vector_print(&membership);
    printf("\n");

    /* Start from existing membership to improve it further */
    igraph_community_leiden(&graph, NULL, NULL, 0.05, 0.01, 1, &membership, &nb_clusters, &quality);

    printf("Iterated Leiden, using CPM (resolution parameter 0.05), quality is %.4f.\n", quality);
    printf("Membership: ");
    igraph_vector_print(&membership);
    printf("\n");

    /* Initialize degree vector to use for optimizing modularity */
    igraph_vector_init(&degree, igraph_vcount(&graph));
    igraph_degree(&graph, &degree, igraph_vss_all(), IGRAPH_ALL, 1);

    /* Perform Leiden algorithm using modularity */
    igraph_community_leiden(&graph, NULL, &degree, 1.0 / (2 * igraph_ecount(&graph)), 0.01, 0, &membership, &nb_clusters, &quality);

    printf("Leiden found %" IGRAPH_PRId " clusters using modularity, quality is %.4f.\n", nb_clusters, quality);
    printf("Membership: ");
    igraph_vector_print(&membership);
    printf("\n");

    igraph_vector_destroy(&degree);
    igraph_vector_destroy(&membership);
    igraph_destroy(&graph);

    return 0;
}


7. Fluid communities

7.1. igraph_community_fluid_communities — Community detection based on fluids interacting on the graph.

int igraph_community_fluid_communities(const igraph_t *graph,
                                       igraph_integer_t no_of_communities,
                                       igraph_vector_t *membership,
                                       igraph_real_t *modularity);

The algorithm is based on the simple idea of several fluids interacting in a non-homogeneous environment (the graph topology), expanding and contracting based on their interaction and density. This function implements the community detection method described in: Parés F, Gasulla DG, et. al. (2018) Fluid Communities: A Competitive, Scalable and Diverse Community Detection Algorithm. In: Complex Networks & Their Applications VI: Proceedings of Complex Networks 2017 (The Sixth International Conference on Complex Networks and Their Applications), Springer, vol 689, p 229.

Arguments: 

graph:

The input graph. The graph must be simple and connected. Empty graphs are not supported as well as single vertex graphs. Edge directions are ignored. Weights are not considered.

no_of_communities:

The number of communities to be found. Must be greater than 0 and fewer than number of vertices in the graph.

membership:

The result vector mapping vertices to the communities they are assigned to.

modularity:

If not a null pointer, then it must be a pointer to a real number. The modularity score of the detected community structure is stored here.

Returns: 

Error code.

Time complexity: O(|E|)

Example 24.8.  File examples/simple/igraph_community_fluid_communities.c

/* -*- mode: C -*-  */
/* vim:set ts=4 sw=4 sts=4 et: */
/*
   IGraph library.
   Copyright (C) 2007-2012  Gabor Csardi <csardi.gabor@gmail.com>
   334 Harvard street, 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>

#include "../../tests/unit/test_utilities.inc"

int main() {
    igraph_t g;
    igraph_integer_t k;
    igraph_vector_t membership;
    igraph_real_t modularity;

    igraph_rng_seed(igraph_rng_default(), 247);

    /* Empty graph */
    igraph_small(&g, 0, IGRAPH_UNDIRECTED, -1);
    igraph_vector_init(&membership, 0);
    igraph_vector_push_back(&membership, 1);
    igraph_community_fluid_communities(&g, 2, &membership, &modularity);
    if (!igraph_is_nan(modularity) || igraph_vector_size(&membership) != 0) {
        return 2;
    }
    igraph_vector_destroy(&membership);
    igraph_destroy(&g);

    /* Graph with one vertex only */
    igraph_small(&g, 1, IGRAPH_UNDIRECTED, -1);
    igraph_vector_init(&membership, 0);
    igraph_community_fluid_communities(&g, 2, &membership, &modularity);
    if (!igraph_is_nan(modularity) || igraph_vector_size(&membership) != 1 || VECTOR(membership)[0] != 0) {
        return 3;
    }
    igraph_vector_destroy(&membership);
    igraph_destroy(&g);

    /* Zachary Karate club -- this is just a quick smoke test */
    igraph_small(&g, 0, IGRAPH_UNDIRECTED,
                 0,  1,  0,  2,  0,  3,  0,  4,  0,  5,
                 0,  6,  0,  7,  0,  8,  0, 10,  0, 11,
                 0, 12,  0, 13,  0, 17,  0, 19,  0, 21,
                 0, 31,  1,  2,  1,  3,  1,  7,  1, 13,
                 1, 17,  1, 19,  1, 21,  1, 30,  2,  3,
                 2,  7,  2,  8,  2,  9,  2, 13,  2, 27,
                 2, 28,  2, 32,  3,  7,  3, 12,  3, 13,
                 4,  6,  4, 10,  5,  6,  5, 10,  5, 16,
                 6, 16,  8, 30,  8, 32,  8, 33,  9, 33,
                 13, 33, 14, 32, 14, 33, 15, 32, 15, 33,
                 18, 32, 18, 33, 19, 33, 20, 32, 20, 33,
                 22, 32, 22, 33, 23, 25, 23, 27, 23, 29,
                 23, 32, 23, 33, 24, 25, 24, 27, 24, 31,
                 25, 31, 26, 29, 26, 33, 27, 33, 28, 31,
                 28, 33, 29, 32, 29, 33, 30, 32, 30, 33,
                 31, 32, 31, 33, 32, 33,
                 -1);

    igraph_vector_init(&membership, 0);
    k = 2;
    igraph_community_fluid_communities(&g, k, &membership,
                                       /*modularity=*/ 0);
    if (!igraph_vector_contains(&membership, 0) || !igraph_vector_contains(&membership, 1)) {
        printf("Resulting graph does not have exactly 2 communities as expected.\n");
        igraph_vector_print(&membership);
        return 1;
    }

    igraph_destroy(&g);
    igraph_vector_destroy(&membership);

    VERIFY_FINALLY_STACK();

    return 0;
}


8. Label propagation

8.1. igraph_community_label_propagation — Community detection based on label propagation.

int igraph_community_label_propagation(const igraph_t *graph,
                                       igraph_vector_t *membership,
                                       const igraph_vector_t *weights,
                                       const igraph_vector_t *initial,
                                       const igraph_vector_bool_t *fixed,
                                       igraph_real_t *modularity);

This function implements the label propagation-based community detection algorithm described by Raghavan, Albert and Kumara. This version extends the original method by the ability to take edge weights into consideration and also by allowing some labels to be fixed.

Weights are taken into account as follows: when the new label of node i is determined, the algorithm iterates over all edges incident on node i and calculate the total weight of edges leading to other nodes with label 0, 1, 2, ..., k - 1 (where k is the number of possible labels). The new label of node i will then be the label whose edges (among the ones incident on node i) have the highest total weight.

Reference:

Raghavan, U.N. and Albert, R. and Kumara, S.: Near linear time algorithm to detect community structures in large-scale networks. Phys Rev E 76, 036106. (2007). https://doi.org/10.1103/PhysRevE.76.036106

Arguments: 

graph:

The input graph, should be undirected to make sense.

membership:

The membership vector, the result is returned here. For each vertex it gives the ID of its community (label).

weights:

The weight vector, it should contain a positive weight for all the edges.

initial:

The initial state. If NULL, every vertex will have a different label at the beginning. Otherwise it must be a vector with an entry for each vertex. Non-negative values denote different labels, negative entries denote vertices without labels. Unlabeled vertices which are not reachable from any labeled ones will remain unlabeled at the end of the label propagation process, and will be labeled in an additional step to avoid returning negative values in membership. In undirected graphs, this happens when entire connected components are unlabeled. Then, each unlabeled component will receive its own separate label. In directed graphs, the outcome of the additional labeling should be considered undefined and may change in the future; please do not rely on it.

fixed:

Boolean vector denoting which labels are fixed. Of course this makes sense only if you provided an initial state, otherwise this element will be ignored. Also note that vertices without labels cannot be fixed. If they are, this vector will be modified to make it consistent with initial.

modularity:

If not a null pointer, then it must be a pointer to a real number. The modularity score of the detected community structure is stored here.

Returns: 

Error code.

Time complexity: O(m+n)

Example 24.9.  File examples/simple/igraph_community_label_propagation.c

/* -*- mode: C -*-  */
/* vim:set ts=4 sw=4 sts=4 et: */
/*
   IGraph library.
   Copyright (C) 2007-2020  The igraph development team <igraph@igraph.org>

   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, see <https://www.gnu.org/licenses/>.
*/

#include <igraph.h>
#include <stdio.h>

int main() {
    igraph_t graph;
    igraph_vector_t membership;
    igraph_real_t modularity;

    igraph_famous(&graph, "Zachary"); /* We use Zachary's karate club network. */

    /* Label propagation is a stochastic method; the result will depend on the random seed. */
    igraph_rng_seed(igraph_rng_default(), 123);

    /* All igraph functions that returns their result in an igraph_vector_t must be given
       an already initialized vector. */
    igraph_vector_init(&membership, 0);
    igraph_community_label_propagation(
                &graph, &membership,
                /* weights= */ NULL, /* initial= */ NULL, /* fixed= */ NULL,
                &modularity);

    printf("%ld communities found; modularity score is %g.\n",
           (long int) (igraph_vector_max(&membership) + 1),
           modularity);

    /* Destroy data structures at the end. */
    igraph_vector_destroy(&membership);
    igraph_destroy(&graph);

    return 0;
}


9. The InfoMAP algorithm

9.1. igraph_community_infomap — Find community structure that minimizes the expected

int igraph_community_infomap(const igraph_t * graph,
                             const igraph_vector_t *e_weights,
                             const igraph_vector_t *v_weights,
                             int nb_trials,
                             igraph_vector_t *membership,
                             igraph_real_t *codelength);

description length of a random walker trajectory. Implementation of the InfoMap community detection algorithm.of Martin Rosvall and Carl T. Bergstrom. See : Visualization of the math and the map generator: www.mapequation.org [2] The original paper: M. Rosvall and C. T. Bergstrom, Maps of information flow reveal community structure in complex networks, PNAS 105, 1118 (2008) [http://dx.doi.org/10.1073/pnas.0706851105 , http://arxiv.org/abs/0707.0609 ] [3] A more detailed paper: M. Rosvall, D. Axelsson, and C. T. Bergstrom, The map equation, Eur. Phys. J. Special Topics 178, 13 (2009). [http://dx.doi.org/10.1140/epjst/e2010-01179-1 , http://arxiv.org/abs/0906.1405 ]

The original C++ implementation of Martin Rosvall is used, see http://www.tp.umu.se/~rosvall/downloads/infomap_undir.tgz . Intergation in igraph has be done by Emmanuel Navarro (who is grateful to * Martin Rosvall and Carl T. Bergstrom for providing this source code.)

Note that the graph must not contain isolated vertices.

If you want to specify a random seed (as in original implementation) you can use igraph_rng_seed().

Arguments: 

graph:

The input graph.

e_weights:

Numeric vector giving the weights of the edges. If it is a NULL pointer then all edges will have equal weights. The weights are expected to be positive.

v_weights:

Numeric vector giving the weights of the vertices. If it is a NULL pointer then all vertices will have equal weights. The weights are expected to be positive.

nb_trials:

The number of attempts to partition the network (can be any integer value equal or larger than 1).

membership:

Pointer to a vector. The membership vector is stored here.

codelength:

Pointer to a real. If not NULL the code length of the partition is stored here.

Returns: 

Error code.

See also: 

Time complexity: TODO.