_laplacian.py 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555
  1. """
  2. Laplacian of a compressed-sparse graph
  3. """
  4. import numpy as np
  5. from scipy.sparse import isspmatrix
  6. from scipy.sparse.linalg import LinearOperator
  7. ###############################################################################
  8. # Graph laplacian
  9. def laplacian(
  10. csgraph,
  11. normed=False,
  12. return_diag=False,
  13. use_out_degree=False,
  14. *,
  15. copy=True,
  16. form="array",
  17. dtype=None,
  18. symmetrized=False,
  19. ):
  20. """
  21. Return the Laplacian of a directed graph.
  22. Parameters
  23. ----------
  24. csgraph : array_like or sparse matrix, 2 dimensions
  25. compressed-sparse graph, with shape (N, N).
  26. normed : bool, optional
  27. If True, then compute symmetrically normalized Laplacian.
  28. Default: False.
  29. return_diag : bool, optional
  30. If True, then also return an array related to vertex degrees.
  31. Default: False.
  32. use_out_degree : bool, optional
  33. If True, then use out-degree instead of in-degree.
  34. This distinction matters only if the graph is asymmetric.
  35. Default: False.
  36. copy: bool, optional
  37. If False, then change `csgraph` in place if possible,
  38. avoiding doubling the memory use.
  39. Default: True, for backward compatibility.
  40. form: 'array', or 'function', or 'lo'
  41. Determines the format of the output Laplacian:
  42. * 'array' is a numpy array;
  43. * 'function' is a pointer to evaluating the Laplacian-vector
  44. or Laplacian-matrix product;
  45. * 'lo' results in the format of the `LinearOperator`.
  46. Choosing 'function' or 'lo' always avoids doubling
  47. the memory use, ignoring `copy` value.
  48. Default: 'array', for backward compatibility.
  49. dtype: None or one of numeric numpy dtypes, optional
  50. The dtype of the output. If ``dtype=None``, the dtype of the
  51. output matches the dtype of the input csgraph, except for
  52. the case ``normed=True`` and integer-like csgraph, where
  53. the output dtype is 'float' allowing accurate normalization,
  54. but dramatically increasing the memory use.
  55. Default: None, for backward compatibility.
  56. symmetrized: bool, optional
  57. If True, then the output Laplacian is symmetric/Hermitian.
  58. The symmetrization is done by ``csgraph + csgraph.T.conj``
  59. without dividing by 2 to preserve integer dtypes if possible
  60. prior to the construction of the Laplacian.
  61. The symmetrization will increase the memory footprint of
  62. sparse matrices unless the sparsity pattern is symmetric or
  63. `form` is 'function' or 'lo'.
  64. Default: False, for backward compatibility.
  65. Returns
  66. -------
  67. lap : ndarray, or sparse matrix, or `LinearOperator`
  68. The N x N Laplacian of csgraph. It will be a NumPy array (dense)
  69. if the input was dense, or a sparse matrix otherwise, or
  70. the format of a function or `LinearOperator` if
  71. `form` equals 'function' or 'lo', respectively.
  72. diag : ndarray, optional
  73. The length-N main diagonal of the Laplacian matrix.
  74. For the normalized Laplacian, this is the array of square roots
  75. of vertex degrees or 1 if the degree is zero.
  76. Notes
  77. -----
  78. The Laplacian matrix of a graph is sometimes referred to as the
  79. "Kirchhoff matrix" or just the "Laplacian", and is useful in many
  80. parts of spectral graph theory.
  81. In particular, the eigen-decomposition of the Laplacian can give
  82. insight into many properties of the graph, e.g.,
  83. is commonly used for spectral data embedding and clustering.
  84. The constructed Laplacian doubles the memory use if ``copy=True`` and
  85. ``form="array"`` which is the default.
  86. Choosing ``copy=False`` has no effect unless ``form="array"``
  87. or the matrix is sparse in the ``coo`` format, or dense array, except
  88. for the integer input with ``normed=True`` that forces the float output.
  89. Sparse input is reformatted into ``coo`` if ``form="array"``,
  90. which is the default.
  91. If the input adjacency matrix is not symmetic, the Laplacian is
  92. also non-symmetric unless ``symmetrized=True`` is used.
  93. Diagonal entries of the input adjacency matrix are ignored and
  94. replaced with zeros for the purpose of normalization where ``normed=True``.
  95. The normalization uses the inverse square roots of row-sums of the input
  96. adjacency matrix, and thus may fail if the row-sums contain
  97. negative or complex with a non-zero imaginary part values.
  98. The normalization is symmetric, making the normalized Laplacian also
  99. symmetric if the input csgraph was symmetric.
  100. References
  101. ----------
  102. .. [1] Laplacian matrix. https://en.wikipedia.org/wiki/Laplacian_matrix
  103. Examples
  104. --------
  105. >>> import numpy as np
  106. >>> from scipy.sparse import csgraph
  107. Our first illustration is the symmetric graph
  108. >>> G = np.arange(4) * np.arange(4)[:, np.newaxis]
  109. >>> G
  110. array([[0, 0, 0, 0],
  111. [0, 1, 2, 3],
  112. [0, 2, 4, 6],
  113. [0, 3, 6, 9]])
  114. and its symmetric Laplacian matrix
  115. >>> csgraph.laplacian(G)
  116. array([[ 0, 0, 0, 0],
  117. [ 0, 5, -2, -3],
  118. [ 0, -2, 8, -6],
  119. [ 0, -3, -6, 9]])
  120. The non-symmetric graph
  121. >>> G = np.arange(9).reshape(3, 3)
  122. >>> G
  123. array([[0, 1, 2],
  124. [3, 4, 5],
  125. [6, 7, 8]])
  126. has different row- and column sums, resulting in two varieties
  127. of the Laplacian matrix, using an in-degree, which is the default
  128. >>> L_in_degree = csgraph.laplacian(G)
  129. >>> L_in_degree
  130. array([[ 9, -1, -2],
  131. [-3, 8, -5],
  132. [-6, -7, 7]])
  133. or alternatively an out-degree
  134. >>> L_out_degree = csgraph.laplacian(G, use_out_degree=True)
  135. >>> L_out_degree
  136. array([[ 3, -1, -2],
  137. [-3, 8, -5],
  138. [-6, -7, 13]])
  139. Constructing a symmetric Laplacian matrix, one can add the two as
  140. >>> L_in_degree + L_out_degree.T
  141. array([[ 12, -4, -8],
  142. [ -4, 16, -12],
  143. [ -8, -12, 20]])
  144. or use the ``symmetrized=True`` option
  145. >>> csgraph.laplacian(G, symmetrized=True)
  146. array([[ 12, -4, -8],
  147. [ -4, 16, -12],
  148. [ -8, -12, 20]])
  149. that is equivalent to symmetrizing the original graph
  150. >>> csgraph.laplacian(G + G.T)
  151. array([[ 12, -4, -8],
  152. [ -4, 16, -12],
  153. [ -8, -12, 20]])
  154. The goal of normalization is to make the non-zero diagonal entries
  155. of the Laplacian matrix to be all unit, also scaling off-diagonal
  156. entries correspondingly. The normalization can be done manually, e.g.,
  157. >>> G = np.array([[0, 1, 1], [1, 0, 1], [1, 1, 0]])
  158. >>> L, d = csgraph.laplacian(G, return_diag=True)
  159. >>> L
  160. array([[ 2, -1, -1],
  161. [-1, 2, -1],
  162. [-1, -1, 2]])
  163. >>> d
  164. array([2, 2, 2])
  165. >>> scaling = np.sqrt(d)
  166. >>> scaling
  167. array([1.41421356, 1.41421356, 1.41421356])
  168. >>> (1/scaling)*L*(1/scaling)
  169. array([[ 1. , -0.5, -0.5],
  170. [-0.5, 1. , -0.5],
  171. [-0.5, -0.5, 1. ]])
  172. Or using ``normed=True`` option
  173. >>> L, d = csgraph.laplacian(G, return_diag=True, normed=True)
  174. >>> L
  175. array([[ 1. , -0.5, -0.5],
  176. [-0.5, 1. , -0.5],
  177. [-0.5, -0.5, 1. ]])
  178. which now instead of the diagonal returns the scaling coefficients
  179. >>> d
  180. array([1.41421356, 1.41421356, 1.41421356])
  181. Zero scaling coefficients are substituted with 1s, where scaling
  182. has thus no effect, e.g.,
  183. >>> G = np.array([[0, 0, 0], [0, 0, 1], [0, 1, 0]])
  184. >>> G
  185. array([[0, 0, 0],
  186. [0, 0, 1],
  187. [0, 1, 0]])
  188. >>> L, d = csgraph.laplacian(G, return_diag=True, normed=True)
  189. >>> L
  190. array([[ 0., -0., -0.],
  191. [-0., 1., -1.],
  192. [-0., -1., 1.]])
  193. >>> d
  194. array([1., 1., 1.])
  195. Only the symmetric normalization is implemented, resulting
  196. in a symmetric Laplacian matrix if and only if its graph is symmetric
  197. and has all non-negative degrees, like in the examples above.
  198. The output Laplacian matrix is by default a dense array or a sparse matrix
  199. inferring its shape, format, and dtype from the input graph matrix:
  200. >>> G = np.array([[0, 1, 1], [1, 0, 1], [1, 1, 0]]).astype(np.float32)
  201. >>> G
  202. array([[0., 1., 1.],
  203. [1., 0., 1.],
  204. [1., 1., 0.]], dtype=float32)
  205. >>> csgraph.laplacian(G)
  206. array([[ 2., -1., -1.],
  207. [-1., 2., -1.],
  208. [-1., -1., 2.]], dtype=float32)
  209. but can alternatively be generated matrix-free as a LinearOperator:
  210. >>> L = csgraph.laplacian(G, form="lo")
  211. >>> L
  212. <3x3 _CustomLinearOperator with dtype=float32>
  213. >>> L(np.eye(3))
  214. array([[ 2., -1., -1.],
  215. [-1., 2., -1.],
  216. [-1., -1., 2.]])
  217. or as a lambda-function:
  218. >>> L = csgraph.laplacian(G, form="function")
  219. >>> L
  220. <function _laplace.<locals>.<lambda> at 0x0000012AE6F5A598>
  221. >>> L(np.eye(3))
  222. array([[ 2., -1., -1.],
  223. [-1., 2., -1.],
  224. [-1., -1., 2.]])
  225. The Laplacian matrix is used for
  226. spectral data clustering and embedding
  227. as well as for spectral graph partitioning.
  228. Our final example illustrates the latter
  229. for a noisy directed linear graph.
  230. >>> from scipy.sparse import diags, random
  231. >>> from scipy.sparse.linalg import lobpcg
  232. Create a directed linear graph with ``N=35`` vertices
  233. using a sparse adjacency matrix ``G``:
  234. >>> N = 35
  235. >>> G = diags(np.ones(N-1), 1, format="csr")
  236. Fix a random seed ``rng`` and add a random sparse noise to the graph ``G``:
  237. >>> rng = np.random.default_rng()
  238. >>> G += 1e-2 * random(N, N, density=0.1, random_state=rng)
  239. Set initial approximations for eigenvectors:
  240. >>> X = rng.random((N, 2))
  241. The constant vector of ones is always a trivial eigenvector
  242. of the non-normalized Laplacian to be filtered out:
  243. >>> Y = np.ones((N, 1))
  244. Alternating (1) the sign of the graph weights allows determining
  245. labels for spectral max- and min- cuts in a single loop.
  246. Since the graph is undirected, the option ``symmetrized=True``
  247. must be used in the construction of the Laplacian.
  248. The option ``normed=True`` cannot be used in (2) for the negative weights
  249. here as the symmetric normalization evaluates square roots.
  250. The option ``form="lo"`` in (2) is matrix-free, i.e., guarantees
  251. a fixed memory footprint and read-only access to the graph.
  252. Calling the eigenvalue solver ``lobpcg`` (3) computes the Fiedler vector
  253. that determines the labels as the signs of its components in (5).
  254. Since the sign in an eigenvector is not deterministic and can flip,
  255. we fix the sign of the first component to be always +1 in (4).
  256. >>> for cut in ["max", "min"]:
  257. ... G = -G # 1.
  258. ... L = csgraph.laplacian(G, symmetrized=True, form="lo") # 2.
  259. ... _, eves = lobpcg(L, X, Y=Y, largest=False, tol=1e-3) # 3.
  260. ... eves *= np.sign(eves[0, 0]) # 4.
  261. ... print(cut + "-cut labels:\\n", 1 * (eves[:, 0]>0)) # 5.
  262. max-cut labels:
  263. [1 0 1 0 1 0 1 0 1 0 1 0 1 0 1 0 1 0 1 0 1 0 1 0 1 0 1 0 1 0 1 0 1 0 1]
  264. min-cut labels:
  265. [1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0]
  266. As anticipated for a (slightly noisy) linear graph,
  267. the max-cut strips all the edges of the graph coloring all
  268. odd vertices into one color and all even vertices into another one,
  269. while the balanced min-cut partitions the graph
  270. in the middle by deleting a single edge.
  271. Both determined partitions are optimal.
  272. """
  273. if csgraph.ndim != 2 or csgraph.shape[0] != csgraph.shape[1]:
  274. raise ValueError('csgraph must be a square matrix or array')
  275. if normed and (
  276. np.issubdtype(csgraph.dtype, np.signedinteger)
  277. or np.issubdtype(csgraph.dtype, np.uint)
  278. ):
  279. csgraph = csgraph.astype(np.float64)
  280. if form == "array":
  281. create_lap = (
  282. _laplacian_sparse if isspmatrix(csgraph) else _laplacian_dense
  283. )
  284. else:
  285. create_lap = (
  286. _laplacian_sparse_flo
  287. if isspmatrix(csgraph)
  288. else _laplacian_dense_flo
  289. )
  290. degree_axis = 1 if use_out_degree else 0
  291. lap, d = create_lap(
  292. csgraph,
  293. normed=normed,
  294. axis=degree_axis,
  295. copy=copy,
  296. form=form,
  297. dtype=dtype,
  298. symmetrized=symmetrized,
  299. )
  300. if return_diag:
  301. return lap, d
  302. return lap
  303. def _setdiag_dense(m, d):
  304. step = len(d) + 1
  305. m.flat[::step] = d
  306. def _laplace(m, d):
  307. return lambda v: v * d[:, np.newaxis] - m @ v
  308. def _laplace_normed(m, d, nd):
  309. laplace = _laplace(m, d)
  310. return lambda v: nd[:, np.newaxis] * laplace(v * nd[:, np.newaxis])
  311. def _laplace_sym(m, d):
  312. return (
  313. lambda v: v * d[:, np.newaxis]
  314. - m @ v
  315. - np.transpose(np.conjugate(np.transpose(np.conjugate(v)) @ m))
  316. )
  317. def _laplace_normed_sym(m, d, nd):
  318. laplace_sym = _laplace_sym(m, d)
  319. return lambda v: nd[:, np.newaxis] * laplace_sym(v * nd[:, np.newaxis])
  320. def _linearoperator(mv, shape, dtype):
  321. return LinearOperator(matvec=mv, matmat=mv, shape=shape, dtype=dtype)
  322. def _laplacian_sparse_flo(graph, normed, axis, copy, form, dtype, symmetrized):
  323. # The keyword argument `copy` is unused and has no effect here.
  324. del copy
  325. if dtype is None:
  326. dtype = graph.dtype
  327. graph_sum = graph.sum(axis=axis).getA1()
  328. graph_diagonal = graph.diagonal()
  329. diag = graph_sum - graph_diagonal
  330. if symmetrized:
  331. graph_sum += graph.sum(axis=1 - axis).getA1()
  332. diag = graph_sum - graph_diagonal - graph_diagonal
  333. if normed:
  334. isolated_node_mask = diag == 0
  335. w = np.where(isolated_node_mask, 1, np.sqrt(diag))
  336. if symmetrized:
  337. md = _laplace_normed_sym(graph, graph_sum, 1.0 / w)
  338. else:
  339. md = _laplace_normed(graph, graph_sum, 1.0 / w)
  340. if form == "function":
  341. return md, w.astype(dtype, copy=False)
  342. elif form == "lo":
  343. m = _linearoperator(md, shape=graph.shape, dtype=dtype)
  344. return m, w.astype(dtype, copy=False)
  345. else:
  346. raise ValueError(f"Invalid form: {form!r}")
  347. else:
  348. if symmetrized:
  349. md = _laplace_sym(graph, graph_sum)
  350. else:
  351. md = _laplace(graph, graph_sum)
  352. if form == "function":
  353. return md, diag.astype(dtype, copy=False)
  354. elif form == "lo":
  355. m = _linearoperator(md, shape=graph.shape, dtype=dtype)
  356. return m, diag.astype(dtype, copy=False)
  357. else:
  358. raise ValueError(f"Invalid form: {form!r}")
  359. def _laplacian_sparse(graph, normed, axis, copy, form, dtype, symmetrized):
  360. # The keyword argument `form` is unused and has no effect here.
  361. del form
  362. if dtype is None:
  363. dtype = graph.dtype
  364. needs_copy = False
  365. if graph.format in ('lil', 'dok'):
  366. m = graph.tocoo()
  367. else:
  368. m = graph
  369. if copy:
  370. needs_copy = True
  371. if symmetrized:
  372. m += m.T.conj()
  373. w = m.sum(axis=axis).getA1() - m.diagonal()
  374. if normed:
  375. m = m.tocoo(copy=needs_copy)
  376. isolated_node_mask = (w == 0)
  377. w = np.where(isolated_node_mask, 1, np.sqrt(w))
  378. m.data /= w[m.row]
  379. m.data /= w[m.col]
  380. m.data *= -1
  381. m.setdiag(1 - isolated_node_mask)
  382. else:
  383. if m.format == 'dia':
  384. m = m.copy()
  385. else:
  386. m = m.tocoo(copy=needs_copy)
  387. m.data *= -1
  388. m.setdiag(w)
  389. return m.astype(dtype, copy=False), w.astype(dtype)
  390. def _laplacian_dense_flo(graph, normed, axis, copy, form, dtype, symmetrized):
  391. if copy:
  392. m = np.array(graph)
  393. else:
  394. m = np.asarray(graph)
  395. if dtype is None:
  396. dtype = m.dtype
  397. graph_sum = m.sum(axis=axis)
  398. graph_diagonal = m.diagonal()
  399. diag = graph_sum - graph_diagonal
  400. if symmetrized:
  401. graph_sum += m.sum(axis=1 - axis)
  402. diag = graph_sum - graph_diagonal - graph_diagonal
  403. if normed:
  404. isolated_node_mask = diag == 0
  405. w = np.where(isolated_node_mask, 1, np.sqrt(diag))
  406. if symmetrized:
  407. md = _laplace_normed_sym(m, graph_sum, 1.0 / w)
  408. else:
  409. md = _laplace_normed(m, graph_sum, 1.0 / w)
  410. if form == "function":
  411. return md, w.astype(dtype, copy=False)
  412. elif form == "lo":
  413. m = _linearoperator(md, shape=graph.shape, dtype=dtype)
  414. return m, w.astype(dtype, copy=False)
  415. else:
  416. raise ValueError(f"Invalid form: {form!r}")
  417. else:
  418. if symmetrized:
  419. md = _laplace_sym(m, graph_sum)
  420. else:
  421. md = _laplace(m, graph_sum)
  422. if form == "function":
  423. return md, diag.astype(dtype, copy=False)
  424. elif form == "lo":
  425. m = _linearoperator(md, shape=graph.shape, dtype=dtype)
  426. return m, diag.astype(dtype, copy=False)
  427. else:
  428. raise ValueError(f"Invalid form: {form!r}")
  429. def _laplacian_dense(graph, normed, axis, copy, form, dtype, symmetrized):
  430. if form != "array":
  431. raise ValueError(f'{form!r} must be "array"')
  432. if dtype is None:
  433. dtype = graph.dtype
  434. if copy:
  435. m = np.array(graph)
  436. else:
  437. m = np.asarray(graph)
  438. if dtype is None:
  439. dtype = m.dtype
  440. if symmetrized:
  441. m += m.T.conj()
  442. np.fill_diagonal(m, 0)
  443. w = m.sum(axis=axis)
  444. if normed:
  445. isolated_node_mask = (w == 0)
  446. w = np.where(isolated_node_mask, 1, np.sqrt(w))
  447. m /= w
  448. m /= w[:, np.newaxis]
  449. m *= -1
  450. _setdiag_dense(m, 1 - isolated_node_mask)
  451. else:
  452. m *= -1
  453. _setdiag_dense(m, w)
  454. return m.astype(dtype, copy=False), w.astype(dtype, copy=False)