File: test_graph_laplacian.py

package info (click to toggle)
python-scipy 1.1.0-7
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 93,828 kB
  • sloc: python: 156,854; ansic: 82,925; fortran: 80,777; cpp: 7,505; makefile: 427; sh: 294
file content (136 lines) | stat: -rw-r--r-- 4,388 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
# Author: Gael Varoquaux <gael.varoquaux@normalesup.org>
#         Jake Vanderplas <vanderplas@astro.washington.edu>
# License: BSD
from __future__ import division, print_function, absolute_import

import numpy as np
from numpy.testing import assert_allclose, assert_array_almost_equal
from pytest import raises as assert_raises
from scipy import sparse

from scipy.sparse import csgraph


def _explicit_laplacian(x, normed=False):
    if sparse.issparse(x):
        x = x.todense()
    x = np.asarray(x)
    y = -1.0 * x
    for j in range(y.shape[0]):
        y[j,j] = x[j,j+1:].sum() + x[j,:j].sum()
    if normed:
        d = np.diag(y).copy()
        d[d == 0] = 1.0
        y /= d[:,None]**.5
        y /= d[None,:]**.5
    return y


def _check_symmetric_graph_laplacian(mat, normed):
    if not hasattr(mat, 'shape'):
        mat = eval(mat, dict(np=np, sparse=sparse))

    if sparse.issparse(mat):
        sp_mat = mat
        mat = sp_mat.todense()
    else:
        sp_mat = sparse.csr_matrix(mat)

    laplacian = csgraph.laplacian(mat, normed=normed)
    n_nodes = mat.shape[0]
    if not normed:
        assert_array_almost_equal(laplacian.sum(axis=0), np.zeros(n_nodes))
    assert_array_almost_equal(laplacian.T, laplacian)
    assert_array_almost_equal(laplacian,
            csgraph.laplacian(sp_mat, normed=normed).todense())

    assert_array_almost_equal(laplacian,
            _explicit_laplacian(mat, normed=normed))


def test_laplacian_value_error():
    for t in int, float, complex:
        for m in ([1, 1],
                  [[[1]]],
                  [[1, 2, 3], [4, 5, 6]],
                  [[1, 2], [3, 4], [5, 5]]):
            A = np.array(m, dtype=t)
            assert_raises(ValueError, csgraph.laplacian, A)


def test_symmetric_graph_laplacian():
    symmetric_mats = ('np.arange(10) * np.arange(10)[:, np.newaxis]',
            'np.ones((7, 7))',
            'np.eye(19)',
            'sparse.diags([1, 1], [-1, 1], shape=(4,4))',
            'sparse.diags([1, 1], [-1, 1], shape=(4,4)).todense()',
            'np.asarray(sparse.diags([1, 1], [-1, 1], shape=(4,4)).todense())',
            'np.vander(np.arange(4)) + np.vander(np.arange(4)).T')
    for mat_str in symmetric_mats:
        for normed in True, False:
            _check_symmetric_graph_laplacian(mat_str, normed)


def _assert_allclose_sparse(a, b, **kwargs):
    # helper function that can deal with sparse matrices
    if sparse.issparse(a):
        a = a.toarray()
    if sparse.issparse(b):
        b = a.toarray()
    assert_allclose(a, b, **kwargs)


def _check_laplacian(A, desired_L, desired_d, normed, use_out_degree):
    for arr_type in np.array, sparse.csr_matrix, sparse.coo_matrix:
        for t in int, float, complex:
            adj = arr_type(A, dtype=t)
            L = csgraph.laplacian(adj, normed=normed, return_diag=False,
                                  use_out_degree=use_out_degree)
            _assert_allclose_sparse(L, desired_L, atol=1e-12)
            L, d = csgraph.laplacian(adj, normed=normed, return_diag=True,
                                  use_out_degree=use_out_degree)
            _assert_allclose_sparse(L, desired_L, atol=1e-12)
            _assert_allclose_sparse(d, desired_d, atol=1e-12)


def test_asymmetric_laplacian():
    # adjacency matrix
    A = [[0, 1, 0],
         [4, 2, 0],
         [0, 0, 0]]

    # Laplacian matrix using out-degree
    L = [[1, -1, 0],
         [-4, 4, 0],
         [0, 0, 0]]
    d = [1, 4, 0]
    _check_laplacian(A, L, d, normed=False, use_out_degree=True)

    # normalized Laplacian matrix using out-degree
    L = [[1, -0.5, 0],
         [-2, 1, 0],
         [0, 0, 0]]
    d = [1, 2, 1]
    _check_laplacian(A, L, d, normed=True, use_out_degree=True)

    # Laplacian matrix using in-degree
    L = [[4, -1, 0],
         [-4, 1, 0],
         [0, 0, 0]]
    d = [4, 1, 0]
    _check_laplacian(A, L, d, normed=False, use_out_degree=False)

    # normalized Laplacian matrix using in-degree
    L = [[1, -0.5, 0],
         [-2, 1, 0],
         [0, 0, 0]]
    d = [2, 1, 1]
    _check_laplacian(A, L, d, normed=True, use_out_degree=False)


def test_sparse_formats():
    for fmt in ('csr', 'csc', 'coo', 'lil', 'dok', 'dia', 'bsr'):
        mat = sparse.diags([1, 1], [-1, 1], shape=(4,4), format=fmt)
        for normed in True, False:
            _check_symmetric_graph_laplacian(mat, normed)