__init__.py 52.4 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) 2006-2013 Tiago de Paula Peixoto <tiago@skewed.de>
7
8
9
10
11
12
13
14
15
16
17
18
#
# 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
Tiago Peixoto's avatar
Tiago Peixoto committed
19
# along with this program.  If not, see <http://www.gnu.org/licenses/>.s
20

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
Tiago Peixoto's avatar
Tiago Peixoto committed
40
   complete_graph
Tiago Peixoto's avatar
Tiago Peixoto committed
41
   circular_graph
42
43
44

Contents
++++++++
45
46
"""

47
48
from __future__ import division, absolute_import, print_function

Tiago Peixoto's avatar
Tiago Peixoto committed
49
from .. dl_import import dl_import
50
dl_import("from . import libgraph_tool_generation")
51

52
from .. import Graph, GraphView, _check_prop_scalar, _prop, _limit_args, _gt_type, _get_rng
Tiago Peixoto's avatar
Tiago Peixoto committed
53
from .. stats import label_parallel_edges, label_self_loops
54
55
import inspect
import types
56
import sys, numpy, numpy.random
57

Tiago Peixoto's avatar
Tiago Peixoto committed
58
__all__ = ["random_graph", "random_rewire", "predecessor_tree", "line_graph",
59
           "graph_union", "triangulation", "lattice", "geometric_graph",
Tiago Peixoto's avatar
Tiago Peixoto committed
60
           "price_network", "complete_graph", "circular_graph"]
61

Tiago Peixoto's avatar
Tiago Peixoto committed
62

63
def random_graph(N, deg_sampler, deg_corr=None, cache_probs=True, directed=True,
64
                 parallel_edges=False, self_loops=False, blockmodel=None,
65
                 block_type="int", degree_block=False,
66
                 random=True, mix_time=10, verbose=False):
Tiago Peixoto's avatar
Tiago Peixoto committed
67
68
69
70
71
72
73
74
75
76
77
78
79
    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.
80

81
82
83
84
85
        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.
86
    deg_corr : function (optional, default: ``None``)
Tiago Peixoto's avatar
Tiago Peixoto committed
87
        A function which gives the degree correlation of the graph. It should be
Tiago Peixoto's avatar
Tiago Peixoto committed
88
89
90
91
92
        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.
93
94
95

        If ``blockmodel != None``, the value passed to the function will be the
        block value of the respective vertices, not the in/out-degree pairs.
96
97
98
99
100
101
102
    cache_probs : bool (optional, default: ``True``)
        If ``True``, the probabilities returned by the ``deg_corr`` parameter
        will be cached internally. This is crucial for good performance, since
        in this case the supplied python function is called only a few times,
        and not at every attempted edge rewire move. However, in the case were
        the different parameter combinations to the probability function is very
        large, the memory requirements to keep the cache may be very large.
103
    directed : bool (optional, default: ``True``)
Tiago Peixoto's avatar
Tiago Peixoto committed
104
        Whether the generated graph should be directed.
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``)
110
111
112
        If supplied, the graph will be sampled from a stochastic blockmodel
        ensemble. If the value is a list or a :class:`~numpy.ndarray`, it must
        have ``len(blockmodel) == N``, and the values will define to which block
113
114
115
116
117
        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
118
119
120
121
122
123
124
125
        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.
126
127
128
129
    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``)
130
131
132
        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"``.
133
134
    verbose : bool (optional, default: ``False``)
        If ``True``, verbose information is displayed.
Tiago Peixoto's avatar
Tiago Peixoto committed
135
136
137

    Returns
    -------
138
    random_graph : :class:`~graph_tool.Graph`
Tiago Peixoto's avatar
Tiago Peixoto committed
139
        The generated graph.
140
141
142
    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
143
144
145
146
147
148
149

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

    Notes
    -----
Tiago Peixoto's avatar
Tiago Peixoto committed
150
151
152
    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
153
154
    the :func:`~graph_tool.generation.random_rewire` function, with the
    ``mix_time`` parameter passed as ``n_iter``.
Tiago Peixoto's avatar
Tiago Peixoto committed
155

156
    The complexity is :math:`O(V + E)` if parallel edges are allowed, and
157
    :math:`O(V + E \times\text{mix-time})` if parallel edges are not allowed.
158
159
160
161
162
163
164
165
166
167
168
169
170


    .. 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
171
172
173
        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
174
175
176

    Examples
    --------
177
178
179
180
181
182
183
184

    .. testcode::
       :hide:

       from numpy.random import randint, random, seed, poisson
       from pylab import *
       seed(43)
       gt.seed_rng(42)
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

    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),
206
207
    ...                     lambda i, k: 1.0 / (1 + abs(i - k)), directed=False,
    ...                     mix_time=100)
Tiago Peixoto's avatar
Tiago Peixoto committed
208
    >>> gt.scalar_assortativity(g, "out")
Tiago Peixoto's avatar
Tiago Peixoto committed
209
    (0.6285094791115295, 0.010745128857935755)
Tiago Peixoto's avatar
Tiago Peixoto committed
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224

    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)
225
    ...
Tiago Peixoto's avatar
Tiago Peixoto committed
226
227
228
229
230
231
232

    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")
233
234
    >>>
    >>> clf()
235
    >>> imshow(hist[0].T, interpolation="nearest", origin="lower")
Tiago Peixoto's avatar
Tiago Peixoto committed
236
237
238
    <...>
    >>> colorbar()
    <...>
239
    >>> xlabel("in-degree")
Tiago Peixoto's avatar
Tiago Peixoto committed
240
    <...>
241
    >>> ylabel("out-degree")
Tiago Peixoto's avatar
Tiago Peixoto committed
242
    <...>
243
    >>> savefig("combined-deg-hist.pdf")
Tiago Peixoto's avatar
Tiago Peixoto committed
244

245
246
247
248
249
    .. testcode::
       :hide:

       savefig("combined-deg-hist.png")

250
    .. figure:: combined-deg-hist.*
Tiago Peixoto's avatar
Tiago Peixoto committed
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
        :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)),
269
270
271
    ...                     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
272
273
274

    Lets plot the average degree correlations to check.

275
    >>> clf()
276
277
    >>> axes([0.1,0.15,0.63,0.8])
    <...>
Tiago Peixoto's avatar
Tiago Peixoto committed
278
    >>> corr = gt.avg_neighbour_corr(g, "in", "in")
279
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
280
    ...         label=r"$\left<\text{in}\right>$ vs in")
281
    <...>
Tiago Peixoto's avatar
Tiago Peixoto committed
282
    >>> corr = gt.avg_neighbour_corr(g, "in", "out")
283
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
284
    ...         label=r"$\left<\text{out}\right>$ vs in")
285
    <...>
Tiago Peixoto's avatar
Tiago Peixoto committed
286
    >>> corr = gt.avg_neighbour_corr(g, "out", "in")
287
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
288
    ...          label=r"$\left<\text{in}\right>$ vs out")
289
    <...>
Tiago Peixoto's avatar
Tiago Peixoto committed
290
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
291
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
292
    ...          label=r"$\left<\text{out}\right>$ vs out")
Tiago Peixoto's avatar
Tiago Peixoto committed
293
    <...>
294
295
296
    >>> legend(bbox_to_anchor=(1.01, 0.5), loc="center left", borderaxespad=0.)
    <...>
    >>> xlabel("Source degree")
Tiago Peixoto's avatar
Tiago Peixoto committed
297
    <...>
298
    >>> ylabel("Average target degree")
Tiago Peixoto's avatar
Tiago Peixoto committed
299
    <...>
300
    >>> savefig("deg-corr-dir.pdf")
Tiago Peixoto's avatar
Tiago Peixoto committed
301

302
303
304
305
306
    .. testcode::
       :hide:

       savefig("deg-corr-dir.png")

307
    .. figure:: deg-corr-dir.*
Tiago Peixoto's avatar
Tiago Peixoto committed
308
309
310
        :align: center

        Average nearest neighbour correlations.
311
312
313
314
315


    **Blockmodels**


316
317
318
    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
319
320
321
322
323
324
325
326
327
328
    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.

Tiago Peixoto's avatar
Tiago Peixoto committed
329
    >>> g, bm = gt.random_graph(5000, lambda: poisson(10), directed=False,
330
331
    ...                         blockmodel=lambda: randint(10), deg_corr=corr,
    ...                         mix_time=500)
332
    >>> gt.graph_draw(g, vertex_fill_color=bm, output="blockmodel.pdf")
333
334
    <...>

335
336
337
338
339
    .. testcode::
       :hide:

       gt.graph_draw(g, vertex_fill_color=bm, output="blockmodel.png")

340
    .. figure:: blockmodel.*
341
342
343
344
345
346
347
348
349
350
        :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
351
       (6): 1087-1092 (1953). :doi:`10.1063/1.1699114`
352
    .. [hastings-monte-carlo-1970] Hastings, W.K. "Monte Carlo Sampling Methods
353
       Using Markov Chains and Their Applications". Biometrika 57 (1): 97-109 (1970).
354
       :doi:`10.1093/biomet/57.1.97`
355
356
357
358
359
360
    .. [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
361
    """
362

363
364
365
366
367
    g = Graph()
    if deg_corr == None:
        uncorrelated = True
    else:
        uncorrelated = False
368
369
370

    if (type(blockmodel) is types.FunctionType or
        type(blockmodel) is types.LambdaType):
371
372
        btype = block_type
        bm = []
373
        if len(inspect.getargspec(blockmodel)[0]) == 0:
374
            for i in range(N):
375
                bm.append(blockmodel())
376
        else:
377
            for i in range(N):
378
379
                bm.append(blockmodel(i))
        blockmodel = bm
Tiago Peixoto's avatar
Tiago Peixoto committed
380
    elif blockmodel is not None:
381
        btype = _gt_type(blockmodel[0])
382
383
384

    if len(inspect.getargspec(deg_sampler)[0]) > 0:
        if blockmodel is not None:
385
            sampler = lambda i: deg_sampler(i, blockmodel[i])
386
        else:
Tiago Peixoto's avatar
Tiago Peixoto committed
387
            sampler = deg_sampler
388
389
390
391
    else:
        sampler = lambda i: deg_sampler()

    libgraph_tool_generation.gen_graph(g._Graph__graph, N, sampler,
392
393
                                       uncorrelated, not parallel_edges,
                                       not self_loops, not directed,
394
                                       _get_rng(), verbose, True)
395
396
    g.set_directed(directed)

397
398
399
400
401
402
403
404
405
406
407
408
    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>"

409
410
411
412
    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():
413
414
415
416
417
418
419
420
                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())
421
422
423
424
425
426
427
428
429
        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
430

431
    if parallel_edges and self_loops and deg_corr is None:
432
        mix_time = 1
Tiago Peixoto's avatar
Tiago Peixoto committed
433
    if random:
434
435
        if deg_corr is not None:
            random_rewire(g, strat="probabilistic", n_iter=mix_time,
Tiago Peixoto's avatar
Tiago Peixoto committed
436
                          parallel_edges=parallel_edges, deg_corr=deg_corr,
437
438
                          cache_probs=cache_probs, self_loops=self_loops,
                          blockmodel=bm, verbose=verbose)
439
440
441
442
        else:
            random_rewire(g, parallel_edges=parallel_edges, n_iter=mix_time,
                          self_loops=self_loops, verbose=verbose)

443
444
445
446
    if bm is None:
        return g
    else:
        return g, bm
447

Tiago Peixoto's avatar
Tiago Peixoto committed
448

449
450
@_limit_args({"strat": ["erdos", "correlated", "uncorrelated",
                        "probabilistic"]})
451
452
def random_rewire(g, strat="uncorrelated", n_iter=1, edge_sweep=True,
                  parallel_edges=False, self_loops=False, deg_corr=None,
453
454
                  cache_probs=True, blockmodel=None, ret_fail=False,
                  verbose=False):
455
    r"""
456
457
458
459
460
461
462
463
    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``.
464
465
466

    Parameters
    ----------
467
    g : :class:`~graph_tool.Graph`
468
        Graph to be shuffled. The graph will be modified.
469
470
471
472
    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
473
        new source and target of each edge both have the same in and out-degree.
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
        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``)
490
491
492
493
494
495
        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,
496
        unless ``strat == "probabilistic"``.
497
498
499

        If ``blockmodel != None``, the value passed to the function will be the
        block value of the respective vertices, not the in/out-degree pairs.
500
501
502
503
504
505
506
    cache_probs : bool (optional, default: ``True``)
        If ``True``, the probabilities returned by the ``deg_corr`` parameter
        will be cached internally. This is crucial for good performance, since
        in this case the supplied python function is called only a few times,
        and not at every attempted edge rewire move. However, in the case were
        the different parameter combinations to the probability function is very
        large, the memory requirements to keep the cache may be very large.
507
508
509
510
    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.
511
512
513
514
515
516
517
518
519
520
521
    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).
522
523
524
525
526
527
528

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

    Notes
    -----
Tiago Peixoto's avatar
Tiago Peixoto committed
529
    This algorithm iterates through all the edges in the network and tries to
530
    swap its target or source with the target or source of another edge.
Tiago Peixoto's avatar
Tiago Peixoto committed
531
532

    .. note::
533

534
535
536
537
538
539
540
541
542
543
        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
544
545
546
        the Markov chain is implemented using the Metropolis-Hastings
        [metropolis-equations-1953]_ [hastings-monte-carlo-1970]_
        acceptance/rejection algorithm.
547

Tiago Peixoto's avatar
Tiago Peixoto committed
548

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

553
554
555
556
557
    Examples
    --------

    Some small graphs for visualization.

558
559
560
561
562
563
564
565
    .. testcode::
       :hide:

       from numpy.random import random, seed
       from pylab import *
       seed(43)
       gt.seed_rng(42)

566
    >>> g, pos = gt.triangulation(random((1000,2)))
567
    >>> pos = gt.arf_layout(g)
568
    >>> gt.graph_draw(g, pos=pos, output="rewire_orig.pdf", output_size=(300, 300))
569
    <...>
570
571
572
573
574
575

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output="rewire_orig.png", output_size=(300, 300))

576
    >>> gt.random_rewire(g, "correlated")
577
    >>> pos = gt.arf_layout(g)
578
    >>> gt.graph_draw(g, pos=pos, output="rewire_corr.pdf", output_size=(300, 300))
579
    <...>
580
581
582
583
584
585

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output="rewire_corr.png", output_size=(300, 300))

586
    >>> gt.random_rewire(g)
587
    >>> pos = gt.arf_layout(g)
588
    >>> gt.graph_draw(g, pos=pos, output="rewire_uncorr.pdf", output_size=(300, 300))
589
    <...>
590
591
592
593
594
595

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output="rewire_uncorr.png", output_size=(300, 300))

596
    >>> gt.random_rewire(g, "erdos")
597
    >>> pos = gt.arf_layout(g)
598
    >>> gt.graph_draw(g, pos=pos, output="rewire_erdos.pdf", output_size=(300, 300))
599
    <...>
600

601
602
603
604
605
    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output="rewire_erdos.png", output_size=(300, 300))

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

608
609
610
611
    .. image:: rewire_orig.*
    .. image:: rewire_corr.*
    .. image:: rewire_uncorr.*
    .. image:: rewire_erdos.*
612

613
614
    **From left to right**: Original graph; Shuffled graph, with degree correlations;
    Shuffled graph, without degree correlations; Shuffled graph, with random degrees.
615

616
    We can try with larger graphs to get better statistics, as follows.
617

618
619
    >>> figure()
    <...>
620
    >>> g = gt.random_graph(30000, lambda: sample_k(20),
621
622
    ...                     lambda i, j: exp(abs(i-j)), directed=False,
    ...                     mix_time=100)
623
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
624
625
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-", label="Original")
    <...>
626
627
    >>> gt.random_rewire(g, "correlated")
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
628
629
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="*", label="Correlated")
    <...>
630
631
    >>> gt.random_rewire(g)
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
632
633
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-", label="Uncorrelated")
    <...>
634
635
    >>> gt.random_rewire(g, "erdos")
    >>> corr = gt.avg_neighbour_corr(g, "out", "out")
636
637
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-", label=r"Erd\H{o}s")
    <...>
638
639
640
641
642
643
    >>> xlabel("$k$")
    <...>
    >>> ylabel(r"$\left<k_{nn}\right>$")
    <...>
    >>> legend(loc="best")
    <...>
644
    >>> savefig("shuffled-stats.pdf")
645

646
647
648
649
650
651
    .. testcode::
       :hide:

       savefig("shuffled-stats.png")


652
    .. figure:: shuffled-stats.*
653
654
655
656
657
658
659
660
661
662
663
        :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)),
664
665
    ...                     lambda a, b: (p.pmf(a[0], b[1]) * p.pmf(a[1], 20 - b[0])),
    ...                     mix_time=100)
666
    >>> figure()
667
668
669
    <...>
    >>> axes([0.1,0.15,0.6,0.8])
    <...>
670
    >>> corr = gt.avg_neighbour_corr(g, "in", "out")
671
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
672
    ...          label=r"$\left<\text{o}\right>$ vs i")
673
    <...>
674
    >>> corr = gt.avg_neighbour_corr(g, "out", "in")
675
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
676
    ...          label=r"$\left<\text{i}\right>$ vs o")
677
    <...>
678
679
    >>> gt.random_rewire(g, "correlated")
    >>> corr = gt.avg_neighbour_corr(g, "in", "out")
680
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
681
    ...          label=r"$\left<\text{o}\right>$ vs i, corr.")
682
    <...>
683
    >>> corr = gt.avg_neighbour_corr(g, "out", "in")
684
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
685
    ...          label=r"$\left<\text{i}\right>$ vs o, corr.")
686
    <...>
687
688
    >>> gt.random_rewire(g, "uncorrelated")
    >>> corr = gt.avg_neighbour_corr(g, "in", "out")
689
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
690
    ...          label=r"$\left<\text{o}\right>$ vs i, uncorr.")
691
    <...>
692
    >>> corr = gt.avg_neighbour_corr(g, "out", "in")
693
    >>> errorbar(corr[2][:-1], corr[0], yerr=corr[1], fmt="o-",
694
    ...          label=r"$\left<\text{i}\right>$ vs o, uncorr.")
695
    <...>
696
697
698
    >>> legend(bbox_to_anchor=(1.01, 0.5), loc="center left", borderaxespad=0.)
    <...>
    >>> xlabel("Source degree")
699
    <...>
700
    >>> ylabel("Average target degree")
701
    <...>
702
    >>> savefig("shuffled-deg-corr-dir.pdf")
703

704
705
706
707
708
    .. testcode::
       :hide:

       savefig("shuffled-deg-corr-dir.png")

709
    .. figure:: shuffled-deg-corr-dir.*
710
711
712
713
714
715
        :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.

716
717
718
719
720
    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
721
       (6): 1087-1092 (1953). :doi:`10.1063/1.1699114`
722
    .. [hastings-monte-carlo-1970] Hastings, W.K. "Monte Carlo Sampling Methods
723
       Using Markov Chains and Their Applications". Biometrika 57 (1): 97-109 (1970).
724
       :doi:`10.1093/biomet/57.1.97`
725
726
727
728
729
730
    .. [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`
731
732

    """
Tiago Peixoto's avatar
Tiago Peixoto committed
733
734
735
736
737
738
739
740
741
742
743
744
745
    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!")

746
    if (deg_corr is not None and not g.is_directed()) and blockmodel is None:
Tiago Peixoto's avatar
Tiago Peixoto committed
747
        corr = lambda i, j: deg_corr(i[1], j[1])
748
749
750
    else:
        corr = deg_corr

751
752
    if strat != "probabilistic":
        g = GraphView(g, reversed=False)
753
754
    elif blockmodel is not None:
        strat = "blockmodel"
755
756
757
    pcount = libgraph_tool_generation.random_rewire(g._Graph__graph, strat,
                                                    n_iter, not edge_sweep,
                                                    self_loops, parallel_edges,
758
                                                    corr, _prop("v", g, blockmodel),
759
                                                    cache_probs,
760
                                                    _get_rng(), verbose)
761
762
    if ret_fail:
        return pcount
Tiago Peixoto's avatar
Tiago Peixoto committed
763

Tiago Peixoto's avatar
Tiago Peixoto committed
764

Tiago Peixoto's avatar
Tiago Peixoto committed
765
def predecessor_tree(g, pred_map):
Tiago Peixoto's avatar
Tiago Peixoto committed
766
    """Return a graph from a list of predecessors given by the ``pred_map`` vertex property."""
Tiago Peixoto's avatar
Tiago Peixoto committed
767
768
769
770
771
772
773

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

Tiago Peixoto's avatar
Tiago Peixoto committed
775

776
def line_graph(g):
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
    """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.

793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819

    Examples
    --------

    >>> g = gt.collection.data["lesmis"]
    >>> lg, vmap = gt.line_graph(g)
    >>> gt.graph_draw(g, pos=g.vp["pos"], output="lesmis.pdf")
    <...>
    >>> pos = gt.graph_draw(lg, output="lesmis-lg.pdf")

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=g.vp["pos"], output="lesmis.png")
       pos = gt.graph_draw(lg, pos=pos, output="lesmis-lg.png")


    .. figure:: lesmis.png
       :align: left

       Coappearances of characters in Victor Hugo's novel "Les Miserables".

    .. figure:: lesmis-lg.png
       :align: right

       Line graph of the coappearance network on the left.

820
821
822
823
    References
    ----------
    .. [line-wiki] http://en.wikipedia.org/wiki/Line_graph
    """
824
825
826
827
828
829
830
831
    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
832

Tiago Peixoto's avatar
Tiago Peixoto committed
833

834
def graph_union(g1, g2, intersection=None, props=None, include=False):
835
836
837
838
839
840
841
842
843
    """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.
844
845
846
847
848
    intersection : :class:`~graph_tool.PropertyMap` (optional, default: ``None``)
       Vertex property map owned by `g1` which maps each of each of its vertices
       to vertex indexes belonging to `g2`. Negative values mean no mapping
       exists, and thus both vertices in `g1` and `g2` will be present in the
       union graph.
849
    props : list of tuples of :class:`~graph_tool.PropertyMap` (optional, default: ``[]``)
850
851
852
853
       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.
854
    include : bool (optional, default: ``False``)
855
856
857
858
859
860
861
862
863
864
       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.
865
866
867
868

    Examples
    --------

869
870
871
872
873
874
875
876
    .. testcode::
       :hide:

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

877
878
879
    >>> g = gt.triangulation(random((300,2)))[0]
    >>> ug = gt.graph_union(g, g)
    >>> uug = gt.graph_union(g, ug)
880
    >>> pos = gt.sfdp_layout(g)
881
    >>> gt.graph_draw(g, pos=pos, output_size=(300,300), output="graph_original.pdf")
882
    <...>
883
884
885
886
887
888
889

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output_size=(300,300), output="graph_original.png")

    >>> pos = gt.sfdp_layout(ug)
890
    >>> gt.graph_draw(ug, pos=pos, output_size=(300,300), output="graph_union.pdf")
891
    <...>
892
893
894
895
896
897
898

    .. testcode::
       :hide:

       gt.graph_draw(ug, pos=pos, output_size=(300,300), output="graph_union.png")

    >>> pos = gt.sfdp_layout(uug)
899
    >>> gt.graph_draw(uug, pos=pos, output_size=(300,300), output="graph_union2.pdf")
900
901
    <...>

902
903
904
905
906
907
    .. testcode::
       :hide:

       gt.graph_draw(uug, pos=pos, output_size=(300,300), output="graph_union2.png")


908
909
910
    .. image:: graph_original.*
    .. image:: graph_union.*
    .. image:: graph_union2.*
911

912
    """
Tiago Peixoto's avatar
Tiago Peixoto committed
913
914
    if props == None:
        props = []
Tiago Peixoto's avatar
Tiago Peixoto committed
915
916
    if not include:
        g1 = Graph(g1)
917
918
919
920
921
922
923
924
    if intersection is None:
        intersection = g1.new_vertex_property("int32_t")
        intersection.a = 0
    else:
        intersection = intersection.copy("int32_t")
        intersection.a[intersection.a >= 0] += 1
        intersection.a[intersection.a < 0] = 0

Tiago Peixoto's avatar
Tiago Peixoto committed
925
926
927
928
929
930
931
932
    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,
933
934
935
936
                                                          g2._Graph__graph,
                                                          _prop("v", g1,
                                                                intersection))
        for p1, p2 in props:
Tiago Peixoto's avatar
Tiago Peixoto committed
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
            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
962

Tiago Peixoto's avatar
Tiago Peixoto committed
963
964

@_limit_args({"type": ["simple", "delaunay"]})
965
def triangulation(points, type="simple", periodic=False):
966
967
968
969
970
971
972
973
    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).
974
    type : string (optional, default: ``'simple'``)
975
        Type of triangulation. May be either 'simple' or 'delaunay'.
976
977
978
    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.
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993

    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
994
    A triangulation [cgal-triang]_ is a division of the convex hull of a point
995
    set into triangles, using only that set as triangle vertices.
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012

    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
    --------
1013
1014
1015
1016
1017
1018
1019
    .. testcode::
       :hide:

       from numpy.random import random, seed
       from pylab import *
       seed(42)
       gt.seed_rng(42)
1020
    >>> points = random((500, 2)) * 4
1021
    >>> g, pos = gt.triangulation(points)
1022
1023
1024
1025
1026
1027
1028
    >>> 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
1029
1030
    >>> gt.graph_draw(g, pos=pos, output_size=(300,300), vertex_fill_color=b[0],
    ...               edge_pen_width=b[1], output="triang.pdf")
1031
    <...>
1032
1033
1034
1035
1036
1037
1038

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output_size=(300,300), vertex_fill_color=b[0],
                     edge_pen_width=b[1], output="triang.png")

1039
    >>> g, pos = gt.triangulation(points, type="delaunay")
1040
1041
1042
1043
1044
1045
    >>> 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
1046
1047
    >>> gt.graph_draw(g, pos=pos, output_size=(300,300), vertex_fill_color=b[0],
    ...               edge_pen_width=b[1], output="triang-delaunay.pdf")
1048
1049
    <...>

1050
1051
1052
1053
1054
1055
1056
    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output_size=(300,300), vertex_fill_color=b[0],
                     edge_pen_width=b[1], output="triang-delaunay.png")


1057
1058
    2D triangulation of random points:

1059
1060
    .. image:: triang.*
    .. image:: triang-delaunay.*
1061

1062
1063
1064
    *Left:* Simple triangulation. *Right:* Delaunay triangulation. The vertex
    colors and the edge thickness correspond to the weighted betweenness
    centrality.
1065
1066
1067

    References
    ----------
Tiago Peixoto's avatar
Tiago Peixoto committed
1068
    .. [cgal-triang] http://www.cgal.org/Manual/last/doc_html/cgal_manual/Triangulation_3/Chapter_main.html
1069
1070
1071

    """

Tiago Peixoto's avatar
Tiago Peixoto committed
1072
    if points.shape[1] not in [2, 3]:
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
        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,
1083
                                           _prop("v", g, pos), type, periodic)
1084
    return g, pos
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094


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.
1095
    periodic : bool (optional, default: ``False``)
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
        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
    --------
1110
1111
1112
1113
1114
    .. testcode::
       :hide:

       gt.seed_rng(42)

1115
    >>> g = gt.lattice([10,10])
1116
1117
    >>> pos = gt.sfdp_layout(g, cooling_step=0.95, epsilon=1e-2)
    >>> gt.graph_draw(g, pos=pos, output_size=(300,300), output="lattice.pdf")
1118
    <...>
1119
1120
1121
1122
1123
1124

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output_size=(300,300), output="lattice.png")

1125
    >>> g = gt.lattice([10,20], periodic=True)
1126
1127
    >>> pos = gt.sfdp_layout(g, cooling_step=0.95, epsilon=1e-2)
    >>> gt.graph_draw(g, pos=pos, output_size=(300,300), output="lattice_periodic.pdf")
1128
    <...>
1129
1130
1131
1132
1133
1134

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output_size=(300,300), output="lattice_periodic.png")

1135
    >>> g = gt.lattice([10,10,10])
1136
1137
    >>> pos = gt.sfdp_layout(g, cooling_step=0.95, epsilon=1e-2)
    >>> gt.graph_draw(g, pos=pos, output_size=(300,300), output="lattice_3d.pdf")
1138
1139
    <...>

1140
1141
1142
1143
1144
1145
    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output_size=(300,300), output="lattice_3d.png")


1146
1147
1148
    .. image:: lattice.*
    .. image:: lattice_periodic.*
    .. image:: lattice_3d.*
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162

    *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

Tiago Peixoto's avatar
Tiago Peixoto committed
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
def complete_graph(N, self_loops=False, directed=False):
    r"""
    Generate complete graph.

    Parameters
    ----------
    N : ``int``
        Number of vertices.
    self_loops : bool (optional, default: ``False``)
        If ``True``, self-loops are included.
    directed : bool (optional, default: ``False``)
        If ``True``, a directed graph is generated.

    Returns
    -------
    complete_graph : :class:`~graph_tool.Graph`
        A complete graph.

    Examples
    --------

    >>> g = gt.complete_graph(30)
    >>> pos = gt.sfdp_layout(g, cooling_step=0.95, epsilon=1e-2)
    >>> gt.graph_draw(g, pos=pos, output_size=(300,300), output="complete.pdf")
    <...>

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output_size=(300,300), output="complete.png")


    .. figure:: complete.*

       A complete graph with :math:`N=30` vertices.

    References
    ----------
    .. [complete] http://en.wikipedia.org/wiki/Complete_graph

    """

    g = Graph(directed=directed)
    libgraph_tool_generation.complete(g._Graph__graph, N, directed, self_loops)
    return g

Tiago Peixoto's avatar
Tiago Peixoto committed
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
def circular_graph(N, k=1, self_loops=False, directed=False):
    r"""
    Generate a circular graph.

    Parameters
    ----------
    N : ``int``
        Number of vertices.
    k : ``int`` (optional, default: ``True``)
        Number of nearest neighbours to be connected.
    self_loops : bool (optional, default: ``False``)
        If ``True``, self-loops are included.
    directed : bool (optional, default: ``False``)
        If ``True``, a directed graph is generated.

    Returns
    -------
    circular_graph : :class:`~graph_tool.Graph`
        A circular graph.

    Examples
    --------

    >>> g = gt.circular_graph(30, 2)
    >>> pos = gt.sfdp_layout(g, cooling_step=0.95)
    >>> gt.graph_draw(g, pos=pos, output_size=(300,300), output="circular.pdf")
    <...>

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output_size=(300,300), output="circular.png")

    .. figure:: circular.*

       A circular graph with :math:`N=30` vertices, and :math:`k=2`.

    """

    g = Graph(directed=directed)
    libgraph_tool_generation.circular(g._Graph__graph, N, k, directed, self_loops)
    return g

1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264

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.
1265
    ranges : list or :class:`~numpy.ndarray` (optional, default: ``None``)
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
        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
    --------
1292
1293
1294
1295
1296
1297
1298
1299
    .. testcode::
       :hide:

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

1300
1301
    >>> points = random((500, 2)) * 4
    >>> g, pos = gt.geometric_graph(points, 0.3)
1302
    >>> gt.graph_draw(g, pos=pos, output_size=(300,300), output="geometric.pdf")
1303
    <...>
1304
1305
1306
1307
1308
1309

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output_size=(300,300), output="geometric.png")

1310
    >>> g, pos = gt.geometric_graph(points, 0.3, [(0,4), (0,4)])
1311
1312
1313
1314
1315
1316
1317
    >>> pos = gt.graph_draw(g, output_size=(300,300), output="geometric_periodic.pdf")

    .. testcode::
       :hide:

       gt.graph_draw(g, pos=pos, output_size=(300,300), output="geometric_periodic.png")

1318

1319
1320
    .. image:: geometric.*
    .. image:: geometric_periodic.*
1321
1322
1323
1324
1325
1326
1327

    *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
1328
       graphs", Phys. Rev. E 66, 016121 (2002), :doi:`10.1103/PhysRevE.66.016121`
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351

    """

    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
1352
1353
1354
1355
1356
1357
1358
1359
1360


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.