algebraicconnectivity.py 20 KB

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