__init__.py 16.2 KB
Newer Older
Tiago Peixoto's avatar
Tiago Peixoto committed
1
#! /usr/bin/env python
2
# -*- coding: utf-8 -*-
Tiago Peixoto's avatar
Tiago Peixoto committed
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>
Tiago Peixoto's avatar
Tiago Peixoto committed
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
#
# 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/>.

"""
``graph_tool.spectral`` - Spectral properties
---------------------------------------------
24
25
26
27
28
29
30
31
32
33

Summary
+++++++

.. autosummary::
   :nosignatures:

   adjacency
   laplacian
   incidence
34
35
   transition
   modularity_matrix
36
37
38

Contents
++++++++
Tiago Peixoto's avatar
Tiago Peixoto committed
39
40
"""

41
42
from __future__ import division, absolute_import, print_function

43
from .. import _degree, _prop, Graph, _limit_args
44
45
from .. stats import label_self_loops
import numpy
Tiago Peixoto's avatar
Tiago Peixoto committed
46
import scipy.sparse
47
import scipy.sparse.linalg
Tiago Peixoto's avatar
Tiago Peixoto committed
48

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

52
__all__ = ["adjacency", "laplacian", "incidence", "transition", "modularity_matrix"]
Tiago Peixoto's avatar
Tiago Peixoto committed
53

54

55
def adjacency(g, weight=None, index=None):
56
57
58
59
    r"""Return the adjacency matrix of the graph.

    Parameters
    ----------
60
    g : :class:`~graph_tool.Graph`
61
        Graph to be used.
62
    weight : :class:`~graph_tool.PropertyMap` (optional, default: True)
63
        Edge property map with the edge weights.
64
65
66
    index : :class:`~graph_tool.PropertyMap` (optional, default: None)
        Vertex property map specifying the row/column indexes. If not provided, the
        internal vertex index is used.
67
68
69

    Returns
    -------
Tiago Peixoto's avatar
Tiago Peixoto committed
70
    a : :class:`~scipy.sparse.csr_matrix`
71
        The (sparse) adjacency matrix.
72
73
74
75
76
77
78
79
80
81

    Notes
    -----
    The adjacency matrix is defined as

    .. math::

        a_{i,j} =
        \begin{cases}
            1 & \text{if } v_i \text{ is adjacent to } v_j, \\
Tiago Peixoto's avatar
Tiago Peixoto committed
82
            2 & \text{if } i = j, \text{ the graph is undirected and there is a self-loop incident in } v_i, \\
83
84
85
            0 & \text{otherwise}
        \end{cases}

Tiago Peixoto's avatar
Tiago Peixoto committed
86
87
    In the case of weighted edges, the entry values are multiplied by the weight
    of the respective edge.
88

89
90
91
    In the case of networks with parallel edges, the entries in the matrix
    become simply the edge multiplicities.

92
93
    Examples
    --------
94
95
    .. testsetup::

Tiago Peixoto's avatar
Tiago Peixoto committed
96
97
       import scipy.linalg
       from pylab import *
98

Tiago Peixoto's avatar
Tiago Peixoto committed
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
    >>> g = gt.collection.data["polblogs"]
    >>> A = gt.adjacency(g)
    >>> ew, ev = scipy.linalg.eig(A.todense())

    >>> figure(figsize=(8, 2))
    <...>
    >>> scatter(real(ew), imag(ew), c=abs(ew))
    <...>
    >>> xlabel(r"$\operatorname{Re}(\lambda)$")
    <...>
    >>> ylabel(r"$\operatorname{Im}(\lambda)$")
    <...>
    >>> tight_layout()
    >>> savefig("adjacency-spectrum.pdf")

    .. testcode::
       :hide:

       savefig("adjacency-spectrum.png")

    .. figure:: adjacency-spectrum.*
        :align: center

        Adjacency matrix spectrum for the political blog network.
123
124
125

    References
    ----------
126
    .. [wikipedia-adjacency] http://en.wikipedia.org/wiki/Adjacency_matrix
127
128
    """

129
130
131
132
133
134
    if index is None:
        if g.get_vertex_filter()[0] != None:
            index = g.new_vertex_property("int64_t")
            index.fa = numpy.arange(g.num_vertices())
        else:
            index = g.vertex_index
Tiago Peixoto's avatar
Tiago Peixoto committed
135

136
137
138
139
    E = g.num_edges() if g.is_directed() else 2 * g.num_edges()
    data = numpy.zeros(E, dtype="double")
    i = numpy.zeros(E, dtype="int32")
    j = numpy.zeros(E, dtype="int32")
Tiago Peixoto's avatar
Tiago Peixoto committed
140

141
142
    libgraph_tool_spectral.adjacency(g._Graph__graph, _prop("v", g, index),
                                     _prop("e", g, weight), data, i, j)
Tiago Peixoto's avatar
Tiago Peixoto committed
143
144
    V = max(g.num_vertices(), max(i.max() + 1, j.max() + 1))
    m = scipy.sparse.coo_matrix((data, (i,j)), shape=(V, V))
145
146
    m = m.tocsr()
    return m
Tiago Peixoto's avatar
Tiago Peixoto committed
147

Tiago Peixoto's avatar
Tiago Peixoto committed
148
149

@_limit_args({"deg": ["total", "in", "out"]})
Tiago Peixoto's avatar
Tiago Peixoto committed
150
def laplacian(g, deg="total", normalized=False, weight=None, index=None):
151
152
153
154
    r"""Return the Laplacian matrix of the graph.

    Parameters
    ----------
155
    g : :class:`~graph_tool.Graph`
156
157
158
        Graph to be used.
    deg : str (optional, default: "total")
        Degree to be used, in case of a directed graph.
Tiago Peixoto's avatar
Tiago Peixoto committed
159
    normalized : bool (optional, default: False)
160
        Whether to compute the normalized Laplacian.
161
    weight : :class:`~graph_tool.PropertyMap` (optional, default: True)
162
        Edge property map with the edge weights.
163
164
165
    index : :class:`~graph_tool.PropertyMap` (optional, default: None)
        Vertex property map specifying the row/column indexes. If not provided, the
        internal vertex index is used.
166
167
168

    Returns
    -------
Tiago Peixoto's avatar
Tiago Peixoto committed
169
    l : :class:`~scipy.sparse.csr_matrix`
170
        The (sparse) Laplacian matrix.
171
172
173

    Notes
    -----
Tiago Peixoto's avatar
Tiago Peixoto committed
174
    The weighted Laplacian matrix is defined as
175
176
177

    .. math::

Tiago Peixoto's avatar
Tiago Peixoto committed
178
        \ell_{ij} =
179
180
        \begin{cases}
        \Gamma(v_i) & \text{if } i = j \\
Tiago Peixoto's avatar
Tiago Peixoto committed
181
        -w_{ij}     & \text{if } i \neq j \text{ and } v_i \text{ is adjacent to } v_j \\
182
183
184
        0           & \text{otherwise}.
        \end{cases}

185
186
    Where :math:`\Gamma(v_i)=\sum_j A_{ij}w_{ij}` is sum of the weights of the
    edges incident on vertex :math:`v_i`. The normalized version is
187
188
189

    .. math::

Tiago Peixoto's avatar
Tiago Peixoto committed
190
        \ell_{ij} =
191
192
        \begin{cases}
        1         & \text{ if } i = j \text{ and } \Gamma(v_i) \neq 0 \\
Tiago Peixoto's avatar
Tiago Peixoto committed
193
        -\frac{w_{ij}}{\sqrt{\Gamma(v_i)\Gamma(v_j)}} & \text{ if } i \neq j \text{ and } v_i \text{ is adjacent to } v_j \\
194
195
196
        0         & \text{otherwise}.
        \end{cases}

Tiago Peixoto's avatar
Tiago Peixoto committed
197
198
199
200
201
    In the case of unweighted edges, it is assumed :math:`w_{ij} = 1`.

    For directed graphs, it is assumed :math:`\Gamma(v_i)=\sum_j A_{ij}w_{ij} +
    \sum_j A_{ji}w_{ji}` if ``deg=="total"``, :math:`\Gamma(v_i)=\sum_j A_{ij}w_{ij}`
    if ``deg=="out"`` or :math:`\Gamma(v_i)=\sum_j A_{ji}w_{ji}` ``deg=="in"``.
202
203
204

    Examples
    --------
Tiago Peixoto's avatar
Tiago Peixoto committed
205

206
207
    .. testsetup::

Tiago Peixoto's avatar
Tiago Peixoto committed
208
209
       import scipy.linalg
       from pylab import *
210

Tiago Peixoto's avatar
Tiago Peixoto committed
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
    >>> g = gt.collection.data["polblogs"]
    >>> L = gt.laplacian(g)
    >>> ew, ev = scipy.linalg.eig(L.todense())

    >>> figure(figsize=(8, 2))
    <...>
    >>> scatter(real(ew), imag(ew), c=abs(ew))
    <...>
    >>> xlabel(r"$\operatorname{Re}(\lambda)$")
    <...>
    >>> ylabel(r"$\operatorname{Im}(\lambda)$")
    <...>
    >>> tight_layout()
    >>> savefig("laplacian-spectrum.pdf")

    .. testcode::
       :hide:

       savefig("laplacian-spectrum.png")

    .. figure:: laplacian-spectrum.*
        :align: center

        Laplacian matrix spectrum for the political blog network.

    >>> L = gt.laplacian(g, normalized=True)
    >>> ew, ev = scipy.linalg.eig(L.todense())

    >>> figure(figsize=(8, 2))
    <...>
    >>> scatter(real(ew), imag(ew), c=abs(ew))
    <...>
    >>> xlabel(r"$\operatorname{Re}(\lambda)$")
    <...>
    >>> ylabel(r"$\operatorname{Im}(\lambda)$")
    <...>
    >>> tight_layout()
    >>> savefig("norm-laplacian-spectrum.pdf")

    .. testcode::
       :hide:

       savefig("norm-laplacian-spectrum.png")

    .. figure:: norm-laplacian-spectrum.*
        :align: center

        Normalized Laplacian matrix spectrum for the political blog network.
259
260
261

    References
    ----------
262
    .. [wikipedia-laplacian] http://en.wikipedia.org/wiki/Laplacian_matrix
263
264
    """

265
266
267
268
269
270
271
    if index is None:
        if g.get_vertex_filter()[0] != None:
            index = g.new_vertex_property("int64_t")
            index.fa = numpy.arange(g.num_vertices())
        else:
            index = g.vertex_index

Tiago Peixoto's avatar
Tiago Peixoto committed
272
    V = g.num_vertices()
273
274
275
276
277
278
279
280
281
282
283
284
    nself = label_self_loops(g, mark_only=True).a.sum()
    E = g.num_edges() - nself
    if not g.is_directed():
        E *= 2
    N = E + g.num_vertices()
    data = numpy.zeros(N, dtype="double")
    i = numpy.zeros(N, dtype="int32")
    j = numpy.zeros(N, dtype="int32")

    if normalized:
        libgraph_tool_spectral.norm_laplacian(g._Graph__graph, _prop("v", g, index),
                                              _prop("e", g, weight), deg, data, i, j)
Tiago Peixoto's avatar
Tiago Peixoto committed
285
    else:
286
287
        libgraph_tool_spectral.laplacian(g._Graph__graph, _prop("v", g, index),
                                         _prop("e", g, weight), deg, data, i, j)
Tiago Peixoto's avatar
Tiago Peixoto committed
288
289
    V = max(g.num_vertices(), max(i.max() + 1, j.max() + 1))
    m = scipy.sparse.coo_matrix((data, (i, j)), shape=(V, V))
290
    m = m.tocsr()
Tiago Peixoto's avatar
Tiago Peixoto committed
291
292
    return m

Tiago Peixoto's avatar
Tiago Peixoto committed
293

294
def incidence(g, vindex=None, eindex=None):
295
296
297
298
    r"""Return the incidence matrix of the graph.

    Parameters
    ----------
299
    g : :class:`~graph_tool.Graph`
300
        Graph to be used.
301
302
303
304
305
306
    vindex : :class:`~graph_tool.PropertyMap` (optional, default: None)
        Vertex property map specifying the row indexes. If not provided, the
        internal vertex index is used.
    eindex : :class:`~graph_tool.PropertyMap` (optional, default: None)
        Edge property map specifying the column indexes. If not provided, the
        internal edge index is used.
307
308
309

    Returns
    -------
Tiago Peixoto's avatar
Tiago Peixoto committed
310
311
    a : :class:`~scipy.sparse.csr_matrix`
        The (sparse) incidence matrix.
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337

    Notes
    -----
    For undirected graphs, the incidence matrix is defined as

    .. math::

        b_{i,j} =
        \begin{cases}
            1 & \text{if vertex } v_i \text{and edge } e_j \text{ are incident}, \\
            0 & \text{otherwise}
        \end{cases}

    For directed graphs, the definition is

    .. math::

        b_{i,j} =
        \begin{cases}
            1 & \text{if edge } e_j \text{ enters vertex } v_i, \\
            -1 & \text{if edge } e_j \text{ leaves vertex } v_i, \\
            0 & \text{otherwise}
        \end{cases}

    Examples
    --------
338
339
340
341
    .. testsetup::

      gt.seed_rng(42)

342
343
    >>> g = gt.random_graph(100, lambda: (2,2))
    >>> m = gt.incidence(g)
344
    >>> print(m.todense())
Tiago Peixoto's avatar
Tiago Peixoto committed
345
    [[-1. -1.  0. ...,  0.  0.  0.]
Tiago Peixoto's avatar
Tiago Peixoto committed
346
     [ 0.  0.  0. ...,  0.  0.  0.]
347
     [ 0.  0.  0. ...,  0.  0.  0.]
Tiago Peixoto's avatar
Tiago Peixoto committed
348
     ..., 
Tiago Peixoto's avatar
Tiago Peixoto committed
349
     [ 0.  0. -1. ...,  0.  0.  0.]
350
     [ 0.  0.  0. ...,  0.  0.  0.]
Tiago Peixoto's avatar
Tiago Peixoto committed
351
     [ 0.  0.  0. ...,  1.  0.  0.]]
352
353
354

    References
    ----------
355
    .. [wikipedia-incidence] http://en.wikipedia.org/wiki/Incidence_matrix
356
357
    """

358
359
360
361
362
363
    if vindex is None:
        if g.get_edge_filter()[0] != None:
            vindex = g.new_vertex_property("int64_t")
            vindex.fa = numpy.arange(g.num_vertices())
        else:
            vindex = g.vertex_index
Tiago Peixoto's avatar
Tiago Peixoto committed
364

365
366
367
368
369
370
    if eindex is None:
        if g.get_edge_filter()[0] != None:
            eindex = g.new_edge_property("int64_t")
            eindex.fa = numpy.arange(g.num_edges())
        else:
            eindex = g.edge_index
Tiago Peixoto's avatar
Tiago Peixoto committed
371
372

    E = g.num_edges()
373
374
375
376
377
378
379
380
    data = numpy.zeros(2 * E, dtype="double")
    i = numpy.zeros(2 * E, dtype="int32")
    j = numpy.zeros(2 * E, dtype="int32")

    libgraph_tool_spectral.incidence(g._Graph__graph, _prop("v", g, vindex),
                                     _prop("e", g, eindex), data, i, j)
    m = scipy.sparse.coo_matrix((data, (i,j)))
    m = m.tocsr()
Tiago Peixoto's avatar
Tiago Peixoto committed
381
    return m
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575

def transition(g, weight=None, index=None):
    r"""Return the transition matrix of the graph.

    Parameters
    ----------
    g : :class:`~graph_tool.Graph`
        Graph to be used.
    weight : :class:`~graph_tool.PropertyMap` (optional, default: True)
        Edge property map with the edge weights.
    index : :class:`~graph_tool.PropertyMap` (optional, default: None)
        Vertex property map specifying the row/column indexes. If not provided, the
        internal vertex index is used.

    Returns
    -------
    T : :class:`~scipy.sparse.csr_matrix`
        The (sparse) transition matrix.

    Notes
    -----
    The transition matrix is defined as

    .. math::

        T_{ij} = \frac{A_{ij}}{k_i}

    where :math:`k_i = \sum_j A_{ij}`, and :math:`A_{ij}` is the adjacency
    matrix.

    In the case of weighted edges, the values of the adjacency matrix are
    multiplied by the edge weights.

    Examples
    --------
    .. testsetup::

       import scipy.linalg
       from pylab import *

    >>> g = gt.collection.data["polblogs"]
    >>> T = gt.transition(g)
    >>> ew, ev = scipy.linalg.eig(T.todense())

    >>> figure(figsize=(8, 2))
    <...>
    >>> scatter(real(ew), imag(ew), c=abs(ew))
    <...>
    >>> xlabel(r"$\operatorname{Re}(\lambda)$")
    <...>
    >>> ylabel(r"$\operatorname{Im}(\lambda)$")
    <...>
    >>> tight_layout()
    >>> savefig("transition-spectrum.pdf")

    .. testcode::
       :hide:

       savefig("transition-spectrum.png")

    .. figure:: transition-spectrum.*
        :align: center

        Transition matrix spectrum for the political blog network.

    References
    ----------
    .. [wikipedia-transition] https://en.wikipedia.org/wiki/Stochastic_matrix
    """

    if index is None:
        if g.get_vertex_filter()[0] != None:
            index = g.new_vertex_property("int64_t")
            index.fa = numpy.arange(g.num_vertices())
        else:
            index = g.vertex_index

    E = g.num_edges() if g.is_directed() else 2 * g.num_edges()
    data = numpy.zeros(E, dtype="double")
    i = numpy.zeros(E, dtype="int32")
    j = numpy.zeros(E, dtype="int32")

    libgraph_tool_spectral.transition(g._Graph__graph, _prop("v", g, index),
                                      _prop("e", g, weight), data, i, j)
    V = max(g.num_vertices(), max(i.max() + 1, j.max() + 1))
    m = scipy.sparse.coo_matrix((data, (i,j)), shape=(V, V))
    m = m.tocsr()
    return m



def modularity_matrix(g, weight=None, index=None):
    r"""Return the modularity matrix of the graph.

    Parameters
    ----------
    g : :class:`~graph_tool.Graph`
        Graph to be used.
    weight : :class:`~graph_tool.PropertyMap` (optional, default: True)
        Edge property map with the edge weights.
    index : :class:`~graph_tool.PropertyMap` (optional, default: None)
        Vertex property map specifying the row/column indexes. If not provided, the
        internal vertex index is used.

    Returns
    -------
    B : :class:`~scipy.sparse.linalg.LinearOperator`
        The (sparse) modularity matrix, represented as a
        :class:`~scipy.sparse.linalg.LinearOperator`.

    Notes
    -----
    The modularity matrix is defined as

    .. math::

        B_{ij} =  A_{ij} - \frac{k^+_i k^-_j}{2E}

    where :math:`k^+_i = \sum_j A_{ij}`, :math:`k^-_i = \sum_j A_{ji}`,
    :math:`2E=\sum_{ij}A_{ij}` and :math:`A_{ij}` is the adjacency matrix.

    In the case of weighted edges, the values of the adjacency matrix are
    multiplied by the edge weights.

    Examples
    --------

    .. testsetup::

       import scipy.linalg
       from pylab import *

    >>> g = gt.collection.data["polblogs"]
    >>> B = gt.modularity_matrix(g)
    >>> B = B * np.identity(B.shape[0])  # transform to a dense matrix
    >>> ew, ev = scipy.linalg.eig(B)

    >>> figure(figsize=(8, 2))
    <...>
    >>> scatter(real(ew), imag(ew), c=abs(ew))
    <...>
    >>> xlabel(r"$\operatorname{Re}(\lambda)$")
    <...>
    >>> ylabel(r"$\operatorname{Im}(\lambda)$")
    <...>
    >>> tight_layout()
    >>> savefig("modularity-spectrum.pdf")

    .. testcode::
       :hide:

       savefig("modularity-spectrum.png")

    .. figure:: modularity-spectrum.*
        :align: center

        Modularity matrix spectrum for the political blog network.

    References
    ----------
    .. [newman-modularity]  M. E. J. Newman, M. Girvan, "Finding and evaluating
       community structure in networks", Phys. Rev. E 69, 026113 (2004).
       :doi:`10.1103/PhysRevE.69.026113`
    """

    A = adjacency(g, weight=weight, index=index)
    if g.is_directed():
        k_in = g.degree_property_map("in", weight=weight).fa
    else:
        k_in = g.degree_property_map("out", weight=weight).fa
    k_out = g.degree_property_map("out", weight=weight).fa

    N = A.shape[0]
    E2 = float(k_out.sum())

    def matvec(x):
        M = x.shape[0]
        if len(x.shape) > 1:
            x = x.reshape(M)
        nx = A * x - k_out * numpy.dot(k_in, x) / E2
        return nx

    def rmatvec(x):
        M = x.shape[0]
        if len(x.shape) > 1:
            x = x.reshape(M)
        nx = A.T * x - k_in * numpy.dot(k_out, x) / E2
        return nx

    B = scipy.sparse.linalg.LinearOperator((g.num_vertices(), g.num_vertices()),
                                           matvec=matvec, rmatvec=rmatvec,
                                           dtype="float")

    return B