blockmodel.py 107 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
from .. spectral import adjacency
29
30
import random
from numpy import *
31
import numpy
32
33
from scipy.optimize import fsolve, fminbound
import scipy.special
34
from collections import defaultdict
35
36
import copy
import heapq
37
38
39
40

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

41
__test__ = False
42

43
44
45
46
47
48
49
50
51
52
53
54
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

55
56
57
58
59
60
61
62
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`
63
        Graph to be modelled.
64
    eweight : :class:`~graph_tool.PropertyMap` (optional, default: ``None``)
65
        Edge multiplicities (for multigraphs or block graphs).
66
    vweight : :class:`~graph_tool.PropertyMap` (optional, default: ``None``)
67
        Vertex multiplicities (for block graphs).
68
69
70
71
72
73
74
75
    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``)
76
77
        Constraint labels on the vertices. If supplied, vertices with different
        label values will not be clustered in the same group.
78
79
80
    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.
81
82
83
84
    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,
85
86
    """

87
88
    _state_ref_count = 0

89
    def __init__(self, g, eweight=None, vweight=None, b=None,
90
91
92
                 B=None, clabel=None, deg_corr=True,
                 max_BE=1000, **kwargs):

93
94
        BlockState._state_ref_count += 1

95
        # initialize weights to unity, if necessary
96
97
        if eweight is None:
            eweight = g.new_edge_property("int")
98
            eweight.fa = 1
99
100
101
102
        elif eweight.value_type() != "int32_t":
            eweight = eweight.copy(value_type="int32_t")
        if vweight is None:
            vweight = g.new_vertex_property("int")
103
            vweight.fa = 1
104
105
        elif vweight.value_type() != "int32_t":
            vweight = vweight.copy(value_type="int32_t")
106
107
108
109
110
111
112
        self.eweight = g.own_property(eweight)
        self.vweight = g.own_property(vweight)

        self.is_weighted = True if self.eweight.fa.max() > 1 else False

        # configure the main graph and block model parameters
        self.g = g
113

114
115
        self.E = int(self.eweight.fa.sum())
        self.N = int(self.vweight.fa.sum())
116
117
118

        self.deg_corr = deg_corr

119
120
121
122
        # ensure we have at most as many blocks as nodes
        if B is not None:
            B = min(B, self.g.num_vertices())

123
        if b is None:
124
            # create a random partition into B blocks.
125
126
            if B is None:
                B = get_max_B(self.N, self.E, directed=g.is_directed())
127
128
            B = min(B, self.g.num_vertices())
            ba = random.randint(0, B, self.g.num_vertices())
129
            ba[:B] = arange(B)        # avoid empty blocks
130
131
            if B < self.g.num_vertices():
                random.shuffle(ba)
132
            b = g.new_vertex_property("int")
133
            b.fa = ba
134
135
            self.b = b
        else:
136
137
138
139
140
141
            # if a partition is available, we will incorporate it.
            if isinstance(b, numpy.ndarray):
                self.b = g.new_vertex_property("int")
                self.b.fa = b
            else:
                self.b = b = g.own_property(b.copy(value_type="int"))
142
            if B is None:
143
144
145
146
                B = int(self.b.fa.max()) + 1

        # if B > self.N:
        #     raise ValueError("B > N!")
147

148
        if self.b.fa.max() >= B:
149
150
151
            raise ValueError("Maximum value of b is larger or equal to B!")

        # Construct block-graph
152
        self.bg = get_block_graph(g, B, self.b, self.vweight, self.eweight)
153
154
155
156
        self.bg.set_fast_edge_removal()

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

158
159
160
161
162
163
164
165
166
        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
167
168
169

        self.vertices = libcommunity.get_vector(B)
        self.vertices.a = arange(B)
170
        self.B = B
171

172
173
174
175
176
177
178
        if clabel is not None:
            if isinstance(clabel, PropertyMap):
                self.clabel = clabel
            else:
                self.clabel = self.g.new_vertex_property("int")
                self.clabel.a = clabel
        else:
179
180
181
182
183
184
185
            self.clabel = self.g.new_vertex_property("int")

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

186
187
        self.overlap = False

188
189
190
191
        # used by mcmc_sweep()
        self.egroups = None
        self.nsampler = None
        self.sweep_vertices = None
192
193
        self.overlap_stats = libcommunity.overlap_stats()
        self.partition_stats = libcommunity.partition_stats()
194

195
196
197
        # computation cache
        libcommunity.init_safelog(int(5 * max(self.E, self.N)))
        libcommunity.init_xlogx(int(5 * max(self.E, self.N)))
198
        libcommunity.init_lgamma(int(3 * max(self.E, self.N)))
199

200
201
202
203
204
205
206
207
    def __del__(self):
        BlockState._state_ref_count -= 1
        if BlockState._state_ref_count == 0:
            libcommunity.clear_safelog()
            libcommunity.clear_xlogx()
            libcommunity.clear_lgamma()


208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
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
301
302
303
304
305
    def __repr__(self):
        return "<BlockState object with %d blocks,%s for graph %s, at 0x%x>" % \
            (self.B, " degree corrected," if self.deg_corr else "", str(self.g),
             id(self))


    def __init_partition_stats(self, empty=True):
        if not empty:
            self.partition_stats = libcommunity.init_partition_stats(self.g._Graph__graph,
                                                                     _prop("v", self.g, self.b),
                                                                     _prop("e", self.g, self.eweight),
                                                                     self.N, self.B)
        else:
            self.partition_stats = libcommunity.partition_stats()



    def copy(self, b=None, B=None, deg_corr=None, clabel=None, overlap=False):
        r"""Copies the block state. The parameters override the state properties, and
         have the same meaning as in the constructor. If ``overlap=True`` an
         instance of :class:`~graph_tool.community.OverlapBlockState` is
         returned."""

        if not overlap:
            state = BlockState(self.g,
                               eweight=self.eweight,
                               vweight=self.vweight,
                               b=self.b.copy() if b is None else b,
                               B=(self.B if b is None else None) if B is None else B,
                               clabel=self.clabel if clabel is None else clabel,
                               deg_corr=self.deg_corr if deg_corr is None else deg_corr,
                               max_BE=self.max_BE)
        else:
            state = OverlapBlockState(self.g,
                                      b=b if b is not None else self.b,
                                      B=(self.B if b is None else None) if B is None else B,
                                      clabel=self.clabel if clabel is None else clabel,
                                      deg_corr=self.deg_corr if deg_corr is None else deg_corr,
                                      max_BE=self.max_BE)

        if not state.__check_clabel():
            b = state.b.a + state.clabel.a * state.B
            continuous_map(b)
            state = state.copy(b=b)

            if __test__:
                assert state.__check_clabel()

        return state


    def __getstate__(self):
        state = dict(g=self.g,
                     eweight=self.eweight,
                     vweight=self.vweight,
                     b=self.b,
                     B=self.B,
                     clabel=self.clabel,
                     deg_corr=self.deg_corr,
                     max_BE=self.max_BE)
        return state

    def __setstate__(self, state):
        self.__init__(**state)
        return state

    def get_block_state(self, b=None, vweight=False, deg_corr=False, overlap=False):
        r"""Returns a :class:`~graph_tool.community.BlockState`` corresponding to the
        block graph. The parameters have the same meaning as the in the constructor."""


        state = BlockState(self.bg, eweight=self.mrs,
                           vweight=self.wr if vweight else None,
                           b=self.bg.vertex_index.copy("int") if b is None else b,
                           clabel=self.get_bclabel(),
                           deg_corr=deg_corr,
                           max_BE=self.max_BE)
        if overlap:
            state = state.copy(overlap=True)
        n_map = self.b.copy()
        return state, n_map

    def get_bclabel(self):
        r"""Returns a :class:`~graph_tool.PropertyMap`` corresponding to constraint
        labels for the block graph."""

        bclabel = self.bg.new_vertex_property("int")
        reverse_map(self.b, bclabel)
        pmap(bclabel, self.clabel)
        return bclabel

    def __check_clabel(self):
        b = self.b.a + self.clabel.a * self.B
        continuous_map(b)
        b2 = self.b.copy()
        continuous_map(b2.a)
        return (b == b2.a).all()

306
307
308
309
    def __get_emat(self):
        if self.emat is None:
            self.__regen_emat()
        return self.emat
310
311

    def __regen_emat(self):
312
313
314
315
        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)
316

317
    def __build_egroups(self, empty=False):
318
319
        self.esrcpos = self.g.new_edge_property("int")
        self.etgtpos = self.g.new_edge_property("int")
320

321
        self.egroups = libcommunity.build_egroups(self.g._Graph__graph,
322
323
324
325
326
                                                  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),
327
                                                  self.is_weighted, empty)
328

329
    def __build_nsampler(self, empty=False):
330
        self.nsampler = libcommunity.init_neighbour_sampler(self.g._Graph__graph,
331
                                                            _prop("e", self.g, self.eweight),
332
                                                            True, empty)
333
334
335
336
337
338

    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()
339
        self.emat = None
340
341
342
343
344
345
346
347
348
349

    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):
350
351
        r"""Returns the edge property map of the block graph which contains the :math:`e_{rs}` matrix entries.
        For undirected graphs, the diagonal values (self-loops) contain :math:`e_{rr}/2`."""
352
353
354
355
356
357
358
359
        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():
360
            return self.mrp, self.mrm
361
362
363
364
365
366
367
        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

368
369
370
    def entropy(self, complete=True, dl=False, partition_dl=False, dense=False,
                multigraph=True, norm=True, dl_ent=False, **kwargs):
        r"""Calculate the entropy associated with the current block partition.
371
372
373
374
375
376
377
378

        Parameters
        ----------
        complete : ``bool`` (optional, default: ``False``)
            If ``True``, the complete entropy will be returned, including constant
            terms not relevant to the block partition.
        dl : ``bool`` (optional, default: ``False``)
            If ``True``, the full description length will be returned.
379
380
381
        partition_dl : ``bool`` (optional, default: ``False``)
            If ``True``, and ``dl == True`` only the partition description
            length will be considered.
382
383
384
        dense : ``bool`` (optional, default: ``False``)
            If ``True``, the "dense" variant of the entropy will be computed.
        multigraph : ``bool`` (optional, default: ``False``)
385
386
387
388
389
390
391
            If ``True``, the multigraph entropy will be used.
        norm : ``bool`` (optional, default: ``True``)
            If ``True``, the entropy will be "normalized" by dividing by the
            number of edges.
        dl_ent : ``bool`` (optional, default: ``False``)
            If ``True``, the description length of the degree sequence will be
            approximated by its entropy.
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417

        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 =
418
        \sum_se_{sr}` are the numbers of out- and in-edges adjacent to block
419
420
        :math:`r`, respectively.

421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
        If ``dense == False`` and ``multigraph == True``, the entropy used will
        be of the "Poisson" model, with the additional term:

        .. math::

            {\mathcal{S}_{cm}^{(d)}} = \mathcal{S}_c^{(d)} + \sum_{i>j} \ln A_{ij}! + \sum_i \ln A_{ii}!!


        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::

            \mathcal{L}_c = \mathcal{L}_t + \sum_r\min\left(\mathcal{L}^{(1)}_r, \mathcal{L}^{(2)}_r\right),

        with

        .. math::

              \mathcal{L}^{(1)}_r &= \ln{\left(\!\!{n_r \choose e_r}\!\!\right)}, \\
              \mathcal{L}^{(2)}_r &= \ln\Xi_r + \ln n_r! - \sum_k \ln n^r_k!,

        and :math:`\ln\Xi_r \simeq 2\sqrt{\zeta(2)e_r}`, where :math:`\zeta(x)`
        is the `Riemann zeta function
        <https://en.wikipedia.org/wiki/Riemann_zeta_function>`_, and
        :math:`n^r_k` is the number of nodes in block :math:`r` with degree
        :math:`k`. For directed graphs we have instead :math:`k \to (k^-, k^+)`,
        and :math:`\ln\Xi_r \to \ln\Xi^+_r + \ln\Xi^-_r` with :math:`\ln\Xi_r^+
        \simeq 2\sqrt{\zeta(2)e^+_r}` and :math:`\ln\Xi_r^- \simeq
        2\sqrt{\zeta(2)e^-_r}`.

        If ``dl_ent=True`` is passed, this will be approximated instead by
455
456
457

        .. math::

458
            \mathcal{L}_c \simeq \mathcal{L}_t - \sum_rn_r\sum_kp^r_k\ln p^r_k,
459

460
        where :math:`p^r_k = n^r_k / n_r`.
461

462
463
        If the "dense" entropies are requested (``dense=True``), they will be
        computed as
464
465
466
467
468
469
470
471
472
473
474
475
476

        .. 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)},

477
478
479
        for multigraphs (i.e. ``multigraph == True``). A dense entropy for the
        degree-corrected model is not available, and if requested will return a
        :exc:`NotImplementedError`.
480

481
482
        If ``complete == False`` constants in the above equations which do not
        depend on the partition of the nodes will be omitted.
483

484
485
        Note that in all cases if ``norm==True`` the value returned corresponds
        to the entropy `per edge`, i.e. :math:`(\mathcal{S}_{t/c}\; [\,+\,\mathcal{L}_{t/c}])/ E`.
486
487
        """

488
489
490
        xi_fast = kwargs.get("xi_fast", False)
        dl_deg_alt = kwargs.get("dl_deg_alt", True)

491
492
493
        E = self.E
        N = self.N

494
495
        if dense:
            if self.deg_corr:
496
                raise NotImplementedError('A degree-corrected "dense" entropy is not yet implemented')
497

498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
            S = libcommunity.entropy_dense(self.bg._Graph__graph,
                                           _prop("e", self.bg, self.mrs),
                                           _prop("v", self.bg, self.wr),
                                           multigraph)
        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 __test__:
                assert not isnan(S) and not isinf(S), "invalid entropy %g (%s) " % (S, str(dict(complete=complete,
                                                                                                random=random, dl=dl,
                                                                                                partition_dl=partition_dl,
                                                                                                dense=dense, multigraph=multigraph,
                                                                                                norm=norm)))
516
517
518
            if complete:
                if self.deg_corr:
                    S -= E
519
520
521
522
                    S += libcommunity.deg_entropy_term(self.g._Graph__graph,
                                                       libcore.any(),
                                                       self.overlap_stats,
                                                       self.N)
523
524
                else:
                    S += E
525

526
527
528
529
530
531
532
533
534
535
                if multigraph:
                    S += libcommunity.entropy_parallel(self.g._Graph__graph,
                                                       _prop("e", self.g, self.eweight))

                if __test__:
                    assert not isnan(S) and not isinf(S), "invalid entropy %g (%s) " % (S, str(dict(complete=complete,
                                                                                                    random=random, dl=dl,
                                                                                                    partition_dl=partition_dl,
                                                                                                    dense=dense, multigraph=multigraph,
                                                                                                    norm=norm)))
536
        if dl:
537
538
539
540
541
542
543
544
545
546
547
548
549
550
            if partition_dl:
                if self.partition_stats.is_enabled():
                    S += self.partition_stats.get_partition_dl()
                else:
                    self.__init_partition_stats(empty=False)
                    S += self.partition_stats.get_partition_dl()
                    self.__init_partition_stats(empty=True)

                if __test__:
                    assert not isnan(S) and not isinf(S), "invalid entropy %g (%s) " % (S, str(dict(complete=complete,
                                                                                                    random=random, dl=dl,
                                                                                                    partition_dl=partition_dl,
                                                                                                    dense=dense, multigraph=multigraph,
                                                                                                    norm=norm)))
551
552
553
            else:
                S += model_entropy(self.B, N, E, directed=self.g.is_directed(), nr=self.wr.a) * E

554
            if self.deg_corr:
555
556
557
558
559
560
                if self.partition_stats.is_enabled():
                    S_seq = self.partition_stats.get_deg_dl(dl_ent, dl_deg_alt, xi_fast)
                else:
                    self.__init_partition_stats(empty=False)
                    S_seq = self.partition_stats.get_deg_dl(dl_ent, dl_deg_alt, xi_fast)
                    self.__init_partition_stats(empty=True)
561

562
                S += S_seq
563

564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
                if __test__:
                    assert not isnan(S_seq) and not isinf(S_seq), "invalid entropy %g (%s) " % (S_seq, str(dict(complete=complete,
                                                                                                                random=random, dl=dl,
                                                                                                                partition_dl=partition_dl,
                                                                                                                dense=dense, multigraph=multigraph,
                                                                                                                norm=norm)))

        if __test__:
            assert not isnan(S) and not isinf(S), "invalid entropy %g (%s) " % (S, str(dict(complete=complete,
                                                                                            random=random, dl=dl,
                                                                                            partition_dl=partition_dl,
                                                                                            dense=dense, multigraph=multigraph,
                                                                                            norm=norm)))

        if norm:
            return S / E
        else:
            return S
582

583
584
585
    def get_matrix(self):
        r"""Returns the block matrix (as a sparse :class:`~scipy.sparse.csr_matrix`),
        which contains the number of edges between each block pair.
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601

        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)
602
           >>> m = state.get_matrix()
603
604
           >>> figure()
           <...>
605
           >>> matshow(m.todense())
606
607
608
609
610
611
612
613
614
615
616
617
618
           <...>
           >>> savefig("bloc_mat.pdf")

        .. testcleanup:: get_matrix

           savefig("bloc_mat.png")

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

           A  5x5 block matrix.

       """
619

620
        return adjacency(self.bg, weight=self.mrs)
621
622


623
def model_entropy(B, N, E, directed=False, nr=None):
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
    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::

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

652

653
654
    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`.
655

656
657
658
659
    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`.
660
661
662
663

    References
    ----------

Tiago Peixoto's avatar
Tiago Peixoto committed
664
665
    .. [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`.
Tiago Peixoto's avatar
Tiago Peixoto committed
666
667
668
    .. [peixoto-hierarchical-2014] Tiago P. Peixoto, "Hierarchical block structures and high-resolution
       model selection in large networks ", Phys. Rev. X 4, 011047 (2014), :doi:`10.1103/PhysRevX.4.011047`,
       :arxiv:`1310.4377`.
669
670
671

    """

672
673
674
675
676
677
    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
678
679
680
681

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

682
def lbinom(n, k):
683
    return scipy.special.gammaln(float(n + 1)) - scipy.special.gammaln(float(n - k + 1)) - scipy.special.gammaln(float(k + 1))
684
685
686
687
688

def partition_entropy(B, N, nr=None):
    if nr is None:
        S = N * log(B) + log1p(-(1 - 1./B) ** N)
    else:
689
        S = lbinom(B + N - 1, N) + scipy.special.gammaln(N + 1) - scipy.special.gammaln(nr + 1).sum()
690
    return S
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710

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
711
712
    .. [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`.
713
714
715
716
717
718
719
720
721
722
723


    """

    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
724
    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.
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757

    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)
758
    2.414413200430159
759
760
761

    References
    ----------
Tiago Peixoto's avatar
Tiago Peixoto committed
762
763
    .. [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`.
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778

    """
    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

779
780
781
782
def mcmc_sweep(state, beta=1., c=1., dl=False, dense=False, multigraph=False,
               node_coherent=False, nmerges=0, nmerge_sweeps=1, merge_map=None,
               coherent_merge=False, sequential=True, parallel=False,
               vertices=None, verbose=False, **kwargs):
783
    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.
784
785
786

    Parameters
    ----------
787
    state : :class:`~graph_tool.community.BlockState` or :class:`~graph_tool.community.OverlapBlockState`
788
        The block state.
789
    beta : ``float`` (optional, default: `1.0`)
790
        The inverse temperature parameter :math:`\beta`.
791
792
793
794
795
    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.
796
797
798
    dl : ``bool`` (optional, default: ``False``)
        If ``True``, the change in the whole description length will be
        considered after each vertex move, not only the entropy.
799
800
801
802
803
    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``.
804
805
806
807
    node_coherent : ``bool`` (optional, default: ``False``)
        If ``True``, and if the ``state`` is an instance of
        :class:`~graph_tool.community.OverlapBlockState`, then all half-edges
        incident on the same node are moved simultaneously.
808
809
810
811
812
    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.
813
814
815
816
817
818
819
820
821
822
    parallel : ``bool`` (optional, default: ``False``)
        If ``True``, the updates are performed in parallel (multiple
        threads).

        .. warning::

            If this is used, the Markov Chain is not guaranteed to be sampled with
            the correct probabilities. This is better used in conjunction with
            ``beta=float('inf')``, where this is not an issue.

Tiago Peixoto's avatar
Tiago Peixoto committed
823
    vertices : ``list of ints`` (optional, default: ``None``)
824
825
        A list of vertices which will be attempted to be moved. If ``None``, all
        vertices will be attempted.
826
827
828
829
830
831
    verbose : ``bool`` (optional, default: ``False``)
        If ``True``, verbose information is displayed.

    Returns
    -------

832
    dS : ``float``
833
834
835
836
837
838
839
840
       The entropy difference (per edge) after a full sweep.
    nmoves : ``int``
       The number of accepted block membership moves.


    Notes
    -----

841
    This algorithm performs a Markov chain Monte Carlo sweep on the network,
842
843
    where the block memberships are randomly moved, and either accepted or
    rejected, so that after sufficiently many sweeps the partitions are sampled
844
    with probability proportional to :math:`e^{-\beta\mathcal{S}_{t/c}}`, where
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
    :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
872
873
    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
874
    block label of a random neighbour of the vertex being moved. See
875
    [peixoto-efficient-2014]_ for more details.
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916

    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
917
       14, no. 1-2 (1992): 5-61. :doi:`10.1016/0378-8733(92)90013-W`
918
919
920
921
922
923
    .. [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
924
925
    .. [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`.
926
927
928
    .. [peixoto-efficient-2014] Tiago P. Peixoto, "Efficient Monte Carlo and greedy
       heuristic for the inference of stochastic block models", Phys. Rev. E 89, 012804 (2014),
       :doi:`10.1103/PhysRevE.89.012804`, :arxiv:`1310.4378`.
929
930
931
    .. [peixoto-model-2014] Tiago P. Peixoto, "Model selection and hypothesis
       testing for large-scale network models with overlapping groups",
       :arxiv:`1409.3059`.
932
933
    """

934
    if state.B == 1:
935
936
        return 0., 0

937
    if vertices is not None:
938
939
940
        vlist = libcommunity.get_vector(len(vertices))
        vlist.a = vertices
        vertices = vlist
941
        #state.sweep_vertices = vertices
942

943
944
945
946
947
    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

948
949
    random_move = c == float("inf")

950
    if random_move or nmerges > 0:
951
952
953
954
        state._BlockState__build_egroups(empty=True)
    elif state.egroups is None:
        state._BlockState__build_egroups(empty=False)

955
    bclabel = state.get_bclabel()
956

957
958
959
960
961
962
963
964
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
    if nmerges == 0:
        nmerge_sweeps = 1
        if state.nsampler is None:
            state._BlockState__build_nsampler(empty=state.overlap)
        nsampler = state.nsampler
        ncavity_sampler = state.nsampler
        merge_map = state.g.new_vertex_property("int")
    else:
        if kwargs.get("unweighted_merge", False):
            emask = state.mrs
        else:
            emask = state.mrs.copy()
            emask.a = emask.a > 0

        nsampler = libcommunity.init_neighbour_sampler(state.bg._Graph__graph,
                                                       _prop("e", state.bg, emask),
                                                       True, False)
        ncavity_sampler = libcommunity.init_neighbour_sampler(state.bg._Graph__graph,
                                                              _prop("e", state.bg, emask),
                                                              False, False)
        beta = float("inf")

        if merge_map is None:
            merge_map = state.g.vertex_index.copy("int")

    if dl and not state.partition_stats.is_enabled():
        if state.overlap:
            state._OverlapBlockState__init_partition_stats(empty=False)
        else:
            state._BlockState__init_partition_stats(empty=False)

    if not dl and state.partition_stats.is_enabled():
        if state.overlap:
            state._OverlapBlockState__init_partition_stats(empty=True)
        else:
            state._BlockState__init_partition_stats(empty=True)

    if __test__:
        assert state._BlockState__check_clabel(), "clabel already invalid!"
        S = state.entropy(dense=dense, multigraph=multigraph, complete=False, dl=dl, dl_deg_alt=False, xi_fast=True)
        assert not (isinf(S) or isnan(S)), "invalid entropy before sweep: %g" % S
998
999

    try:
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
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
        if not state.overlap:
            dS, nmoves = libcommunity.move_sweep(state.g._Graph__graph,
                                                 state.bg._Graph__graph,
                                                 state._BlockState__get_emat(),
                                                 nsampler, ncavity_sampler,
                                                 _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, 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,
                                                 parallel, random_move,
                                                 c, state.is_weighted,
                                                 nmerges, nmerge_sweeps,
                                                 _prop("v", state.g, merge_map),
                                                 state.partition_stats,
                                                 verbose, _get_rng())
        else:
            dS, nmoves = libcommunity.move_sweep_overlap(state.g._Graph__graph,
                                                         state.bg._Graph__graph,
                                                         state._BlockState__get_emat(),
                                                         nsampler,
                                                         ncavity_sampler,
                                                         _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, bclabel),
                                                         state.sweep_vertices,
                                                         state.deg_corr, dense, multigraph,
                                                         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, parallel,
                                                         random_move, float(c),
                                                         ((nmerges == 0 and node_coherent) or
                                                          (nmerges > 0 and coherent_merge)),
                                                         state.is_weighted,
                                                         nmerges, nmerge_sweeps,
                                                         _prop("v", state.g, merge_map),
                                                         state.overlap_stats,
                                                         state.partition_stats,
                                                         verbose, _get_rng())

1057
1058
1059
    finally:
        if random_move:
            state.egroups = None
1060
1061
1062
        if nmerges > 0:
            state.nsampler = None
            state.egroups = None
1063

1064
1065
1066
1067
1068
1069
1070
1071
    if __test__:
        assert state._BlockState__check_clabel(), "clabel invalidated!"
        assert not (isinf(dS) or isnan(dS)), "invalid after sweep: %g" % dS
        S2 = state.entropy(dense=dense, multigraph=multigraph, complete=False, dl=dl, dl_deg_alt=False, xi_fast=True)
        c_dS = S2 - S
        if not abs(dS / state.E - c_dS) < 1e-6:
            print(dS / state.E, c_dS, nmoves, state.overlap, dense, multigraph, state.deg_corr, state.is_weighted, node_coherent)
        assert abs(dS / state.E - c_dS) < 1e-6, "invalid delta S (%g, %g)" % (dS / state.E, c_dS)
1072

1073
    return dS / state.E, nmoves
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
def pmap(prop, value_map):
    """Maps all the values of `prop` to the values given by `value_map`, which
    is indexed by the values of `prop`."""
    if isinstance(prop, PropertyMap):
        prop = prop.a
    if isinstance(value_map, PropertyMap):
        value_map = value_map.a
    if prop.max() >= len(value_map):
        raise ValueError("value map is not large enough!")
    libcommunity.vector_map(prop, value_map)

def reverse_map(prop, value_map):
    """Modify `value_map` such that the positions indexed by the values in `prop`
    correspond to their index in `prop`."""
    if isinstance(prop, PropertyMap):
        prop = prop.a
    if isinstance(value_map, PropertyMap):
        value_map = value_map.a
    if prop.max() >= len(value_map):
        raise ValueError("value map is not large enough!")
    libcommunity.vector_rmap(prop, value_map)

def continuous_map(prop):
    """Remap the values of ``prop`` in the continuous range :math:`[0, N-1]`."""
    if isinstance(prop, PropertyMap):
        prop = prop.a
    if prop.max() < len(prop):
        rmap = -ones(len(prop), dtype=prop.dtype)
        libcommunity.vector_map(prop, rmap)
    else:
        libcommunity.vector_continuous_map(prop)

def greedy_shrink(state, B, **kwargs):
1109
1110
    if B > state.B:
        raise ValueError("Cannot shrink to a larger size!")
1111

1112
1113
1114
1115
1116
    kwargs = kwargs.copy()
    if kwargs.get("nmerge_sweeps", None) is None:
        kwargs["nmerge_sweeps"] = max((2 * state.g.num_edges()) // state.g.num_vertices(), 1)

    verbose = kwargs.get("verbose", False)
1117

1118
1119
    orig_state = state
    state = state.copy(B=state.B)
1120

1121
1122
1123
1124
    # merge according to indirect neighbourhood; we put all group-nodes in their
    # own groups, and merge/move them until the desired size is reached
    curr_B = (state.wr.a > 0).sum()
    assert curr_B >= B, "shrinking to a larger size ?! (%d, %d)" % (curr_B, B)
1125

1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
    random = kwargs.get("random_move", False)
    old_state = state
    if not state.overlap:
        state, n_map = state.get_block_state(vweight=True, deg_corr=state.deg_corr)
        unilevel_minimize(state, **kwargs)
    merge_map = state.g.vertex_index.copy("int")

    unweighted = False
    kwargs["c"] = 0 if not random else float("inf")
    kwargs["dl"] = False
    while curr_B > B:
        dS, nmoves = mcmc_sweep(state, beta=float("inf"),
                                nmerges=curr_B - B,
                                merge_map=merge_map,
                                unweighted_merge=unweighted,
                                **kwargs)

        #assert nmoves == curr_B - (state.wr.a > 0).sum()
        curr_B = (state.wr.a > 0).sum()
1145

1146
        if verbose:
1147
1148
1149
1150
            print("merging, B=%d" % curr_B, "left:", curr_B - B,
                  "(%g, %d%s%s)" % (dS, nmoves, ", random" if random else "",
                                    ", unweighted" if unweighted else ""))

1151
        if nmoves == 0:
1152
1153
1154
1155
1156
1157
1158
            if not unweighted:
                unweighted = True
            else:
                kwargs["c"] = float("inf")
                random = True

    #assert curr_B == (state.wr.a > 0).sum()
1159

1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
    if not state.overlap:
        unilevel_minimize(state, **kwargs)  # block level moves
        pmap(merge_map, state.b)
        pmap(n_map, merge_map)
        continuous_map(n_map)
        state = orig_state.copy(b=n_map, B=B)
    else:
        pmap(merge_map, state.b)
        continuous_map(merge_map)
        state = orig_state.copy(b=merge_map, B=B)
1170
1171


1172
1173
1174
1175
1176
1177
1178
    if __test__:
        assert state._BlockState__check_clabel(), "clabel already invalidated!"
        assert curr_B == (state.wr.a > 0).sum()
        curr_B = (state.wr.a > 0).sum()
        assert state.B == curr_B
        assert state.B == B

1179
1180
1181
1182
1183
1184
    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.
1185
    This object can be saved to disk via the :mod:`pickle` interface."""
1186
1187
1188
1189

    def __init__(self):
        self.b_cache = {}
        self.checkpoint_state = defaultdict(dict)
1190
        self.init = True
1191
1192

    def clear(self):
1193
        r"""Clear state."""
1194
1195
1196
        self.b_cache.clear()
        self.checkpoint_state.clear()

1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
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
1257
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
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
def unilevel_minimize(state, nsweeps=10, adaptive_sweeps=True, epsilon=0,
                      anneal=(1., 1.), greedy=True, c=0., dl=False, dense=False,
                      multigraph=True, sequential=True, parallel=False,
                      verbose=False, **kwargs):
    kwargs = kwargs.copy()
    kwargs.update(dict(c=c, dl=dl, dense=dense, multigraph=multigraph,
                       sequential=sequential, parallel=parallel))

    t_dS, t_nmoves = 0, 0

    S = state.entropy()

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

        for i in range(ntotal):
            if i < niter:
                continue
            if i < nsweeps and not greedy:
                beta = anneal[0]
            else:
                beta = float("inf")

            delta, nmoves = mcmc_sweep(state, beta=beta, **kwargs)

            S += delta
            t_dS += delta
            t_nmoves += nmoves
            niter += 1
    else:
        # adaptive mode
        min_dl = S
        max_dl = S
        count = 0
        bump = False
        beta =  anneal[0]
        last_min = min_dl
        greedy_step = greedy
        total_nmoves = 0

        if verbose and not greedy:
            print("Performing sweeps for beta = %g, B=%d (N=%d)..." % \
                   (beta, state.B, 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
                    count = 0
                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, state.B, state.g.num_vertices()))

            delta, nmoves = mcmc_sweep(state, beta=beta, **kwargs)

            if state.overlap and beta == float("inf"):
                ds, nm = mcmc_sweep(state, beta=beta, node_coherent=True, **kwargs)
                delta += ds
                nmoves += nm

            S += delta
            niter += 1
            total_nmoves += nmoves

            t_dS += delta
            t_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

        if verbose:
            if not greedy_step:
                print("... performed %d sweeps with %d vertex moves" % (niter, total_nmoves))
            print("Performing sweeps for beta = ∞, B=%d (N=%d)..." % \
                  (state.B, state.g.num_vertices()))

        if not greedy_step:
            min_dl = S
            count = 0

        niter = 0
        total_nmoves = 0
        while count <= nsweeps:
            delta, nmoves = mcmc_sweep(state, beta=float("inf"), **kwargs)

            if state.overlap:
                ds, nm = mcmc_sweep(state, beta=float("inf"), node_coherent=True, **kwargs)
                delta += ds
                nmoves += nm

            S += delta
            niter += 1
            total_nmoves += nmoves

            t_dS += delta
            t_nmoves += nmoves

            if abs(delta) > eps and nmoves / state.g.num_vertices() > epsilon:
                min_dl = S
                count = 0
            else:
                count += 1

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

        bi = state.b

    return t_dS, t_nmoves

1327
1328
1329

def multilevel_minimize(state, B, nsweeps=10, adaptive_sweeps=True, epsilon=0,
                        anneal=(1., 1.), r=2., nmerge_sweeps=10, greedy=True,
1330
1331
1332
                        c=0., dl=False, dense=False, multigraph=True,
                        sequential=True, parallel=False, checkpoint=None,
                        minimize_state=None, verbose=False, **kwargs):
Tiago Peixoto's avatar
Tiago Peixoto committed
1333
    r"""Performs an agglomerative heuristic, which progressively merges blocks together (while allowing individual node moves) to achieve a good partition in ``B`` blocks.
1334
1335
1336

    Parameters
    ----------
1337
    state : :class:`~graph_tool.community.BlockState` or :class:`~graph_tool.community.OverlapBlockState`
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
        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.
1366
    c : ``float`` (optional, default: ``0.0``)
1367
1368
1369
1370
        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.
1371
1372
1373
    dl : ``bool`` (optional, default: ``False``)
        If ``True``, the change in the whole description length will be
        considered after each vertex move, not only the entropy.
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
    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.
1384
1385
1386
    parallel : ``bool`` (optional, default: ``False``)
        If ``True``, the updates are performed in parallel (multiple
        threads).
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
    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
1431
    the :func:`mcmc_sweep` moves, at different scales. See [peixoto-efficient-2014]_
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
    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
    ----------

1465
1466
1467
    .. [peixoto-efficient-2014] Tiago P. Peixoto, "Efficient Monte Carlo and greedy
       heuristic for the inference of stochastic block models", Phys. Rev. E 89, 012804 (2014),
       :doi:`10.1103/PhysRevE.89.012804`, :arxiv:`1310.4378`.
1468
1469
    """

1470
1471
    nonoverlap_compare = kwargs.get("nonoverlap_compare", False)

1472
1473
1474
1475
1476
    if minimize_state is None:
        minimize_state = MinimizeState()
    b_cache = minimize_state.b_cache
    checkpoint_state = minimize_state.checkpoint_state

1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
    nkwargs = dict(nsweeps=nsweeps, epsilon=epsilon, c=c,
                   dl=dl, dense=dense, multigraph=multigraph,
                   nmerge_sweeps=nmerge_sweeps,
                   sequential=sequential, parallel=parallel)
    kwargs = copy.copy(kwargs)
    kwargs.update(nkwargs)

    orig_state = state

    if __test__:
        assert state._BlockState__check_clabel(), "orig clabel already invalidated!"

    # some simple boundary conditions
1490
    if B == 1:
1491
1492
        if state.clabel.fa.max() > 0:
            raise ValueError("Cannot shrink to B = 1 without invalidating constraints")
1493
        bi = state.g.new_vertex_property("int")
1494
        state = state.copy(b=bi)
1495
1496
        return state
    if B == state.g.num_vertices():
1497
1498
        bi = state.g.vertex_index.copy("int")
        state = state.copy(b=bi)
1499
1500
1501
1502
        return state

    Bi = state.B
    while True:
1503
        # keep reducing B by a factor of "r", until desired size is reached
1504
1505
1506
1507
        Bi = max(int(round(Bi / r)), B)
        if Bi == state.B and Bi > B:
            Bi -= 1

1508
        # check cache for previous results
1509
        if b_cache is not None and Bi in b_cache:
1510
1511
1512
1513
1514
            if __test__:
                assert (state.clabel.a == b_cache[Bi][1].clabel.a).all(), "wrong clabel in cache"
                assert state._BlockState__check_clabel(), "clabel already invalidated before cache"
                assert b_cache[Bi][1]._BlockState__check_clabel(), "clabel already invalidated after cache"
            state = b_cache[Bi][1].copy()
1515

1516
1517
1518
1519
        if __test__:
            assert state._BlockState__check_clabel(), "clabel already invalidated!"

        # if necessary, shrink state
1520
1521
1522
1523
        if Bi < state.B:
            if verbose:
                print("Shrinking:", state.B, "->", Bi)

1524
            state = greedy_shrink(state, B=Bi, verbose=verbose, **kwargs)
1525

1526
1527
            if __test__:
                assert state._BlockState__check_clabel(), "clabel invalidated after shrink"
1528

1529
        dS, nmoves = unilevel_minimize(state, checkpoint=checkpoint, verbose=verbose, **kwargs)
1530

1531
1532
        if __test__:
            assert state._BlockState__check_clabel(), "clabel invalidated after unilevel minimize!"
1533

1534
        if state.overlap and state.deg_corr and nonoverlap_compare:
1535
            if verbose:
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
                print("Attempting nonoverlapping minimize...")
            nstate = state.copy(b=state.get_nonoverlap_blocks(), overlap=False)
            assert nstate.B <= nstate.N
            nstate = multilevel_minimize(nstate, B=Bi, verbose=verbose, **kwargs)
            nstate = nstate.copy(overlap=True, clabel=state.clabel.a)
            unilevel_minimize(nstate, **kwargs)

            if nstate.B > Bi:
                nstate = multilevel_minimize(nstate, B=Bi, verbose=verbose,
                                             nonoverlap_compare=False, **kwargs)

            if nstate.entropy(dense=dense, multigraph=multigraph) < state.entropy(dense=dense, multigraph=multigraph):
                if verbose:
                    print("Nonoverlapping minimize improved.")
                state = nstate

                if __test__:
                    assert state._BlockState__check_clabel(), "clabel invalidated after nonoverlap compare!"

        if Bi == B:
            break
1557

1558
    return state
1559

1560
1561
def get_b_dl(state, dense, multigraph, nested_dl, complete=False,
             nested_overlap=False, dl_ent=False):
1562
    if not nested_dl:
1563
1564
        dl = state.entropy(dense=dense, multigraph=multigraph, dl=True,
                           complete=complete, dl_ent=dl_ent)