algebraicconnectivity.py 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659
  1. """
  2. Algebraic connectivity and Fiedler vectors of undirected graphs.
  3. """
  4. from functools import partial
  5. import networkx as nx
  6. from networkx.utils import (
  7. not_implemented_for,
  8. np_random_state,
  9. reverse_cuthill_mckee_ordering,
  10. )
  11. __all__ = [
  12. "algebraic_connectivity",
  13. "fiedler_vector",
  14. "spectral_ordering",
  15. "spectral_bisection",
  16. ]
  17. class _PCGSolver:
  18. """Preconditioned conjugate gradient method.
  19. To solve Ax = b:
  20. M = A.diagonal() # or some other preconditioner
  21. solver = _PCGSolver(lambda x: A * x, lambda x: M * x)
  22. x = solver.solve(b)
  23. The inputs A and M are functions which compute
  24. matrix multiplication on the argument.
  25. A - multiply by the matrix A in Ax=b
  26. M - multiply by M, the preconditioner surrogate for A
  27. Warning: There is no limit on number of iterations.
  28. """
  29. def __init__(self, A, M):
  30. self._A = A
  31. self._M = M
  32. def solve(self, B, tol):
  33. import numpy as np
  34. # Densifying step - can this be kept sparse?
  35. B = np.asarray(B)
  36. X = np.ndarray(B.shape, order="F")
  37. for j in range(B.shape[1]):
  38. X[:, j] = self._solve(B[:, j], tol)
  39. return X
  40. def _solve(self, b, tol):
  41. import numpy as np
  42. import scipy as sp
  43. import scipy.linalg.blas # call as sp.linalg.blas
  44. A = self._A
  45. M = self._M
  46. tol *= sp.linalg.blas.dasum(b)
  47. # Initialize.
  48. x = np.zeros(b.shape)
  49. r = b.copy()
  50. z = M(r)
  51. rz = sp.linalg.blas.ddot(r, z)
  52. p = z.copy()
  53. # Iterate.
  54. while True:
  55. Ap = A(p)
  56. alpha = rz / sp.linalg.blas.ddot(p, Ap)
  57. x = sp.linalg.blas.daxpy(p, x, a=alpha)
  58. r = sp.linalg.blas.daxpy(Ap, r, a=-alpha)
  59. if sp.linalg.blas.dasum(r) < tol:
  60. return x
  61. z = M(r)
  62. beta = sp.linalg.blas.ddot(r, z)
  63. beta, rz = beta / rz, beta
  64. p = sp.linalg.blas.daxpy(p, z, a=beta)
  65. class _LUSolver:
  66. """LU factorization.
  67. To solve Ax = b:
  68. solver = _LUSolver(A)
  69. x = solver.solve(b)
  70. optional argument `tol` on solve method is ignored but included
  71. to match _PCGsolver API.
  72. """
  73. def __init__(self, A):
  74. import scipy as sp
  75. import scipy.sparse.linalg # call as sp.sparse.linalg
  76. self._LU = sp.sparse.linalg.splu(
  77. A,
  78. permc_spec="MMD_AT_PLUS_A",
  79. diag_pivot_thresh=0.0,
  80. options={"Equil": True, "SymmetricMode": True},
  81. )
  82. def solve(self, B, tol=None):
  83. import numpy as np
  84. B = np.asarray(B)
  85. X = np.ndarray(B.shape, order="F")
  86. for j in range(B.shape[1]):
  87. X[:, j] = self._LU.solve(B[:, j])
  88. return X
  89. def _preprocess_graph(G, weight):
  90. """Compute edge weights and eliminate zero-weight edges."""
  91. if G.is_directed():
  92. H = nx.MultiGraph()
  93. H.add_nodes_from(G)
  94. H.add_weighted_edges_from(
  95. ((u, v, e.get(weight, 1.0)) for u, v, e in G.edges(data=True) if u != v),
  96. weight=weight,
  97. )
  98. G = H
  99. if not G.is_multigraph():
  100. edges = (
  101. (u, v, abs(e.get(weight, 1.0))) for u, v, e in G.edges(data=True) if u != v
  102. )
  103. else:
  104. edges = (
  105. (u, v, sum(abs(e.get(weight, 1.0)) for e in G[u][v].values()))
  106. for u, v in G.edges()
  107. if u != v
  108. )
  109. H = nx.Graph()
  110. H.add_nodes_from(G)
  111. H.add_weighted_edges_from((u, v, e) for u, v, e in edges if e != 0)
  112. return H
  113. def _rcm_estimate(G, nodelist):
  114. """Estimate the Fiedler vector using the reverse Cuthill-McKee ordering."""
  115. import numpy as np
  116. G = G.subgraph(nodelist)
  117. order = reverse_cuthill_mckee_ordering(G)
  118. n = len(nodelist)
  119. index = dict(zip(nodelist, range(n)))
  120. x = np.ndarray(n, dtype=float)
  121. for i, u in enumerate(order):
  122. x[index[u]] = i
  123. x -= (n - 1) / 2.0
  124. return x
  125. def _tracemin_fiedler(L, X, normalized, tol, method):
  126. """Compute the Fiedler vector of L using the TraceMIN-Fiedler algorithm.
  127. The Fiedler vector of a connected undirected graph is the eigenvector
  128. corresponding to the second smallest eigenvalue of the Laplacian matrix
  129. of the graph. This function starts with the Laplacian L, not the Graph.
  130. Parameters
  131. ----------
  132. L : Laplacian of a possibly weighted or normalized, but undirected graph
  133. X : Initial guess for a solution. Usually a matrix of random numbers.
  134. This function allows more than one column in X to identify more than
  135. one eigenvector if desired.
  136. normalized : bool
  137. Whether the normalized Laplacian matrix is used.
  138. tol : float
  139. Tolerance of relative residual in eigenvalue computation.
  140. Warning: There is no limit on number of iterations.
  141. method : string
  142. Should be 'tracemin_pcg' or 'tracemin_lu'.
  143. Otherwise exception is raised.
  144. Returns
  145. -------
  146. sigma, X : Two NumPy arrays of floats.
  147. The lowest eigenvalues and corresponding eigenvectors of L.
  148. The size of input X determines the size of these outputs.
  149. As this is for Fiedler vectors, the zero eigenvalue (and
  150. constant eigenvector) are avoided.
  151. """
  152. import numpy as np
  153. import scipy as sp
  154. import scipy.linalg # call as sp.linalg
  155. import scipy.linalg.blas # call as sp.linalg.blas
  156. import scipy.sparse # call as sp.sparse
  157. n = X.shape[0]
  158. if normalized:
  159. # Form the normalized Laplacian matrix and determine the eigenvector of
  160. # its nullspace.
  161. e = np.sqrt(L.diagonal())
  162. # TODO: rm csr_array wrapper when spdiags array creation becomes available
  163. D = sp.sparse.csr_array(sp.sparse.spdiags(1 / e, 0, n, n, format="csr"))
  164. L = D @ L @ D
  165. e *= 1.0 / np.linalg.norm(e, 2)
  166. if normalized:
  167. def project(X):
  168. """Make X orthogonal to the nullspace of L."""
  169. X = np.asarray(X)
  170. for j in range(X.shape[1]):
  171. X[:, j] -= (X[:, j] @ e) * e
  172. else:
  173. def project(X):
  174. """Make X orthogonal to the nullspace of L."""
  175. X = np.asarray(X)
  176. for j in range(X.shape[1]):
  177. X[:, j] -= X[:, j].sum() / n
  178. if method == "tracemin_pcg":
  179. D = L.diagonal().astype(float)
  180. solver = _PCGSolver(lambda x: L @ x, lambda x: D * x)
  181. elif method == "tracemin_lu":
  182. # Convert A to CSC to suppress SparseEfficiencyWarning.
  183. A = sp.sparse.csc_array(L, dtype=float, copy=True)
  184. # Force A to be nonsingular. Since A is the Laplacian matrix of a
  185. # connected graph, its rank deficiency is one, and thus one diagonal
  186. # element needs to modified. Changing to infinity forces a zero in the
  187. # corresponding element in the solution.
  188. i = (A.indptr[1:] - A.indptr[:-1]).argmax()
  189. A[i, i] = float("inf")
  190. solver = _LUSolver(A)
  191. else:
  192. raise nx.NetworkXError(f"Unknown linear system solver: {method}")
  193. # Initialize.
  194. Lnorm = abs(L).sum(axis=1).flatten().max()
  195. project(X)
  196. W = np.ndarray(X.shape, order="F")
  197. while True:
  198. # Orthonormalize X.
  199. X = np.linalg.qr(X)[0]
  200. # Compute iteration matrix H.
  201. W[:, :] = L @ X
  202. H = X.T @ W
  203. sigma, Y = sp.linalg.eigh(H, overwrite_a=True)
  204. # Compute the Ritz vectors.
  205. X = X @ Y
  206. # Test for convergence exploiting the fact that L * X == W * Y.
  207. res = sp.linalg.blas.dasum(W @ Y[:, 0] - sigma[0] * X[:, 0]) / Lnorm
  208. if res < tol:
  209. break
  210. # Compute X = L \ X / (X' * (L \ X)).
  211. # L \ X can have an arbitrary projection on the nullspace of L,
  212. # which will be eliminated.
  213. W[:, :] = solver.solve(X, tol)
  214. X = (sp.linalg.inv(W.T @ X) @ W.T).T # Preserves Fortran storage order.
  215. project(X)
  216. return sigma, np.asarray(X)
  217. def _get_fiedler_func(method):
  218. """Returns a function that solves the Fiedler eigenvalue problem."""
  219. import numpy as np
  220. if method == "tracemin": # old style keyword <v2.1
  221. method = "tracemin_pcg"
  222. if method in ("tracemin_pcg", "tracemin_lu"):
  223. def find_fiedler(L, x, normalized, tol, seed):
  224. q = 1 if method == "tracemin_pcg" else min(4, L.shape[0] - 1)
  225. X = np.asarray(seed.normal(size=(q, L.shape[0]))).T
  226. sigma, X = _tracemin_fiedler(L, X, normalized, tol, method)
  227. return sigma[0], X[:, 0]
  228. elif method == "lanczos" or method == "lobpcg":
  229. def find_fiedler(L, x, normalized, tol, seed):
  230. import scipy as sp
  231. import scipy.sparse # call as sp.sparse
  232. import scipy.sparse.linalg # call as sp.sparse.linalg
  233. L = sp.sparse.csc_array(L, dtype=float)
  234. n = L.shape[0]
  235. if normalized:
  236. # TODO: rm csc_array wrapping when spdiags array becomes available
  237. D = sp.sparse.csc_array(
  238. sp.sparse.spdiags(
  239. 1.0 / np.sqrt(L.diagonal()), [0], n, n, format="csc"
  240. )
  241. )
  242. L = D @ L @ D
  243. if method == "lanczos" or n < 10:
  244. # Avoid LOBPCG when n < 10 due to
  245. # https://github.com/scipy/scipy/issues/3592
  246. # https://github.com/scipy/scipy/pull/3594
  247. sigma, X = sp.sparse.linalg.eigsh(
  248. L, 2, which="SM", tol=tol, return_eigenvectors=True
  249. )
  250. return sigma[1], X[:, 1]
  251. else:
  252. X = np.asarray(np.atleast_2d(x).T)
  253. # TODO: rm csr_array wrapping when spdiags array becomes available
  254. M = sp.sparse.csr_array(sp.sparse.spdiags(1.0 / L.diagonal(), 0, n, n))
  255. Y = np.ones(n)
  256. if normalized:
  257. Y /= D.diagonal()
  258. sigma, X = sp.sparse.linalg.lobpcg(
  259. L, X, M=M, Y=np.atleast_2d(Y).T, tol=tol, maxiter=n, largest=False
  260. )
  261. return sigma[0], X[:, 0]
  262. else:
  263. raise nx.NetworkXError(f"unknown method {method!r}.")
  264. return find_fiedler
  265. @np_random_state(5)
  266. @not_implemented_for("directed")
  267. def algebraic_connectivity(
  268. G, weight="weight", normalized=False, tol=1e-8, method="tracemin_pcg", seed=None
  269. ):
  270. r"""Returns the algebraic connectivity of an undirected graph.
  271. The algebraic connectivity of a connected undirected graph is the second
  272. smallest eigenvalue of its Laplacian matrix.
  273. Parameters
  274. ----------
  275. G : NetworkX graph
  276. An undirected graph.
  277. weight : object, optional (default: None)
  278. The data key used to determine the weight of each edge. If None, then
  279. each edge has unit weight.
  280. normalized : bool, optional (default: False)
  281. Whether the normalized Laplacian matrix is used.
  282. tol : float, optional (default: 1e-8)
  283. Tolerance of relative residual in eigenvalue computation.
  284. method : string, optional (default: 'tracemin_pcg')
  285. Method of eigenvalue computation. It must be one of the tracemin
  286. options shown below (TraceMIN), 'lanczos' (Lanczos iteration)
  287. or 'lobpcg' (LOBPCG).
  288. The TraceMIN algorithm uses a linear system solver. The following
  289. values allow specifying the solver to be used.
  290. =============== ========================================
  291. Value Solver
  292. =============== ========================================
  293. 'tracemin_pcg' Preconditioned conjugate gradient method
  294. 'tracemin_lu' LU factorization
  295. =============== ========================================
  296. seed : integer, random_state, or None (default)
  297. Indicator of random number generation state.
  298. See :ref:`Randomness<randomness>`.
  299. Returns
  300. -------
  301. algebraic_connectivity : float
  302. Algebraic connectivity.
  303. Raises
  304. ------
  305. NetworkXNotImplemented
  306. If G is directed.
  307. NetworkXError
  308. If G has less than two nodes.
  309. Notes
  310. -----
  311. Edge weights are interpreted by their absolute values. For MultiGraph's,
  312. weights of parallel edges are summed. Zero-weighted edges are ignored.
  313. See Also
  314. --------
  315. laplacian_matrix
  316. Examples
  317. --------
  318. For undirected graphs algebraic connectivity can tell us if a graph is connected or not
  319. `G` is connected iff ``algebraic_connectivity(G) > 0``:
  320. >>> G = nx.complete_graph(5)
  321. >>> nx.algebraic_connectivity(G) > 0
  322. True
  323. >>> G.add_node(10) # G is no longer connected
  324. >>> nx.algebraic_connectivity(G) > 0
  325. False
  326. """
  327. if len(G) < 2:
  328. raise nx.NetworkXError("graph has less than two nodes.")
  329. G = _preprocess_graph(G, weight)
  330. if not nx.is_connected(G):
  331. return 0.0
  332. L = nx.laplacian_matrix(G)
  333. if L.shape[0] == 2:
  334. return 2.0 * L[0, 0] if not normalized else 2.0
  335. find_fiedler = _get_fiedler_func(method)
  336. x = None if method != "lobpcg" else _rcm_estimate(G, G)
  337. sigma, fiedler = find_fiedler(L, x, normalized, tol, seed)
  338. return sigma
  339. @np_random_state(5)
  340. @not_implemented_for("directed")
  341. def fiedler_vector(
  342. G, weight="weight", normalized=False, tol=1e-8, method="tracemin_pcg", seed=None
  343. ):
  344. """Returns the Fiedler vector of a connected undirected graph.
  345. The Fiedler vector of a connected undirected graph is the eigenvector
  346. corresponding to the second smallest eigenvalue of the Laplacian matrix
  347. of the graph.
  348. Parameters
  349. ----------
  350. G : NetworkX graph
  351. An undirected graph.
  352. weight : object, optional (default: None)
  353. The data key used to determine the weight of each edge. If None, then
  354. each edge has unit weight.
  355. normalized : bool, optional (default: False)
  356. Whether the normalized Laplacian matrix is used.
  357. tol : float, optional (default: 1e-8)
  358. Tolerance of relative residual in eigenvalue computation.
  359. method : string, optional (default: 'tracemin_pcg')
  360. Method of eigenvalue computation. It must be one of the tracemin
  361. options shown below (TraceMIN), 'lanczos' (Lanczos iteration)
  362. or 'lobpcg' (LOBPCG).
  363. The TraceMIN algorithm uses a linear system solver. The following
  364. values allow specifying the solver to be used.
  365. =============== ========================================
  366. Value Solver
  367. =============== ========================================
  368. 'tracemin_pcg' Preconditioned conjugate gradient method
  369. 'tracemin_lu' LU factorization
  370. =============== ========================================
  371. seed : integer, random_state, or None (default)
  372. Indicator of random number generation state.
  373. See :ref:`Randomness<randomness>`.
  374. Returns
  375. -------
  376. fiedler_vector : NumPy array of floats.
  377. Fiedler vector.
  378. Raises
  379. ------
  380. NetworkXNotImplemented
  381. If G is directed.
  382. NetworkXError
  383. If G has less than two nodes or is not connected.
  384. Notes
  385. -----
  386. Edge weights are interpreted by their absolute values. For MultiGraph's,
  387. weights of parallel edges are summed. Zero-weighted edges are ignored.
  388. See Also
  389. --------
  390. laplacian_matrix
  391. Examples
  392. --------
  393. Given a connected graph the signs of the values in the Fiedler vector can be
  394. used to partition the graph into two components.
  395. >>> G = nx.barbell_graph(5, 0)
  396. >>> nx.fiedler_vector(G, normalized=True, seed=1)
  397. array([-0.32864129, -0.32864129, -0.32864129, -0.32864129, -0.26072899,
  398. 0.26072899, 0.32864129, 0.32864129, 0.32864129, 0.32864129])
  399. The connected components are the two 5-node cliques of the barbell graph.
  400. """
  401. import numpy as np
  402. if len(G) < 2:
  403. raise nx.NetworkXError("graph has less than two nodes.")
  404. G = _preprocess_graph(G, weight)
  405. if not nx.is_connected(G):
  406. raise nx.NetworkXError("graph is not connected.")
  407. if len(G) == 2:
  408. return np.array([1.0, -1.0])
  409. find_fiedler = _get_fiedler_func(method)
  410. L = nx.laplacian_matrix(G)
  411. x = None if method != "lobpcg" else _rcm_estimate(G, G)
  412. sigma, fiedler = find_fiedler(L, x, normalized, tol, seed)
  413. return fiedler
  414. @np_random_state(5)
  415. def spectral_ordering(
  416. G, weight="weight", normalized=False, tol=1e-8, method="tracemin_pcg", seed=None
  417. ):
  418. """Compute the spectral_ordering of a graph.
  419. The spectral ordering of a graph is an ordering of its nodes where nodes
  420. in the same weakly connected components appear contiguous and ordered by
  421. their corresponding elements in the Fiedler vector of the component.
  422. Parameters
  423. ----------
  424. G : NetworkX graph
  425. A graph.
  426. weight : object, optional (default: None)
  427. The data key used to determine the weight of each edge. If None, then
  428. each edge has unit weight.
  429. normalized : bool, optional (default: False)
  430. Whether the normalized Laplacian matrix is used.
  431. tol : float, optional (default: 1e-8)
  432. Tolerance of relative residual in eigenvalue computation.
  433. method : string, optional (default: 'tracemin_pcg')
  434. Method of eigenvalue computation. It must be one of the tracemin
  435. options shown below (TraceMIN), 'lanczos' (Lanczos iteration)
  436. or 'lobpcg' (LOBPCG).
  437. The TraceMIN algorithm uses a linear system solver. The following
  438. values allow specifying the solver to be used.
  439. =============== ========================================
  440. Value Solver
  441. =============== ========================================
  442. 'tracemin_pcg' Preconditioned conjugate gradient method
  443. 'tracemin_lu' LU factorization
  444. =============== ========================================
  445. seed : integer, random_state, or None (default)
  446. Indicator of random number generation state.
  447. See :ref:`Randomness<randomness>`.
  448. Returns
  449. -------
  450. spectral_ordering : NumPy array of floats.
  451. Spectral ordering of nodes.
  452. Raises
  453. ------
  454. NetworkXError
  455. If G is empty.
  456. Notes
  457. -----
  458. Edge weights are interpreted by their absolute values. For MultiGraph's,
  459. weights of parallel edges are summed. Zero-weighted edges are ignored.
  460. See Also
  461. --------
  462. laplacian_matrix
  463. """
  464. if len(G) == 0:
  465. raise nx.NetworkXError("graph is empty.")
  466. G = _preprocess_graph(G, weight)
  467. find_fiedler = _get_fiedler_func(method)
  468. order = []
  469. for component in nx.connected_components(G):
  470. size = len(component)
  471. if size > 2:
  472. L = nx.laplacian_matrix(G, component)
  473. x = None if method != "lobpcg" else _rcm_estimate(G, component)
  474. sigma, fiedler = find_fiedler(L, x, normalized, tol, seed)
  475. sort_info = zip(fiedler, range(size), component)
  476. order.extend(u for x, c, u in sorted(sort_info))
  477. else:
  478. order.extend(component)
  479. return order
  480. def spectral_bisection(
  481. G, weight="weight", normalized=False, tol=1e-8, method="tracemin_pcg", seed=None
  482. ):
  483. """Bisect the graph using the Fiedler vector.
  484. This method uses the Fiedler vector to bisect a graph.
  485. The partition is defined by the nodes which are associated with
  486. either positive or negative values in the vector.
  487. Parameters
  488. ----------
  489. G : NetworkX Graph
  490. weight : str, optional (default: weight)
  491. The data key used to determine the weight of each edge. If None, then
  492. each edge has unit weight.
  493. normalized : bool, optional (default: False)
  494. Whether the normalized Laplacian matrix is used.
  495. tol : float, optional (default: 1e-8)
  496. Tolerance of relative residual in eigenvalue computation.
  497. method : string, optional (default: 'tracemin_pcg')
  498. Method of eigenvalue computation. It must be one of the tracemin
  499. options shown below (TraceMIN), 'lanczos' (Lanczos iteration)
  500. or 'lobpcg' (LOBPCG).
  501. The TraceMIN algorithm uses a linear system solver. The following
  502. values allow specifying the solver to be used.
  503. =============== ========================================
  504. Value Solver
  505. =============== ========================================
  506. 'tracemin_pcg' Preconditioned conjugate gradient method
  507. 'tracemin_lu' LU factorization
  508. =============== ========================================
  509. seed : integer, random_state, or None (default)
  510. Indicator of random number generation state.
  511. See :ref:`Randomness<randomness>`.
  512. Returns
  513. --------
  514. bisection : tuple of sets
  515. Sets with the bisection of nodes
  516. Examples
  517. --------
  518. >>> G = nx.barbell_graph(3, 0)
  519. >>> nx.spectral_bisection(G)
  520. ({0, 1, 2}, {3, 4, 5})
  521. References
  522. ----------
  523. .. [1] M. E. J Newman 'Networks: An Introduction', pages 364-370
  524. Oxford University Press 2011.
  525. """
  526. import numpy as np
  527. v = nx.fiedler_vector(G, weight, normalized, tol, method, seed)
  528. nodes = np.array(list(G))
  529. pos_vals = v >= 0
  530. return set(nodes[~pos_vals]), set(nodes[pos_vals])