_laplacian.py 18 KB

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