__init__.py 44.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
                 parallel_edges=False, self_loops=False, blockmodel=None,
61
                 block_type="int", degree_block=False,
62
                 random=True, mix_time=10, verbose=False):
Tiago Peixoto's avatar
Tiago Peixoto committed
63
64
65
66
67
68
69
70
71
72
73
74
75
    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.
76

77
78
79
80
81
82
        Optionally, you can also pass a function which receives one or two
        arguments: If ``blockmodel == None``, the single argument passed will
        be the index of the vertex which will receive the degree.
        If ``blockmodel != None``, the first value passed will be the vertex
        index, and the second will be the block value of the vertex.
        
83
    deg_corr : function (optional, default: ``None``)
Tiago Peixoto's avatar
Tiago Peixoto committed
84
        A function which gives the degree correlation of the graph. It should be
Tiago Peixoto's avatar
Tiago Peixoto committed
85
86
87
88
89
        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.
90
91
92
93

        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
94
        Whether the generated graph should be directed.
95
96
97
98
99
100
101
102
103
104
105
106
107
    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
108
109
110
111
112
113
114
115
        a type compatible with the ``block_type`` parameter.
    block_type : string (optional, default: ``"int"``)
        Value type of block labels. Valid only if ``blockmodel != None``.
    degree_block : bool (optional, default: ``False``)
        If ``True``, the degree of each vertex will be appended to block labels
        when constructing the blockmodel, such that the resulting block type
        will be a pair :math:`(r, k)`, where :math:`r` is the original block
        label.
116
117
118
119
    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``)
120
121
122
        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"``.
123
124
    verbose : bool (optional, default: ``False``)
        If ``True``, verbose information is displayed.
Tiago Peixoto's avatar
Tiago Peixoto committed
125
126
127

    Returns
    -------
128
    random_graph : :class:`~graph_tool.Graph`
Tiago Peixoto's avatar
Tiago Peixoto committed
129
        The generated graph.
130
131
132
    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
133
134
135
136
137
138
139

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

    Notes
    -----
Tiago Peixoto's avatar
Tiago Peixoto committed
140
141
142
    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
143
144
    the :func:`~graph_tool.generation.random_rewire` function, with the
    ``mix_time`` parameter passed as ``n_iter``.
Tiago Peixoto's avatar
Tiago Peixoto committed
145

146
    The complexity is :math:`O(V + E)` if parallel edges are allowed, and
147
    :math:`O(V + E \times\text{mix-time})` if parallel edges are not allowed.
148
149
150
151
152
153
154
155
156
157
158
159
160


    .. 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
161
162
163
        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
164
165
166
167
168

    Examples
    --------
    >>> from numpy.random import randint, random, seed, poisson
    >>> from pylab import *
169
    >>> seed(43)
Tiago Peixoto's avatar
Tiago Peixoto committed
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190

    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),
191
192
    ...                     lambda i, k: 1.0 / (1 + abs(i - k)), directed=False,
    ...                     mix_time=100)
Tiago Peixoto's avatar
Tiago Peixoto committed
193
    >>> gt.scalar_assortativity(g, "out")
Tiago Peixoto's avatar
Tiago Peixoto committed
194
    (0.6435658697163692, 0.010420519538259333)
Tiago Peixoto's avatar
Tiago Peixoto committed
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221

    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()
    <...>
222
    >>> xlabel("in-degree")
Tiago Peixoto's avatar
Tiago Peixoto committed
223
    <...>
224
    >>> ylabel("out-degree")
Tiago Peixoto's avatar
Tiago Peixoto committed
225
    <...>
226
    >>> savefig("combined-deg-hist.pdf")
Tiago Peixoto's avatar
Tiago Peixoto committed
227

228
    .. figure:: combined-deg-hist.*
Tiago Peixoto's avatar
Tiago Peixoto committed
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
        :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)),
247
248
249
    ...                     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
250
251
252

    Lets plot the average degree correlations to check.

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

281
    .. figure:: deg-corr-dir.*
Tiago Peixoto's avatar
Tiago Peixoto committed
282
283
284
        :align: center

        Average nearest neighbour correlations.
285
286
287
288
289


    **Blockmodels**


290
291
292
    The following example shows how a stochastic blockmodel
    [holland-stochastic-1983]_ [karrer-stochastic-2011]_ can be generated. We
    will consider a system of 10 blocks, which form communities. The connection
293
294
295
296
297
298
299
300
301
302
303
304
305
    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)
306
    >>> gt.graph_draw(g, vcolor=bm, layout="sfdp", output="blockmodel.pdf")
307
308
    <...>

309
    .. figure:: blockmodel.*
310
311
312
313
314
315
316
317
318
319
320
321
322
323
        :align: center

        Simple blockmodel with 10 blocks.


    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`
324
325
326
327
328
329
    .. [holland-stochastic-1983] Paul W. Holland, Kathryn Blackmond Laskey, and
       Samuel Leinhardt, "Stochastic blockmodels: First steps," Social Networks
       5, no. 2: 109-13 (1983) :doi:`10.1016/0378-8733(83)90021-7`
    .. [karrer-stochastic-2011] Brian Karrer and M. E. J. Newman, "Stochastic
       blockmodels and community structure in networks," Physical Review E 83,
       no. 1: 016107 (2011) :doi:`10.1103/PhysRevE.83.016107` :arxiv:`1008.3926`
Tiago Peixoto's avatar
Tiago Peixoto committed
330
    """
331

332
    seed = numpy.random.randint(0, sys.maxint)
333
334
335
336
337
    g = Graph()
    if deg_corr == None:
        uncorrelated = True
    else:
        uncorrelated = False
338
339
340

    if (type(blockmodel) is types.FunctionType or
        type(blockmodel) is types.LambdaType):
341
342
        btype = block_type
        bm = []
343
344
        if len(inspect.getargspec(blockmodel)[0]) == 0:
            for i in xrange(N):
345
                bm.append(blockmodel())
346
347
        else:
            for i in xrange(N):
348
349
                bm.append(blockmodel(i))
        blockmodel = bm
Tiago Peixoto's avatar
Tiago Peixoto committed
350
    elif blockmodel is not None:
351
        btype = _gt_type(blockmodel[0])
352
353
354

    if len(inspect.getargspec(deg_sampler)[0]) > 0:
        if blockmodel is not None:
355
            sampler = lambda i: deg_sampler(i, blockmodel[i])
356
        else:
Tiago Peixoto's avatar
Tiago Peixoto committed
357
            sampler = deg_sampler
358
359
360
361
    else:
        sampler = lambda i: deg_sampler()

    libgraph_tool_generation.gen_graph(g._Graph__graph, N, sampler,
362
363
364
                                       uncorrelated, not parallel_edges,
                                       not self_loops, not directed,
                                       seed, verbose, True)
365
366
    g.set_directed(directed)

367
368
369
370
371
372
373
374
375
376
377
378
    if degree_block:
        if btype in ["object", "string"] or "vector" in btype:
            btype = "object"
        elif btype in ["int", "int32_t", "bool"]:
            btype = "vector<int32_t>"
        elif btype in ["long", "int64_t"]:
            btype = "vector<int64_t>"
        elif btype in ["double"]:
            btype = "vector<double>"
        elif btype in ["long double"]:
            btype = "vector<long double>"

379
380
381
382
    if blockmodel is not None:
        bm = g.new_vertex_property(btype)
        if btype in ["object", "string"] or "vector" in btype:
            for v in g.vertices():
383
384
385
386
387
388
389
390
                if not degree_block:
                    bm[v] = blockmodel[int(v)]
                else:
                    if g.is_directed():
                        bm[v] = (blockmodel[int(v)], v.in_degree(),
                                 v.out_degree())
                    else:
                        bm[v] = (blockmodel[int(v)], v.out_degree())
391
392
393
394
395
396
397
398
399
        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
400

401
    if parallel_edges and self_loops and deg_corr is None:
402
        mix_time = 1
Tiago Peixoto's avatar
Tiago Peixoto committed
403
    if random:
404
405
        if deg_corr is not None:
            random_rewire(g, strat="probabilistic", n_iter=mix_time,
Tiago Peixoto's avatar
Tiago Peixoto committed
406
                          parallel_edges=parallel_edges, deg_corr=deg_corr,
407
408
                          self_loops=self_loops, blockmodel=bm,
                          verbose=verbose)
409
410
411
412
        else:
            random_rewire(g, parallel_edges=parallel_edges, n_iter=mix_time,
                          self_loops=self_loops, verbose=verbose)

413
414
415
416
    if bm is None:
        return g
    else:
        return g, bm
417

Tiago Peixoto's avatar
Tiago Peixoto committed
418

419
420
@_limit_args({"strat": ["erdos", "correlated", "uncorrelated",
                        "probabilistic"]})
421
422
def random_rewire(g, strat="uncorrelated", n_iter=1, edge_sweep=True,
                  parallel_edges=False, self_loops=False, deg_corr=None,
423
                  blockmodel=None, ret_fail=False, verbose=False):
424
    r"""
425
426
427
428
429
430
431
432
    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``.
433
434
435

    Parameters
    ----------
436
    g : :class:`~graph_tool.Graph`
437
        Graph to be shuffled. The graph will be modified.
438
439
440
441
    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
442
        new source and target of each edge both have the same in and out-degree.
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
        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``)
459
460
461
462
463
464
        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,
465
        unless ``strat == "probabilistic"``.
466
467
468
469
470
471
472

        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.
473
474
475
476
477
478
479
480
481
482
483
    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).
484
485
486
487
488
489
490

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

    Notes
    -----
Tiago Peixoto's avatar
Tiago Peixoto committed
491
    This algorithm iterates through all the edges in the network and tries to
492
    swap its target or source with the target or source of another edge.
Tiago Peixoto's avatar
Tiago Peixoto committed
493
494

    .. note::
495

496
497
498
499
500
501
502
503
504
505
        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
506
507
508
        the Markov chain is implemented using the Metropolis-Hastings
        [metropolis-equations-1953]_ [hastings-monte-carlo-1970]_
        acceptance/rejection algorithm.
509

Tiago Peixoto's avatar
Tiago Peixoto committed
510

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

515
516
517
518
519
    Examples
    --------

    Some small graphs for visualization.

520
    >>> from numpy.random import random, seed
521
    >>> from pylab import *
522
    >>> seed(43)
523
    >>> g, pos = gt.triangulation(random((1000,2)))
524
    >>> gt.graph_draw(g, layout="arf", output="rewire_orig.pdf", size=(6,6))
525
    <...>
526
    >>> gt.random_rewire(g, "correlated")
527
    >>> gt.graph_draw(g, layout="arf", output="rewire_corr.pdf", size=(6,6))
528
    <...>
529
    >>> gt.random_rewire(g)
530
    >>> gt.graph_draw(g, layout="arf", output="rewire_uncorr.pdf", size=(6,6))
531
    <...>
532
    >>> gt.random_rewire(g, "erdos")
533
    >>> gt.graph_draw(g, layout="arf", output="rewire_erdos.pdf", size=(6,6))
534
    <...>
535

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

538
539
540
541
    .. image:: rewire_orig.*
    .. image:: rewire_corr.*
    .. image:: rewire_uncorr.*
    .. image:: rewire_erdos.*
542

543
544
    *From left to right:* Original graph --- Shuffled graph, with degree
    correlations --- Shuffled graph, without degree correlations --- Shuffled graph,
545
    with random degrees.
546
547
548

    We can try some larger graphs to get better statistics.

549
550
    >>> figure()
    <...>
551
    >>> g = gt.random_graph(30000, lambda: sample_k(20),
552
553
    ...                     lambda i, j: exp(abs(i-j)), directed=False,
    ...                     mix_time=100)
554
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
555
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-", label="original")
556
557
558
    (...)
    >>> gt.random_rewire(g, "correlated")
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
559
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="*", label="correlated")
560
561
562
    (...)
    >>> gt.random_rewire(g)
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
563
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-", label="uncorrelated")
564
    (...)
565
566
    >>> gt.random_rewire(g, "erdos")
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
567
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-", label="Erdos")
568
    (...)
569
570
571
572
573
574
    >>> xlabel("$k$")
    <...>
    >>> ylabel(r"$\left<k_{nn}\right>$")
    <...>
    >>> legend(loc="best")
    <...>
575
    >>> savefig("shuffled-stats.pdf")
576

577
    .. figure:: shuffled-stats.*
578
579
580
581
582
583
584
585
586
587
588
        :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)),
589
590
    ...                     lambda a, b: (p.pmf(a[0], b[1]) * p.pmf(a[1], 20 - b[0])),
    ...                     mix_time=100)
591
592
593
594
    >>> figure(figsize=(6,3))
    <...>
    >>> axes([0.1,0.15,0.6,0.8])
    <...>
595
    >>> corr = gt.avg_neighbour_corr(g, "in", "out")
596
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
597
    ...          label=r"$\left<\text{o}\right>$ vs i")
598
599
    (...)
    >>> corr = gt.avg_neighbour_corr(g, "out", "in")
600
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
601
    ...          label=r"$\left<\text{i}\right>$ vs o")
602
603
604
    (...)
    >>> gt.random_rewire(g, "correlated")
    >>> corr = gt.avg_neighbour_corr(g, "in", "out")
605
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
606
    ...          label=r"$\left<\text{o}\right>$ vs i, corr.")
607
608
    (...)
    >>> corr = gt.avg_neighbour_corr(g, "out", "in")
609
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
610
    ...          label=r"$\left<\text{i}\right>$ vs o, corr.")
611
612
613
    (...)
    >>> gt.random_rewire(g, "uncorrelated")
    >>> corr = gt.avg_neighbour_corr(g, "in", "out")
614
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
615
    ...          label=r"$\left<\text{o}\right>$ vs i, uncorr.")
616
617
    (...)
    >>> corr = gt.avg_neighbour_corr(g, "out", "in")
618
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
619
    ...          label=r"$\left<\text{i}\right>$ vs o, uncorr.")
620
    (...)
621
    >>> legend(loc=(1.05,0.45))
622
623
624
625
626
    <...>
    >>> xlabel("source degree")
    <...>
    >>> ylabel("average target degree")
    <...>
627
    >>> savefig("shuffled-deg-corr-dir.pdf")
628

629
    .. figure:: shuffled-deg-corr-dir.*
630
631
632
633
634
635
        :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.

636
637
638
639
640
641
642
643
644
    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`
645
646
647
648
649
650
    .. [holland-stochastic-1983] Paul W. Holland, Kathryn Blackmond Laskey, and
       Samuel Leinhardt, "Stochastic blockmodels: First steps," Social Networks
       5, no. 2: 109-13 (1983) :doi:`10.1016/0378-8733(83)90021-7`
    .. [karrer-stochastic-2011] Brian Karrer and M. E. J. Newman, "Stochastic
       blockmodels and community structure in networks," Physical Review E 83,
       no. 1: 016107 (2011) :doi:`10.1103/PhysRevE.83.016107` :arxiv:`1008.3926`
651
652

    """
653
    seed = numpy.random.randint(0, sys.maxint)
654

Tiago Peixoto's avatar
Tiago Peixoto committed
655
656
657
658
659
660
661
662
663
664
665
666
667
    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!")

668
    if (deg_corr is not None and not g.is_directed()) and blockmodel is None:
Tiago Peixoto's avatar
Tiago Peixoto committed
669
        corr = lambda i, j: deg_corr(i[1], j[1])
670
671
672
    else:
        corr = deg_corr

673
674
    if strat != "probabilistic":
        g = GraphView(g, reversed=False)
675
676
    elif blockmodel is not None:
        strat = "blockmodel"
677
678
679
    pcount = libgraph_tool_generation.random_rewire(g._Graph__graph, strat,
                                                    n_iter, not edge_sweep,
                                                    self_loops, parallel_edges,
680
681
                                                    corr, _prop("v", g, blockmodel),
                                                    seed, verbose)
682
683
    if ret_fail:
        return pcount
Tiago Peixoto's avatar
Tiago Peixoto committed
684

Tiago Peixoto's avatar
Tiago Peixoto committed
685

Tiago Peixoto's avatar
Tiago Peixoto committed
686
def predecessor_tree(g, pred_map):
Tiago Peixoto's avatar
Tiago Peixoto committed
687
    """Return a graph from a list of predecessors given by the ``pred_map`` vertex property."""
Tiago Peixoto's avatar
Tiago Peixoto committed
688
689
690
691
692
693
694

    _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
695

Tiago Peixoto's avatar
Tiago Peixoto committed
696

697
def line_graph(g):
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
    """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
    """
718
719
720
721
722
723
724
725
    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
726

Tiago Peixoto's avatar
Tiago Peixoto committed
727
728

def graph_union(g1, g2, props=None, include=False):
729
730
731
732
733
734
735
736
737
    """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.
738
    props : list of tuples of :class:`~graph_tool.PropertyMap` (optional, default: ``[]``)
739
740
741
742
       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.
743
    include : bool (optional, default: ``False``)
744
745
746
747
748
749
750
751
752
753
       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.
754
755
756
757
758
759
760
761
762

    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)
763
    >>> gt.graph_draw(g, layout="arf", size=(8,8), output="graph_original.pdf")
764
    <...>
765
    >>> gt.graph_draw(ug, layout="arf", size=(8,8), output="graph_union.pdf")
766
    <...>
767
    >>> gt.graph_draw(uug, layout="arf", size=(8,8), output="graph_union2.pdf")
768
769
    <...>

770
771
772
    .. image:: graph_original.*
    .. image:: graph_union.*
    .. image:: graph_union2.*
773

774
    """
Tiago Peixoto's avatar
Tiago Peixoto committed
775
776
    if props == None:
        props = []
Tiago Peixoto's avatar
Tiago Peixoto committed
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
    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
815

Tiago Peixoto's avatar
Tiago Peixoto committed
816
817

@_limit_args({"type": ["simple", "delaunay"]})
818
def triangulation(points, type="simple", periodic=False):
819
820
821
822
823
824
825
826
    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).
827
    type : string (optional, default: ``'simple'``)
828
        Type of triangulation. May be either 'simple' or 'delaunay'.
829
830
831
    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.
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846

    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
847
    A triangulation [cgal-triang]_ is a division of the convex hull of a point
848
    set into triangles, using only that set as triangle vertices.
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867

    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)
868
    >>> points = random((500, 2)) * 4
869
    >>> g, pos = gt.triangulation(points)
870
871
872
873
874
875
876
877
    >>> 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],
878
    ...               eprops={"penwidth":b[1]}, output="triang.pdf")
879
880
    <...>
    >>> g, pos = gt.triangulation(points, type="delaunay")
881
882
883
884
885
886
887
    >>> 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],
888
    ...               eprops={"penwidth":b[1]}, output="triang-delaunay.pdf")
889
890
891
892
    <...>

    2D triangulation of random points:

893
894
    .. image:: triang.*
    .. image:: triang-delaunay.*
895

896
897
898
    *Left:* Simple triangulation. *Right:* Delaunay triangulation. The vertex
    colors and the edge thickness correspond to the weighted betweenness
    centrality.
899
900
901

    References
    ----------
Tiago Peixoto's avatar
Tiago Peixoto committed
902
    .. [cgal-triang] http://www.cgal.org/Manual/last/doc_html/cgal_manual/Triangulation_3/Chapter_main.html
903
904
905

    """

Tiago Peixoto's avatar
Tiago Peixoto committed
906
    if points.shape[1] not in [2, 3]:
907
908
909
910
911
912
913
914
915
916
        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,
917
                                           _prop("v", g, pos), type, periodic)
918
    return g, pos
919
920
921
922
923
924
925
926
927
928


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.
929
    periodic : bool (optional, default: ``False``)
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
        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])
945
    >>> gt.graph_draw(g, size=(8,8), output="lattice.pdf")
946
947
    <...>
    >>> g = gt.lattice([10,20], periodic=True)
948
    >>> gt.graph_draw(g, size=(8,8), output="lattice_periodic.pdf")
949
950
    <...>
    >>> g = gt.lattice([10,10,10])
951
    >>> gt.graph_draw(g, size=(8,8), output="lattice_3d.pdf")
952
953
    <...>

954
955
956
    .. image:: lattice.*
    .. image:: lattice_periodic.*
    .. image:: lattice_3d.*
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

    *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.
984
    ranges : list or :class:`~numpy.ndarray` (optional, default: ``None``)
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
        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)
1015
    >>> gt.graph_draw(g, pos=pos, pin=True, size=(8,8), output="geometric.pdf")
1016
1017
    <...>
    >>> g, pos = gt.geometric_graph(points, 0.3, [(0,4), (0,4)])
1018
    >>> gt.graph_draw(g, size=(8,8), output="geometric_periodic.pdf")
1019
1020
    <...>

1021
1022
    .. image:: geometric.*
    .. image:: geometric_periodic.*
1023
1024
1025
1026
1027
1028
1029

    *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
1030
       graphs", Phys. Rev. E 66, 016121 (2002), :doi:`10.1103/PhysRevE.66.016121`
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053

    """

    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:
        raise ValueError("points list must be a two-dimensional array!")
    if ranges is not None:
        periodic = True
        if type(ranges) != numpy.ndarray:
            ranges = numpy.array(ranges, dtype="float")
        else:
            ranges = array(ranges, dtype="float")
    else:
        periodic = False
        ranges = ()

    libgraph_tool_generation.geometric(g._Graph__graph, points, float(radius),
                                       ranges, periodic,
                                       _prop("v", g, pos))
    return g, pos
1054
1055
1056
1057
1058
1059
1060
1061
1062


def price_network(N, m=1, c=None, gamma=1, directed=True, seed_graph=None):
    r"""A generalized version of Price's -- or Barabási-Albert if undirected -- preferential attachment network model.

    Parameters
    ----------
    N : int
        Size of the network.
1063
    m : int (optional, default: ``1``)
1064
        Out-degree of newly added vertices.
1065
    c : float (optional, default: ``1 if directed == True else 0``)
1066
1067
        Constant factor added to the probability of a vertex receiving an edge
        (see notes below).
1068
    gamma : float (optional, default: ``1``)
1069
        Preferential attachment power (see notes below).
1070
    directed : bool (optional, default: ``True``)
1071
1072
        If ``True``, a Price network is generated. If ``False``, a
        Barabási-Albert network is generated.
1073
    seed_graph : :class:`~graph_tool.Graph` (optional, default: ``None``)
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
        If provided, this graph will be used as the starting point of the
        algorithm.

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

    Notes
    -----

    The (generalized) [price]_ network is either a directed or undirected graph
    (the latter is called a Barabási-Albert network), generated dynamically by
    at each step adding a new vertex, and connecting it to :math:`m` other
1088
    vertices, chosen with probability :math:`\pi` defined as:
1089
1090
1091

    .. math::

1092
        \pi \propto k^\gamma + c
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110

    where :math:`k` is the in-degree of the vertex (or simply the degree in the
    undirected case). If :math:`\gamma=1`, the tail of resulting in-degree
    distribution of the directed case is given by

    .. math::

        P_{k_\text{in}} \sim k_\text{in}^{-(2 + c/m)},

    or for the undirected case

    .. math::

        P_{k} \sim k^{-(3 + c/m)}.

    However, if :math:`\gamma \ne 1`, the in-degree distribution is not
    scale-free (see [dorogovtsev-evolution]_ for details).

1111
1112
1113
1114
1115
1116
1117
    Note that if `seed_graph` is not given, the algorithm will *always* start
    with one node if :math:`c > 0`, or with two nodes with a link between them
    otherwise. If :math:`m > 1`, the degree of the newly added vertices will be
    vary dynamically as :math:`m'(t) = \min(m, N(t))`, where :math:`N(t)` is the
    number of vertices added so far. If this behaviour is undesired, a proper
    seed graph with :math:`N \ge m` vertices must be provided.

1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
    This algorithm runs in :math:`O(N\log N)` time.

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

    Examples
    --------
    >>> from numpy.random import seed, random
    >>> seed(42)
    >>> g = gt.price_network(100000)
    >>> gt.graph_draw(g, layout="sfdp", size=(12,12), vcolor=g.vertex_index,
1133
    ...               output="price-network.pdf")
1134
1135
1136
    <...>
    >>> g = gt.price_network(100000, c=0.1)
    >>> gt.graph_draw(g, layout="sfdp", size=(12,12), vcolor=g.vertex_index,
1137
    ...               output="price-network-broader.pdf")
1138
1139
    <...>

1140
1141
    .. image:: price-network.*
    .. image:: price-network-broader.*
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152

    Price networks with :math:`N=10^5` nodes. **Left:** :math:`c=1`, **Right:**
    :math:`c=0.1`. The colors represent the order in which vertices were
    added.

    References
    ----------

    .. [yule] Yule, G. U. "A Mathematical Theory of Evolution, based on the
       Conclusions of Dr. J. C. Willis, F.R.S.". Philosophical Transactions of
       the Royal Society of London, Ser. B 213: 21–87, 1925,
Tiago Peixoto's avatar
Tiago Peixoto committed
1153
       :doi:`10.1098/rstb.1925.0002`
1154
1155
1156
    .. [price] Derek De Solla Price, "A general theory of bibliometric and other
       cumulative advantage processes", Journal of the American Society for
       Information Science, Volume 27, Issue 5, pages 292–306, September 1976,
Tiago Peixoto's avatar
Tiago Peixoto committed
1157
       :doi:`10.1002/asi.4630270505`
1158
    .. [barabasi-albert] Barabási, A.-L., and Albert, R., "Emergence of
Tiago Peixoto's avatar
Tiago Peixoto committed
1159
1160
       scaling in random networks", Science, 286, 509, 1999,
       :doi:`10.1126/science.286.5439.509`
1161
1162
    .. [dorogovtsev-evolution] S. N. Dorogovtsev and J. F. F. Mendes, "Evolution
       of networks", Advances in Physics, 2002, Vol. 51, No. 4, 1079-1187,
Tiago Peixoto's avatar
Tiago Peixoto committed
1163
       :doi:`10.1080/00018730110112519`
1164
1165
1166
1167
1168
1169
    """

    if c is None:
        c = 1 if directed else 0

    if seed_graph is None:
1170
1171
1172
        g = Graph(directed=directed)
        if c > 0:
            g.add_vertex()
1173
        else:
1174
1175
            g.add_vertex(2)
            g.add_edge(g.vertex(1), g.vertex(0))
1176
1177
1178
1179
1180
1181
        N -= g.num_vertices()
    else:
        g = seed_graph
    seed = numpy.random.randint(0, sys.maxint)
    libgraph_tool_generation.price(g._Graph__graph, N, gamma, c, m, seed)
    return g