blockmodel.py 89.6 KB
Newer Older
1
2
3
4
5
#! /usr/bin/env python
# -*- coding: utf-8 -*-
#
# graph_tool -- a general graph manipulation python module
#
Tiago Peixoto's avatar
Tiago Peixoto committed
6
# Copyright (C) 2006-2014 Tiago de Paula Peixoto <tiago@skewed.de>
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
#
# 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 3 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 <http://www.gnu.org/licenses/>.

from __future__ import division, absolute_import, print_function
import sys
if sys.version_info < (3,):
    range = xrange

26
27
from .. import _degree, _prop, Graph, GraphView, libcore, _get_rng, PropertyMap
from .. stats import label_self_loops
28
29
import random
from numpy import *
30
import numpy
31
32
from scipy.optimize import fsolve, fminbound
import scipy.special
33
from collections import defaultdict
34
35
import copy
import heapq
36
37
38
39
40

from .. dl_import import dl_import
dl_import("from . import libgraph_tool_community as libcommunity")


41
42
43
44
45
46
47
48
49
50
51
52
def get_block_graph(g, B, b, vcount, ecount):
    cg, br, vcount, ecount = condensation_graph(g, b,
                                                vweight=vcount,
                                                eweight=ecount,
                                                self_loops=True)[:4]
    cg.vp["count"] = vcount
    cg.ep["count"] = ecount
    cg = Graph(cg, vorder=br)

    cg.add_vertex(B - cg.num_vertices())
    return cg

53
54
55
56
57
58
59
60
class BlockState(object):
    r"""This class encapsulates the block state of a given graph.

    This must be instantiated and used by functions such as :func:`mcmc_sweep`.

    Parameters
    ----------
    g : :class:`~graph_tool.Graph`
61
        Graph to be modelled.
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
    eweight : :class:`~graph_tool.PropertyMap` (optional, default: ``None``)
        Edge weights (i.e. multiplicity).
    vweight : :class:`~graph_tool.PropertyMap` (optional, default: ``None``)
        Vertex weights (i.e. multiplicity).
    b : :class:`~graph_tool.PropertyMap` (optional, default: ``None``)
        Initial block labels on the vertices. If not supplied, it will be
        randomly sampled.
    B : ``int`` (optional, default: ``None``)
        Number of blocks. If not supplied it will be either obtained from the
        parameter ``b``, or set to the maximum possible value according to the
        minimum description length.
    clabel : :class:`~graph_tool.PropertyMap` (optional, default: ``None``)
        This parameter provides a constraint label, such that vertices with
        different labels will not be allowed to belong to the same block. If not given,
        all labels will be assumed to be the same.
    deg_corr : ``bool`` (optional, default: ``True``)
        If ``True``, the degree-corrected version of the blockmodel ensemble will
        be assumed, otherwise the traditional variant will be used.
80
81
82
83
    max_BE : ``int`` (optional, default: ``1000``)
        If the number of blocks exceeds this number, a sparse representation of
        the block graph is used, which is slightly less efficient, but uses less
        memory,
84
85
    """

86
87
    def __init__(self, g, eweight=None, vweight=None, b=None,
                 B=None, clabel=None, deg_corr=True, max_BE=1000):
88
89
90
91
92
93
94
95
96
97
98
99
100
101
        self.g = g
        if eweight is None:
            eweight = g.new_edge_property("int")
            eweight.a = 1
        elif eweight.value_type() != "int32_t":
            eweight = eweight.copy(value_type="int32_t")
        if vweight is None:
            vweight = g.new_vertex_property("int")
            vweight.a = 1
        elif vweight.value_type() != "int32_t":
            vweight = vweight.copy(value_type="int32_t")
        self.eweight = eweight
        self.vweight = vweight

102
103
        self.E = int(self.eweight.fa.sum())
        self.N = int(self.vweight.fa.sum())
104
105
106
107
108
109

        self.deg_corr = deg_corr

        if b is None:
            if B is None:
                B = get_max_B(self.N, self.E, directed=g.is_directed())
110
111
112
            ba = random.randint(0, B, g.num_vertices())
            ba[:B] = arange(B)        # avoid empty blocks
            random.shuffle(ba)
113
            b = g.new_vertex_property("int")
114
            b.fa = ba
115
116
117
            self.b = b
        else:
            if B is None:
118
                B = int(b.fa.max()) + 1
119
120
            self.b = b = b.copy(value_type="int32_t")

121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
        if b.fa.max() >= B:
            raise ValueError("Maximum value of b is larger or equal to B!")

        # Construct block-graph
        self.bg = get_block_graph(g, B, b, vweight, eweight)
        self.bg.set_fast_edge_removal()

        self.mrs = self.bg.ep["count"]
        self.wr = self.bg.vp["count"]
        del self.bg.ep["count"]
        del self.bg.vp["count"]

        self.mrp = self.bg.degree_property_map("out", weight=self.mrs)

        if g.is_directed():
            self.mrm = self.bg.degree_property_map("in", weight=self.mrs)
        else:
            self.mrm = self.mrp
139
140
141

        self.vertices = libcommunity.get_vector(B)
        self.vertices.a = arange(B)
142
        self.B = B
143

144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
        self.clabel = clabel
        if self.clabel is None:
            self.clabel = self.g.new_vertex_property("int")

        self.bclabel = self.bg.new_vertex_property("int")
        libcommunity.vector_rmap(self.b.a, self.bclabel.a)
        libcommunity.vector_map(self.bclabel.a, self.clabel.a)

        self.emat = None
        if max_BE is None:
            max_BE = 1000
        self.max_BE = max_BE

        # used by mcmc_sweep()
        self.egroups = None
        self.nsampler = None
        self.sweep_vertices = None

        # used by merge_sweep()
        self.bnsampler = None
164

165
166
167
        libcommunity.init_safelog(int(2 * max(self.E, self.N)))
        libcommunity.init_xlogx(int(2 * max(self.E, self.N)))
        libcommunity.init_lgamma(int(3 * max(self.E, self.N)))
168

169
170
171
172
    def __get_emat(self):
        if self.emat is None:
            self.__regen_emat()
        return self.emat
173
174

    def __regen_emat(self):
175
176
177
178
        if self.B <= self.max_BE:
            self.emat = libcommunity.create_emat(self.bg._Graph__graph)
        else:
            self.emat = libcommunity.create_ehash(self.bg._Graph__graph)
179

180
    def __build_egroups(self, empty=False):
181
182
183
        self.esrcpos = self.g.new_edge_property("int")
        self.etgtpos = self.g.new_edge_property("int")
        self.is_weighted = True if self.eweight.fa.max() > 1 else False
184
        self.egroups = libcommunity.build_egroups(self.g._Graph__graph,
185
186
187
188
189
                                                  self.bg._Graph__graph,
                                                  _prop("v", self.g, self.b),
                                                  _prop("e", self.g, self.eweight),
                                                  _prop("e", self.g, self.esrcpos),
                                                  _prop("e", self.g, self.etgtpos),
190
                                                  self.is_weighted, empty)
191
192
193
194
195
196
197
198
199
200
201
202
203

    def __build_nsampler(self):
        self.nsampler = libcommunity.init_neighbour_sampler(self.g._Graph__graph,
                                                            _prop("e", self.g, self.eweight))
    def __build_bnsampler(self):
        self.bnsampler = libcommunity.init_neighbour_sampler(self.bg._Graph__graph,
                                                             _prop("e", self.bg, self.mrs))

    def __cleanup_bg(self):
        emask = self.bg.new_edge_property("bool")
        emask.a = self.mrs.a[:len(emask.a)] > 0
        self.bg.set_edge_filter(emask)
        self.bg.purge_edges()
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222

    def get_blocks(self):
        r"""Returns the property map which contains the block labels for each vertex."""
        return self.b

    def get_bg(self):
        r"""Returns the block graph."""
        return self.bg

    def get_ers(self):
        r"""Returns the edge property map of the block graph which contains the :math:`e_{rs}` matrix entries."""
        return self.mrs

    def get_er(self):
        r"""Returns the vertex property map of the block graph which contains the number
        :math:`e_r` of half-edges incident on block :math:`r`. If the graph is
        directed, a pair of property maps is returned, with the number of
        out-edges :math:`e^+_r` and in-edges :math:`e^-_r`, respectively."""
        if self.bg.is_directed():
223
            return self.mrp, self.mrm
224
225
226
227
228
229
230
231
232
233
234
235
236
237
        else:
            return self.mrp

    def get_nr(self):
        r"""Returns the vertex property map of the block graph which contains the block sizes :math:`n_r`."""
        return self.wr

    def get_eweight(self):
        r"""Returns the block edge counts associated with the block matrix
        :math:`e_{rs}`. For directed graphs it is identical to :math:`e_{rs}`,
        but for undirected graphs it is identical except for the diagonal, which
        is :math:`e_{rr}/2`."""
        eweight = self.mrs.copy()
        if not self.g.is_directed():
238
239
            sl = label_self_loops(self.bg, mark_only=True)
            eweight.a[sl.a > 0] /= 2
240
241
        return eweight

242
243
    def entropy(self, complete=False, random=False, dl=False, dense=False,
                multigraph=False):
244
245
246
247
248
249
250
251
252
253
254
255
        r"""Calculate the entropy per edge associated with the current block partition.

        Parameters
        ----------
        complete : ``bool`` (optional, default: ``False``)
            If ``True``, the complete entropy will be returned, including constant
            terms not relevant to the block partition.
        random : ``bool`` (optional, default: ``False``)
            If ``True``, the entropy entropy corresponding to an equivalent random
            graph (i.e. no block partition) will be returned.
        dl : ``bool`` (optional, default: ``False``)
            If ``True``, the full description length will be returned.
256
257
258
259
260
        dense : ``bool`` (optional, default: ``False``)
            If ``True``, the "dense" variant of the entropy will be computed.
        multigraph : ``bool`` (optional, default: ``False``)
            If ``True``, the multigraph entropy will be used. Only has an effect
            if ``dense == True``.
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300

        Notes
        -----

        For the traditional blockmodel (``deg_corr == False``), the entropy is
        given by

        .. math::

          \mathcal{S}_t &\cong E - \frac{1}{2} \sum_{rs}e_{rs}\ln\left(\frac{e_{rs}}{n_rn_s}\right), \\
          \mathcal{S}^d_t &\cong E - \sum_{rs}e_{rs}\ln\left(\frac{e_{rs}}{n_rn_s}\right),

        for undirected and directed graphs, respectively, where :math:`e_{rs}`
        is the number of edges from block :math:`r` to :math:`s` (or the number
        of half-edges for the undirected case when :math:`r=s`), and :math:`n_r`
        is the number of vertices in block :math:`r` .


        For the degree-corrected variant with "hard" degree constraints the
        equivalent expressions are

        .. math::

            \mathcal{S}_c &\cong -E -\sum_kN_k\ln k! - \frac{1}{2} \sum_{rs}e_{rs}\ln\left(\frac{e_{rs}}{e_re_s}\right), \\
            \mathcal{S}^d_c &\cong -E -\sum_{k^+}N_{k^+}\ln k^+!  -\sum_{k^-}N_{k^-}\ln k^-! - \sum_{rs}e_{rs}\ln\left(\frac{e_{rs}}{e^+_re^-_s}\right),

        where :math:`e_r = \sum_se_{rs}` is the number of half-edges incident on
        block :math:`r`, and :math:`e^+_r = \sum_se_{rs}` and :math:`e^-_r =
        \sum_se_{sr}` are the number of out- and in-edges adjacent to block
        :math:`r`, respectively.

        If ``complete == False`` only the last term of the equations above will
        be returned. If ``random == True`` it will be assumed that :math:`B=1`
        despite the actual :math:`e_{rs}` matrix.  If ``dl == True``, the
        description length :math:`\mathcal{L}_t` of the model will be returned
        as well, as described in :func:`model_entropy`. Note that for the
        degree-corrected version the description length is

        .. math::

301
            \mathcal{L}_c = \mathcal{L}_t - \sum_rn_r\sum_kp^r_k\ln p^r_k,
302

303
304
        where :math:`p^r_k` is the fraction of nodes in block $r$ with degree :math:`k`. For directed
        graphs we have instead :math:`k \to (k^-, k^+)`.
305

306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
        If the "dense" entropies are requested, they will be computed as

        .. math::

            \mathcal{S}_t  &= \sum_{r>s} \ln{\textstyle {n_rn_s \choose e_{rs}}} + \sum_r \ln{\textstyle {{n_r\choose 2}\choose e_{rr}/2}}\\
            \mathcal{S}^d_t  &= \sum_{rs} \ln{\textstyle {n_rn_s \choose e_{rs}}},

        for simple graphs, and

        .. math::

            \mathcal{S}_m  &= \sum_{r>s} \ln{\textstyle \left(\!\!{n_rn_s \choose e_{rs}}\!\!\right)} + \sum_r \ln{\textstyle \left(\!\!{\left(\!{n_r\choose 2}\!\right)\choose e_{rr}/2}\!\!\right)}\\
            \mathcal{S}^d_m  &= \sum_{rs} \ln{\textstyle \left(\!\!{n_rn_s \choose e_{rs}}\!\!\right)},

        for multigraphs (i.e. ``multigraph == True``).

        Note that in all cases the value returned corresponds to the entropy `per edge`,
323
324
325
326
327
328
329
        i.e. :math:`(\mathcal{S}_{t/c}\; [\,+\, \mathcal{L}_{t/c}])/ E`.

        """

        E = self.E
        N = self.N

330
331
332
333
334
335
336
337
338
        if dense:
            if random:
                bg = get_block_graph(self.bg, 1,
                                     self.bg.new_vertex_property("int"),
                                     self.wr, self.mrs)
                S = libcommunity.entropy_dense(bg._Graph__graph,
                                               _prop("e", bg, bg.ep["count"]),
                                               _prop("v", bg, bg.vp["count"]),
                                               multigraph)
339
            else:
340
341
342
343
                S = libcommunity.entropy_dense(self.bg._Graph__graph,
                                               _prop("e", self.bg, self.mrs),
                                               _prop("v", self.bg, self.wr),
                                               multigraph)
344
        else:
345
346
347
348
349
350
351
352
            if self.deg_corr:
                if self.g.is_directed():
                    S_rand = E * log(E)
                else:
                    S_rand = E * log(2 * E)
            else:
                ak = E / float(N) if self.g.is_directed() else  2 * E / float(N)
                S_rand = E * log (N / ak)
353

354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
            if random:
                S = S_rand
            else:
                S = libcommunity.entropy(self.bg._Graph__graph,
                                         _prop("e", self.bg, self.mrs),
                                         _prop("v", self.bg, self.mrp),
                                         _prop("v", self.bg, self.mrm),
                                         _prop("v", self.bg, self.wr),
                                         self.deg_corr)

            if complete:
                if self.deg_corr:
                    S -= E
                    for v in self.g.vertices():
                        S -= scipy.special.gammaln(v.out_degree() + 1)
                        if self.g.is_directed():
                            S -= scipy.special.gammaln(v.in_degree() + 1)
                else:
                    S += E
            else:
                S -= S_rand
375

376
377
378
379
380
381
        if dl:
            if random:
                S += model_entropy(1, N, E, directed=self.g.is_directed()) * E
            else:
                S += model_entropy(self.B, N, E, directed=self.g.is_directed(), nr=self.wr.a) * E

382
            if self.deg_corr:
383
384
385
                S_seq = libcommunity.deg_entropy(self.g._Graph__graph,
                                                 _prop("v", self.g, self.b),
                                                 self.B)
386
387
388
389
390
                S += S_seq

        return S / E

    def remove_vertex(self, v):
391
        r"""Remove vertex ``v`` from its current block."""
392
393
394
395
396
397
398
399
        libcommunity.remove_vertex(self.g._Graph__graph,
                                   self.bg._Graph__graph,
                                   int(v),
                                   _prop("e", self.bg, self.mrs),
                                   _prop("v", self.bg, self.mrp),
                                   _prop("v", self.bg, self.mrm),
                                   _prop("v", self.bg, self.wr),
                                   _prop("v", self.g, self.b))
400
401
402
        self.egroups = None
        self.nb_list = None
        self.nb_count = None
403
404
405


    def add_vertex(self, v, r):
406
        r"""Add vertex ``v`` to block ``r``."""
407
408
409
410
411
412
413
414
        libcommunity.add_vertex(v.get_graph()._Graph__graph,
                                self.bg._Graph__graph,
                                int(v), int(r),
                                _prop("e", self.bg, self.mrs),
                                _prop("v", self.bg, self.mrp),
                                _prop("v", self.bg, self.mrm),
                                _prop("v", self.bg, self.wr),
                                _prop("v", self.g, self.b))
415
416
417
        self.egroups = None
        self.nb_list = None
        self.nb_count = None
418
419

    def move_vertex(self, v, nr):
420
        r"""Move vertex ``v`` to block ``nr``, and return the entropy difference."""
421
422
        dS = libcommunity.move_vertex(self.g._Graph__graph,
                                      self.bg._Graph__graph,
423
                                      self.__get_emat(),
424
425
426
427
428
429
430
431
432
                                      int(v), int(nr),
                                      _prop("e", self.bg, self.mrs),
                                      _prop("v", self.bg, self.mrp),
                                      _prop("v", self.bg, self.mrm),
                                      _prop("v", self.bg, self.wr),
                                      _prop("v", self.g, self.b),
                                      self.deg_corr,
                                      _prop("e", self.bg, self.eweight),
                                      _prop("v", self.bg, self.vweight))
433
434
435
        self.egroups = None
        self.nb_list = None
        self.nb_count = None
436
437
        return dS / float(self.E)

438
    def get_matrix(self, reorder=False, niter=0, ret_order=False):
439
440
        r"""Returns the block matrix, which contains the number of edges between
        each block pair.
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486

        Parameters
        ----------
        reorder : ``bool`` (optional, default: ``False``)
            If ``True``, the matrix is reordered so that blocks which are
            'similar' are close together.
        niter : ``int`` (optional, default: `0`)
            Number of iterations performed to obtain the best ordering. If
            ``niter == 0`` it will automatically determined. Only has effect
            if ``reorder == True``.
        ret_order : ``bool`` (optional, default: ``False``)
            If ``True``, the vertex ordering is returned. Only has effect if
            ``reorder == True``.

        Examples
        --------

        .. testsetup:: get_matrix

           gt.seed_rng(42)
           np.random.seed(42)
           from pylab import *

        .. doctest:: get_matrix

           >>> g = gt.collection.data["polbooks"]
           >>> state = gt.BlockState(g, B=5, deg_corr=True)
           >>> for i in range(1000):
           ...     ds, nmoves = gt.mcmc_sweep(state)
           >>> m = state.get_matrix(reorder=True)
           >>> figure()
           <...>
           >>> matshow(m)
           <...>
           >>> savefig("bloc_mat.pdf")

        .. testcleanup:: get_matrix

           savefig("bloc_mat.png")

        .. figure:: bloc_mat.*
           :align: center

           A  5x5 block matrix.

       """
487
        B = self.B
488
489
490
491
492
493
        vmap = {}
        for r in range(len(self.vertices)):
            vmap[self.vertices[r]] = r

        if reorder:
            if niter == 0:
494
                niter = 10
495
496
497
498
499

            states = []

            label = None
            states = [self]
500
            Bi = self.B // 2
501
502

            while Bi > 1:
503
504
505
506
507

                state = BlockState(states[-1].bg,
                                   b=states[-1].bg.vertex_index.copy("int"),
                                   B=states[-1].bg.num_vertices(),
                                   clabel=states[-1].bclabel,
508
                                   vweight=states[-1].wr,
509
510
511
512
513
514
                                   eweight=states[-1].mrs,
                                   deg_corr=states[-1].deg_corr,
                                   max_BE=states[-1].max_BE)

                state = greedy_shrink(state, B=Bi, nsweeps=niter,
                                      epsilon=1e-3, c=0,
515
                                      nmerge_sweeps=niter)
516
517

                for i in range(niter):
518
                    mcmc_sweep(state, c=0, beta=float("inf"))
519
520
521

                states.append(state)

522
                Bi //= 2
523

524
                if Bi < self.bclabel.a.max() + 1:
525
526
                    break

527
            vorder = list(range(len(states[-1].vertices)))
528
            for state in reversed(states[1:]):
529
                norder = [[] for i in range(state.B)]
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
                for v in state.g.vertices():
                    pos = vorder.index(state.b[v])
                    norder[pos].append(int(v))
                vorder = [item for sublist in norder for item in sublist]
        else:
            vorder = self.vertices

        order_map = zeros(B, dtype="int")
        for i, v in enumerate(vorder):
            order_map[vmap[v]] = i

        m = zeros((B, B))
        rmap = {}
        for e in self.bg.edges():
            r, s = vmap[int(e.source())], vmap[int(e.target())]
            r = order_map[r]
            s = order_map[s]
            rmap[r] = int(e.source())
            rmap[s] = int(e.target())
            m[r, s] = self.mrs[e]
            if not self.bg.is_directed():
                m[s, r] = m[r, s]

553
554
555
556
        for r in range(B):
            if r not in rmap:
                rmap[r] = r

557
558
559
560
561
562
        if ret_order:
            return m, rmap
        else:
            return m


563
def model_entropy(B, N, E, directed=False, nr=None):
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
    r"""Computes the amount of information necessary for the parameters of the traditional blockmodel ensemble, for ``B`` blocks, ``N`` vertices, ``E`` edges, and either a directed or undirected graph.

    A traditional blockmodel is defined as a set of :math:`N` vertices which can
    belong to one of :math:`B` blocks, and the matrix :math:`e_{rs}` describes
    the number of edges from block :math:`r` to :math:`s` (or twice that number
    if :math:`r=s` and the graph is undirected).

    For an undirected graph, the number of distinct :math:`e_{rs}` matrices is given by,

    .. math::

       \Omega_m = \left(\!\!{\left(\!{B \choose 2}\!\right) \choose E}\!\!\right)

    and for a directed graph,

    .. math::
       \Omega_m = \left(\!\!{B^2 \choose E}\!\!\right)


    where :math:`\left(\!{n \choose k}\!\right) = {n+k-1\choose k}` is the
    number of :math:`k` combinations with repetitions from a set of size :math:`n`.

    The total information necessary to describe the model is then,

    .. math::

590
591
       \mathcal{L}_t = \ln\Omega_m + \ln\left(\!\!{B \choose N}\!\!\right) + \ln N! - \sum_r \ln n_r!,

592

593
594
    where the remaining term is the information necessary to describe the
    block partition, where :math:`n_r` is the number of nodes in block :math:`r`.
595

596
597
598
599
    If ``nr`` is ``None``, it is assumed :math:`n_r=N/B`.

    The value returned corresponds to the information per edge, i.e.
    :math:`\mathcal{L}_t/E`.
600
601
602
603

    References
    ----------

Tiago Peixoto's avatar
Tiago Peixoto committed
604
605
606
    .. [peixoto-parsimonious-2013] Tiago P. Peixoto, "Parsimonious module inference in large networks",
       Phys. Rev. Lett. 110, 148701 (2013), :doi:`10.1103/PhysRevLett.110.148701`, :arxiv:`1212.4794`.
    .. [peixoto-hierarchical-2013] Tiago P. Peixoto, "Hierarchical block structures and high-resolution
607
       model selection in large networks ", :arxiv:`1310.4377`.
608
609
610

    """

611
612
613
614
615
616
    if directed:
        x = (B * B);
    else:
        x = (B * (B + 1)) / 2;
    L = lbinom(x + E - 1, E) + partition_entropy(B, N, nr)
    return L / E
617
618
619
620

def Sdl(B, S, N, E, directed=False):
    return S + model_entropy(B, N, E, directed)

621
622
623
624
625
626
627
def lbinom(n, k):
    return scipy.special.gammaln(n + 1) - scipy.special.gammaln(n - k + 1) - scipy.special.gammaln(k + 1)

def partition_entropy(B, N, nr=None):
    if nr is None:
        S = N * log(B) + log1p(-(1 - 1./B) ** N)
    else:
628
        S = lbinom(B + N - 1, N) + scipy.special.gammaln(N + 1) - scipy.special.gammaln(nr + 1).sum()
629
    return S
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649

def get_max_B(N, E, directed=False):
    r"""Return the maximum detectable number of blocks, obtained by minimizing:

    .. math::

        \mathcal{L}_t(B, N, E) - E\ln B

    where :math:`\mathcal{L}_t(B, N, E)` is the information necessary to
    describe a traditional blockmodel with `B` blocks, `N` nodes and `E`
    edges (see :func:`model_entropy`).

    Examples
    --------

    >>> gt.get_max_B(N=1e6, E=5e6)
    1572

    References
    ----------
Tiago Peixoto's avatar
Tiago Peixoto committed
650
651
    .. [peixoto-parsimonious-2013] Tiago P. Peixoto, "Parsimonious module inference in large networks",
       Phys. Rev. Lett. 110, 148701 (2013), :doi:`10.1103/PhysRevLett.110.148701`, :arxiv:`1212.4794`.
652
653
654
655
656
657
658
659
660
661
662


    """

    B = fminbound(lambda B: Sdl(B, -log(B), N, E, directed), 1, E,
                  xtol=1e-6, maxfun=1500, disp=0)
    if isnan(B):
        B = 1
    return max(int(ceil(B)), 2)

def get_akc(B, I, N=float("inf"), directed=False):
Tiago Peixoto's avatar
Tiago Peixoto committed
663
    r"""Return the minimum value of the average degree of the network, so that some block structure with :math:`B` blocks can be detected, according to the minimum description length criterion.
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696

    This is obtained by solving

    .. math::

       \Sigma_b = \mathcal{L}_t(B, N, E) - E\mathcal{I}_{t/c} = 0,

    where :math:`\mathcal{L}_{t}` is the necessary information to describe the
    blockmodel parameters (see :func:`model_entropy`), and
    :math:`\mathcal{I}_{t/c}` characterizes the planted block structure, and is
    given by

    .. math::

        \mathcal{I}_t &= \sum_{rs}m_{rs}\ln\left(\frac{m_{rs}}{w_rw_s}\right),\\
        \mathcal{I}_c &= \sum_{rs}m_{rs}\ln\left(\frac{m_{rs}}{m_rm_s}\right),

    where :math:`m_{rs} = e_{rs}/2E` (or :math:`m_{rs} = e_{rs}/E` for directed
    graphs) and :math:`w_r=n_r/N`. We note that :math:`\mathcal{I}_{t/c}\in[0,
    \ln B]`. In the case where :math:`E \gg B^2`, this simplifies to

    .. math::

       \left<k\right>_c &= \frac{2\ln B}{\mathcal{I}_{t/c}},\\
       \left<k^{-/+}\right>_c &= \frac{\ln B}{\mathcal{I}_{t/c}},

    for undirected and directed graphs, respectively. This limit is assumed if
    ``N == inf``.

    Examples
    --------

    >>> gt.get_akc(10, log(10) / 100, N=100)
697
    2.414413200430159
698
699
700

    References
    ----------
Tiago Peixoto's avatar
Tiago Peixoto committed
701
702
    .. [peixoto-parsimonious-2013] Tiago P. Peixoto, "Parsimonious module inference in large networks",
       Phys. Rev. Lett. 110, 148701 (2013), :doi:`10.1103/PhysRevLett.110.148701`, :arxiv:`1212.4794`.
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717

    """
    if N != float("inf"):
        if directed:
            get_dl = lambda ak: model_entropy(B, N, N * ak, directed) - N * ak * I
        else:
            get_dl = lambda ak: model_entropy(B, N, N * ak / 2., directed) - N * ak * I / 2.
        ak = fsolve(lambda ak: get_dl(ak), 10)
        ak = float(ak)
    else:
        ak = 2 * log(B) / S
        if directed:
            ak /= 2
    return ak

718
719
720
721
def mcmc_sweep(state, beta=1., random_move=False, c=1., dense=False,
               multigraph=False, sequential=True, vertices=None,
               verbose=False):
    r"""Performs a Markov chain Monte Carlo sweep on the network, to sample the block partition according to a probability :math:`\propto e^{-\beta \mathcal{S}_{t/c}}`, where :math:`\mathcal{S}_{t/c}` is the blockmodel entropy.
722
723
724
725
726

    Parameters
    ----------
    state : :class:`~graph_tool.community.BlockState`
        The block state.
727
    beta : ``float`` (optional, default: `1.0`)
728
        The inverse temperature parameter :math:`\beta`.
729
730
731
732
733
    random_move : ``bool`` (optional, default: ``False``)
        If ``True``, the proposed moves will attempt to place the vertices in
        fully randomly-chosen blocks. If ``False``, the proposed moves will be
        chosen with a probability depending on the membership of the neighbours
        and the currently-inferred block structure.
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
    c : ``float`` (optional, default: ``1.0``)
        This parameter specifies how often fully random moves are attempted,
        instead of more likely moves based on the inferred block partition.
        For ``c == 0``, no fully random moves are attempted, and for ``c == inf``
        they are always attempted.
    dense : ``bool`` (optional, default: ``False``)
        If ``True``, the "dense" variant of the entropy will be computed.
    multigraph : ``bool`` (optional, default: ``False``)
        If ``True``, the multigraph entropy will be used. Only has an effect
        if ``dense == True``.
    sequential : ``bool`` (optional, default: ``True``)
        If ``True``, the move attempts on the vertices are done in sequential
        random order. Otherwise a total of `N` moves attempts are made, where
        `N` is the number of vertices, where each vertex can be selected with
        equal probability.
Tiago Peixoto's avatar
Tiago Peixoto committed
749
    vertices : ``list of ints`` (optional, default: ``None``)
750
751
        A list of vertices which will be attempted to be moved. If ``None``, all
        vertices will be attempted.
752
753
754
755
756
757
    verbose : ``bool`` (optional, default: ``False``)
        If ``True``, verbose information is displayed.

    Returns
    -------

758
    dS : ``float``
759
760
761
762
763
764
765
766
       The entropy difference (per edge) after a full sweep.
    nmoves : ``int``
       The number of accepted block membership moves.


    Notes
    -----

767
    This algorithm performs a Markov chain Monte Carlo sweep on the network,
768
769
    where the block memberships are randomly moved, and either accepted or
    rejected, so that after sufficiently many sweeps the partitions are sampled
770
    with probability proportional to :math:`e^{-\beta\mathcal{S}_{t/c}}`, where
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
    :math:`\mathcal{S}_{t/c}` is the blockmodel entropy, given by

    .. math::

      \mathcal{S}_t &\cong - \frac{1}{2} \sum_{rs}e_{rs}\ln\left(\frac{e_{rs}}{n_rn_s}\right), \\
      \mathcal{S}^d_t &\cong - \sum_{rs}e_{rs}\ln\left(\frac{e_{rs}}{n_rn_s}\right),

    for undirected and directed traditional blockmodels (``deg_corr == False``),
    respectively, where :math:`e_{rs}` is the number of edges from block
    :math:`r` to :math:`s` (or the number of half-edges for the undirected case
    when :math:`r=s`), and :math:`n_r` is the number of vertices in block
    :math:`r`, and constant terms which are independent of the block partition
    were dropped (see :meth:`BlockState.entropy` for the complete entropy). For
    the degree-corrected variant with "hard" degree constraints the equivalent
    expressions are

    .. math::

       \mathcal{S}_c &\cong  - \frac{1}{2} \sum_{rs}e_{rs}\ln\left(\frac{e_{rs}}{e_re_s}\right), \\
       \mathcal{S}^d_c &\cong - \sum_{rs}e_{rs}\ln\left(\frac{e_{rs}}{e^+_re^-_s}\right),

    where :math:`e_r = \sum_se_{rs}` is the number of half-edges incident on
    block :math:`r`, and :math:`e^+_r = \sum_se_{rs}` and :math:`e^-_r =
    \sum_se_{sr}` are the number of out- and in-edges adjacent to block
    :math:`r`, respectively.

    The Monte Carlo algorithm employed attempts to improve the mixing time of
798
799
    the Markov chain by proposing membership moves :math:`r\to s` with
    probability :math:`p(r\to s|t) \propto e_{ts} + c`, where :math:`t` is the
800
    block label of a random neighbour of the vertex being moved. See
Tiago Peixoto's avatar
Tiago Peixoto committed
801
    [peixoto-efficient-2013]_ for more details.
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842

    This algorithm has a complexity of :math:`O(E)`, where :math:`E` is the
    number of edges in the network.

    Examples
    --------
    .. testsetup:: mcmc

       gt.seed_rng(42)
       np.random.seed(42)

    .. doctest:: mcmc

       >>> g = gt.collection.data["polbooks"]
       >>> state = gt.BlockState(g, B=3, deg_corr=True)
       >>> pv = None
       >>> for i in range(1000):        # remove part of the transient
       ...     ds, nmoves = gt.mcmc_sweep(state)
       >>> for i in range(1000):
       ...     ds, nmoves = gt.mcmc_sweep(state)
       ...     pv = gt.collect_vertex_marginals(state, pv)
       >>> gt.graph_draw(g, pos=g.vp["pos"], vertex_shape="pie", vertex_pie_fractions=pv, output="polbooks_blocks_soft.pdf")
       <...>

    .. testcleanup:: mcmc

       gt.graph_draw(g, pos=g.vp["pos"], vertex_shape="pie", vertex_pie_fractions=pv, output="polbooks_blocks_soft.png")

    .. figure:: polbooks_blocks_soft.*
       :align: center

       "Soft" block partition of a political books network with :math:`B=3`.

     References
    ----------

    .. [holland-stochastic-1983] Paul W. Holland, Kathryn Blackmond Laskey,
       Samuel Leinhardt, "Stochastic blockmodels: First steps",
       Carnegie-Mellon University, Pittsburgh, PA 15213, U.S.A., :doi:`10.1016/0378-8733(83)90021-7`
    .. [faust-blockmodels-1992] Katherine Faust, and Stanley
       Wasserman. "Blockmodels: Interpretation and Evaluation." Social Networks
843
       14, no. 1-2 (1992): 5-61. :doi:`10.1016/0378-8733(92)90013-W`
844
845
846
847
848
849
    .. [karrer-stochastic-2011] Brian Karrer, and M. E. J. Newman. "Stochastic
       Blockmodels and Community Structure in Networks." Physical Review E 83,
       no. 1 (2011): 016107. :doi:`10.1103/PhysRevE.83.016107`.
    .. [peixoto-entropy-2012] Tiago P. Peixoto "Entropy of Stochastic Blockmodel
       Ensembles." Physical Review E 85, no. 5 (2012): 056122. :doi:`10.1103/PhysRevE.85.056122`,
       :arxiv:`1112.6028`.
Tiago Peixoto's avatar
Tiago Peixoto committed
850
851
852
    .. [peixoto-parsimonious-2013] Tiago P. Peixoto, "Parsimonious module inference in large networks",
       Phys. Rev. Lett. 110, 148701 (2013), :doi:`10.1103/PhysRevLett.110.148701`, :arxiv:`1212.4794`.
    .. [peixoto-efficient-2013] Tiago P. Peixoto, "Efficient Monte Carlo and greedy
853
       heuristic for the inference of stochastic block models", :arxiv:`1310.4378`.
854
855
    """

856
    if state.B == 1:
857
858
        return 0., 0

859
    if vertices is not None:
860
861
862
        vlist = libcommunity.get_vector(len(vertices))
        vlist.a = vertices
        vertices = vlist
863
        state.sweep_vertices = vertices
864

865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
    if state.sweep_vertices is None:
        vertices = libcommunity.get_vector(state.g.num_vertices())
        vertices.a = state.g.vertex_index.copy("int").fa
        state.sweep_vertices = vertices

    if random_move:
        state._BlockState__build_egroups(empty=True)
    elif state.egroups is None:
        state._BlockState__build_egroups(empty=False)

    if state.nsampler is None:
        state._BlockState__build_nsampler()

    state.bnsampler = None

    try:
        dS, nmoves = libcommunity.move_sweep(state.g._Graph__graph,
                                             state.bg._Graph__graph,
                                             state._BlockState__get_emat(),
                                             state.nsampler,
                                             _prop("e", state.bg, state.mrs),
                                             _prop("v", state.bg, state.mrp),
                                             _prop("v", state.bg, state.mrm),
                                             _prop("v", state.bg, state.wr),
                                             _prop("v", state.g, state.b),
                                             _prop("v", state.bg, state.bclabel),
                                             state.sweep_vertices,
                                             state.deg_corr, dense, multigraph,
                                             _prop("e", state.g, state.eweight),
                                             _prop("v", state.g, state.vweight),
                                             state.egroups,
                                             _prop("e", state.g, state.esrcpos),
                                             _prop("e", state.g, state.etgtpos),
                                             float(beta), sequential, random_move,
899
900
                                             c, state.is_weighted, verbose,
                                             _get_rng())
901
902
903
    finally:
        if random_move:
            state.egroups = None
904
905
906
    return dS / state.E, nmoves


907
908
def merge_sweep(state, bm, nmerges, nsweeps=10, dense=False, multigraph=False,
                random_moves=False, verbose=False):
909

910
911
    if state.B == 1:
        return 0., 0
912

913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
    if state.bnsampler is None:
        state._BlockState__build_bnsampler()

    state.egroups = None
    state.nsampler = None

    dS, nmoves = libcommunity.merge_sweep(state.bg._Graph__graph,
                                          state._BlockState__get_emat(),
                                          state.bnsampler,
                                          _prop("e", state.bg, state.mrs),
                                          _prop("v", state.bg, state.mrp),
                                          _prop("v", state.bg, state.mrm),
                                          _prop("v", state.bg, state.wr),
                                          _prop("v", state.bg, bm),
                                          _prop("v", state.bg, state.bclabel),
                                          state.deg_corr, dense, multigraph,
                                          nsweeps, nmerges, random_moves,
                                          verbose, _get_rng())
931

932
    return dS / state.E, nmoves
933

934

935
936
937
938
939
940
def greedy_shrink(state, B, nsweeps=10, adaptive_sweeps=True, nmerge_sweeps=None,
                  epsilon=0, r=2, greedy=True, anneal=(1, 1), c=1, dense=False,
                  multigraph=False, random_move=False, verbose=False,
                  sequential=True):
    if B > state.B:
        raise ValueError("Cannot shrink to a larger size!")
941

942
943
    if nmerge_sweeps is None:
        nmerge_sweeps = max((2 * state.g.num_edges()) // state.g.num_vertices(), 1)
944

945
    nmerged = 0
946

947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
    state = BlockState(state.g, b=state.b, B=state.B,
                       clabel=state.clabel, vweight=state.vweight,
                       eweight=state.eweight, deg_corr=state.deg_corr,
                       max_BE=state.max_BE)

    cg = state.bg.copy()
    cg_vweight = cg.own_property(state.wr.copy())
    cg_eweight = cg.own_property(state.mrs.copy())
    cg_clabel = cg.own_property(state.bclabel.copy())

    # merge according to indirect neighbourhood
    bm = state.bg.vertex_index.copy("int")
    random = random_move
    while nmerged < state.B - B:
        dS, nmoves = merge_sweep(state, bm, nmerges=state.B - B - nmerged,
                                 nsweeps=nmerge_sweeps, random_moves=random)
        nmerged += nmoves
964
        if verbose:
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
            print("merging", dS, nmoves, nmerged)
        if nmoves == 0:
            random = True
            if verbose:
                print("can't merge... switching to random")

    # Merged block-level state
    bmap = -ones(len(bm.a), dtype=bm.a.dtype)
    libcommunity.vector_map(bm.a, bmap)

    bm = cg.own_property(bm)
    bg_state = BlockState(cg, b=bm, B=B, clabel=cg_clabel,
                          vweight=cg_vweight, eweight=cg_eweight,
                          deg_corr=state.deg_corr, max_BE=state.max_BE)

    if bg_state.g.num_vertices() != state.g.num_vertices() and nsweeps > 0:
        # Perform block-level moves
        if verbose:
            print("Performing block-level moves...")
        multilevel_minimize(bg_state, B=B, nsweeps=nsweeps,
                            adaptive_sweeps=adaptive_sweeps,
                            epsilon=epsilon, r=r, greedy=greedy,
                            anneal=anneal, c=c, dense=dense,
                            multigraph=multigraph, random_move=random_move,
                            sequential=sequential, verbose=verbose)

    bm = bg_state.b
    libcommunity.vector_map(state.b.a, bm.a)

    state = BlockState(state.g, b=state.b, B=B, clabel=state.clabel,
                       vweight=state.vweight, eweight=state.eweight,
                       deg_corr=state.deg_corr, max_BE=state.max_BE)
    return state


class MinimizeState(object):
    r"""This object stores information regarding the current entropy minimization
    state, so that the algorithms can resume previously started runs.
    This object can be saved to disk via the  :mod:`pickle` interface."""

    def __init__(self):
        self.b_cache = {}
        self.checkpoint_state = defaultdict(dict)

    def clear(self):
        self.b_cache.clear()
        self.checkpoint_state.clear()


def multilevel_minimize(state, B, nsweeps=10, adaptive_sweeps=True, epsilon=0,
                        anneal=(1., 1.), r=2., nmerge_sweeps=10, greedy=True,
                        random_move=False, c=1., dense=False, multigraph=False,
                        sequential=True, checkpoint=None,
                        minimize_state=None, verbose=False):
Tiago Peixoto's avatar
Tiago Peixoto committed
1019
    r"""Performs an agglomerative heuristic, which progressively merges blocks together (while allowing individual node moves) to achieve a good partition in ``B`` blocks.
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115

    Parameters
    ----------
    state : :class:`~graph_tool.community.BlockState`
        The block state.
    B : ``int``
        The desired number of blocks.
    nsweeps : ``int`` (optional, default: ``10``)
        The number of sweeps done after each merge step to reach the local
        minimum.
    adaptive_sweeps : ``bool`` (optional, default: ``True``)
        If ``True``, the number of sweeps necessary for the local minimum will
        be estimated to be enough so that no more than ``epsilon * N`` nodes
        changes their states in the last ``nsweeps`` sweeps.
    epsilon : ``float`` (optional, default: ``0``)
        Converge criterion for ``adaptive_sweeps``.
    anneal : pair of ``floats`` (optional, default: ``(1., 1.)``)
        The first value specifies the starting value for  ``beta`` of the MCMC
        steps, and the second value is the factor which is multiplied to ``beta``
        after each estimated equilibration (according to ``nsweeps`` and
        ``adaptive_sweeps``).
    r : ``float`` (optional, default: ``2.``)
        Agglomeration ratio for the merging steps. Each merge step will attempt
        to find the best partition into :math:`B_{i-1} / r` blocks, where
        :math:`B_{i-1}` is the number of blocks in the last step.
    nmerge_sweeps : `int` (optional, default: `10`)
        The number of merge sweeps done, where in each sweep a better merge
        candidate is searched for every block.
    greedy : ``bool`` (optional, default: ``True``)
        If ``True``, the value of ``beta`` of the MCMC steps are kept at
        infinity for all steps. Otherwise they change according to the ``anneal``
        parameter.
    random_move : ``bool`` (optional, default: ``False``)
        If ``True``, the proposed moves will attempt to place the vertices in
        fully randomly-chosen blocks. If ``False``, the proposed moves will be
        chosen with a probability depending on the membership of the neighbours
        and the currently-inferred block structure.
    c : ``float`` (optional, default: ``1.0``)
        This parameter specifies how often fully random moves are attempted,
        instead of more likely moves based on the inferred block partition.
        For ``c == 0``, no fully random moves are attempted, and for ``c == inf``
        they are always attempted.
    dense : ``bool`` (optional, default: ``False``)
        If ``True``, the "dense" variant of the entropy will be computed.
    multigraph : ``bool`` (optional, default: ``False``)
        If ``True``, the multigraph entropy will be used. Only has an effect
        if ``dense == True``.
    sequential : ``bool`` (optional, default: ``True``)
        If ``True``, the move attempts on the vertices are done in sequential
        random order. Otherwise a total of `N` moves attempts are made, where
        `N` is the number of vertices, where each vertex can be selected with
        equal probability.
    vertices: ``list of ints`` (optional, default: ``None``)
        A list of vertices which will be attempted to be moved. If ``None``, all
        vertices will be attempted.
    checkpoint : function (optional, default: ``None``)
        If provided, this function will be called after each call to
        :func:`mcmc_sweep`. This can be used to store the current state, so it
        can be continued later. The function must have the following signature:

        .. code-block:: python

            def checkpoint(state, S, delta, nmoves, minimize_state):
                ...

        where `state` is either a :class:`~graph_tool.community.BlockState`
        instance or ``None``, `S` is the current entropy value, `delta` is
        the entropy difference in the last MCMC sweep, and `nmoves` is the
        number of accepted block membership moves. The ``minimize_state``
        argument is a :class:`MinimizeState` instance which specifies the current
        state of the algorithm, which can be stored via :mod:`pickle`, and
        supplied via the ``minimize_state`` option below to continue from an
        interrupted run.

        This function will also be called when the MCMC has finished for the
        current value of :math:`B`, in which case ``state == None``, and the
        remaining parameters will be zero, except the last.
    minimize_state : :class:`MinimizeState` (optional, default: ``None``)
        If provided, this will specify an exact point of execution from which
        the algorithm will continue. The expected object is a :class:`MinimizeState`
        instance which will be passed to the callback of the ``checkpoint``
        option above, and  can be stored by :mod:`pickle`.
    verbose : ``bool`` (optional, default: ``False``)
        If ``True``, verbose information is displayed.

    Returns
    -------

    state : :class:`~graph_tool.community.BlockState`
        The new :class:`~graph_tool.community.BlockState` with ``B`` blocks.

    Notes
    -----

    This algorithm performs an agglomerative heuristic on the current block state,
    where blocks are progressively merged together, using repeated applications of
Tiago Peixoto's avatar
Tiago Peixoto committed
1116
    the :func:`mcmc_sweep` moves, at different scales. See [peixoto-efficient-2013]_
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
    for more details.

    This algorithm has a complexity of :math:`O(N\ln^2 N)`, where :math:`N` is the
    number of nodes in the network.

    Examples
    --------
    .. testsetup:: multilevel_minimize

       gt.seed_rng(42)
       np.random.seed(42)

    .. doctest:: multilevel_minimize

       >>> g = gt.collection.data["polblogs"]
       >>> g = gt.GraphView(g, vfilt=gt.label_largest_component(gt.GraphView(g, directed=False)))
       >>> state = gt.BlockState(g, B=g.num_vertices(), deg_corr=True)
       >>> state = gt.multilevel_minimize(state, B=2)
       >>> gt.graph_draw(g, pos=g.vp["pos"], vertex_fill_color=state.get_blocks(), output="polblogs_agg.pdf")
       <...>

    .. testcleanup:: multilevel_minimize

       gt.graph_draw(g, pos=g.vp["pos"], vertex_fill_color=state.get_blocks(), output="polblogs_agg.png")

    .. figure:: polblogs_agg.*
       :align: center

       Block partition of a political blogs network with :math:`B=2`.

     References
    ----------

Tiago Peixoto's avatar
Tiago Peixoto committed
1150
    .. [peixoto-efficient-2013] Tiago P. Peixoto, "Efficient Monte Carlo and greedy
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
       heuristic for the inference of stochastic block models", :arxiv:`1310.4378`.
    """

    if minimize_state is None:
        minimize_state = MinimizeState()
    b_cache = minimize_state.b_cache
    checkpoint_state = minimize_state.checkpoint_state

    # some trivial boundary conditions
    if B == 1:
        bi = state.g.new_vertex_property("int")
        state = BlockState(state.g, vweight=state.vweight, eweight=state.eweight,
                           b=bi, clabel=state.clabel, deg_corr=state.deg_corr,
                           max_BE=state.max_BE)
        return state
    if B == state.g.num_vertices():
        bi = state.g.new_vertex_property("int")
        bi.fa = range(state.g.num_vertices())
        state = BlockState(state.g, vweight=state.vweight, eweight=state.eweight,
                           B=state.g.num_vertices(), b=bi,
                           clabel=state.clabel, deg_corr=state.deg_corr,
                           max_BE=state.max_BE)
        return state

    Bi = state.B
    while True:
        Bi = max(int(round(Bi / r)), B)
        if Bi == state.B and Bi > B:
            Bi -= 1

        if b_cache is not None and Bi in b_cache:
            bi = state.g.new_vertex_property("int")
            bi.fa = b_cache[Bi][1]
            state = BlockState(state.g, B=Bi, b=bi,
                               vweight=state.vweight, eweight=state.eweight,
                               clabel=state.clabel, deg_corr=state.deg_corr,
                               max_BE=state.max_BE)

        if Bi < state.B:
            if verbose:
                print("Shrinking:", state.B, "->", Bi)
            state = greedy_shrink(state, B=Bi, nsweeps=nsweeps, epsilon=epsilon, c=c,
                                  dense=dense, multigraph=multigraph,
                                  nmerge_sweeps=nmerge_sweeps, sequential=sequential,
                                  verbose=verbose)

        if "S" in checkpoint_state[Bi]:
            S = checkpoint_state[Bi]["S"]
            niter = checkpoint_state[Bi]["niter"]
        else:
            S = state.entropy(dl=True)
            checkpoint_state[Bi]["S"] = S
            niter = 0
            checkpoint_state[Bi]["niter"] = niter

        if b_cache is not None and Bi not in b_cache:
            b_cache[Bi] = [float("inf"), array(state.b.fa), None]

        if not adaptive_sweeps:
            ntotal = nsweeps if greedy else 2 * nsweeps
            if verbose:
                print("Performing %d sweeps for B=%d..." % (ntotal, Bi))

            for i in range(ntotal):
                if i < niter:
                    continue
                if i < nsweeps and not greedy:
                    beta = anneal[0]
1219
                else:
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
                    beta = float("inf")
                delta, nmoves = mcmc_sweep(state, beta=beta, c=c,
                                           dense=dense, multigraph=multigraph,
                                           sequential=sequential,
                                           random_move=random_move)
                S += delta
                niter += 1
                checkpoint_state[Bi]["S"] = S
                checkpoint_state[Bi]["niter"] = niter
                if b_cache is not None:
                    b_cache[Bi][1] = array(state.b.fa)
                if checkpoint is not None:
                    checkpoint(state, S, delta, nmoves, minimize_state)
        else:
            # adaptive mode
            min_dl = checkpoint_state[Bi].get("min_dl", S)
            max_dl = checkpoint_state[Bi].get("max_dl", S)
            count = checkpoint_state[Bi].get("count", 0)
            bump = checkpoint_state[Bi].get("bump", False)
            beta =  checkpoint_state[Bi].get("beta", anneal[0])
            last_min = checkpoint_state[Bi].get("last_min", min_dl)
            greedy_step = checkpoint_state[Bi].get("greedy_step", greedy)
            total_nmoves = checkpoint_state[Bi].get("total_nmoves", 0)

            if verbose and not greedy:
                print("Performing sweeps for beta = %g, B=%d (N=%d)..." % \
                       (beta, Bi, state.g.num_vertices()))

            eps = 1e-8
            niter = 0
            while True:
                if greedy_step:
                    break
                if count > nsweeps:
                    if not bump:
                        min_dl = max_dl = S
                        bump = True
1257
                        count = 0
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
                    else:
                        if anneal[1] <= 1 or min_dl == last_min:
                            break
                        else:
                            beta *= anneal[1]
                            count = 0
                            last_min = min_dl
                            if verbose:
                                print("Performing sweeps for beta = %g, B=%d (N=%d)..." % \
                                       (beta, Bi, state.g.num_vertices()))

                delta, nmoves = mcmc_sweep(state, beta=beta, c=c,
                                           dense=dense, multigraph=multigraph,
                                           sequential=sequential,
                                           random_move=random_move)
                S += delta
                niter += 1
                total_nmoves += nmoves

                if S > max_dl + eps:
                    max_dl = S
                    count = 0
                elif S < min_dl - eps:
                    min_dl = S
                    count = 0
                else:
                    count += 1

                checkpoint_state[B]["S"] = S
                checkpoint_state[B]["niter"] = niter
                checkpoint_state[B]["min_dl"] = min_dl
                checkpoint_state[B]["max_dl"] = max_dl
                checkpoint_state[B]["count"] = count
                checkpoint_state[B]["bump"] = bump
                checkpoint_state[B]["total_nmoves"] = total_nmoves

                if b_cache is not None:
                    b_cache[Bi][0] = float("inf")
                    b_cache[Bi][1] = array(state.b.fa)
                if checkpoint is not None:
                    checkpoint(state, S, delta, nmoves, minimize_state)
1299

1300
1301
1302
1303
1304
            if verbose:
                if not greedy_step:
                    print("... performed %d sweeps with %d vertex moves" % (niter, total_nmoves))
                print(u"Performing sweeps for beta = ∞, B=%d (N=%d)..." % \
                       (Bi, state.g.num_vertices()))
1305

1306
1307
            if not greedy_step:
                checkpoint_state[Bi]["greedy_step"] = True
1308
1309
1310
                min_dl = S
                count = 0

1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
            niter = 0
            total_nmoves = 0
            while count <= nsweeps:
                delta, nmoves = mcmc_sweep(state, beta=float("inf"), c=c,
                                           dense=dense, multigraph=multigraph,
                                           sequential=sequential,
                                           random_move=random_move)
                S += delta
                niter += 1
                total_nmoves += nmoves

                # if verbose:
                #     print("Moved:", delta, nmoves,
                #           nmoves / state.g.num_vertices(),
                #           epsilon, count)

                #if nmoves > epsilon * state.g.num_vertices():
                if abs(delta) > eps and nmoves / state.g.num_vertices() > epsilon:
                    min_dl = S
                    count = 0
                else:
                    count += 1
                checkpoint_state[Bi]["S"] = S
                checkpoint_state[Bi]["min_dl"] = min_dl
                checkpoint_state[Bi]["count"] = count
                checkpoint_state[B]["total_nmoves"] = total_nmoves
                if b_cache is not None:
                    b_cache[Bi][1] = array(state.b.fa)
                if checkpoint is not None:
                    checkpoint(state, S, delta, nmoves, minimize_state)
1341

1342
1343
            if verbose:
                print("... performed %d sweeps with %d vertex moves" % (niter, total_nmoves))
1344

1345
1346
1347
            bi = state.b
            if Bi == B:
                break
1348

1349
    return state
1350

1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
def get_state_dl(state, dense, nested_dl, clabel=None):
    if not nested_dl:
        dl = state.entropy(dense=dense, multigraph=dense, dl=True)
    else:
        dl = state.entropy(dense=dense, multigraph=dense, dl=False) + \
             partition_entropy(B=state.B, N=state.N, nr=state.wr.a) / state.E
        if clabel is None:
            bclabel = state.bclabel
        else:
            bclabel = state.bg.new_vertex_property("int")
            libcommunity.vector_rmap(state.b.a, bclabel.a)
            libcommunity.vector_map(bclabel.a, clabel.a)

        bstate = BlockState(state.bg, b=bclabel, eweight=state.mrs,
                            deg_corr=False)
        dl += bstate.entropy(dl=False, dense=True, multigraph=True) + \
              partition_entropy(B=bstate.B, N=bstate.N, nr=bstate.wr.a) / state.E
    return dl


def get_b_dl(g, vweight, eweight, B, nsweeps, adaptive_sweeps, c, random_move,
             sequential, shrink, r, anneal, greedy, epsilon, nmerge_sweeps, clabel,
             deg_corr, dense, sparse_heuristic, checkpoint, minimize_state,
             max_BE, nested_dl,  verbose):
    bs = minimize_state.b_cache
    checkpoint_state = minimize_state.checkpoint_state
    previous = None
1378
    if B in bs and checkpoint_state[B].get("done", False):
1379
1380
1381
        # A previous finished result is available. Use that and keep going.
        if verbose:
            print("(using previous finished result for B=%d)" % B)
1382
        return bs[B][0]
1383
    elif B in bs:
1384
        # A previous unfinished result is available. Use that as the starting point.
1385
        if verbose:
1386
1387
1388
1389
1390
1391
            print("(starting from previous result for B=%d)" % B)
        b = g.new_vertex_property("int")
        b.fa = bs[B][1]
        state = BlockState(g, b=b, B=B, vweight=vweight, eweight=eweight,
                           clabel=clabel, deg_corr=deg_corr, max_BE=max_BE)
        previous = bs[B]
1392
    else:
1393
        # No previous result is available.
1394
        bs_keys = [k for k in bs.keys() if type(k) != str]
1395
        B_sup = max(max(bs_keys), B) if len(bs_keys) > 0 else B
1396
1397
1398
        for Bi in bs_keys:
            if Bi > B and Bi < B_sup:
                B_sup = Bi
1399
1400
1401
1402
1403
1404
1405
1406
        if B_sup == B or not shrink:
            # Start from scratch.
            bi = g.new_vertex_property("int")
            bi.fa = range(g.num_vertices())
            state = BlockState(g, B=g.num_vertices(), b=bi,
                               vweight=vweight,
                               eweight=eweight, clabel=clabel,
                               deg_corr=deg_corr, max_BE=max_BE)
1407
        else:
1408
            # Start from result with B_sup > B, and shrink it.
1409
            if verbose:
1410
1411
1412
                print("(shrinking from B=%d to B=%d)" % (B_sup, B))
            b = g.new_vertex_property("int")
            b.fa = bs[B_sup][1]
1413

1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
            if B > 1:
                state = BlockState(g, B=B_sup, b=b, vweight=vweight, eweight=eweight,
                                   clabel=clabel, deg_corr=deg_corr,
                                   max_BE=max_BE)
            else:
                bi = g.new_vertex_property("int")
                bi.fa = range(g.num_vertices())
                state = BlockState(g, B=g.num_vertices(), b=bi,
                                   vweight=vweight,
                                   eweight=eweight, clabel=clabel,
                                   deg_corr=deg_corr, max_BE=max_BE)

    # perform the actual minimization
    state = multilevel_minimize(state, B, nsweeps=nsweeps,
                                adaptive_sweeps=adaptive_sweeps,
                                epsilon=epsilon, r=r, greedy=greedy,
                                nmerge_sweeps=nmerge_sweeps, anneal=anneal,
                                c=c, random_move=random_move,
                                dense=dense and not sparse_heuristic,
                                multigraph=dense,
                                sequential=sequential,
                                minimize_state=minimize_state,
                                checkpoint=checkpoint,
                                verbose=verbose)
    dl = get_state_dl(state, dense, nested_dl)

    if previous is None or dl < previous[0]:
        # the current result improved the previous one
        bs[B] = [dl, array(state.b.fa)]
        if verbose:
            print("(using new result for B=%d with L=%g)" % (B, dl))
    else:
        # the previous result is better than the current one
        if verbose:
            print("(kept old result for B=%d with L=%g [vs L=%g])" % (B, previous[0], dl))
        dl = previous[0]
1450
    checkpoint_state[B]["done"] = True
1451
    assert(not isinf(dl))
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
    return dl

def fibo(n):
    phi = (1 + sqrt(5)) / 2
    return int(round(phi ** n / sqrt(5)))

def fibo_n_floor(x):
    phi = (1 + sqrt(5)) / 2
    n = floor(log(x * sqrt(5) + 0.5) / log(phi))
    return int(n)

def get_mid(a, b):
    n = fibo_n_floor(b - a)
    return b - fibo(n - 1)

def is_fibo(x):
    return fibo(fibo_n_floor(x)) == x

1470
1471
1472
1473
1474
1475
1476
1477
def minimize_blockmodel_dl(g, eweight=None, vweight=None, deg_corr=True, dense=False,
                           sparse_heuristic=False, random_move=False, c=0, nsweeps=100,
                           adaptive_sweeps=True, epsilon=0., anneal=(1., 1.),
                           greedy_cooling=True, sequential=True, r=2,
                           nmerge_sweeps=10, max_B=None, min_B=1, mid_B=None,
                           clabel=None, checkpoint=None, minimize_state=None,
                           exhaustive=False, max_BE=None, nested_dl=False,
                           verbose=False):
1478
1479
1480
1481
1482
1483
1484
1485
    r"""Find the block partition of an unspecified size which minimizes the description
    length of the network, according to the stochastic blockmodel ensemble which
    best describes it.

    Parameters
    ----------
    g : :class:`~graph_tool.Graph`
        Graph being used.
1486
1487
1488
1489
    eweight : :class:`~graph_tool.PropertyMap` (optional, default: ``None``)
        Edge weights (i.e. multiplicity).
    vweight : :class:`~graph_tool.PropertyMap` (optional, default: ``None``)
        Vertex weights (i.e. multiplicity).
1490
1491
1492
    deg_corr : ``bool`` (optional, default: ``True``)
        If ``True``, the degree-corrected version of the blockmodel ensemble will
        be assumed, otherwise the traditional variant will be used.
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
    dense : ``bool`` (optional, default: ``False``)
        If ``True``, the "dense" variant of the entropy will be computed.
    sparse_heuristic : ``bool`` (optional, default: ``False``)
        If ``True``, the sparse entropy will be used to find the best partition,
        but the dense entropy will be used to compare different partitions. This
        has an effect only if ``dense == True``.
    random_move : ``bool`` (optional, default: ``False``)
        If ``True``, the proposed moves will attempt to place the vertices in
        fully randomly-chosen blocks. If ``False``, the proposed moves will be
        chosen with a probability depending on the membership of the neighbours
        and the currently-inferred block structure.
    c : ``float`` (optional, default: ``1.0``)
        This parameter specifies how often fully random moves are attempted,
        instead of more likely moves based on the inferred block partition.
        For ``c == 0``, no fully random moves are attempted, and for ``c == inf``
        they are always attempted.
    nsweeps : ``int`` (optional, default: ``10``)
        The number of sweeps done after each merge step to reach the local
        minimum.
    adaptive_sweeps : ``bool`` (optional, default: ``True``)
        If ``True``, the number of sweeps necessary for the local minimum will
        be estimated to be enough so that no more than ``epsilon * N`` nodes
        changes their states in the last ``nsweeps`` sweeps.
    epsilon : ``float`` (optional, default: ``0``)
        Converge criterion for ``adaptive_sweeps``.
    anneal : pair of ``floats`` (optional, default: ``(1., 1.)``)
        The first value specifies the starting value for  ``beta`` of the MCMC
        steps, and the second value is the factor which is multiplied to ``beta``
        after each estimated equilibration (according to ``nsweeps`` and
        ``adaptive_sweeps``).
1523
    greedy_cooling : ``bool`` (optional, default: ``True``)
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
        If ``True``, the value of ``beta`` of the MCMC steps are kept at
        infinity for all steps. Otherwise they change according to the ``anneal``
        parameter.
    sequential : ``bool`` (optional, default: ``True``)
        If ``True``, the move attempts on the vertices are done in sequential
        random order. Otherwise a total of `N` moves attempts are made, where
        `N` is the number of vertices, where each vertex can be selected with
        equal probability.
    r : ``float`` (optional, default: ``2.``)
        Agglomeration ratio for the merging steps. Each merge step will attempt
        to find the best partition into :math:`B_{i-1} / r` blocks, where
        :math:`B_{i-1}` is the number of blocks in the last step.
    nmerge_sweeps : `int` (optional, default: `10`)
        The number of merge sweeps done, where in each sweep a better merge
        candidate is searched for every block.
1539
1540
1541
1542
1543
1544
1545
1546
    max_B : ``int`` (optional, default: ``None``)
        Maximum number of blocks tried. If not supplied, it will be
        automatically determined.
    min_B : ``int`` (optional, default: `1`)
        Minimum number of blocks tried.
    mid_B : ``int`` (optional, default: ``None``)
        Middle of the range which brackets the minimum. If not supplied, will be
        automatically determined.
1547
1548
1549
    clabel : :class:`~graph_tool.PropertyMap` (optional, default: ``None``)
        Constraint labels on the vertices, such that vertices with different
        labels cannot belong to the same block.
1550
1551
1552
1553
1554
1555
1556
    checkpoint : function (optional, default: ``None``)
        If provided, this function will be called after each call to
        :func:`mcmc_sweep`. This can be used to store the current state, so it
        can be continued later. The function must have the following signature:

        .. code-block:: python

1557
            def checkpoint(state, L, delta, nmoves, minimize_state):
1558
1559
1560
1561
1562
                ...

        where `state` is either a :class:`~graph_tool.community.BlockState`
        instance or ``None``, `L` is the current description length, `delta` is
        the entropy difference in the last MCMC sweep, and `nmoves` is the
1563
1564
1565
1566
1567
        number of accepted block membership moves. The ``minimize_state``
        argument is a :class:`~graph_tool.community.MinimizeState` instance
        which specifies the current state of the algorithm, which can be stored
        via :mod:`pickle`, and supplied via the ``minimize_state`` option below
        to continue from an interrupted run.
1568
1569
1570

        This function will also be called when the MCMC has finished for the
        current value of :math:`B`, in which case ``state == None``, and the
1571
        remaining parameters will be zero, except the last.
1572
    minimize_state : :class:`~graph_tool.community.MinimizeState` (optional, default: ``None``)
1573
        If provided, this will specify an exact point of execution from which
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
        the algorithm will continue. The expected object is a
        :class:`~graph_tool.community.MinimizeState`
        instance which will be passed to the callback of the ``checkpoint``
        option above, and  can be stored by :mod:`pickle`.
    exhaustive : ``bool`` (optional, default: ``False``)
        If ``True``, the best value of ``B`` will be found by testing all possible
        values, instead of performing a bisection search.
    max_BE : ``int`` (optional, default: ``1000``)
        If the number of blocks exceeds this number, a sparse representation of
        the block graph is used, which is slightly less efficient, but uses less
        memory,
1585
1586
1587
1588
1589
1590
1591
    verbose : ``bool`` (optional, default: ``False``)
        If ``True``, verbose information is displayed.

    Returns
    -------
    b : :class:`~graph_tool.PropertyMap`
       Vertex property map with the best block partition.
1592
1593
    min_dl : ``float``
       Minimum value of the description length (in `nats <http://en.wikipedia.org/wiki/Nat_%28information%29>`_ per edge).
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614

    Notes
    -----

    This algorithm attempts to find a block partition of an unspecified size
    which minimizes the description length of the network,

    .. math::

       \Sigma_{t/c} = \mathcal{S}_{t/c} + \mathcal{L}_{t/c},

    where :math:`\mathcal{S}_{t/c}` is the blockmodel entropy (as described in
    the docstring of :func:`mcmc_sweep` and :meth:`BlockState.entropy`) and
    :math:`\mathcal{L}_{t/c}` is the information necessary to describe the model
    (as described in the docstring of :func:`model_entropy` and
    :meth:`BlockState.entropy`).

    The algorithm works by minimizing the entropy :math:`\mathcal{S}_{t/c}` for
    specific values of :math:`B` via :func:`mcmc_sweep` (with :math:`\beta = 1`
    and :math:`\beta\to\infty`), and minimizing :math:`\Sigma_{t/c}` via an
    one-dimensional Fibonacci search on :math:`B`. See
Tiago Peixoto's avatar
Tiago Peixoto committed
1615
    [peixoto-parsimonious-2013]_ for more details.
1616

1617
1618
    This algorithm has a complexity of :math:`O(\tau N\ln^2 B_{\text{max}})`,
    where :math:`N` is the number of nodes in the network, :math:`\tau` is the
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
    mixing time of the MCMC, and :math:`B_{\text{max}}` is the maximum number of
    blocks considered. If :math:`B_{\text{max}}` is not supplied, it is computed
    as :math:`\sim\sqrt{E}` via :func:`get_max_B`, in which case the complexity
    becomes :math:`O(\tau E\ln E)`.


    Examples
    --------
    .. testsetup:: mdl

       gt.seed_rng(42)
       np.random.seed(42)

    .. doctest:: mdl

       >>> g = gt.collection.data["polbooks"]
1635
       >>> b, mdl = gt.minimize_blockmodel_dl(g)
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
       >>> gt.graph_draw(g, pos=g.vp["pos"], vertex_fill_color=b, vertex_shape=b, output="polbooks_blocks_mdl.pdf")
       <...>

    .. testcleanup:: mdl

       gt.graph_draw(g, pos=g.vp["pos"], vertex_fill_color=b, vertex_shape=b, output="polbooks_blocks_mdl.png")

    .. figure:: polbooks_blocks_mdl.*
       :align: center

       Block partition of a political books network, which minimizes the description
1647
       length of the network according to the degree-corrected stochastic blockmodel.
1648
1649
1650
1651
1652
1653
1654
1655
1656

    References
    ----------

    .. [holland-stochastic-1983] Paul W. Holland, Kathryn Blackmond Laskey,
       Samuel Leinhardt, "Stochastic blockmodels: First steps",
       Carnegie-Mellon University, Pittsburgh, PA 15213, U.S.A., :doi:`10.1016/0378-8733(83)90021-7`
    .. [faust-blockmodels-1992] Katherine Faust, and Stanley
       Wasserman. "Blockmodels: Interpretation and Evaluation." Social Networks
1657
       14, no. 1-2 (1992): 5-61. :doi:`10.1016/0378-8733(92)90013-W`
1658
1659
1660
1661
1662
1663
    .. [karrer-stochastic-2011] Brian Karrer, and M. E. J. Newman. "Stochastic
       Blockmodels and Community Structure in Networks." Physical Review E 83,
       no. 1 (2011): 016107. :doi:`10.1103/PhysRevE.83.016107`.
    .. [peixoto-entropy-2012] Tiago P. Peixoto "Entropy of Stochastic Blockmodel
       Ensembles." Physical Review E 85, no. 5 (2012): 056122. :doi:`10.1103/PhysRevE.85.056122`,
       :arxiv:`1112.6028`.
Tiago Peixoto's avatar
Tiago Peixoto committed
1664
1665
1666
    .. [peixoto-parsimonious-2013] Tiago P. Peixoto, "Parsimonious module inference in large networks",
       Phys. Rev. Lett. 110, 148701 (2013), :doi:`10.1103/PhysRevLett.110.148701`, :arxiv:`1212.4794`.
    .. [peixoto-efficient-2013] Tiago P. Peixoto, "Efficient Monte Carlo and greedy
1667
       heuristic for the inference of stochastic block models", :arxiv:`1310.4378`.
1668
1669
1670
1671

    """

    if max_B is None:
1672
1673
1674
1675
        if dense:
            max_B = max(g.num_vertices(), 1)
        else:
            max_B = get_max_B(g.num_vertices(), g.num_edges(), g.is_directed())
1676
1677
1678
1679
1680
1681
1682
1683
1684
        if verbose:
            print("max_B:", max_B)
    if min_B is None:
        min_B = 1

    if mid_B is None:
        mid_B = get_mid(min_B, max_B)

    greedy = greedy_cooling
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723