__init__.py 42.2 KB
Newer Older
1
#! /usr/bin/env python
2
# -*- coding: utf-8 -*-
3
#
4
5
# graph_tool -- a general graph manipulation python module
#
Tiago Peixoto's avatar
Tiago Peixoto committed
6
# Copyright (C) 2007-2011 Tiago de Paula Peixoto <tiago@skewed.de>
7
8
9
10
11
12
13
14
15
16
17
18
19
20
#
# 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/>.

21
"""
22
``graph_tool.generation`` - Random graph generation
23
---------------------------------------------------
24
25
26
27
28
29
30
31
32
33
34
35

Summary
+++++++

.. autosummary::
   :nosignatures:

   random_graph
   random_rewire
   predecessor_tree
   line_graph
   graph_union
36
   triangulation
37
38
   lattice
   geometric_graph
39
   price_network
40
41
42

Contents
++++++++
43
44
"""

Tiago Peixoto's avatar
Tiago Peixoto committed
45
46
from .. dl_import import dl_import
dl_import("import libgraph_tool_generation")
47

48
from .. import Graph, GraphView, _check_prop_scalar, _prop, _limit_args, _gt_type
Tiago Peixoto's avatar
Tiago Peixoto committed
49
from .. stats import label_parallel_edges, label_self_loops
50
51
import inspect
import types
52
import sys, numpy, numpy.random
53

Tiago Peixoto's avatar
Tiago Peixoto committed
54
__all__ = ["random_graph", "random_rewire", "predecessor_tree", "line_graph",
55
56
           "graph_union", "triangulation", "lattice", "geometric_graph",
           "price_network"]
57

Tiago Peixoto's avatar
Tiago Peixoto committed
58

59
def random_graph(N, deg_sampler, deg_corr=None, directed=True,
60
61
                 parallel_edges=False, self_loops=False, blockmodel=None,
                 random=True, mix_time=10, verbose=False):
Tiago Peixoto's avatar
Tiago Peixoto committed
62
63
64
65
66
67
68
69
70
71
72
73
74
    r"""
    Generate a random graph, with a given degree distribution and correlation.

    Parameters
    ----------
    N : int
        Number of vertices in the graph.
    deg_sampler : function
        A degree sampler function which is called without arguments, and returns
        a tuple of ints representing the in and out-degree of a given vertex (or
        a single int for undirected graphs, representing the out-degree). This
        function is called once per vertex, but may be called more times, if the
        degree sequence cannot be used to build a graph.
75
76
77
78
79
80

        Optionally, you can also pass a function which receives one argument. In
        this case the argument passed will be the index of the vertex which will
        receive the degree. If ``blockmodel != None``, the value passed will be
        the block value of the respective vertex instead.
    deg_corr : function (optional, default: ``None``)
Tiago Peixoto's avatar
Tiago Peixoto committed
81
        A function which gives the degree correlation of the graph. It should be
Tiago Peixoto's avatar
Tiago Peixoto committed
82
83
84
85
86
        callable with two parameters: the in,out-degree pair of the source
        vertex an edge, and the in,out-degree pair of the target of the same
        edge (for undirected graphs, both parameters are single values). The
        function should return a number proportional to the probability of such
        an edge existing in the generated graph.
87
88
89
90

        If ``blockmodel != None``, the value passed to the function will be the
        block value of the respective vertices, not the in/out-degree pairs.
    directed : bool (optional, default: ``True``)
Tiago Peixoto's avatar
Tiago Peixoto committed
91
        Whether the generated graph should be directed.
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
    parallel_edges : bool (optional, default: ``False``)
        If ``True``, parallel edges are allowed.
    self_loops : bool (optional, default: ``False``)
        If ``True``, self-loops are allowed.
    blockmodel : list or :class:`~numpy.ndarray` or function (optional, default: ``None``)
        If supplied, the graph will be sampled from a blockmodel ensemble. If
        the value is a list or a :class:`~numpy.ndarray`, it must have
        ``len(block_model) == N``, and the values will define to which block
        each vertex belongs.

        If this value is a function, it will be used to sample the block
        types. It must be callable either with no arguments or with a single
        argument which will be the vertex index. In either case it must return
        an integer.
    random : bool (optional, default: ``True``)
        If ``True``, the returned graph is randomized. Otherwise a deterministic
        placement of the edges will be used.
    mix_time : int (optional, default: ``10``)
110
111
112
        Number of edge sweeps to perform in order to mix the graph. This value
        is ignored if ``parallel_edges == self_loops == True`` and
        ``strat != "probabilistic"``.
113
114
    verbose : bool (optional, default: ``False``)
        If ``True``, verbose information is displayed.
Tiago Peixoto's avatar
Tiago Peixoto committed
115
116
117

    Returns
    -------
118
    random_graph : :class:`~graph_tool.Graph`
Tiago Peixoto's avatar
Tiago Peixoto committed
119
        The generated graph.
120
121
122
    blocks : :class:`~graph_tool.PropertyMap`
        A vertex property map with the block values. This is only returned if
        ``blockmodel != None``.
Tiago Peixoto's avatar
Tiago Peixoto committed
123
124
125
126
127
128
129

    See Also
    --------
    random_rewire: in place graph shuffling

    Notes
    -----
Tiago Peixoto's avatar
Tiago Peixoto committed
130
131
132
    The algorithm makes sure the degree sequence is graphical (i.e. realizable)
    and keeps re-sampling the degrees if is not. With a valid degree sequence,
    the edges are placed deterministically, and later the graph is shuffled with
133
134
    the :func:`~graph_tool.generation.random_rewire` function, with the
    ``mix_time`` parameter passed as ``n_iter``.
Tiago Peixoto's avatar
Tiago Peixoto committed
135

136
    The complexity is :math:`O(V + E)` if parallel edges are allowed, and
137
    :math:`O(V + E \times\text{mix-time})` if parallel edges are not allowed.
138
139
140
141
142
143
144
145
146
147
148
149
150


    .. note ::

        If ``parallel_edges == False`` this algorithm only guarantees that the
        returned graph will be a random sample from the desired ensemble if
        ``mix_time`` is sufficiently large. The algorithm implements an
        efficient Markov chain based on edge swaps, with a mixing time which
        depends on the degree distribution and correlations desired. If degree
        correlations are provided, the mixing time tends to be larger.

        If ``strat == "probabilistic"``, the Markov chain still needs to be
        mixed, even if parallel edges and self-loops are allowed. In this case
151
152
153
        the Markov chain is implemented using the Metropolis-Hastings
        [metropolis-equations-1953]_ [hastings-monte-carlo-1970]_
        acceptance/rejection algorithm.
Tiago Peixoto's avatar
Tiago Peixoto committed
154
155
156
157
158

    Examples
    --------
    >>> from numpy.random import randint, random, seed, poisson
    >>> from pylab import *
159
    >>> seed(43)
Tiago Peixoto's avatar
Tiago Peixoto committed
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180

    This is a degree sampler which uses rejection sampling to sample from the
    distribution :math:`P(k)\propto 1/k`, up to a maximum.

    >>> def sample_k(max):
    ...     accept = False
    ...     while not accept:
    ...         k = randint(1,max+1)
    ...         accept = random() < 1.0/k
    ...     return k
    ...

    The following generates a random undirected graph with degree distribution
    :math:`P(k)\propto 1/k` (with k_max=40) and an *assortative* degree
    correlation of the form:

    .. math::

        P(i,k) \propto \frac{1}{1+|i-k|}

    >>> g = gt.random_graph(1000, lambda: sample_k(40),
181
182
    ...                     lambda i, k: 1.0 / (1 + abs(i - k)), directed=False,
    ...                     mix_time=100)
Tiago Peixoto's avatar
Tiago Peixoto committed
183
    >>> gt.scalar_assortativity(g, "out")
Tiago Peixoto's avatar
Tiago Peixoto committed
184
    (0.6435658697163692, 0.010420519538259333)
Tiago Peixoto's avatar
Tiago Peixoto committed
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211

    The following samples an in,out-degree pair from the joint distribution:

    .. math::

        p(j,k) = \frac{1}{2}\frac{e^{-m_1}m_1^j}{j!}\frac{e^{-m_1}m_1^k}{k!} +
                 \frac{1}{2}\frac{e^{-m_2}m_2^j}{j!}\frac{e^{-m_2}m_2^k}{k!}

    with :math:`m_1 = 4` and :math:`m_2 = 20`.

    >>> def deg_sample():
    ...    if random() > 0.5:
    ...        return poisson(4), poisson(4)
    ...    else:
    ...        return poisson(20), poisson(20)
    ...

    The following generates a random directed graph with this distribution, and
    plots the combined degree correlation.

    >>> g = gt.random_graph(20000, deg_sample)
    >>>
    >>> hist = gt.combined_corr_hist(g, "in", "out")
    >>> imshow(hist[0], interpolation="nearest")
    <...>
    >>> colorbar()
    <...>
212
    >>> xlabel("in-degree")
Tiago Peixoto's avatar
Tiago Peixoto committed
213
    <...>
214
    >>> ylabel("out-degree")
Tiago Peixoto's avatar
Tiago Peixoto committed
215
    <...>
216
    >>> savefig("combined-deg-hist.pdf")
Tiago Peixoto's avatar
Tiago Peixoto committed
217

218
    .. figure:: combined-deg-hist.*
Tiago Peixoto's avatar
Tiago Peixoto committed
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
        :align: center

        Combined degree histogram.

    A correlated directed graph can be build as follows. Consider the following
    degree correlation:

    .. math::

         P(j',k'|j,k)=\frac{e^{-k}k^{j'}}{j'!}
         \frac{e^{-(20-j)}(20-j)^{k'}}{k'!}

    i.e., the in->out correlation is "disassortative", the out->in correlation
    is "assortative", and everything else is uncorrelated.
    We will use a flat degree distribution in the range [1,20).

    >>> p = scipy.stats.poisson
    >>> g = gt.random_graph(20000, lambda: (sample_k(19), sample_k(19)),
237
238
239
    ...                     lambda a,b: (p.pmf(a[0], b[1]) *
    ...                                  p.pmf(a[1], 20 - b[0])),
    ...                     mix_time=100)
Tiago Peixoto's avatar
Tiago Peixoto committed
240
241
242

    Lets plot the average degree correlations to check.

243
244
245
246
    >>> figure(figsize=(6,3))
    <...>
    >>> axes([0.1,0.15,0.63,0.8])
    <...>
Tiago Peixoto's avatar
Tiago Peixoto committed
247
    >>> corr = gt.avg_neighbour_corr(g, "in", "in")
248
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
249
    ...         label=r"$\left<\text{in}\right>$ vs in")
Tiago Peixoto's avatar
Tiago Peixoto committed
250
251
    (...)
    >>> corr = gt.avg_neighbour_corr(g, "in", "out")
252
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
253
    ...         label=r"$\left<\text{out}\right>$ vs in")
Tiago Peixoto's avatar
Tiago Peixoto committed
254
255
    (...)
    >>> corr = gt.avg_neighbour_corr(g, "out", "in")
256
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
257
    ...          label=r"$\left<\text{in}\right>$ vs out")
Tiago Peixoto's avatar
Tiago Peixoto committed
258
259
    (...)
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
260
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
261
    ...          label=r"$\left<\text{out}\right>$ vs out")
Tiago Peixoto's avatar
Tiago Peixoto committed
262
    (...)
263
    >>> legend(loc=(1.05,0.5))
Tiago Peixoto's avatar
Tiago Peixoto committed
264
265
266
267
268
    <...>
    >>> xlabel("source degree")
    <...>
    >>> ylabel("average target degree")
    <...>
269
    >>> savefig("deg-corr-dir.pdf")
Tiago Peixoto's avatar
Tiago Peixoto committed
270

271
    .. figure:: deg-corr-dir.*
Tiago Peixoto's avatar
Tiago Peixoto committed
272
273
274
        :align: center

        Average nearest neighbour correlations.
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294


    **Blockmodels**


    The following example shows how a blockmodel [#]_ can be generated. We will
    consider a system of 10 blocks, which form communities. The connection
    probability will be given by

    >>> def corr(a, b):
    ...    if a == b:
    ...        return 0.999
    ...    else:
    ...        return 0.001

    The blockmodel can be generated as follows.

    >>> g, bm = gt.random_graph(1000, lambda: poisson(10), directed=False,
    ...                         blockmodel=lambda: randint(10), deg_corr=corr,
    ...                         mix_time=500)
295
    >>> gt.graph_draw(g, vcolor=bm, layout="sfdp", output="blockmodel.pdf")
296
297
    <...>

298
    .. figure:: blockmodel.*
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
        :align: center

        Simple blockmodel with 10 blocks.


    .. [#] Blockmodels are equivalent to the hidden variable model
       [boguna-correlated-2003]_.


    References
    ----------
    .. [metropolis-equations-1953]  Metropolis, N.; Rosenbluth, A.W.;
       Rosenbluth, M.N.; Teller, A.H.; Teller, E. "Equations of State
       Calculations by Fast Computing Machines". Journal of Chemical Physics 21
       (6): 1087–1092 (1953). :doi:`10.1063/1.1699114`
    .. [hastings-monte-carlo-1970] Hastings, W.K. "Monte Carlo Sampling Methods
       Using Markov Chains and Their Applications". Biometrika 57 (1): 97–109 (1970).
       :doi:`10.1093/biomet/57.1.97`
    .. [boguna-correlated-2003] M. Boguñá and R. Pastor-Satorras, "Class of
       correlated random networks with hidden variables" Physical Review E
       68, 036112 (2003) :doi:`10.1103/PhysRevE.68.036112`
Tiago Peixoto's avatar
Tiago Peixoto committed
320
    """
321

322
    seed = numpy.random.randint(0, sys.maxint)
323
324
325
326
327
    g = Graph()
    if deg_corr == None:
        uncorrelated = True
    else:
        uncorrelated = False
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343

    if (type(blockmodel) is types.FunctionType or
        type(blockmodel) is types.LambdaType):
        bm = numpy.zeros(N, dtype="int")
        if len(inspect.getargspec(blockmodel)[0]) == 0:
            for i in xrange(N):
                bm[i] = blockmodel()
        else:
            for i in xrange(N):
                bm[i] = blockmodel(i)
        blockmodel = numpy.array(bm)

    if len(inspect.getargspec(deg_sampler)[0]) > 0:
        if blockmodel is not None:
            sampler = lambda i: deg_sampler(blockmodel[i])
        else:
Tiago Peixoto's avatar
Tiago Peixoto committed
344
            sampler = deg_sampler
345
346
347
348
    else:
        sampler = lambda i: deg_sampler()

    libgraph_tool_generation.gen_graph(g._Graph__graph, N, sampler,
349
350
351
                                       uncorrelated, not parallel_edges,
                                       not self_loops, not directed,
                                       seed, verbose, True)
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
    g.set_directed(directed)

    if blockmodel is not None:
        btype = _gt_type(blockmodel[0])
        bm = g.new_vertex_property(btype)
        if btype in ["object", "string"] or "vector" in btype:
            for v in g.vertices():
                bm[v] = blockmodel[int(v)]
        else:
            try:
                bm.a = blockmodel
            except ValueError:
                bm = g.new_vertex_property("object")
                for v in g.vertices():
                    bm[v] = blockmodel[int(v)]
    else:
        bm = None
369

370
    if parallel_edges and self_loops and deg_corr is None:
371
        mix_time = 1
Tiago Peixoto's avatar
Tiago Peixoto committed
372
    if random:
373
374
        if deg_corr is not None:
            random_rewire(g, strat="probabilistic", n_iter=mix_time,
Tiago Peixoto's avatar
Tiago Peixoto committed
375
                          parallel_edges=parallel_edges, deg_corr=deg_corr,
376
377
                          self_loops=self_loops, blockmodel=bm,
                          verbose=verbose)
378
379
380
381
        else:
            random_rewire(g, parallel_edges=parallel_edges, n_iter=mix_time,
                          self_loops=self_loops, verbose=verbose)

382
383
384
385
    if bm is None:
        return g
    else:
        return g, bm
386

Tiago Peixoto's avatar
Tiago Peixoto committed
387

Tiago Peixoto's avatar
Tiago Peixoto committed
388
@_limit_args({"strat": ["erdos", "correlated", "uncorrelated", "probabilistic"]})
389
390
def random_rewire(g, strat="uncorrelated", n_iter=1, edge_sweep=True,
                  parallel_edges=False, self_loops=False, deg_corr=None,
391
                  blockmodel=None, ret_fail=False, verbose=False):
392
    r"""
393
394
395
396
397
398
399
400
    Shuffle the graph in-place.

    If ``strat != "erdos"``, the degrees (either in or out) of each vertex are
    always the same, but otherwise the edges are randomly placed. If
    ``strat == "correlated"``, the degree correlations are also maintained: The
    new source and target of each edge both have the same in and out-degree. If
    ``strat == "probabilistic"``, then edges are rewired according to the degree
    correlation given by the parameter ``deg_corr``.
401
402
403

    Parameters
    ----------
404
    g : :class:`~graph_tool.Graph`
405
        Graph to be shuffled. The graph will be modified.
406
407
408
409
    strat : string (optional, default: ``"uncorrelated"``)
        If ``strat == "erdos"``, the resulting graph will be entirely random. If
        ``strat == "uncorrelated"`` only the degrees of the vertices will be
        maintained, nothing else. If ``strat == "correlated"``, additionally the
410
        new source and target of each edge both have the same in and out-degree.
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
        If ``strat == "probabilistic"``, than edges are rewired according to the
        degree correlation given by the parameter ``deg_corr``.
    n_iter : int (optional, default: ``1``)
        Number of iterations. If ``edge_sweep == True``, each iteration
        corresponds to an entire "sweep" over all edges. Otherwise this
        corresponds to the total number of edges which are randomly chosen for a
        swap attempt (which may repeat).
    edge_sweep : bool (optional, default: ``True``)
        If ``True``, each iteration will perform an entire "sweep" over the
        edges, where each edge is visited once in random order, and a edge swap
        is attempted.
    parallel : bool (optional, default: ``False``)
        If ``True``, parallel edges are allowed.
    self_loops : bool (optional, default: ``False``)
        If ``True``, self-loops are allowed.
    deg_corr : function (optional, default: ``None``)
427
428
429
430
431
432
        A function which gives the degree correlation of the graph. It should be
        callable with two parameters: the in,out-degree pair of the source
        vertex an edge, and the in,out-degree pair of the target of the same
        edge (for undirected graphs, both parameters are single values). The
        function should return a number proportional to the probability of such
        an edge existing in the generated graph. This parameter is ignored,
433
        unless ``strat == "probabilistic"``.
434
435
436
437
438
439
440

        If ``blockmodel != None``, the value passed to the function will be the
        block value of the respective vertices, not the in/out-degree pairs.
    blockmodel : :class:`~graph_tool.PropertyMap` (optional, default: ``None``)
        If supplied, the graph will be rewired to conform to a blockmodel
        ensemble. The value must be a vertex property map which defines the
        block of each vertex.
441
442
443
444
445
446
447
448
449
450
451
    ret_fail : bool (optional, default: ``False``)
        If ``True``, the number of failed edge moves (due to parallel edges or
        self-loops) is returned.
    verbose : bool (optional, default: ``False``)
        If ``True``, verbose information is displayed.


    Returns
    -------
    fail_count : int
        Number of failed edge moves (due to parallel edges or self-loops).
452
453
454
455
456
457
458

    See Also
    --------
    random_graph: random graph generation

    Notes
    -----
Tiago Peixoto's avatar
Tiago Peixoto committed
459
    This algorithm iterates through all the edges in the network and tries to
460
    swap its target or source with the target or source of another edge.
Tiago Peixoto's avatar
Tiago Peixoto committed
461
462

    .. note::
463

464
465
466
467
468
469
470
471
472
473
        If ``parallel_edges = False``, parallel edges are not placed during
        rewiring. In this case, the returned graph will be a uncorrelated sample
        from the desired ensemble only if ``n_iter`` is sufficiently large. The
        algorithm implements an efficient Markov chain based on edge swaps, with
        a mixing time which depends on the degree distribution and correlations
        desired. If degree probabilistic correlations are provided, the mixing
        time tends to be larger.

        If ``strat == "probabilistic"``, the Markov chain still needs to be
        mixed, even if parallel edges and self-loops are allowed. In this case
474
475
476
        the Markov chain is implemented using the Metropolis-Hastings
        [metropolis-equations-1953]_ [hastings-monte-carlo-1970]_
        acceptance/rejection algorithm.
477

Tiago Peixoto's avatar
Tiago Peixoto committed
478

479
    Each edge is tentatively swapped once per iteration, so the overall
480
481
    complexity is :math:`O(V + E \times \text{n-iter})`. If ``edge_sweep ==
    False``, the complexity becomes :math:`O(V + E + \text{n-iter})`.
482

483
484
485
486
487
    Examples
    --------

    Some small graphs for visualization.

488
    >>> from numpy.random import random, seed
489
    >>> from pylab import *
490
    >>> seed(43)
491
    >>> g, pos = gt.triangulation(random((1000,2)))
492
    >>> gt.graph_draw(g, layout="arf", output="rewire_orig.pdf", size=(6,6))
493
    <...>
494
    >>> gt.random_rewire(g, "correlated")
495
    >>> gt.graph_draw(g, layout="arf", output="rewire_corr.pdf", size=(6,6))
496
    <...>
497
    >>> gt.random_rewire(g)
498
    >>> gt.graph_draw(g, layout="arf", output="rewire_uncorr.pdf", size=(6,6))
499
    <...>
500
    >>> gt.random_rewire(g, "erdos")
501
    >>> gt.graph_draw(g, layout="arf", output="rewire_erdos.pdf", size=(6,6))
502
    <...>
503

504
    Some `ridiculograms <http://www.youtube.com/watch?v=YS-asmU3p_4>`_ :
505

506
507
508
509
    .. image:: rewire_orig.*
    .. image:: rewire_corr.*
    .. image:: rewire_uncorr.*
    .. image:: rewire_erdos.*
510

511
512
    *From left to right:* Original graph --- Shuffled graph, with degree
    correlations --- Shuffled graph, without degree correlations --- Shuffled graph,
513
    with random degrees.
514
515
516

    We can try some larger graphs to get better statistics.

517
518
    >>> figure()
    <...>
519
    >>> g = gt.random_graph(30000, lambda: sample_k(20),
520
521
    ...                     lambda i, j: exp(abs(i-j)), directed=False,
    ...                     mix_time=100)
522
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
523
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-", label="original")
524
525
526
    (...)
    >>> gt.random_rewire(g, "correlated")
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
527
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="*", label="correlated")
528
529
530
    (...)
    >>> gt.random_rewire(g)
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
531
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-", label="uncorrelated")
532
    (...)
533
534
    >>> gt.random_rewire(g, "erdos")
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
535
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-", label="Erdos")
536
    (...)
537
538
539
540
541
542
    >>> xlabel("$k$")
    <...>
    >>> ylabel(r"$\left<k_{nn}\right>$")
    <...>
    >>> legend(loc="best")
    <...>
543
    >>> savefig("shuffled-stats.pdf")
544

545
    .. figure:: shuffled-stats.*
546
547
548
549
550
551
552
553
554
555
556
        :align: center

        Average degree correlations for the different shuffled and non-shuffled
        graphs. The shuffled graph with correlations displays exactly the same
        correlation as the original graph.

    Now let's do it for a directed graph. See
    :func:`~graph_tool.generation.random_graph` for more details.

    >>> p = scipy.stats.poisson
    >>> g = gt.random_graph(20000, lambda: (sample_k(19), sample_k(19)),
557
558
    ...                     lambda a, b: (p.pmf(a[0], b[1]) * p.pmf(a[1], 20 - b[0])),
    ...                     mix_time=100)
559
560
561
562
    >>> figure(figsize=(6,3))
    <...>
    >>> axes([0.1,0.15,0.6,0.8])
    <...>
563
    >>> corr = gt.avg_neighbour_corr(g, "in", "out")
564
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
565
    ...          label=r"$\left<\text{o}\right>$ vs i")
566
567
    (...)
    >>> corr = gt.avg_neighbour_corr(g, "out", "in")
568
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
569
    ...          label=r"$\left<\text{i}\right>$ vs o")
570
571
572
    (...)
    >>> gt.random_rewire(g, "correlated")
    >>> corr = gt.avg_neighbour_corr(g, "in", "out")
573
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
574
    ...          label=r"$\left<\text{o}\right>$ vs i, corr.")
575
576
    (...)
    >>> corr = gt.avg_neighbour_corr(g, "out", "in")
577
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
578
    ...          label=r"$\left<\text{i}\right>$ vs o, corr.")
579
580
581
    (...)
    >>> gt.random_rewire(g, "uncorrelated")
    >>> corr = gt.avg_neighbour_corr(g, "in", "out")
582
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
583
    ...          label=r"$\left<\text{o}\right>$ vs i, uncorr.")
584
585
    (...)
    >>> corr = gt.avg_neighbour_corr(g, "out", "in")
586
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
587
    ...          label=r"$\left<\text{i}\right>$ vs o, uncorr.")
588
    (...)
589
    >>> legend(loc=(1.05,0.45))
590
591
592
593
594
    <...>
    >>> xlabel("source degree")
    <...>
    >>> ylabel("average target degree")
    <...>
595
    >>> savefig("shuffled-deg-corr-dir.pdf")
596

597
    .. figure:: shuffled-deg-corr-dir.*
598
599
600
601
602
603
        :align: center

        Average degree correlations for the different shuffled and non-shuffled
        directed graphs. The shuffled graph with correlations displays exactly
        the same correlation as the original graph.

604
605
606
607
608
609
610
611
612
613
614
    References
    ----------
    .. [metropolis-equations-1953]  Metropolis, N.; Rosenbluth, A.W.;
       Rosenbluth, M.N.; Teller, A.H.; Teller, E. "Equations of State
       Calculations by Fast Computing Machines". Journal of Chemical Physics 21
       (6): 1087–1092 (1953). :doi:`10.1063/1.1699114`
    .. [hastings-monte-carlo-1970] Hastings, W.K. "Monte Carlo Sampling Methods
       Using Markov Chains and Their Applications". Biometrika 57 (1): 97–109 (1970).
       :doi:`10.1093/biomet/57.1.97`

    """
615
    seed = numpy.random.randint(0, sys.maxint)
616

Tiago Peixoto's avatar
Tiago Peixoto committed
617
618
619
620
621
622
623
624
625
626
627
628
629
    if not parallel_edges:
        p = label_parallel_edges(g)
        if p.a.max() != 0:
            raise ValueError("Parallel edge detected. Can't rewire " +
                             "graph without parallel edges if it " +
                             "already contains parallel edges!")
    if not self_loops:
        l = label_self_loops(g)
        if l.a.max() != 0:
            raise ValueError("Self-loop detected. Can't rewire graph " +
                             "without self-loops if it already contains" +
                             " self-loops!")

630
    if (deg_corr is not None and not g.is_directed()) and blockmodel is None:
Tiago Peixoto's avatar
Tiago Peixoto committed
631
        corr = lambda i, j: deg_corr(i[1], j[1])
632
633
634
    else:
        corr = deg_corr

635
636
    if strat != "probabilistic":
        g = GraphView(g, reversed=False)
637
638
    elif blockmodel is not None:
        strat = "blockmodel"
639
640
641
    pcount = libgraph_tool_generation.random_rewire(g._Graph__graph, strat,
                                                    n_iter, not edge_sweep,
                                                    self_loops, parallel_edges,
642
643
                                                    corr, _prop("v", g, blockmodel),
                                                    seed, verbose)
644
645
    if ret_fail:
        return pcount
Tiago Peixoto's avatar
Tiago Peixoto committed
646

Tiago Peixoto's avatar
Tiago Peixoto committed
647

Tiago Peixoto's avatar
Tiago Peixoto committed
648
def predecessor_tree(g, pred_map):
Tiago Peixoto's avatar
Tiago Peixoto committed
649
    """Return a graph from a list of predecessors given by the ``pred_map`` vertex property."""
Tiago Peixoto's avatar
Tiago Peixoto committed
650
651
652
653
654
655
656

    _check_prop_scalar(pred_map, "pred_map")
    pg = Graph()
    libgraph_tool_generation.predecessor_graph(g._Graph__graph,
                                               pg._Graph__graph,
                                               _prop("v", g, pred_map))
    return pg
657

Tiago Peixoto's avatar
Tiago Peixoto committed
658

659
def line_graph(g):
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
    """Return the line graph of the given graph `g`.

    Notes
    -----
    Given an undirected graph G, its line graph L(G) is a graph such that

        * each vertex of L(G) represents an edge of G; and
        * two vertices of L(G) are adjacent if and only if their corresponding
          edges share a common endpoint ("are adjacent") in G.

    For a directed graph, the second criterion becomes:

       * Two vertices representing directed edges from u to v and from w to x in
         G are connected by an edge from uv to wx in the line digraph when v =
         w.

    References
    ----------
    .. [line-wiki] http://en.wikipedia.org/wiki/Line_graph
    """
680
681
682
683
684
685
686
687
    lg = Graph(directed=g.is_directed())

    vertex_map = lg.new_vertex_property("int64_t")

    libgraph_tool_generation.line_graph(g._Graph__graph,
                                        lg._Graph__graph,
                                        _prop("v", lg, vertex_map))
    return lg, vertex_map
Tiago Peixoto's avatar
Tiago Peixoto committed
688

Tiago Peixoto's avatar
Tiago Peixoto committed
689
690

def graph_union(g1, g2, props=None, include=False):
691
692
693
694
695
696
697
698
699
    """Return the union of graphs g1 and g2, composed of all edges and vertices
    of g1 and g2, without overlap.

    Parameters
    ----------
    g1 : :class:`~graph_tool.Graph`
       First graph in the union.
    g2 : :class:`~graph_tool.Graph`
       Second graph in the union.
700
    props : list of tuples of :class:`~graph_tool.PropertyMap` (optional, default: ``[]``)
701
702
703
704
       Each element in this list must be a tuple of two PropertyMap objects. The
       first element must be a property of `g1`, and the second of `g2`. The
       values of the property maps are propagated into the union graph, and
       returned.
705
    include : bool (optional, default: ``False``)
706
707
708
709
710
711
712
713
714
715
       If true, graph `g2` is inserted into `g1` which is modified. If false, a
       new graph is created, and both graphs remain unmodified.

    Returns
    -------
    ug : :class:`~graph_tool.Graph`
        The union graph
    props : list of :class:`~graph_tool.PropertyMap` objects
        List of propagated properties.  This is only returned if `props` is not
        empty.
716
717
718
719
720
721
722
723
724

    Examples
    --------

    >>> from numpy.random import random, seed
    >>> seed(42)
    >>> g = gt.triangulation(random((300,2)))[0]
    >>> ug = gt.graph_union(g, g)
    >>> uug = gt.graph_union(g, ug)
725
    >>> gt.graph_draw(g, layout="arf", size=(8,8), output="graph_original.pdf")
726
    <...>
727
    >>> gt.graph_draw(ug, layout="arf", size=(8,8), output="graph_union.pdf")
728
    <...>
729
    >>> gt.graph_draw(uug, layout="arf", size=(8,8), output="graph_union2.pdf")
730
731
    <...>

732
733
734
    .. image:: graph_original.*
    .. image:: graph_union.*
    .. image:: graph_union2.*
735

736
    """
Tiago Peixoto's avatar
Tiago Peixoto committed
737
738
    if props == None:
        props = []
Tiago Peixoto's avatar
Tiago Peixoto committed
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
    if not include:
        g1 = Graph(g1)
    g1.stash_filter(directed=True)
    g1.set_directed(True)
    g2.stash_filter(directed=True)
    g2.set_directed(True)
    n_props = []

    try:
        vmap, emap = libgraph_tool_generation.graph_union(g1._Graph__graph,
                                                          g2._Graph__graph)
        for p in props:
            p1, p2 = p
            if not include:
                p1 = g1.copy_property(p1)
            if p2.value_type() != p1.value_type():
                p2 = g2.copy_property(p2, value_type=p1.value_type())
            if p1.key_type() == 'v':
                libgraph_tool_generation.\
                      vertex_property_union(g1._Graph__graph, g2._Graph__graph,
                                            vmap, emap,
                                            _prop(p1.key_type(), g1, p1),
                                            _prop(p2.key_type(), g2, p2))
            else:
                libgraph_tool_generation.\
                      edge_property_union(g1._Graph__graph, g2._Graph__graph,
                                          vmap, emap,
                                          _prop(p1.key_type(), g1, p1),
                                          _prop(p2.key_type(), g2, p2))
            n_props.append(p1)
    finally:
        g1.pop_filter(directed=True)
        g2.pop_filter(directed=True)

    if len(n_props) > 0:
        return g1, n_props
    else:
        return g1
777

Tiago Peixoto's avatar
Tiago Peixoto committed
778
779

@_limit_args({"type": ["simple", "delaunay"]})
780
def triangulation(points, type="simple", periodic=False):
781
782
783
784
785
786
787
788
    r"""
    Generate a 2D or 3D triangulation graph from a given point set.

    Parameters
    ----------
    points : :class:`~numpy.ndarray`
        Point set for the triangulation. It may be either a N x d array, where N
        is the number of points, and d is the space dimension (either 2 or 3).
789
    type : string (optional, default: ``'simple'``)
790
        Type of triangulation. May be either 'simple' or 'delaunay'.
791
792
793
    periodic : bool (optional, default: ``False``)
        If ``True``, periodic boundary conditions will be used. This is
        parameter is valid only for type="delaunay", and is otherwise ignored.
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808

    Returns
    -------
    triangulation_graph : :class:`~graph_tool.Graph`
        The generated graph.
    pos : :class:`~graph_tool.PropertyMap`
        Vertex property map with the Cartesian coordinates.

    See Also
    --------
    random_graph: random graph generation

    Notes
    -----

Tiago Peixoto's avatar
Tiago Peixoto committed
809
    A triangulation [cgal-triang]_ is a division of the convex hull of a point
810
    set into triangles, using only that set as triangle vertices.
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829

    In simple triangulations (`type="simple"`), the insertion of a point is done
    by locating a face that contains the point, and splitting this face into
    three new faces (the order of insertion is therefore important). If the
    point falls outside the convex hull, the triangulation is restored by
    flips. Apart from the location, insertion takes a time O(1). This bound is
    only an amortized bound for points located outside the convex hull.

    Delaunay triangulations (`type="delaunay"`) have the specific empty sphere
    property, that is, the circumscribing sphere of each cell of such a
    triangulation does not contain any other vertex of the triangulation in its
    interior. These triangulations are uniquely defined except in degenerate
    cases where five points are co-spherical. Note however that the CGAL
    implementation computes a unique triangulation even in these cases.

    Examples
    --------
    >>> from numpy.random import seed, random
    >>> seed(42)
830
    >>> points = random((500, 2)) * 4
831
    >>> g, pos = gt.triangulation(points)
832
833
834
835
836
837
838
839
    >>> weight = g.new_edge_property("double") # Edge weights corresponding to
    ...                                        # Euclidean distances
    >>> for e in g.edges():
    ...    weight[e] = sqrt(sum((array(pos[e.source()]) -
    ...                          array(pos[e.target()]))**2))
    >>> b = gt.betweenness(g, weight=weight)
    >>> b[1].a *= 100
    >>> gt.graph_draw(g, pos=pos, pin=True, size=(8,8), vsize=0.07, vcolor=b[0],
840
    ...               eprops={"penwidth":b[1]}, output="triang.pdf")
841
842
    <...>
    >>> g, pos = gt.triangulation(points, type="delaunay")
843
844
845
846
847
848
849
    >>> weight = g.new_edge_property("double")
    >>> for e in g.edges():
    ...    weight[e] = sqrt(sum((array(pos[e.source()]) -
    ...                          array(pos[e.target()]))**2))
    >>> b = gt.betweenness(g, weight=weight)
    >>> b[1].a *= 120
    >>> gt.graph_draw(g, pos=pos, pin=True, size=(8,8), vsize=0.07, vcolor=b[0],
850
    ...               eprops={"penwidth":b[1]}, output="triang-delaunay.pdf")
851
852
853
854
    <...>

    2D triangulation of random points:

855
856
    .. image:: triang.*
    .. image:: triang-delaunay.*
857

858
859
860
    *Left:* Simple triangulation. *Right:* Delaunay triangulation. The vertex
    colors and the edge thickness correspond to the weighted betweenness
    centrality.
861
862
863

    References
    ----------
Tiago Peixoto's avatar
Tiago Peixoto committed
864
    .. [cgal-triang] http://www.cgal.org/Manual/last/doc_html/cgal_manual/Triangulation_3/Chapter_main.html
865
866
867

    """

Tiago Peixoto's avatar
Tiago Peixoto committed
868
    if points.shape[1] not in [2, 3]:
869
870
871
872
873
874
875
876
877
878
        raise ValueError("points array must have shape N x d, with d either 2 or 3.")
    # copy points to ensure continuity and correct data type
    points = numpy.array(points, dtype='float64')
    if points.shape[1] == 2:
        npoints = numpy.zeros((points.shape[0], 3))
        npoints[:,:2] = points
        points = npoints
    g = Graph(directed=False)
    pos = g.new_vertex_property("vector<double>")
    libgraph_tool_generation.triangulation(g._Graph__graph, points,
879
                                           _prop("v", g, pos), type, periodic)
880
    return g, pos
881
882
883
884
885
886
887
888
889
890


def lattice(shape, periodic=False):
    r"""
    Generate a N-dimensional square lattice.

    Parameters
    ----------
    shape : list or :class:`~numpy.ndarray`
        List of sizes in each dimension.
891
    periodic : bool (optional, default: ``False``)
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
        If ``True``, periodic boundary conditions will be used.

    Returns
    -------
    lattice_graph : :class:`~graph_tool.Graph`
        The generated graph.

    See Also
    --------
    triangulation: 2D or 3D triangulation
    random_graph: random graph generation

    Examples
    --------
    >>> g = gt.lattice([10,10])
907
    >>> gt.graph_draw(g, size=(8,8), output="lattice.pdf")
908
909
    <...>
    >>> g = gt.lattice([10,20], periodic=True)
910
    >>> gt.graph_draw(g, size=(8,8), output="lattice_periodic.pdf")
911
912
    <...>
    >>> g = gt.lattice([10,10,10])
913
    >>> gt.graph_draw(g, size=(8,8), output="lattice_3d.pdf")
914
915
    <...>

916
917
918
    .. image:: lattice.*
    .. image:: lattice_periodic.*
    .. image:: lattice_3d.*
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

    *Left:* 10x10 2D lattice. *Middle:* 10x20 2D periodic lattice (torus).
    *Right:* 10x10x10 3D lattice.

    References
    ----------
    .. [lattice] http://en.wikipedia.org/wiki/Square_lattice

    """

    g = Graph(directed=False)
    libgraph_tool_generation.lattice(g._Graph__graph, shape, periodic)
    return g


def geometric_graph(points, radius, ranges=None):
    r"""
    Generate a geometric network form a set of N-dimensional points.

    Parameters
    ----------
    points : list or :class:`~numpy.ndarray`
        List of points. This must be a two-dimensional array, where the rows are
        coordinates in a N-dimensional space.
    radius : float
        Pairs of points with an euclidean distance lower than this parameters
        will be connected.
946
    ranges : list or :class:`~numpy.ndarray` (optional, default: ``None``)
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
        If provided, periodic boundary conditions will be assumed, and the
        values of this parameter it will be used as the ranges in all
        dimensions. It must be a two-dimensional array, where each row will
        cointain the lower and upper bound of each dimension.

    Returns
    -------
    geometric_graph : :class:`~graph_tool.Graph`
        The generated graph.
    pos : :class:`~graph_tool.PropertyMap`
        A vertex property map with the position of each vertex.

    Notes
    -----
    A geometric graph [geometric-graph]_ is generated by connecting points
    embedded in a N-dimensional euclidean space which are at a distance equal to
    or smaller than a given radius.

    See Also
    --------
    triangulation: 2D or 3D triangulation
    random_graph: random graph generation
    lattice : N-dimensional square lattice

    Examples
    --------
    >>> from numpy.random import seed, random
    >>> seed(42)
    >>> points = random((500, 2)) * 4
    >>> g, pos = gt.geometric_graph(points, 0.3)
977
    >>> gt.graph_draw(g, pos=pos, pin=True, size=(8,8), output="geometric.pdf")
978
979
    <...>
    >>> g, pos = gt.geometric_graph(points, 0.3, [(0,4), (0,4)])
980
    >>> gt.graph_draw(g, size=(8,8), output="geometric_periodic.pdf")
981
982
    <...>

983
984
    .. image:: geometric.*
    .. image:: geometric_periodic.*
985
986
987
988
989
990
991

    *Left:* Geometric network with random points. *Right:* Same network, but
     with periodic boundary conditions.

    References
    ----------
    .. [geometric-graph] Jesper Dall and Michael Christensen, "Random geometric
Tiago Peixoto's avatar
Tiago Peixoto committed
992
       graphs", Phys. Rev. E 66, 016121 (2002), :doi:`10.1103/PhysRevE.66.016121`
993
994
995
996
997
998
999
1000

    """

    g = Graph(directed=False)
    pos = g.new_vertex_property("vector<double>")
    if type(points) != numpy.ndarray:
        points = numpy.array(points)
    if len(points.shape) < 2: