graph_blockmodel.hh 55.9 KB
Newer Older
1 2
// graph-tool -- a general graph modification and manipulation thingy
//
Tiago Peixoto's avatar
Tiago Peixoto committed
3
// Copyright (C) 2006-2014 Tiago de Paula Peixoto <tiago@skewed.de>
4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22
//
// 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/>.

#ifndef GRAPH_BLOCKMODEL_HH
#define GRAPH_BLOCKMODEL_HH

#include <cmath>
#include <iostream>
23
#include <queue>
24

25
#include "config.h"
Tiago Peixoto's avatar
Tiago Peixoto committed
26 27 28
#include <unordered_set>
#include <unordered_map>
#include <tuple>
29

30
#ifdef HAVE_SPARSEHASH
31 32
#include SPARSEHASH_INCLUDE(dense_hash_set)
#include SPARSEHASH_INCLUDE(dense_hash_map)
33
#endif
34

35
#include "../generation/sampler.hh"
36
#include "../generation/dynamic_sampler.hh"
37

38 39 40
namespace graph_tool
{

41
#ifdef HAVE_SPARSEHASH
42
using google::dense_hash_set;
43
using google::dense_hash_map;
44
#endif
45 46 47 48 49 50 51

using namespace boost;

// ====================
// Entropy calculation
// ====================

52 53 54 55 56 57
// Repeated computation of x*log(x) and log(x) actually adds up to a lot of
// time. A significant speedup can be made by caching pre-computed values. This
// is doable since the values of mrse are bounded in [0, 2E], where E is the
// total number of edges in the network. Here we simply grow the cache as
// needed.

58 59 60 61 62 63 64 65 66 67
template <class Type>
__attribute__((always_inline))
inline double safelog(Type x)
{
    if (x == 0)
        return 0;
    return log(x);
}

__attribute__((always_inline))
68 69 70 71 72 73 74
inline double safelog(size_t x);

__attribute__((always_inline))
inline double xlogx(size_t x);

__attribute__((always_inline))
inline double lgamma_fast(size_t x);
75

76 77 78 79

//
// Sparse entropy
//
80 81

// "edge" term of the entropy
82
template <class Graph>
83
__attribute__((always_inline))
84
inline double eterm(size_t r, size_t s, size_t mrs, const Graph& g)
85
{
86 87
    if (!is_directed::apply<Graph>::type::value && r == s)
        mrs *= 2;
88

89
    double val = xlogx(mrs);
90 91 92 93 94 95 96 97

    if (is_directed::apply<Graph>::type::value || r != s)
        return -val;
    else
        return -val / 2;
}

// "vertex" term of the entropy
98 99
template <class Graph, class Vertex>
inline double vterm(Vertex v, size_t mrp, size_t mrm, size_t wr, bool deg_corr,
100
                    Graph&)
101 102 103 104 105 106 107
{
    double one = 0.5;

    if (is_directed::apply<Graph>::type::value)
        one = 1;

    if (deg_corr)
108
        return one * (xlogx(mrm) + xlogx(mrp));
109
    else
110
        return one * (mrm * safelog(wr) + mrp * safelog(wr));
111 112 113 114 115 116 117 118 119 120 121
}

struct entropy
{
    template <class Graph, class Eprop, class Vprop>
    void operator()(Eprop mrs, Vprop mrp, Vprop mrm, Vprop wr, bool deg_corr,
                    Graph& g, double& S) const
    {
        S = 0;
        typename graph_traits<Graph>::edge_iterator e, e_end;
        for (tie(e, e_end) = edges(g); e != e_end; ++e)
122
            S += eterm(source(*e, g), target(*e, g), mrs[*e], g);
123 124
        typename graph_traits<Graph>::vertex_iterator v, v_end;
        for (tie(v, v_end) = vertices(g); v != v_end; ++v)
125
            S += vterm(*v, mrp[*v], mrm[*v], wr[*v], deg_corr, g);
126 127 128 129
    }
};


130 131 132
//
// Dense entropy
//
133

134
double lbinom(double N, double k)
135
{
136
    return lgamma(N + 1) - lgamma(N - k + 1) - lgamma(k + 1);
137 138
}

139
double lbinom_fast(int N, int k)
140
{
141 142
    return lgamma_fast(N + 1) - lgamma_fast(N - k + 1) - lgamma_fast(k + 1);
}
143 144


145 146 147 148 149 150 151 152
// "edge" term of the entropy
template <class Graph>
__attribute__((always_inline))
inline double eterm_dense(size_t r, size_t s, int mrs, int wr_r,
                          int wr_s, bool multigraph, const Graph& g)
{
    int nrns;
    int ers = mrs;
153

154 155
    if (ers == 0)
        return 0.;
156

157 158 159
    if (r != s || is_directed::apply<Graph>::type::value)
    {
        nrns = wr_r * wr_s;
160
    }
161
    else
162
    {
163 164 165 166
        if (multigraph)
            nrns = (wr_r * (wr_r + 1)) / 2;
        else
            nrns = (wr_r * (wr_r - 1)) / 2;
167 168
    }

169 170 171 172 173 174 175
    double S;
    if (multigraph)
        S = lbinom_fast(nrns + ers - 1, ers);
    else
        S = lbinom_fast(nrns, ers);
    return S;
}
176

177 178
struct entropy_dense
{
179
    template <class Graph, class Eprop, class Vprop>
180
    void operator()(Eprop mrs, Vprop wr, bool multigraph, Graph& g, double& S) const
181
    {
182 183 184
        S = 0;
        typename graph_traits<Graph>::edge_iterator e, e_end;
        for (tie(e, e_end) = edges(g); e != e_end; ++e)
185
        {
186 187 188 189
            typename graph_traits<Graph>::vertex_descriptor r, s;
            r = source(*e, g);
            s = target(*e, g);
            S += eterm_dense(r, s, mrs[*e], wr[r], wr[s], multigraph, g);
190 191 192 193
        }
    }
};

194

195 196 197 198
// ===============================
// Block moves
// ===============================

199 200
// this structure speeds up the access to the edges between given blocks,
// since we're using an adjacency list to store the block structure (the emat_t
201 202 203 204 205 206 207 208 209 210 211 212 213 214
// is simply a corresponding adjacency matrix)
struct get_emat_t
{
    template <class Graph>
    struct apply
    {
        typedef multi_array<pair<typename graph_traits<Graph>::edge_descriptor, bool>, 2> type;
    };
};


struct create_emat
{
    template <class Graph>
215
    void operator()(Graph& g, boost::any& oemap) const
216 217
    {
        typedef typename get_emat_t::apply<Graph>::type emat_t;
218
        size_t B = num_vertices(g);
219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237
        emat_t emat(boost::extents[B][B]);

        for (size_t i = 0; i < B; ++i)
            for (size_t j = 0; j < B; ++j)
                emat[i][j].second = false;
        typename graph_traits<Graph>::edge_iterator e, e_end;
        for (tie(e, e_end) = edges(g); e != e_end; ++e)
        {
            if (source(*e, g) >= B || target(*e, g) >= B)
                throw GraphException("incorrect number of blocks when creating emat!");
            emat[source(*e, g)][target(*e, g)] = make_pair(*e, true);
            if (!is_directed::apply<Graph>::type::value)
                emat[target(*e, g)][source(*e, g)] = make_pair(*e, true);
        }

        oemap = emat;
    }
};

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 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383
template <class Graph>
inline __attribute__((always_inline))
pair<typename graph_traits<Graph>::edge_descriptor, bool>
get_me(typename graph_traits<Graph>::vertex_descriptor r,
       typename graph_traits<Graph>::vertex_descriptor s,
       const typename get_emat_t::apply<Graph>::type& emat, const Graph& g)
{
    return emat[r][s];
}

template <class Graph>
inline __attribute__((always_inline))
void
put_me(typename graph_traits<Graph>::vertex_descriptor r,
       typename graph_traits<Graph>::vertex_descriptor s,
       const typename graph_traits<Graph>::edge_descriptor& e,
       typename get_emat_t::apply<Graph>::type& emat, const Graph& g)
{
    emat[r][s] = make_pair(e, true);
    if (!is_directed::apply<Graph>::type::value && r != s)
        emat[s][r] = make_pair(e, true);
}

template <class Graph>
inline __attribute__((always_inline))
void
remove_me(typename graph_traits<Graph>::vertex_descriptor r,
          typename graph_traits<Graph>::vertex_descriptor s,
          const typename graph_traits<Graph>::edge_descriptor& e,
          typename get_emat_t::apply<Graph>::type& emat, Graph& g,
          bool delete_edge=true)
{
    if (!delete_edge)
    {
        emat[r][s].second = false;
        if (!is_directed::apply<Graph>::type::value && r != s)
            emat[s][r].second = false;
    }
}


// this structure speeds up the access to the edges between given blocks, since
// we're using an adjacency list to store the block structure (this is like
// emat_t above, but takes less space and is slower)
struct get_ehash_t
{
    template <class Graph>
    struct apply
    {
        typedef typename graph_traits<Graph>::vertex_descriptor vertex_t;
        typedef typename graph_traits<Graph>::edge_descriptor edge_t;
#ifdef HAVE_SPARSEHASH
        typedef dense_hash_map<vertex_t, edge_t> map_t;
#else
        typedef unordered_map<vertex_t, edge_t> map_t;
#endif
        typedef vector<map_t> type;
    };
};


template<class Graph>
inline __attribute__((always_inline))
pair<typename graph_traits<Graph>::edge_descriptor, bool>
get_me(typename graph_traits<Graph>::vertex_descriptor r,
       typename graph_traits<Graph>::vertex_descriptor s,
       const typename get_ehash_t::apply<Graph>::type& ehash, const Graph& bg)
{
    const typename get_ehash_t::apply<Graph>::map_t& map = ehash[r];
    typeof(map.begin()) iter = map.find(s);
    if (iter == map.end())
        return (make_pair(typename graph_traits<Graph>::edge_descriptor(), false));
    return make_pair(iter->second, true);
}

template<class Graph>
inline __attribute__((always_inline))
void
put_me(typename graph_traits<Graph>::vertex_descriptor r,
       typename graph_traits<Graph>::vertex_descriptor s,
       const typename graph_traits<Graph>::edge_descriptor& e,
       typename get_ehash_t::apply<Graph>::type& ehash,
       const Graph& bg)
{
    ehash[r][s] = e;
    if (!is_directed::apply<Graph>::type::value)
        ehash[s][r] = e;
}

template<class Graph>
inline __attribute__((always_inline))
void
remove_me(typename graph_traits<Graph>::vertex_descriptor r,
          typename graph_traits<Graph>::vertex_descriptor s,
          const typename graph_traits<Graph>::edge_descriptor& e,
          typename get_ehash_t::apply<Graph>::type& ehash, Graph& bg,
          bool delete_edge=true)
{
    ehash[r].erase(s);
    if (!is_directed::apply<Graph>::type::value)
        ehash[s].erase(r);
    if (delete_edge)
        remove_edge(e, bg);
}

struct create_ehash
{
    template <class Graph>
    void operator()(Graph& g, boost::any& oemap) const
    {
        typedef typename get_ehash_t::apply<Graph>::type emat_t;
        typedef typename get_ehash_t::apply<Graph>::map_t map_t;
        typedef typename graph_traits<Graph>::vertex_descriptor vertex_t;

        emat_t emat(num_vertices(g));

#ifdef HAVE_SPARSEHASH
        typename graph_traits<Graph>::vertex_iterator v, v_end;
        for (tie(v, v_end) = vertices(g); v != v_end; ++v)
        {
            emat[*v].set_empty_key(numeric_limits<vertex_t>::max());
            emat[*v].set_deleted_key(num_vertices(g));
        }
#endif

        typename graph_traits<Graph>::edge_iterator e, e_end;
        for (tie(e, e_end) = edges(g); e != e_end; ++e)
            put_me(source(*e, g), target(*e, g), *e, emat, g);

        oemap = emat;
    }
};

template <class Vertex, class Eprop, class Emat, class BGraph>
__attribute__((always_inline))
inline size_t get_mrs(Vertex r, Vertex s, const Eprop& mrs, Emat& emat,
                      BGraph& bg)
{
    const pair<typename graph_traits<BGraph>::edge_descriptor, bool> me =
        get_me(r, s, emat, bg);
    if (me.second)
        return mrs[me.first];
    else
        return 0;
}

384 385 386 387 388 389 390 391 392 393
// remove a vertex from its current block
template <class Graph, class BGraph, class Eprop, class Vprop, class EWprop,
          class VWprop, class EMat>
void remove_vertex(size_t v, Eprop& mrs, Vprop& mrp, Vprop& mrm, Vprop& wr,
                   Vprop& b, const EWprop& eweight, const VWprop& vweight,
                   Graph& g, BGraph& bg, EMat& emat)
{
    typedef typename graph_traits<Graph>::vertex_descriptor vertex_t;
    vertex_t r = b[v];

394
    bool self_count = false;
395
    typename graph_traits<Graph>::out_edge_iterator ei, ei_end;
396
    for (tie(ei, ei_end) = out_edges(v, g); ei != ei_end; ++ei)
397
    {
398 399
        typename graph_traits<Graph>::edge_descriptor e = *ei;
        vertex_t u = target(e, g);
400 401 402 403 404 405
        if (u == v && !is_directed::apply<Graph>::type::value)
        {
            if (self_count)
                continue;
            self_count = true;
        }
406 407
        vertex_t s = b[u];

408 409 410
        // if (!emat[r][s].second)
        //     throw GraphException("no edge? " + lexical_cast<string>(r) +
        //                          " " + lexical_cast<string>(s));
411

412 413
        typename graph_traits<BGraph>::edge_descriptor me =
            get_me(r, s, emat, bg).first;
414

415
        size_t ew = eweight[e];
416 417 418
        mrs[me] -= ew;

        mrp[r] -= ew;
419 420 421 422
        mrm[s] -= ew;

        if (mrs[me] == 0)
            remove_me(r, s, me, emat, bg);
423 424 425 426 427 428 429 430
    }

    if (is_directed::apply<Graph>::type::value)
    {
        typename in_edge_iteratorS<Graph>::type ie, ie_end;
        for (tie(ie, ie_end) = in_edge_iteratorS<Graph>::get_edges(vertex(v, g), g);
             ie != ie_end; ++ie)
        {
431 432
            typename graph_traits<Graph>::edge_descriptor e = *ie;
            vertex_t u = source(e, g);
433 434 435 436
            if (u == v)
                continue;
            vertex_t s = b[u];

437 438 439
            // if (!emat[s][r].second)
            //     throw GraphException("no edge? " + lexical_cast<string>(s) +
            //                          " " + lexical_cast<string>(r));
440

441 442
            typename graph_traits<BGraph>::edge_descriptor me =
                get_me(s, r, emat, bg).first;
443

444
            size_t ew = eweight[e];
445 446 447 448
            mrs[me] -= ew;

            mrp[s] -= ew;
            mrm[r] -= ew;
449 450 451

            if (mrs[me] == 0)
                remove_me(s, r, me, emat, bg);
452 453 454
        }
    }

455
    wr[r] -= vweight[v];
456 457 458 459 460
}

// add a vertex to block rr
template <class Graph, class BGraph, class Eprop, class Vprop, class EWprop,
          class VWprop,class EMat>
461
void add_vertex(size_t v, size_t r, Eprop& mrs, Vprop& mrp, Vprop& mrm,
462 463 464 465 466
                Vprop& wr, Vprop& b, const EWprop& eweight,
                const VWprop& vweight, Graph& g, BGraph& bg, EMat& emat)
{
    typedef typename graph_traits<Graph>::vertex_descriptor vertex_t;

467
    bool self_count = false;
468 469
    typename graph_traits<Graph>::out_edge_iterator ei, ei_end;
    for (tie(ei, ei_end) = out_edges(vertex(v, g), g); ei != ei_end; ++ei)
470
    {
471 472
        typename graph_traits<Graph>::edge_descriptor e = *ei;
        vertex_t u = target(e, g);
473 474 475 476 477 478
        if (u == v && !is_directed::apply<Graph>::type::value )
        {
            if (self_count)
                continue;
            self_count = true;
        }
479 480 481 482 483 484 485 486 487
        vertex_t s;

        if (u != v)
            s = b[u];
        else
            s = r;

        typename graph_traits<BGraph>::edge_descriptor me;

488 489
        pair<typename graph_traits<BGraph>::edge_descriptor, bool> mep =
                get_me(r, s, emat, bg);
490 491 492

        if (!mep.second)
        {
493 494 495
            mep = add_edge(r, s, bg);
            put_me(r, s, mep.first, emat, bg);
            mrs[mep.first] = 0;
496
        }
497
        me = mep.first;
498

499
        size_t ew = eweight[e];
500 501 502

        mrs[me] += ew;

503 504
        mrp[r] += ew;
        mrm[s] += ew;
505 506 507 508 509 510
    }

    typename in_edge_iteratorS<Graph>::type ie, ie_end;
    for (tie(ie, ie_end) = in_edge_iteratorS<Graph>::get_edges(vertex(v, g), g);
         ie != ie_end; ++ie)
    {
511 512
        typename graph_traits<Graph>::edge_descriptor e = *ie;
        vertex_t u = source(e, g);
513 514 515 516 517 518
        if (u == v)
            continue;

        vertex_t s = b[u];

        typename graph_traits<BGraph>::edge_descriptor me;
519 520 521
        pair<typename graph_traits<BGraph>::edge_descriptor, bool> mep =
                get_me(s, r, emat, bg);

522 523 524

        if (!mep.second)
        {
525 526 527
            mep = add_edge(s, r, bg);
            put_me(s, r, mep.first, emat, bg);
            mrs[mep.first] = 0;
528
        }
529
        me = mep.first;
530

531
        size_t ew = eweight[e];
532 533 534 535 536 537 538

        mrs[me] += ew;

        mrp[s] += ew;
        mrm[r] += ew;
    }

539
    wr[r] += vweight[v];
540 541 542
    b[v] = r;
}

543 544 545 546 547 548 549
// move a vertex from its current block to block nr
template <class Graph, class BGraph, class Eprop, class Vprop, class EWprop,
          class VWprop, class EMat>
void move_vertex(size_t v, size_t nr, Eprop& mrs, Vprop& mrp, Vprop& mrm,
                 Vprop& wr, Vprop& b, bool deg_corr, const EWprop& eweight,
                 const VWprop& vweight, Graph& g, BGraph& bg, EMat& emat)

550
{
551 552
    remove_vertex(v, mrs, mrp, mrm, wr, b, eweight, vweight, g, bg, emat);
    add_vertex(v, nr, mrs, mrp, mrm, wr, b, eweight, vweight, g, bg, emat);
553 554
}

555 556 557 558

template <class Type1, class Type2, class Graph>
__attribute__((always_inline))
inline pair<Type1,Type2> make_ordered_pair(const Type1& v1, const Type2& v2, const Graph&)
559
{
560
    if (!is_directed::apply<Graph>::type::value)
561
    {
562 563 564 565
        if (v1 < v2)
            return make_pair(v1, v2);
        else
            return make_pair(v2, v1);
566
    }
567
    return make_pair(v1, v2);
568 569
}

570 571
template <class Graph>
class EntrySet
572
{
573 574
public:
    EntrySet(size_t B)
575
    {
576 577 578 579 580 581 582 583 584
        _null = numeric_limits<size_t>::max();
        _r_field_t.resize(B, _null);
        _nr_field_t.resize(B, _null);

        if (is_directed::apply<Graph>::type::value)
        {
            _r_field_s.resize(B, _null);
            _nr_field_s.resize(B, _null);
        }
585 586
    }

587
    void SetMove(size_t r, size_t nr)
588
    {
589
        _rnr = make_pair(r, nr);
590 591
    }

592
    void InsertDeltaTarget(size_t r, size_t s, int delta)
593
    {
594 595 596 597 598 599
        if (!is_directed::apply<Graph>::type::value &&
            (s == _rnr.first || s == _rnr.second) && s < r)
        {
            InsertDeltaTarget(s, r, delta);
            return;
        }
600

601 602 603 604 605 606
        if (is_directed::apply<Graph>::type::value &&
            (s == _rnr.first || s == _rnr.second) && s < r)
        {
            InsertDeltaSource(r, s, delta);
            return;
        }
607

608 609 610 611 612 613 614 615 616 617 618 619
        vector<size_t>& field = (_rnr.first == r) ? _r_field_t : _nr_field_t;
        if (field[s] == _null)
        {
            field[s] = _entries.size();
            _entries.push_back(make_pair(r, s));
            _delta.push_back(delta);
        }
        else
        {
            _delta[field[s]] += delta;
        }
    }
620

621 622 623 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 650 651 652 653 654 655 656 657 658 659 660 661 662 663 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 697 698 699 700 701 702 703 704 705 706 707 708 709
    void InsertDeltaSource(size_t s, size_t r, int delta)
    {
        if (s == r)
        {
            InsertDeltaTarget(r, s, delta);
            return;
        }

        if ((s == _rnr.first || s == _rnr.second) && s < r)
        {
            InsertDeltaTarget(s, r, delta);
            return;
        }

        vector<size_t>& field = (_rnr.first == r) ? _r_field_s : _nr_field_s;
        if (field[s] == _null)
        {
            field[s] = _entries.size();
            _entries.push_back(make_pair(s, r));
            _delta.push_back(delta);
        }
        else
        {
            _delta[field[s]] += delta;
        }
    }

    int GetDelta(size_t t, size_t s)
    {
        if (is_directed::apply<Graph>::type::value)
        {
            if (t == _rnr.first || t == _rnr.second)
                return GetDeltaTarget(t, s);
            if (s == _rnr.first || s == _rnr.second)
                return GetDeltaSource(t, s);
        }
        else
        {
            if (t == _rnr.first || t == _rnr.second)
                return GetDeltaTarget(t, s);
            if (s == _rnr.first || s == _rnr.second)
                return GetDeltaTarget(s, t);
        }
        return 0;
    }

    int GetDeltaTarget(size_t r, size_t s)
    {
        vector<size_t>& field = (_rnr.first == r) ? _r_field_t : _nr_field_t;
        if (field[s] == _null)
        {
            return 0;
        }
        else
        {
            return _delta[field[s]];
        }
    }

    int GetDeltaSource(size_t s, size_t r)
    {
        vector<size_t>& field = (_rnr.first == r) ? _r_field_s : _nr_field_s;
        if (field[s] == _null)
        {
            return 0;
        }
        else
        {
            return _delta[field[s]];
        }
    }

    void Clear()
    {
        size_t r, s;
        for (size_t i = 0; i < _entries.size(); ++i)
        {
            tie(r, s) = _entries[i];
            _r_field_t[r] = _nr_field_t[r] = _null;
            _r_field_t[s] = _nr_field_t[s] = _null;
            if (is_directed::apply<Graph>::type::value)
            {
                _r_field_s[r] = _nr_field_s[r] = _null;
                _r_field_s[s] = _nr_field_s[s] = _null;
            }
        }
        _entries.clear();
        _delta.clear();
    }
710

711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729
    vector<pair<size_t, size_t> >& GetEntries() { return _entries; }
    vector<int>& GetDelta() { return _delta; }

private:
    pair<size_t, size_t> _rnr;
    size_t _null;
    vector<size_t> _r_field_t;
    vector<size_t> _nr_field_t;
    vector<size_t> _r_field_s;
    vector<size_t> _nr_field_s;
    vector<pair<size_t, size_t> > _entries;
    vector<int> _delta;
};

// obtain the necessary entries in the e_rs matrix which need to be modified
// after the move
template <class Graph, class BGraph, class Vertex, class Vprop, class Eprop>
void move_entries(Vertex v, Vertex nr, Vprop b, Eprop eweights, Graph& g,
                  BGraph& bg, EntrySet<Graph>& m_entries)
730 731
{
    typedef typename graph_traits<Graph>::vertex_descriptor vertex_t;
732
    vertex_t r = b[v];
733

734
    m_entries.SetMove(r, nr);
735

736 737 738 739 740 741 742 743 744 745 746 747 748 749
    bool self_count = false;
    typename graph_traits<Graph>::out_edge_iterator e, e_end;
    for (tie(e, e_end) = out_edges(v, g); e != e_end; ++e)
    {
        vertex_t u = target(*e, g);
        if (u == v && !is_directed::apply<Graph>::type::value)
        {
            if (self_count)
                continue;
            self_count = true;
        }
        vertex_t s = b[u];
        int ew = eweights[*e];
        //assert(ew > 0);
750

751 752 753 754 755 756 757 758 759 760 761 762 763
        m_entries.InsertDeltaTarget(r, s, -ew);

        //insert_m_entry(r,  s, -ew, m_entries, m_entries_set, g);
        if (u == v)
            s = nr;
        m_entries.InsertDeltaTarget(nr, s, +ew);

        //insert_m_entry(nr, s, +ew, m_entries, m_entries_set, g);
    }

    typename in_edge_iteratorS<Graph>::type ie, ie_end;
    for (tie(ie, ie_end) = in_edge_iteratorS<Graph>::get_edges(v, g);
         ie != ie_end; ++ie)
764
    {
765 766 767 768 769
        vertex_t u = source(*ie, g);
        if (u == v)
            continue;
        vertex_t s = b[u];
        int ew = eweights[*ie];
770

771 772
        m_entries.InsertDeltaSource(s,  r, -ew);
        m_entries.InsertDeltaSource(s, nr, +ew);
773

774 775
        // insert_m_entry(s,  r, -ew, m_entries, m_entries_set, g);
        // insert_m_entry(s, nr, +ew, m_entries, m_entries_set, g);
776 777
    }

778 779 780 781 782 783 784 785 786 787
    //clear_entries_set(m_entries, m_entries_set);
}

// obtain the entropy difference given a set of entries in the e_rs matrix
template <class Graph, class Eprop, class BGraph, class EMat>
double entries_dS(EntrySet<Graph>& m_entries, Eprop& mrs, BGraph& bg, EMat& emat)
{
    typedef typename graph_traits<BGraph>::vertex_descriptor vertex_t;
    vector<pair<size_t, size_t> >& entries = m_entries.GetEntries();
    vector<int>& delta = m_entries.GetDelta();
788

789 790
    double dS = 0;
    for (size_t i = 0; i < entries.size(); ++i)
791
    {
792 793 794 795 796 797 798
        vertex_t er = entries[i].first;
        vertex_t es = entries[i].second;
        int d = delta[i];

        int ers = get_mrs(er, es, mrs, emat, bg);
        //assert(ers + delta >= 0);
        dS += eterm(er, es, ers + d, bg) - eterm(er, es, ers, bg);
799
    }
800
    return dS;
801
}
802

803 804 805 806 807 808 809
// compute the entropy difference of a virtual move of vertex r to block nr
template <class Graph, class BGraph, class Eprop, class Vprop, class EWprop,
          class VWprop, class EMat>
double virtual_move(size_t v, size_t nr, Eprop& mrs, Vprop& mrp, Vprop& mrm,
                    Vprop& wr, Vprop& b, bool deg_corr, const EWprop& eweight,
                    const VWprop& vweight, Graph& g, BGraph& bg, EMat& emat,
                    EntrySet<Graph>& m_entries)
810 811

{
812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831
    typedef typename graph_traits<Graph>::vertex_descriptor vertex_t;
    vertex_t r = b[v];

    if (r == nr)
        return 0.;

    m_entries.Clear();
    move_entries(v, nr, b, eweight, g, bg, m_entries);
    double dS = entries_dS(m_entries, mrs, bg, emat);
    int kout = out_degreeS()(v, g, eweight);
    int kin = kout;
    if (is_directed::apply<Graph>::type::value)
        kin = in_degreeS()(v, g, eweight);
    //assert(mrm[r]  - kin >= 0);
    //assert(mrp[r]  - kout >= 0);
    dS += vterm(r,  mrp[r]  - kout, mrm[r]  - kin, wr[r]  - vweight[v], deg_corr, bg);
    dS += vterm(nr, mrp[nr] + kout, mrm[nr] + kin, wr[nr] + vweight[v], deg_corr, bg);
    dS -= vterm(r,  mrp[r]        , mrm[r]       , wr[r]              , deg_corr, bg);
    dS -= vterm(nr, mrp[nr]       , mrm[nr]      , wr[nr]             , deg_corr, bg);
    return dS;
832 833
}

834 835 836 837 838 839 840 841 842

// compute the entropy difference of a virtual move of vertex r to block nr
template <class Graph, class BGraph, class Eprop, class Vprop, class EWprop,
          class VWprop, class EMat>
double virtual_move_dense(size_t v, size_t nr, Eprop& mrs, Vprop& mrp,
                          Vprop& mrm, Vprop& wr, Vprop& b, bool deg_corr,
                          const EWprop& eweight, const VWprop& vweight,
                          Graph& g, BGraph& bg, EMat& emat, bool multigraph)

843
{
844 845 846 847
    typedef typename graph_traits<Graph>::vertex_descriptor vertex_t;
    vertex_t r = b[v];

    if (r == nr)
848
        return 0;
849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 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 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954

    if (deg_corr)
        throw GraphException("Dense entropy for degree corrected model not implemented!");

    vector<int> deltap(num_vertices(bg), 0);
    int deltal = 0;
    bool self_count = false;
    typename graph_traits<Graph>::out_edge_iterator e, e_end;
    for (tie(e, e_end) = out_edges(v, g); e != e_end; ++e)
    {
        vertex_t u = target(*e, g);
        vertex_t s = b[u];
        if (u == v)
        {
            if (self_count)
                continue;
            self_count = true;
            deltal += eweight[*e];
        }
        else
        {
            deltap[s] += eweight[*e];
        }
    }

    vector<int> deltam(num_vertices(bg), 0);
    typename in_edge_iteratorS<Graph>::type ie, ie_end;
    for (tie(ie, ie_end) = in_edge_iteratorS<Graph>::get_edges(v, g);
         ie != ie_end; ++ie)
    {
        vertex_t u = source(*ie, g);
        if (u == v)
            continue;
        vertex_t s = b[u];
        deltam[s] += eweight[*ie];
    }

    double Si = 0, Sf = 0;
    for (vertex_t s = 0; s < num_vertices(bg); ++s)
    {
        int ers = get_mrs(r, s, mrs, emat, bg);
        int enrs = get_mrs(nr, s, mrs, emat, bg);

        if (!is_directed::apply<Graph>::type::value)
        {
            if (s != nr && s != r)
            {
                Si += eterm_dense(r,  s, ers,              wr[r],               wr[s], multigraph, bg);
                Sf += eterm_dense(r,  s, ers - deltap[s],  wr[r] - vweight[v],  wr[s], multigraph, bg);
                Si += eterm_dense(nr, s, enrs,             wr[nr],              wr[s], multigraph, bg);
                Sf += eterm_dense(nr, s, enrs + deltap[s], wr[nr] + vweight[v], wr[s], multigraph, bg);
            }

            if (s == r)
            {
                Si += eterm_dense(r, r, ers,                      wr[r],              wr[r],              multigraph, bg);
                Sf += eterm_dense(r, r, ers - deltap[r] - deltal, wr[r] - vweight[v], wr[r] - vweight[v], multigraph, bg);
            }

            if (s == nr)
            {
                Si += eterm_dense(nr, nr, enrs,                       wr[nr],              wr[nr],              multigraph, bg);
                Sf += eterm_dense(nr, nr, enrs + deltap[nr] + deltal, wr[nr] + vweight[v], wr[nr] + vweight[v], multigraph, bg);

                Si += eterm_dense(r, nr, ers,                          wr[r],              wr[nr], multigraph, bg);
                Sf += eterm_dense(r, nr, ers - deltap[nr] + deltap[r], wr[r] - vweight[v], wr[nr] + vweight[v], multigraph, bg);
            }
        }
        else
        {
            int esr = get_mrs(s, r, mrs, emat, bg);
            int esnr = get_mrs(s, nr, mrs, emat, bg);

            if (s != nr && s != r)
            {
                Si += eterm_dense(r, s, ers            , wr[r]             , wr[s]             , multigraph, bg);
                Sf += eterm_dense(r, s, ers - deltap[s], wr[r] - vweight[v], wr[s]             , multigraph, bg);
                Si += eterm_dense(s, r, esr            , wr[s]             , wr[r]             , multigraph, bg);
                Sf += eterm_dense(s, r, esr - deltam[s], wr[s]             , wr[r] - vweight[v], multigraph, bg);

                Si += eterm_dense(nr, s, enrs            , wr[nr]             , wr[s]              , multigraph, bg);
                Sf += eterm_dense(nr, s, enrs + deltap[s], wr[nr] + vweight[v], wr[s]              , multigraph, bg);
                Si += eterm_dense(s, nr, esnr            , wr[s]              , wr[nr]             , multigraph, bg);
                Sf += eterm_dense(s, nr, esnr + deltam[s], wr[s]              , wr[nr] + vweight[v], multigraph, bg);
            }

            if(s == r)
            {
                Si += eterm_dense(r, r, ers                                  , wr[r]             , wr[r]             , multigraph, bg);
                Sf += eterm_dense(r, r, ers - deltap[r]  - deltam[r] - deltal, wr[r] - vweight[v], wr[r] - vweight[v], multigraph, bg);

                Si += eterm_dense(r, nr, esnr                         , wr[r]             , wr[nr]             , multigraph, bg);
                Sf += eterm_dense(r, nr, esnr - deltap[nr] + deltam[r], wr[r] - vweight[v], wr[nr] + vweight[v], multigraph, bg);
            }

            if(s == nr)
            {
                Si += eterm_dense(nr, nr, esnr                                   , wr[nr]             , wr[nr]             , multigraph, bg);
                Sf += eterm_dense(nr, nr, esnr + deltap[nr] + deltam[nr] + deltal, wr[nr] + vweight[v], wr[nr] + vweight[v], multigraph, bg);

                Si += eterm_dense(nr, r, esr                         , wr[nr]             , wr[r]             , multigraph, bg);
                Sf += eterm_dense(nr, r, esr + deltap[r] - deltam[nr], wr[nr] + vweight[v], wr[r] - vweight[v], multigraph, bg);
            }
        }
    }
    return Sf - Si;
955 956
}

957

958 959 960
// ====================================
// Construct and manage half-edge lists
// ====================================
961 962 963

//the following guarantees a stable (source, target) ordering even for
//undirected graphs
964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983
template <class Edge, class Graph>
inline typename graph_traits<Graph>::vertex_descriptor
get_source(const Edge& e, const Graph &g)
{
    typename graph_traits<Graph>::vertex_descriptor u, v;
    u = source(e, g);
    v = target(e, g);
    return min(u, v);
}

template <class Edge, class Graph>
inline typename graph_traits<Graph>::vertex_descriptor
get_target(const Edge& e, const Graph &g)
{
    typename graph_traits<Graph>::vertex_descriptor u, v;
    u = source(e, g);
    v = target(e, g);
    return max(u, v);
}

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 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 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 1116 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 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172
struct egroups_manage
{

    template <class Graph, class Weighted>
    struct get_sampler
    {
        typedef typename mpl::if_<Weighted,
                                  DynamicSampler<std::tuple<typename graph_traits<Graph>::edge_descriptor, bool>>,
                                  vector<std::tuple<typename graph_traits<Graph>::edge_descriptor, bool>>>::type type;
    };

    template <class Eprop, class Vprop, class VEprop, class Graph, class VertexIndex>
    static void build(Vprop b, boost::any& oegroups, VEprop esrcpos,
                      VEprop etgtpos, Eprop eweight, Graph& g,
                      VertexIndex vertex_index, bool weighted)
    {
        if (weighted)
        {
            typedef typename get_sampler<Graph, mpl::true_>::type sampler_t;
            typedef typename property_map_type::apply<sampler_t,
                                                      VertexIndex>::type vemap_t;
            vemap_t egroups_checked(vertex_index);
            oegroups = egroups_checked;
            build_dispatch(b, egroups_checked.get_unchecked(num_vertices(g)),
                           esrcpos, etgtpos, eweight, g, vertex_index,
                           mpl::true_());
        }
        else
        {
            typedef typename get_sampler<Graph, mpl::false_>::type sampler_t;
            typedef typename property_map_type::apply<sampler_t,
                                                      VertexIndex>::type vemap_t;
            vemap_t egroups_checked(vertex_index);
            oegroups = egroups_checked;
            build_dispatch(b, egroups_checked.get_unchecked(num_vertices(g)),
                           esrcpos, etgtpos, eweight, g, vertex_index,
                           mpl::true_());
        }
    }

    template <class Eprop, class Vprop, class VEprop, class Graph, class VertexIndex, class Egroups>
    static void build_dispatch(Vprop b, Egroups egroups, VEprop esrcpos,
                               VEprop etgtpos, Eprop eweight, Graph& g,
                               VertexIndex vertex_index, mpl::true_)
    {

        typename graph_traits<Graph>::edge_iterator e, e_end;
        for (tie(e, e_end) = edges(g); e != e_end; ++e)
        {
            size_t r = b[get_source(*e, g)];
            auto& r_elist = egroups[r];
            esrcpos[*e] = insert_edge(std::make_tuple(*e, true), r_elist, eweight[*e]);

            size_t s = b[get_target(*e, g)];
            auto& s_elist = egroups[s];
            etgtpos[*e] = insert_edge(std::make_tuple(*e, false), s_elist, eweight[*e]);
        }
    }

    template <class Edge, class EV>
    static size_t insert_edge(const Edge& e, EV& elist, size_t weight)
    {
        elist.push_back(e);
        return elist.size() - 1;
    }

    template <class Edge>
    static size_t insert_edge(const Edge& e, DynamicSampler<Edge>& elist,
                              size_t weight)
    {
        return elist.insert(e, weight);
    }


    template <class Edge, class Epos>
    static void remove_edge(size_t pos, Epos& esrcpos, Epos& etgtpos,
                            DynamicSampler<Edge>& elist)
    {
        elist.remove(pos);
    }

    template <class Edge, class Epos>
    static void remove_edge(size_t pos, Epos& esrcpos, Epos& etgtpos,
                            vector<Edge>& elist)
    {
        if (get<1>(elist.back()))
            esrcpos[get<0>(elist.back())] = pos;
        else
            etgtpos[get<0>(elist.back())] = pos;
        elist[pos] = elist.back();
        elist.pop_back();
    }

    template <class Vertex, class Graph, class EVprop, class Eprop, class VEprop>
    static void remove_egroups(Vertex v, Vertex r, Eprop eweight,
                               EVprop egroups, VEprop esrcpos, VEprop etgtpos,
                               Graph& g)
    {
        typedef Vertex vertex_t;
        bool self_count = false;

        //update the half-edge lists
        typename all_edges_iteratorS<Graph>::type e, e_end;
        for (tie(e, e_end) = all_edges_iteratorS<Graph>::get_edges(v, g);
             e != e_end; ++e)
        {
            vertex_t src = get_source(*e, g);
            vertex_t tgt = get_target(*e, g);

            bool is_src = (src == v);

            // self-loops will appear twice
            if (src == tgt)
            {
                is_src = !self_count;
                self_count = true;
            }

            auto& elist = egroups[r];
            size_t pos = (is_src) ? esrcpos[*e] : etgtpos[*e];
            remove_edge(pos, esrcpos, etgtpos, elist);
        }
    }

    template <class Vertex, class Graph, class EVprop, class Eprop, class VEprop>
    static void add_egroups(Vertex v, Vertex s, Eprop eweight, EVprop egroups,
                            VEprop esrcpos, VEprop etgtpos, Graph& g)
    {
        typedef Vertex vertex_t;
        bool self_count = false;

        //update the half-edge lists
        typename all_edges_iteratorS<Graph>::type e, e_end;
        for (tie(e, e_end) = all_edges_iteratorS<Graph>::get_edges(v, g);
             e != e_end; ++e)
        {
            vertex_t src = get_source(*e, g);
            vertex_t tgt = get_target(*e, g);

            bool is_src = (src == v);

            // self-loops will appear twice
            if (src == tgt)
            {
                is_src = !self_count;
                self_count = true;
            }

            auto& elist = egroups[s];
            typedef typename tuple_element<0, typename property_traits<EVprop>::value_type::value_type>::type e_type;
            size_t pos = insert_edge(std::make_tuple(e_type(*e), is_src),
                                     elist, size_t(eweight[*e]));
            if (is_src)
                esrcpos[*e] = pos;
            else
                etgtpos[*e] = pos;
        }
    }

    template <class Vertex, class Graph, class EVprop, class Eprop, class VEprop>
    static void update_egroups(Vertex v, Vertex r, Vertex s, Eprop eweight,
                               EVprop egroups, VEprop esrcpos, VEprop etgtpos,
                               Graph& g)
    {
        remove_egroups(v, r, eweight, egroups, esrcpos, etgtpos, g);
        add_egroups(v, s, eweight, egroups, esrcpos, etgtpos, g);
    }

    template <class Edge, class RNG>
    static typename std::tuple_element<0, Edge>::type
    sample_edge(DynamicSampler<Edge>& elist, RNG& rng)
    {
        return get<0>(elist.sample(rng));
    }

    template <class Edge, class RNG>
    static typename std::tuple_element<0, Edge>::type
    sample_edge(vector<Edge>& elist, RNG& rng)
    {
        std::uniform_int_distribution<size_t> urand(0, elist.size() - 1);
        size_t ur = urand(rng);
        return get<0>(elist[ur]);
    }
};

//============
// Main loops
//============

1173
//computes the move proposal probability
1174 1175
template <class Vertex, class Graph, class Vprop, class Eprop, class Emat,
          class BGraph>
1176
inline double
1177 1178 1179
get_move_prob(Vertex v, Vertex r, Vertex s, double c, Vprop b, Eprop mrs,
              Vprop mrp, Vprop mrm, Emat& emat, Eprop eweight, Graph& g,
              BGraph& bg, EntrySet<Graph>& m_entries, bool reverse)
1180 1181
{
    typedef typename graph_traits<Graph>::vertex_descriptor vertex_t;
1182
    size_t B = num_vertices(bg);
1183
    double p = 0;
1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194
    size_t w = 0;

    int kout = 0, kin = 0;
    if (reverse)
    {
        kout = out_degreeS()(v, g, eweight);
        kin = kout;
        if (is_directed::apply<Graph>::type::value)
            kin = in_degreeS()(v, g, eweight);
    }

1195 1196 1197 1198
    typename all_edges_iteratorS<Graph>::type e, e_end;
    for (tie(e, e_end) = all_edges_iteratorS<Graph>::get_edges(v, g);
         e != e_end; ++e)
    {
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
        vertex_t u = target(*e, g);
        if (is_directed::apply<Graph>::type::value && u == v)
            u = source(*e, g);
        vertex_t t = b[u];
        if (u == v)
            t = r;
        size_t ew = eweight[*e];
        w += ew;

        int mts = get_mrs(t, s, mrs, emat, bg);
        int mtp = mrp[t];
        int mst = mts;
        int mtm = mtp;

        if (is_directed::apply<Graph>::type::value)
        {
            mst = get_mrs(s, t, mrs, emat, bg);
            mtm = mrm[t];
        }

        if (reverse)
        {
            int dts = m_entries.GetDelta(t, s);
            int dst = dts;
            if (is_directed::apply<Graph>::type::value)
                dst = m_entries.GetDelta(s, t);

            mts += dts;
            mst += dst;

            if (t == s)
            {
                mtp -= kout;
                mtm -= kin;
            }

            if (t == r)
            {
                mtp += kout;
                mtm += kin;
            }
        }

1242
        if (is_directed::apply<Graph>::type::value)
1243 1244 1245
        {
            p += ew * ((mts + mst + c) / (mtp + mtm + c * B));
        }
1246
        else
1247 1248 1249 1250 1251
        {
            if (t == s)
                mts *= 2;
            p += ew * (mts + c) / (mtp + c * B);
        }
1252
    }
1253
    return p / w;
1254 1255
}

1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270

//A single Monte Carlo Markov chain sweep
template <class Graph, class BGraph, class EMprop, class Eprop, class Vprop,
          class EMat, class EVprop, class VEprop, class SamplerMap, class RNG>
void move_sweep(EMprop mrs, Vprop mrp, Vprop mrm, Vprop wr, Vprop b,
                Vprop clabel, vector<int>& vlist, bool deg_corr, bool dense,
                bool multigraph, double beta, Eprop eweight, Vprop vweight,
                EVprop egroups, VEprop esrcpos, VEprop etgtpos, Graph& g,
                BGraph& bg, EMat& emat, SamplerMap neighbour_sampler,
                bool sequential, bool random_move, double c, bool verbose,
                RNG& rng, double& S, size_t& nmoves)
{
    typedef typename graph_traits<Graph>::vertex_descriptor vertex_t;

    size_t B = num_vertices(bg);
1271

Tiago Peixoto's avatar
Tiago Peixoto committed
1272 1273
    typedef std::uniform_real_distribution<> rdist_t;
    auto rand_real = std::bind(rdist_t(), std::ref(rng));
1274

1275 1276
    // it is useful to shuffle the vertex order even in the parallel case, so
    // threads become more balanced
1277
    std::shuffle(vlist.begin(), vlist.end(), rng);
1278

Tiago Peixoto's avatar
Tiago Peixoto committed
1279
    std::uniform_int_distribution<size_t> s_rand(0, B - 1);
1280

1281 1282
    nmoves = 0;
    S = 0;
1283

1284 1285
    EntrySet<Graph> m_entries(B);

1286 1287 1288
    int i = 0, N = vlist.size();
    for (i = 0; i < N; ++i)
    {
1289 1290 1291 1292 1293 1294 1295 1296 1297 1298
        vertex_t v;
        if (sequential)
        {
            v = vertex(vlist[i], g);
        }
        else
        {
            std::uniform_int_distribution<size_t> v_rand(0, N - 1);
            v = vertex(vlist[v_rand(rng)], g);
        }
1299 1300 1301

        vertex_t r = b[v];

1302 1303 1304 1305
        // blocks can't become empty
        if (wr[r] == 1)
            continue;

1306
        // attempt random block