Algorithms#
Module to compute connected components on topological domains.
- toponetx.algorithms.components.connected_component_subcomplexes(domain: ComplexTypeVar, return_singletons: bool = True) Generator[ComplexTypeVar, None, None] [source]#
Compute connected component subcomplexes with s=1.
Same as
s_component_subcomplexes()
with s=1.- Parameters:
- domainCellComplex or CombinaorialComplex or ColoredHyperGraph
The domain for which to compute the the connected subcomplexes.
- return_singletonsbool, optional
When True, returns singletons connected components.
- Yields:
- Complex
The subcomplexes obtained as the resstriction of the input complex on its connected components.
See also
s_component_subcomplexes
Method implemented in the library to get a generator for s-component subcomplexes.
Examples
>>> CC = tnx.CellComplex() >>> CC.add_cell([2, 3, 4], rank=2) >>> CC.add_cell([5, 6, 7], rank=2) >>> list(tnx.connected_component_subcomplexes(CC)) >>> CC.add_cell([4, 5], rank=1) >>> list(tnx.connected_component_subcomplexes(CC))
- toponetx.algorithms.components.connected_components(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph, cells: Literal[True] = False, return_singletons: bool = True) Generator[set[tuple[Hashable, ...]], None, None] [source]#
- toponetx.algorithms.components.connected_components(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph, cells: Literal[False], return_singletons: bool = True) Generator[set[Hashable], None, None]
- toponetx.algorithms.components.connected_components(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph, cells: bool, return_singletons: bool = True) Generator[set[Hashable] | set[tuple[Hashable, ...]], None, None]
Compute s-connected components with s=1.
Same as s_connected_component` with s=1, but nodes returned.
- Parameters:
- domainComplex
Supported complexes are cell/combintorial and hypegraphs.
- cellsbool, default=False
If True will return cell components, if False will return node components.
- return_singletonsbool, default=True
When True, returns singletons connected components.
- Yields:
- set[Hashable] | set[tuple[Hashable, …]]
Yields subcomplexes generated by the cells (or nodes) in the cell(node) components of complex.
See also
s_connected_components
Method implemented in the library for s-connected-components.
Examples
>>> CC = tnx.CellComplex() >>> CC.add_cell([2, 3, 4], rank=2) >>> CC.add_cell([5, 6, 7], rank=2) >>> list(tnx.connected_components(CC, cells=False)) >>> CC.add_cell([4, 5], rank=1) >>> list(tnx.CC.connected_components(CC, cells=False))
- toponetx.algorithms.components.s_component_subcomplexes(domain: ComplexTypeVar, s: int = 1, cells: bool = True, return_singletons: bool = False) Generator[ComplexTypeVar, None, None] [source]#
Return a generator for the induced subcomplexes of s_connected components.
Removes singletons unless return_singletons is set to True.
- Parameters:
- domainCellComplex or CombinatorialComplex or ColoredHyperGraph
The domain for which to compute the the s-connected subcomplexes.
- sint, default=1
The number of intersections between pairwise consecutive cells.
- cellsbool, default=True
Determines if cell or node components are desired. Returns subcomplexes equal to the cell complex restricted to each set of nodes(cells) in the s-connected components or s-cell-connected components.
- return_singletonsbool, default=False
When True, returns singletons connected components.
- Yields:
- Complex
Returns subcomplexes generated by the cells (or nodes) in the s-cell(node) components of complex.
Examples
>>> CC = tnx.CellComplex() >>> CC.add_cell([2, 3, 4], rank=2) >>> CC.add_cell([5, 6, 7], rank=2) >>> list(tnx.s_component_subcomplexes(CC, 1, cells=False)) >>> CCC = CC.to_combinatorial_complex() >>> list(tnx.s_component_subcomplexes(CCC, s=1, cells=False)) >>> CHG = CC.to_colored_hypergraph() >>> list(tnx.s_component_subcomplexes(CHG, s=1, cells=False)) >>> CC.add_cell([4, 5], rank=1) >>> list(tnx.s_component_subcomplexes(CC, s=1, cells=False))
- toponetx.algorithms.components.s_connected_components(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph, s: int, cells: Literal[True] = True, return_singletons: bool = False) Generator[set[tuple[Hashable, ...]], None, None] [source]#
- toponetx.algorithms.components.s_connected_components(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph, s: int, cells: Literal[False], return_singletons: bool = False) Generator[set[Hashable], None, None]
- toponetx.algorithms.components.s_connected_components(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph, s: int, cells: bool, return_singletons: bool = False) Generator[set[Hashable] | set[tuple[Hashable, ...]], None, None]
Return generator for the s-connected components.
- Parameters:
- domainCellComplex or CombinatorialComplex or ColoredHyperGraph
The domain on which to compute the s-connected components.
- sint
The number of intersections between pairwise consecutive cells.
- cellsbool, default=True
If True will return cell components, if False will return node components.
- return_singletonsbool, default=False
When True, returns singleton connected components.
- Yields:
- set[Hashable] or set[tuple[Hashable, …]]
Returns sets of uids of the cells (or nodes) in the s-cells(node) components of the complex.
Notes
If cells=True, this method returns the s-cell-connected components as lists of lists of cell uids. An s-cell-component has the property that for any two cells e1 and e2 there is a sequence of cells starting with e1 and ending with e2 such that pairwise adjacent cells in the sequence intersect in at least s nodes. If s=1 these are the path components of the cell complex.
If cells=False this method returns s-node-connected components. A list of sets of uids of the nodes which are s-walk connected. Two nodes v1 and v2 are s-walk-connected if there is a sequence of nodes starting with v1 and ending with v2 such that pairwise adjacent nodes in the sequence share s cells. If s=1 these are the path components of the cell complex .
Examples
>>> CC = tnx.CellComplex() >>> CC.add_cell([2, 3, 4], rank=2) >>> CC.add_cell([5, 6, 7], rank=2) >>> list(tnx.s_connected_components(CC, s=1, cells=False)) [{2, 3, 4}, {5, 6, 7}] >>> list(tnx.s_connected_components(CC, s=1, cells=True)) [{(2, 3), (2, 3, 4), (2, 4), (3, 4)}, {(5, 6), (5, 6, 7), (5, 7), (6, 7)}] >>> CHG = CC.to_colored_hypergraph() >>> list(tnx.s_connected_components(CHG, s=1, cells=False)) >>> CC.add_cell([4, 5], rank=1) >>> list(tnx.s_connected_components(CC, s=1, cells=False)) [{2, 3, 4, 5, 6, 7}] >>> CCC = CC.to_combinatorial_complex() >>> list(tnx.s_connected_components(CCC, s=1, cells=False))
Module to distance measures on topological domains.
- toponetx.algorithms.distance_measures.cell_diameter(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph, s: int = 1) int [source]#
Return the length of the longest shortest s-walk between cells.
- Parameters:
- domainCellComplex or CombinatorialComplex or ColoredHyperGraph
Supported complexes are cell/combintorial and hypegraphs.
- sint, default=1
The number of intersections between pairwise consecutive cells.
- Returns:
- int
Returns the length of the longest shortest s-walk between cells.
- Raises:
- RuntimeError
If cell complex is not s-cell-connected
Notes
Two cells are s-coadjacent if they share s nodes. Two nodes e_start and e_end are s-walk connected if there is a sequence of cells (one or two dimensional) e_start, e_1, e_2, … e_n-1, e_end such that consecutive cells are s-coadjacent. If the cell complex is not connected, an error will be raised.
Examples
>>> CC = tnx.CellComplex() >>> CC.add_cell([2, 3, 4], rank=2) >>> CC.add_cell([5, 6, 7], rank=2) >>> CC.add_cell([2, 5], rank=1) >>> tnx.cell_diameter(CC) >>> CCC = CC.to_combinatorial_complex() >>> tnx.cell_diameter(CCC) >>> CHG = CC.to_colored_hypergraph() >>> tnx.cell_diameter(CHG)
- toponetx.algorithms.distance_measures.cell_diameters(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph, s: int = 1) tuple[list[int], list[set[int]]] [source]#
Return the cell diameters of the s_cell_connected component subgraphs.
- Parameters:
- domainCellComplex or CombinatorialComplex or ColoredHyperGraph
Supported complexes are cell/combintorial and hypegraphs.
- sint, optional
The number of intersections between pairwise consecutive cells.
- Returns:
- list of diameterslist
List of cell_diameters for s-cell component subcomplexes in the cell complex.
- list of componentlist
List of the cell uids in the s-cell component subcomplexes.
Examples
>>> CC = CellComplex() >>> CC.add_cell([2, 3, 4], rank=2) >>> CC.add_cell([5, 6, 7], rank=2) >>> tnx.cell_diameters(CC) >>> CCC = CC.to_combinatorial_complex() >>> tnx.cell_diameters(CCC) >>> CHG = CC.to_colored_hypergraph() >>> tnx.cell_diameters(CHG)
- toponetx.algorithms.distance_measures.diameter(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph) int [source]#
Return length of the longest shortest s-walk between nodes.
- Parameters:
- domainCellComplex or CombinatorialComplex or ColoredHyperGraph
Supported complexes are cell/combintorial and hypegraphs.
- Returns:
- int
The diameter of the longest shortest s-walk between nodes.
- Raises:
- RuntimeError
If the cell complex is not s-cell-connected
Notes
Two nodes are s-adjacent if they share s cells. Two nodes v_start and v_end are s-walk connected if there is a sequence of nodes v_start, v_1, v_2, … v_n-1, v_end such that consecutive nodes are s-adjacent. If the cell complex is not connected, an error will be raised.
Examples
>>> CC = tnx.CellComplex() >>> CC.add_cell([2, 3, 4], rank=2) >>> CC.add_cell([5, 6, 7], rank=2) >>> CC.add_cell([2, 5], rank=2) >>> tnx.diameter(CC) >>> CCC = CC.to_combinatorial_complex() >>> tnx.diameter(CCC) >>> CHG = CC.to_colored_hypergraph() >>> tnx.diameter(CHG)
- toponetx.algorithms.distance_measures.node_diameters(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph) tuple[list[int], list[set[Hashable]]] [source]#
Return the node diameters of the connected components in cell complex.
- Parameters:
- domainCellComplex or CombinaorialComplex or ColoredHyperGraph
The complex to be used to generate the node diameters for.
- Returns:
- diameterslist
List of the diameters of the s-components.
- componentslist
List of the s-component nodes.
Examples
>>> CC = tnx.CellComplex() >>> CC.add_cell([2, 3, 4], rank=2) >>> CC.add_cell([5, 6, 7], rank=2) >>> tnx.node_diameters(CC) >>> CCC = CC.to_combinatorial_complex() >>> tnx.node_diameters(CCC) >>> CHG = CC.to_colored_hypergraph() >>> tnx.node_diameters(CHG)
Module to compute distance between nodes or cells on topological domains.
- toponetx.algorithms.distance.cell_distance(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph, source: Iterable | HyperEdge | Cell, target: Iterable | HyperEdge | Cell, s: int = 1) int [source]#
Return the shortest s-walk distance between two cells in the cell complex.
- Parameters:
- domainCellComplex or CombinatorialComplex or ColoredHyperGraph
The domain on which to compute the s-walk distance between source and target cells.
- sourceIterable or HyperEdge or Cell
An Iterable representing a cell in the input complex cell complex.
- targetIterable or HyperEdge or Cell
An Iterable representing a cell in the input complex cell complex.
- sint
The number of intersections between pairwise consecutive cells.
- Returns:
- int
Shortest s-walk cell distance between source and target. A shortest s-walk is computed as a sequence of cells, the s-walk distance is the number of cells in the sequence minus 1.
- Raises:
- TopoNetXNoPath
If no s-walk path exists between source and target cells.
See also
distance
Method implemented in the library that returns shortest s-walk distance between two nodes in the cell complex.
Notes
The s-distance is the shortest s-walk length between the cells. An s-walk between cells is a sequence of cells such that consecutive pairwise cells intersect in at least s nodes. The length of the shortest s-walk is 1 less than the number of cells in the path sequence.
Uses the networkx shortest_path_length method on the graph generated by the s-cell_adjacency matrix.
Examples
>>> CC = tnx.CellComplex() >>> CC.add_cell([2, 3, 4], rank=2) >>> CC.add_cell([5, 6, 7], rank=2) >>> CC.add_cell([5, 2], rank=1) >>> tnx.cell_distance(CC, [2, 3], [6, 7]) >>> CHG = CC.to_colored_hypergraph() >>> tnx.cell_distance(CHG, (frozenset({2, 3}), 0), (frozenset({6, 7}), 0)) >>> CCC = CC.to_combinatorial_complex() >>> tnx.cell_distance(CCC, frozenset({2, 3}), frozenset({6, 7}))
- toponetx.algorithms.distance.distance(domain: CellComplex | CombinatorialComplex | ColoredHyperGraph, source: Hashable, target: Hashable, s: int = 1) int [source]#
Return shortest s-walk distance between two nodes in the cell complex.
- Parameters:
- domainCellComplex or CombinaorialComplex or ColoredHyperGraph
The domain on which to compute the s-walk distance between source and target.
- sourceHashable
A node in the input complex.
- targetHashable
A node in the input complex.
- sint, optional
The number of intersections between pairwise consecutive cells.
- Returns:
- int
Shortest s-walk distance between two nodes in the cell complex.
- Raises:
- TopoNetXNoPath
If no s-walk path exists between source and target nodes.
See also
cell_distance
Method implemented in the library that returns the shortest s-walk distance between two cells in the cell complex.
Notes
The s-distance is the shortest s-walk length between the nodes. An s-walk between nodes is a sequence of nodes that pairwise share at least s cells. The length of the shortest s-walk is 1 less than the number of nodes in the path sequence.
Uses the networkx shortest_path_length method on the graph generated by the s-adjacency matrix.
Examples
>>> CC = tnx.CellComplex() >>> CC.add_cell([2, 3, 4], rank=2) >>> CC.add_cell([5, 6, 7], rank=2) >>> list(tnx.node_diameters(CC)) >>> CCC = CC.to_combinatorial_complex() >>> list(tnx.node_diameters(CCC)) >>> CHG = CC.to_colored_hypergraph() >>> list(tnx.node_diameters(CHG))
Module to compute spectra.
- toponetx.algorithms.spectrum.cell_complex_adjacency_spectrum(CC: CellComplex, rank: int)[source]#
Return eigenvalues of the adjacency matrix of the cell complex.
- Parameters:
- CCtoponetx.classes.CellComplex
The cell complex for which to compute the spectrum.
- rankint
Rank of the cells to take the adjacency from: - 0-cells are nodes - 1-cells are edges - 2-cells are polygons Currently, no cells of rank > 2 are supported.
- Returns:
- numpy.ndarray
Eigenvalues.
Examples
>>> CC = tnx.CellComplex() >>> CC.add_cell([1, 2, 3, 4], rank=2) >>> CC.add_cell([2, 3, 4, 5], rank=2) >>> CC.add_cell([5, 6, 7, 8], rank=2) >>> tnx.cell_complex_adjacency_spectrum(CC, 1)
- toponetx.algorithms.spectrum.cell_complex_hodge_laplacian_spectrum(CC: CellComplex, rank: int, weight: str | None = None) ndarray [source]#
Return eigenvalues of the Laplacian of G.
- Parameters:
- CCtoponetx.classes.CellComplex
The cell complex for which to compute the spectrum.
- rankint
Rank of the cells to compute the Hodge Laplacian for.
- weightstr, optional
If None, then each cell has weight 1.
- Returns:
- numpy.ndarray
Eigenvalues.
Examples
>>> CC = tnx.CellComplex() >>> CC.add_cell([1, 2, 3, 4], rank=2) >>> CC.add_cell([2, 3, 4, 5], rank=2) >>> CC.add_cell([5, 6, 7, 8], rank=2) >>> tnx.cell_complex_hodge_laplacian_spectrum(CC, 1)
- toponetx.algorithms.spectrum.combinatorial_complex_adjacency_spectrum(CCC: CombinatorialComplex, rank: int, via_rank: int) ndarray [source]#
Return eigenvalues of the adjacency matrix of the combinatorial complex.
- Parameters:
- CCCtoponetx.classes.CombinatorialComplex
The combinatorial complex for which to compute the spectrum.
- rank, via_rankint
Rank of the cells to compute the adjacency spectrum for.
- Returns:
- numpy.ndarray
Eigenvalues.
Examples
>>> CCC = tnx.CombinatorialComplex(cells=[[1, 2, 3], [2, 3], [0]], ranks=[2, 1, 0]) >>> tnx.combinatorial_complex_adjacency_spectrum(CCC, 0, 2)
- toponetx.algorithms.spectrum.hodge_laplacian_eigenvectors(hodge_laplacian, n_components: int) tuple[ndarray, ndarray] [source]#
Compute the first k eigenvectors of the hodge laplacian matrix.
- Parameters:
- hodge_laplacianscipy sparse matrix
Hodge laplacian.
- n_componentsint
Number of eigenvectors that should be computed.
- Returns:
- numpy.ndarray
First k eigevals and eigenvec associated with the hodge laplacian matrix.
Examples
>>> SC = tnx.SimplicialComplex([[1, 2, 3], [2, 3, 5], [0, 1]]) >>> L1 = SC.hodge_laplacian_matrix(1) >>> vals, vecs = tnx.hodge_laplacian_eigenvectors(L1, 2)
- toponetx.algorithms.spectrum.laplacian_beltrami_eigenvectors(SC: SimplicialComplex, mode: str = 'fem') tuple[ndarray, ndarray] [source]#
Compute the first k eigenvectors of the laplacian beltrami matrix.
- Parameters:
- SCtoponetx.classes.SimplicialComplex
The simplicial complex for which to compute the beltrami eigenvectors.
- mode{“fem”, “sphara”}, default=”fem”
Mode of the spharapy basis.
- Returns:
- k_eigenvectorsnumpy.ndarray
First k Eigenvectors associated with the hodge laplacian matrix.
- k_eigenvalsnumpy.ndarray
First k Eigenvalues associated with the hodge laplacian matrix.
- Raises:
- RuntimeError
If package spharapy is not installed.
Examples
>>> SC = tnx.datasets.stanford_bunny("simplicial") >>> eigenvectors, eigenvalues = tnx.laplacian_beltrami_eigenvectors(SC)
- toponetx.algorithms.spectrum.laplacian_spectrum(matrix) ndarray [source]#
Return eigenvalues of the Laplacian matrix.
- Parameters:
- matrixscipy sparse matrix
The Laplacian matrix.
- Returns:
- numpy.ndarray
Eigenvalues.
- toponetx.algorithms.spectrum.set_hodge_laplacian_eigenvector_attrs(SC: SimplicialComplex, dim: int, n_components: int, laplacian_type: Literal['up', 'down', 'hodge'] = 'hodge', normalized: bool = True) None [source]#
Set the hodge laplacian eigenvectors as simplex attributes.
- Parameters:
- SCSimplicialComplex
The simplicial complex for which to compute the hodge laplacian eigenvectors.
- dimint
Dimension of the hodge laplacian to be computed.
- n_componentsint
The number of eigenvectors to be computed.
- laplacian_type{“up”, “down”, “hodge”}, default=”hodge”
Type of hodge matrix to be computed.
- normalizedbool, default=True
Normalize the eigenvector or not.
Examples
>>> SC = tnx.SimplicialComplex([[1, 2, 3], [2, 3, 5], [0, 1]]) >>> tnx.set_hodge_laplacian_eigenvector_attrs(SC, 1, 2, "down") >>> SC.get_simplex_attributes("0.th_eigen", 1)
- toponetx.algorithms.spectrum.simplicial_complex_adjacency_spectrum(SC: SimplicialComplex, rank: int, weight: str | None = None) ndarray [source]#
Return eigenvalues of the Laplacian of G.
- Parameters:
- SCtoponetx.classes.SimplicialComplex
The simplicial complex for which to compute the spectrum.
- rankint
Rank of the adjacency matrix.
- weightstr, optional
If None, then each cell has weight 1.
- Returns:
- numpy.ndarray
Eigenvalues.
- toponetx.algorithms.spectrum.simplicial_complex_hodge_laplacian_spectrum(SC: SimplicialComplex, rank: int, weight: str = 'weight') ndarray [source]#
Return eigenvalues of the Laplacian of G.
- Parameters:
- SCtoponetx.classes.SimplicialComplex
The simplicial complex for which to compute the spectrum.
- rankint
Rank of the Hodge Laplacian.
- weightstr or None, optional (default=’weight’)
If None, then each cell has weight 1.
- Returns:
- numpy.ndarray
Eigenvalues.
Examples
>>> SC = tnx.SimplicialComplex([[1, 2, 3], [2, 3, 5], [0, 1]]) >>> spectrum = tnx.simplicial_complex_hodge_laplacian_spectrum(SC, 1)