toponetx.SimplicialComplex#

class toponetx.SimplicialComplex(simplices=None, **kwargs)[source]#

Bases: Complex

Class representing a simplicial complex.

Class for construction boundary operators, Hodge Laplacians, higher order (co)adjacency operators from a collection of simplices.

A simplicial complex is a topological space of a specific kind, constructed by “gluing together” points, line segments, triangles, and their higher-dimensional counterparts. It is a generalization of the notion of a triangle in a triangulated surface, or a tetrahedron in a tetrahedralized 3-dimensional manifold. Simplicial complexes are the basic objects of study in combinatorial topology.

For example, a triangle is a simplicial complex because it is a collection of three points that are connected to each other in a specific way. Similarly, a tetrahedron is a simplicial complex because it is a collection of four points that are connected to each other in a specific way. These simplices can be thought of as the “building blocks” of a simplicial complex, and the complex itself is constructed by combining these building blocks in a specific way. For example, a 2-dimensional simplicial complex could be a collection of triangles that are connected to each other to form a surface, while a 3-dimensional simplicial complex could be a collection of tetrahedra that are connected to each other to form a solid object.

The SimplicialComplex class is a class for representing simplicial complexes, which are a type of topological space constructed by “gluing together” points, line segments, triangles, and higher-dimensional counterparts. The class provides methods for computing boundary operators, Hodge Laplacians, and higher-order adjacency operators on the simplicial complex. It also allows for compatibility with NetworkX and the GUDHI library.

Features:

  1. The SimplicialComplex class allows for the dynamic construction of simplicial complexes,

    enabling users to add or remove simplices from the complex after its initial creation.

  2. The class provides methods for computing boundary operators, Hodge Laplacians,

    and higher-order adjacency operators on the simplicial complex.

  3. The class is compatible with the gudhi library, allowing users to leverage the powerful

    algorithms and data structures provided by this package.

  4. The class supports the attachment of arbitrary attributes and data to simplices,

    enabling users to store and manipulate additional information about these objects.

  5. The class has robust error handling and input validation, ensuring reliable and easy use of the class.

Parameters:
simplicesiterable, optional

Iterable of maximal simplices that define the simplicial complex.

**kwargskeyword arguments, optional

Attributes to add to the complex as key=value pairs.

Attributes:
complexdict

A dictionary that can be used to store additional information about the complex.

Methods

add_elements_from_nx_graph(G)

Add elements from a networkx graph to this simplicial complex.

add_node(node, **kwargs)

Add node to simplicial complex.

add_simplex(simplex, **kwargs)

Add simplex to simplicial complex.

add_simplices_from(simplices)

Add simplices from iterable to simplicial complex.

adjacency_matrix(rank[, signed, weight, index])

Compute the adjacency matrix of the simplicial complex.

clone()

Return a copy of the simplicial complex.

coadjacency_matrix(rank[, signed, weight, index])

Compute the coadjacency matrix of the simplicial complex.

coincidence_matrix(rank[, signed, weight, index])

Compute coincidence matrix of the simplicial complex.

dirac_operator_matrix([signed, weight, index])

Compute dirac operator matrix matrix.

down_laplacian_matrix(rank[, signed, ...])

Compute the down Laplacian matrix of the simplicial complex.

from_gudhi(tree)

Import from gudhi.

from_nx(G)

Convert from netwrokx graph.

from_spharpy(mesh)

Import from sharpy.

from_trimesh(mesh)

Import from trimesh.

get_all_maximal_simplices()

Get all maximal simplices of this simplicial complex.

get_boundaries(simplices[, min_dim, max_dim])

Get boundaries of the given simplices.

get_cofaces(simplex, codimension)

Get cofaces of simplex.

get_edges_from_matrix(matrix)

Get edges from matrix.

get_maximal_simplices_of_simplex(simplex)

Get maximal simplices of simplex.

get_node_attributes(name)

Get node attributes from combinatorial complex.

get_simplex_attributes(name[, rank])

Get node attributes from simplical complex.

get_star(simplex)

Get star.

graph_skeleton()

Return the graph-skeleton of this simplicial complex.

hodge_laplacian_matrix(rank[, signed, ...])

Compute hodge-laplacian matrix for the simplicial complex.

incidence_matrix(rank[, signed, weight, index])

Compute incidence matrix of the simplicial complex.

is_connected()

Check if the simplicial complex is connected.

is_maximal(simplex)

Check if simplex is maximal.

is_triangular_mesh()

Check if the simplicial complex is a triangular mesh.

laplace_beltrami_operator([mode])

Compute a laplacian matrix for a triangular mesh.

load_mesh(file_path[, process])

Load a mesh.

normalized_laplacian_matrix(rank[, weight])

Return the normalized hodge Laplacian matrix of simplicial complex .

remove_maximal_simplex(simplex)

Remove maximal simplex from simplicial complex.

remove_nodes(node_set)

Remove the given nodes from the simplicial complex.

restrict_to_nodes(node_set)

Construct a new simplicial complex by restricting the simplices.

restrict_to_simplices(cell_set)

Construct a simplicial complex using a subset of the simplices.

set_simplex_attributes(values[, name])

Set simplex attributes.

simplicial_closure_of_hypergraph(H)

Compute the simplicial complex closure of a hypergraph.

skeleton(rank)

Compute skeleton.

to_cell_complex()

Convert a simplicial complex to a cell complex.

to_combinatorial_complex()

Convert a simplicial complex to a combinatorial complex.

to_hasse_graph()

Create the hasse graph corresponding to this simplicial complex.

to_hypergraph()

Convert a simplicial complex to a hypergraph.

to_spharapy([vertex_position_name])

Convert to sharapy.

to_trimesh([vertex_position_name])

Convert simplicial complex to trimesh object.

up_laplacian_matrix(rank[, signed, weight, ...])

Compute the up Laplacian matrix of the simplicial complex.

Notes

A simplicial complex is determined by its maximal simplices, simplices that are not contained in any other simplices. If a maximal simplex is inserted, all faces of this simplex will be inserted automatically.

Examples

Define a simplicial complex using a set of maximal simplices:

>>> SC = tnx.SimplicialComplex([[1, 2, 3], [2, 3, 5], [0, 1]])

TopoNetX is also compatible with NetworkX, allowing users to create a simplicial complex from a NetworkX graph. Existing node and edge attributes are copied to the simplicial complex:

>>> G = nx.Graph()
>>> G.add_edge(0, 1, weight=4)
>>> G.add_edges_from([(0, 3), (0, 4), (1, 4)])
>>> SC = tnx.SimplicialComplex(simplices=G)
>>> SC.add_simplex([1, 2, 3])
>>> SC.simplices
SimplexView([(0,), (1,), (3,), (4,), (2,), (0, 1), (0, 3), (0, 4), (1, 4), (1, 2), (1, 3), (2, 3), (1, 2, 3)])
>>> SC[(0, 1)]["weight"]
4
__init__(simplices=None, **kwargs) None[source]#

Initialize a new instance of the Complex class.

Parameters:
**kwargskeyword arguments, optional

Attributes to add to the complex as key=value pairs.

Methods

__init__([simplices])

Initialize a new instance of the Complex class.

add_elements_from_nx_graph(G)

Add elements from a networkx graph to this simplicial complex.

add_node(node, **kwargs)

Add node to simplicial complex.

add_simplex(simplex, **kwargs)

Add simplex to simplicial complex.

add_simplices_from(simplices)

Add simplices from iterable to simplicial complex.

adjacency_matrix(rank[, signed, weight, index])

Compute the adjacency matrix of the simplicial complex.

clone()

Return a copy of the simplicial complex.

coadjacency_matrix(rank[, signed, weight, index])

Compute the coadjacency matrix of the simplicial complex.

coincidence_matrix(rank[, signed, weight, index])

Compute coincidence matrix of the simplicial complex.

dirac_operator_matrix([signed, weight, index])

Compute dirac operator matrix matrix.

down_laplacian_matrix(rank[, signed, ...])

Compute the down Laplacian matrix of the simplicial complex.

from_gudhi(tree)

Import from gudhi.

from_nx(G)

Convert from netwrokx graph.

from_spharpy(mesh)

Import from sharpy.

from_trimesh(mesh)

Import from trimesh.

get_all_maximal_simplices()

Get all maximal simplices of this simplicial complex.

get_boundaries(simplices[, min_dim, max_dim])

Get boundaries of the given simplices.

get_cofaces(simplex, codimension)

Get cofaces of simplex.

get_edges_from_matrix(matrix)

Get edges from matrix.

get_maximal_simplices_of_simplex(simplex)

Get maximal simplices of simplex.

get_node_attributes(name)

Get node attributes from combinatorial complex.

get_simplex_attributes(name[, rank])

Get node attributes from simplical complex.

get_star(simplex)

Get star.

graph_skeleton()

Return the graph-skeleton of this simplicial complex.

hodge_laplacian_matrix(rank[, signed, ...])

Compute hodge-laplacian matrix for the simplicial complex.

incidence_matrix(rank[, signed, weight, index])

Compute incidence matrix of the simplicial complex.

is_connected()

Check if the simplicial complex is connected.

is_maximal(simplex)

Check if simplex is maximal.

is_triangular_mesh()

Check if the simplicial complex is a triangular mesh.

laplace_beltrami_operator([mode])

Compute a laplacian matrix for a triangular mesh.

load_mesh(file_path[, process])

Load a mesh.

normalized_laplacian_matrix(rank[, weight])

Return the normalized hodge Laplacian matrix of simplicial complex .

remove_maximal_simplex(simplex)

Remove maximal simplex from simplicial complex.

remove_nodes(node_set)

Remove the given nodes from the simplicial complex.

restrict_to_nodes(node_set)

Construct a new simplicial complex by restricting the simplices.

restrict_to_simplices(cell_set)

Construct a simplicial complex using a subset of the simplices.

set_simplex_attributes(values[, name])

Set simplex attributes.

simplicial_closure_of_hypergraph(H)

Compute the simplicial complex closure of a hypergraph.

skeleton(rank)

Compute skeleton.

to_cell_complex()

Convert a simplicial complex to a cell complex.

to_combinatorial_complex()

Convert a simplicial complex to a combinatorial complex.

to_hasse_graph()

Create the hasse graph corresponding to this simplicial complex.

to_hypergraph()

Convert a simplicial complex to a hypergraph.

to_spharapy([vertex_position_name])

Convert to sharapy.

to_trimesh([vertex_position_name])

Convert simplicial complex to trimesh object.

up_laplacian_matrix(rank[, signed, weight, ...])

Compute the up Laplacian matrix of the simplicial complex.

Attributes

dim

Dimension of the simplicial complex.

maxdim

Maximum dimension of the simplicial complex.

nodes

Return the list of nodes in the simplicial complex.

shape

Shape of simplicial complex.

simplices

Set of all simplices in the simplicial complex.

complex

add_elements_from_nx_graph(G: Graph) None[source]#

Add elements from a networkx graph to this simplicial complex.

Parameters:
Gnetworkx.Graph

A networkx graph instance.

add_node(node: Hashable, **kwargs) None[source]#

Add node to simplicial complex.

Parameters:
nodeHashable or Simplex

The node to be added to the simplicial complex.

**kwargskeyword arguments, optional

Additional attributes to be associated with the node.

add_simplex(simplex: Collection, **kwargs) None[source]#

Add simplex to simplicial complex.

In case sub-simplices are missing, they are added without attributes to the simplicial complex to fulfill the simplex requirements.

If the given simplex is already part of the simplicial complex, its attributes will be updated by the provided ones.

Parameters:
simplexCollection

The simplex to be added to the simplicial complex.

If a Simplex object is given, its attributes will be copied to the simplicial complex. kwargs take precedence over the attributes of the Simplex object.

**kwargskeyword arguments, optional

Additional attributes to be associated with the simplex.

add_simplices_from(simplices) None[source]#

Add simplices from iterable to simplicial complex.

Parameters:
simplicesiterable

Iterable of simplices to be added to the simplicial complex.

adjacency_matrix(rank, signed: bool = False, weight: str | None = None, index: bool = False)[source]#

Compute the adjacency matrix of the simplicial complex.

Parameters:
rankint

Rank of the adjacency matrix.

signedbool

Whether to return the signed or unsigned adjacency matrix.

weightstr, optional

The name of the simplex attribute to use as weights for the hodge laplacian matrix. If None, the matrix is binary.

indexbool, default=False

Whether to return the indices of the rows and columns of the adjacency matrix.

Returns:
indicesdict

Dictionary mapping each row and column of the coadjacency matrix to a simplex. Only returned if index is True.

adjacency_matrixscipy.sparse.csr.csr_matrix

The adjacency matrix of rank rank of this simplicial complex.

Examples

>>> SC = tnx.SimplicialComplex()
>>> SC.add_simplex([1, 2, 3, 4])
>>> SC.add_simplex([1, 2, 4])
>>> SC.add_simplex([3, 4, 8])
>>> adj1 = SC.adjacency_matrix(1)
clone() Self[source]#

Return a copy of the simplicial complex.

The clone method by default returns an independent shallow copy of the simplicial complex. Use Python’s copy.deepcopy for new containers.

Returns:
SimplicialComplex

A shallow copy of this simplicial complex.

coadjacency_matrix(rank: int, signed: bool = False, weight=None, index: bool = False)[source]#

Compute the coadjacency matrix of the simplicial complex.

Parameters:
rankint

Rank of the coadjacency matrix.

signedbool

Whether to return the signed or unsigned coadjacency matrix.

weightstr, optional

The name of the simplex attribute to use as weights for the hodge laplacian matrix. If None, the matrix is binary.

indexbool, default=False

Whether to return the indices of the rows and columns of the coadjacency matrix.

Returns:
indiceslist

List identifying the rows and columns of the coadjacency matrix with the simplices of the simplicial complex. Only returned if index is True.

coadjacency_matrixscipy.sparse.csr.csr_matrix

The coadjacency matrix of this simplicial complex.

coincidence_matrix(rank, signed: bool = True, weight=None, index: bool = False)[source]#

Compute coincidence matrix of the simplicial complex.

This is also called the coboundary matrix.

Parameters:
rankint

For which rank of simplices to compute the coincidence matrix.

signedbool, default=True

Whether to return the signed or unsigned coincidence matrix.

weightstr, optional

The name of the simplex attribute to use as weights for the coincidence matrix. If None, the matrix is binary.

indexbool, default=False

Whether to return the indices of the rows and columns of the coincidence matrix.

Returns:
row_indices, col_indicesdict

Dictionary assigning each row and column of the coincidence matrix to a simplex. Only returned if index is True.

coincidence_matrixscipy.sparse.csr.csr_matrix

The coincidence matrix.

property dim: int#

Dimension of the simplicial complex.

This is the highest dimension of any simplex in the complex.

Returns:
int

The dimension of the simplicial complex.

dirac_operator_matrix(signed: bool = True, weight: str | None = None, index: bool = False)[source]#

Compute dirac operator matrix matrix.

Parameters:
signedbool, default=False

Whether the returned dirac matrix should be signed (i.e., respect orientations) or unsigned.

weightstr, optional

The name of the simplex attribute to use as weights for the dirac operator matrix. If None, the matrix is binary.

indexbool, default=False

Whether to return the indices of the rows and columns of the dirac operator matrix.

Returns:
row_indices, col_indicesdict

List identifying rows and columns of the dirac operator matrix. Only returned if index is True.

scipy.sparse.csr.csc_matrix

The dirac operator matrix of this simplicial complex.

Examples

>>> SC = tnx.SimplicialComplex()
>>> SC.add_simplex([1, 2, 3, 4])
>>> SC.add_simplex([1, 2, 4])
>>> SC.add_simplex([3, 4, 8])
>>> SC.dirac_operator_matrix()
down_laplacian_matrix(rank, signed: bool = True, weight=None, index: bool = False)[source]#

Compute the down Laplacian matrix of the simplicial complex.

Parameters:
rankint

Rank of the down Laplacian matrix.

signedbool

Whether to return the signed or unsigned down laplacian.

weightstr, optional

The name of the simplex attribute to use as weights for the hodge laplacian matrix. If None, the matrix is binary.

indexbool, default=False

Whether to return the indices of the rows and columns of the laplacian.

Returns:
indicesdict

Dictionary assigning each row and column of the laplacian matrix to a simplex. Only returned if index is True.

down_laplacianscipy.sparse.csr.csr_matrix

The down laplacian matrix of this simplicial complex.

classmethod from_gudhi(tree: None) Self[source]#

Import from gudhi.

Parameters:
treegudhi.SimplexTree

The input gudhi simplex tree.

Returns:
SimplicialComplex

The resulting SimplicialComplex.

Examples

>>> from gudhi import SimplexTree
>>> tree = SimplexTree()
>>> _ = tree.insert([1, 2, 3, 5])
>>> SC = tnx.SimplicialComplex.from_gudhi(tree)
classmethod from_nx(G: Graph) Self[source]#

Convert from netwrokx graph.

Parameters:
Gnx.Graph

The graph to convert to a simplicial complex.

Returns:
SimplicialComplex

The simplicial complex.

Examples

>>> G = nx.Graph()  # or DiGraph, MultiGraph, MultiDiGraph, etc
>>> G.add_edge(1, 2, weight=2)
>>> G.add_edge(3, 4, weight=4)
>>> SC = tnx.SimplicialComplex.from_nx(G)
>>> SC[(1, 2)]["weight"]
2
classmethod from_spharpy(mesh) Self[source]#

Import from sharpy.

Parameters:
meshspharapy.trimesh.TriMesh

The input spharapy object.

Returns:
SimplicialComplex

The resulting SimplicialComplex.

Examples

>>> import spharapy.trimesh as tm
>>> import spharapy.spharabasis as sb
>>> import spharapy.datasets as sd
>>> mesh = tm.TriMesh(
...     [[0, 1, 2]], [[1.0, 0.0, 0.0], [0.0, 2.0, 0.0], [0.0, 0.0, 3.0]]
... )
>>> SC = tnx.SimplicialComplex.from_spharpy(mesh)
classmethod from_trimesh(mesh) Self[source]#

Import from trimesh.

Parameters:
meshtrimesh.Trimesh

The input trimesh object.

Returns:
SimplicialComplex

The resulting SimplicialComplex.

Examples

>>> import trimesh
>>> mesh = trimesh.Trimesh(
...     vertices=[[0, 0, 0], [0, 0, 1], [0, 1, 0]],
...     faces=[[0, 1, 2]],
...     process=False,
... )
>>> SC = tnx.SimplicialComplex.from_trimesh(mesh)
>>> print(SC.nodes)
>>> print(SC.simplices)
>>> SC[(0)]["position"]
get_all_maximal_simplices()[source]#

Get all maximal simplices of this simplicial complex.

A simplex is maximal if it is not a face of any other simplex in the complex.

Returns:
list of tuples

List of maximal simplices in this simplicial complex.

Examples

>>> SC = tnx.SimplicialComplex([(1, 2), (1, 2, 3), (1, 2, 4), (2, 5)])
>>> SC.get_all_maximal_simplices()
[(2, 5), (1, 2, 3), (1, 2, 4)]
static get_boundaries(simplices: Iterable, min_dim=None, max_dim=None) set[frozenset][source]#

Get boundaries of the given simplices.

Parameters:
simplicesIterable

Iterable of simplices for which to compute the boundaries.

min_dimint, optional

Constrain the max dimension of faces.

max_dimint, optional

Constrain the max dimension of faces.

Returns:
set of frozensets

Set of simplices that are boundaries of the given simplices. If min_dim or max_dim are given, only simplices with dimension in the given range are returned.

get_cofaces(simplex: Iterable[Hashable], codimension: int) list[frozenset][source]#

Get cofaces of simplex.

Parameters:
simplexlist, tuple or simplex

The simplex to get the cofaces of.

codimensionint

The codimension. If codimension = 0, all cofaces are returned.

Returns:
list of tuples

The cofaces of the given simplex.

static get_edges_from_matrix(matrix)[source]#

Get edges from matrix.

Parameters:
matrixnumpy or scipy array

The matrix to get the edges from.

Returns:
list of int

List of indices where the operator is not zero.

Notes

Most operaters (e.g. adjacencies/(co)boundary maps) that describe connectivity of the simplicial complex can be described as a G whose nodes are the simplices used to construct the operator and whose edges correspond to the entries in the matrix where the operator is not zero.

This property implies that many computations on simplicial complexes can be reduced to G computations.

get_maximal_simplices_of_simplex(simplex: Iterable[Hashable]) set[frozenset][source]#

Get maximal simplices of simplex.

Parameters:
simplexIterable

The simplex for which to compute the maximal simplices.

Returns:
set of frozensets

Set of maximal simplices of the given simplex.

get_node_attributes(name: str) dict[Hashable, Any][source]#

Get node attributes from combinatorial complex.

Parameters:
namestr

Attribute name.

Returns:
dict

Dictionary mapping each node to the value of the given attribute name.

Examples

>>> SC = tnx.SimplicialComplex()
>>> SC.add_simplex([1, 2, 3, 4])
>>> SC.add_simplex([1, 2, 4])
>>> SC.add_simplex([3, 4, 8])
>>> SC.set_simplex_attributes({1: "red", 2: "blue", 3: "black"}, name="color")
>>> SC.get_node_attributes("color")
{1: 'red', 2: 'blue', 3: 'black'}
get_simplex_attributes(name: str, rank: int | None = None) dict[tuple[Hashable, ...], Any][source]#

Get node attributes from simplical complex.

Parameters:
namestr

Attribute name.

rankint, optional

Restrict the returned attribute values to simplices of a specific rank.

Returns:
dict[tuple[Hashable, …], Any]

Dictionary mapping each simplex to the value of the given attribute name.

Examples

>>> SC = tnx.SimplicialComplex()
>>> SC.add_simplex([1, 2, 3, 4])
>>> SC.add_simplex([1, 2, 4])
>>> SC.add_simplex([3, 4, 8])
>>> d = {(1, 2): "red", (2, 3): "blue", (3, 4): "black"}
>>> SC.set_simplex_attributes(d, name="color")
>>> SC.get_simplex_attributes("color")
{frozenset({1, 2}): 'red', frozenset({2, 3}): 'blue', frozenset({3, 4}): 'black'}
get_star(simplex) list[frozenset][source]#

Get star.

Parameters:
simplexlist, tuple or simplex

The simplex represented by a list of its nodes.

Returns:
list[frozenset]

The star of the given simplex.

Notes

This function is equivalent to get_cofaces(simplex, 0).

graph_skeleton() Graph[source]#

Return the graph-skeleton of this simplicial complex.

The graph-skeleton consists of the 0 and 1-simplices (i.e., nodes and edges) of the simplicial complex.

Returns:
nx.Graph

The graph-skeleton of this simplicial complex.

hodge_laplacian_matrix(rank: int, signed: bool = True, weight=None, index: bool = False)[source]#

Compute hodge-laplacian matrix for the simplicial complex.

Parameters:
rankint

Dimension of the Laplacian matrix.

signedbool, default=True

Whether to return the signed or unsigned hodge laplacian.

weightstr, optional

The name of the simplex attribute to use as weights for the hodge laplacian matrix. If None, the matrix is binary.

indexbool, default=False

Indicates whether to return the indices that define the hodge laplacian matrix.

Returns:
indexlist

List assigning each row and column of the laplacian matrix to a simplex. Only available when index is True.

laplacianscipy.sparse.csr.csr_matrix

The hodge laplacian matrix of rank rank.

Examples

>>> SC = tnx.SimplicialComplex()
>>> SC.add_simplex([1, 2, 3, 4])
>>> SC.add_simplex([1, 2, 4])
>>> SC.add_simplex([3, 4, 8])
>>> L1 = SC.hodge_laplacian_matrix(1)
incidence_matrix(rank, signed: bool = True, weight: str | None = None, index: bool = False) csr_matrix | tuple[dict, dict, csr_matrix][source]#

Compute incidence matrix of the simplicial complex.

Getting the matrix that correpodnds to the boundary matrix of the input SC.

Parameters:
rankint

For which rank of simplices to compute the incidence matrix.

signedbool, default=True

Whether to return the signed or unsigned incidence matrix.

weightstr, optional

The name of the simplex attribute to use as weights for the incidence matrix. If None, the matrix is binary.

indexbool, default=False

Whether to return the indices of the rows and columns of the incidence matrix.

Returns:
row_indices, col_indicesdict

Dictionary assigning each row and column of the incidence matrix to a simplex. Only returned if index is True.

incidence_matrixscipy.sparse.csr.csr_matrix

The incidence matrix.

Raises:
ValueError

If rank is negative or larger than the dimension of the simplicial complex.

Examples

>>> SC = tnx.SimplicialComplex()
>>> SC.add_simplex([1, 2, 3, 4])
>>> SC.add_simplex([1, 2, 4])
>>> SC.add_simplex([3, 4, 8])
>>> B1 = SC.incidence_matrix(1)
>>> B2 = SC.incidence_matrix(2)
is_connected() bool[source]#

Check if the simplicial complex is connected.

Returns:
bool

True if the simplicial complex is connected, False otherwise.

Notes

A simplicial complex is connected iff its 1-skeleton is connected.

is_maximal(simplex: Iterable) bool[source]#

Check if simplex is maximal.

Parameters:
simplexIterable

The Simplex to check.

Returns:
bool

True if the given simplex is maximal in this simplicial complex, False otherwise.

Raises:
ValueError

If simplex is not in the simplicial complex.

Examples

>>> SC = tnx.SimplicialComplex([[1, 2, 3]])
>>> SC.is_maximal([1, 2, 3])
True
>>> SC.is_maximal([1, 2])
False
is_triangular_mesh() bool[source]#

Check if the simplicial complex is a triangular mesh.

Returns:
bool

True if the simplicial complex can be converted to a triangular mesh, False otherwise.

laplace_beltrami_operator(mode: str = 'inv_euclidean')[source]#

Compute a laplacian matrix for a triangular mesh.

The method creates a laplacian matrix for a triangular mesh using different weighting function.

Parameters:
mode{‘unit’, ‘inv_euclidean’, ‘half_cotangent’}, optional

The methods for determining the edge weights. Using the option ‘unit’ all edges of the mesh are weighted by unit weighting function, the result is an adjacency matrix. The option ‘inv_euclidean’ results in edge weights corresponding to the inverse Euclidean distance of the edge lengths. The option ‘half_cotangent’ uses the half of the cotangent of the two angles opposed to an edge as weighting function. the default weighting function is ‘inv_euclidean’.

Returns:
numpy.ndarray, shape (n_vertices, n_vertices)

Matrix, which contains the discrete laplace operator for data defined at the vertices of a triangular mesh. The number of vertices of the triangular mesh is n_points.

classmethod load_mesh(file_path, process: bool = False) Self[source]#

Load a mesh.

Parameters:
file_pathstr or pathlib.Path

The source of the data to be loaded.

processbool

Whether trimesh should try to process the mesh before loading it.

Returns:
SimplicialComplex

The output simplicial complex stores the same structure stored in the mesh input file.

Notes

mesh files supported : obj, off, glb

Examples

>>> SC = tnx.SimplicialComplex.load_mesh("C:/temp/stanford-bunny.obj")
>>> SC.nodes
property maxdim: int#

Maximum dimension of the simplicial complex.

This is the highest dimension of any simplex in the complex.

Returns:
int

The maximum dimension of the simplicial complex.

property nodes#

Return the list of nodes in the simplicial complex.

Returns:
NodeView

A NodeView object representing the nodes of the simplicial complex.

normalized_laplacian_matrix(rank: int, weight: str | None = None)[source]#

Return the normalized hodge Laplacian matrix of simplicial complex .

The normalized hodge Laplacian is the matrix

\[N_d = D^{-1/2} L_d D^{-1/2}\]

where L is the simplicial complex Laplacian and D is the diagonal matrix of simplices of rank d.

Parameters:
rankint

Rank of the hodge laplacian matrix.

weightstr, optional

The name of the simplex attribute to use as weights for the hodge laplacian matrix. If None, the matrix is binary.

Returns:
Scipy sparse matrix

The normalized hodge Laplacian matrix.

Examples

>>> SC = tnx.SimplicialComplex([[1, 2, 3], [2, 3, 5], [0, 1]])
>>> SC.normalized_laplacian_matrix(1)
remove_maximal_simplex(simplex: Collection) None[source]#

Remove maximal simplex from simplicial complex.

Parameters:
simplexCollection

The simplex to be removed from the simplicial complex.

Raises:
KeyError

If simplex is not in simplicial complex.

ValueError

If simplex is not maximal.

Examples

>>> SC = tnx.SimplicialComplex()
>>> SC.add_simplex((1, 2, 3, 4), weight=1)
>>> SC.add_simplex((1, 2, 3, 4, 5))
>>> SC.remove_maximal_simplex((1, 2, 3, 4, 5))
remove_nodes(node_set: Iterable[Hashable]) None[source]#

Remove the given nodes from the simplicial complex.

Any simplices that become invalid due to the removal of nodes are also removed.

Parameters:
node_setIterable

The nodes to be removed from the simplicial complex.

Examples

>>> SC = tnx.SimplicialComplex([(1, 2), (2, 3), (3, 4)])
>>> SC.remove_nodes([1])
>>> SC.simplices
SimplexView([(2,), (3,), (4,), (2, 3), (3, 4)])
restrict_to_nodes(node_set)[source]#

Construct a new simplicial complex by restricting the simplices.

The simplices are restricted to the nodes referenced by node_set.

Parameters:
node_setiterable of hashables

A subset of nodes of the simplicial complex to restrict to.

Returns:
SimplicialComplex

A new simplicial complex restricted to the nodes in node_set.

Examples

>>> c1 = tnx.Simplex((1, 2, 3))
>>> c2 = tnx.Simplex((1, 2, 4))
>>> c3 = tnx.Simplex((1, 2, 5))
>>> SC = tnx.SimplicialComplex([c1, c2, c3])
>>> new_complex = SC.restrict_to_nodes([1, 2, 3, 4])
>>> new_complex.simplices
SimplexView([(1,), (2,), (3,), (4,), (1, 2), (1, 3), (1, 4), (2, 3), (2, 4), (1, 2, 3), (1, 2, 4)])
restrict_to_simplices(cell_set) Self[source]#

Construct a simplicial complex using a subset of the simplices.

Parameters:
cell_setiterable of hashables or simplices

A subset of elements of the simplicial complex that should be in the new simplicial complex.

Returns:
SimplicialComplex

New simplicial complex restricted to the simplices in cell_set.

Examples

>>> c1 = tnx.Simplex((1, 2, 3))
>>> c2 = tnx.Simplex((1, 2, 4))
>>> c3 = tnx.Simplex((1, 2, 5))
>>> SC = tnx.SimplicialComplex([c1, c2, c3])
>>> SC1 = SC.restrict_to_simplices([c1, (2, 4)])
>>> SC1.simplices
SimplexView([(1,), (2,), (3,), (4,), (1, 2), (1, 3), (2, 3), (2, 4), (1, 2, 3)])
set_simplex_attributes(values, name: str | None = None) None[source]#

Set simplex attributes.

Parameters:
valuesdict

Either provide a mapping from simplices to values or a dict of dicts. In the former case, the attribute name for each simplex is set to the corresponding value. In the latter case, the entire dictionary is used to update the attributes of the simplices.

namestr, optional

The name of the attribute to set.

Notes

If the dict contains simplices that are not in self.simplices, they are silently ignored.

Examples

After computing some property of the simplex of a simplicial complex, you may want to assign a simplex attribute to store the value of that property for each simplex:

>>> SC = tnx.SimplicialComplex()
>>> SC.add_simplex([1, 2, 3, 4])
>>> SC.add_simplex([1, 2, 4])
>>> SC.add_simplex([3, 4, 8])
>>> d = {(1, 2, 3): "red", (1, 2, 4): "blue"}
>>> SC.set_simplex_attributes(d, name="color")
>>> SC[(1, 2, 3)]["color"]
'red'

If you provide a dictionary of dictionaries as the second argument, the entire dictionary will be used to update simplex attributes:

>>> SC = tnx.SimplicialComplex()
>>> SC.add_simplex([1, 3, 4])
>>> SC.add_simplex([1, 2, 3])
>>> SC.add_simplex([1, 2, 4])
>>> d = {
...     (1, 3, 4): {"color": "red", "attr2": 1},
...     (1, 2, 4): {"color": "blue", "attr2": 3},
... }
>>> SC.set_simplex_attributes(d)
>>> SC[(1, 3, 4)]["color"]
'red'
property shape: tuple[int, ...]#

Shape of simplicial complex.

(number of simplices[i], for i in range(0,dim(Sc)) )

Returns:
tuple of ints

This gives the number of cells in each rank.

property simplices: SimplexView#

Set of all simplices in the simplicial complex.

Returns:
SimplexView

A SimplexView object representing the set of all simplices in the simplicial complex.

classmethod simplicial_closure_of_hypergraph(H) Self[source]#

Compute the simplicial complex closure of a hypergraph.

Parameters:
Hhyernetx hypergraph

The hypergraph to compute the simplicial complex closure of.

Returns:
SimplicialComplex

Simplicial complex closure of the hypergraph.

Examples

>>> import hypernetx as hnx
>>> hg = hnx.Hypergraph([[1, 2, 3, 4], [1, 2, 3]], static=True)
>>> sc = tnx.SimplicialComplex.simplicial_closure_of_hypergraph(hg)
>>> sc.simplices
SimplexView([(1,), (2,), (3,), (4,), (1, 2), (1, 3), (1, 4), (2, 3), (2, 4), (3, 4), (1, 2, 3), (1, 2, 4), (1, 3, 4), (2, 3, 4), (1, 2, 3, 4)])
skeleton(rank: int) list[tuple[Hashable, ...]][source]#

Compute skeleton.

Parameters:
rankint

The rank of the skeleton to compute.

Returns:
list of tuples

Simplices of rank rank in the simplicial complex.

to_cell_complex()[source]#

Convert a simplicial complex to a cell complex.

Returns:
toponetx.classes.CellComplex

The cell complex corresponding to this simplicial complex.

Examples

>>> c1 = tnx.Simplex((1, 2, 3))
>>> c2 = tnx.Simplex((1, 2, 4))
>>> c3 = tnx.Simplex((2, 5))
>>> SC = tnx.SimplicialComplex([c1, c2, c3])
>>> SC.to_cell_complex()
to_combinatorial_complex()[source]#

Convert a simplicial complex to a combinatorial complex.

Returns:
toponetx.classes.CombinatorialComplex

The combinatorial complex equivalent to this simplicial complex.

Examples

>>> c1 = tnx.Simplex((1, 2, 3))
>>> c2 = tnx.Simplex((1, 2, 3))
>>> c3 = tnx.Simplex((1, 2, 4))
>>> SC = tnx.SimplicialComplex([c1, c2, c3])
>>> CCC = SC.to_combinatorial_complex()
to_hasse_graph() DiGraph[source]#

Create the hasse graph corresponding to this simplicial complex.

Returns:
nx.DiGraph

A NetworkX Digraph representing the Hasse graph of the Simplicial Complex.

Examples

>>> SC = tnx.SimplicialComplex()
>>> SC.add_simplex([0, 1, 2])
>>> G = SC.to_hasse_graph()
to_hypergraph() None[source]#

Convert a simplicial complex to a hypergraph.

Returns:
Hypergraph

The hypergraph corresponding to this simplicial complex.

Examples

>>> c1 = tnx.Simplex((1, 2, 3))
>>> c2 = tnx.Simplex((1, 2, 4))
>>> c3 = tnx.Simplex((2, 5))
>>> SC = tnx.SimplicialComplex([c1, c2, c3])
>>> SC.to_hypergraph()
Hypergraph({'e0': [1, 2], 'e1': [1, 3], 'e2': [1, 4], 'e3': [2, 3], 'e4': [2, 4], 'e5': [2, 5], 'e6': [1, 2, 3], 'e7': [1, 2, 4]},name=)
to_spharapy(vertex_position_name: str = 'position')[source]#

Convert to sharapy.

Parameters:
vertex_position_namestr, default=”position”

The simplex attribute name that contains the vertex positions.

Returns:
spharapy.trimesh.TriMesh

The spharapy object corresponding to this simplicial complex.

Raises:
RuntimeError

If package spharapy is not installed.

Examples

>>> import spharapy.trimesh as tm
>>> import spharapy.spharabasis as sb
>>> import spharapy.datasets as sd
>>> mesh = tm.TriMesh([[0, 1, 2]], [[0, 0, 0], [0, 0, 1], [0, 1, 0]])
>>> SC = tnx.SimplicialComplex.from_spharpy(mesh)
>>> mesh2 = SC.to_spharapy()
>>> mesh2.vertlist == mesh.vertlist
>>> mesh2.trilist == mesh.trilist
to_trimesh(vertex_position_name: str = 'position')[source]#

Convert simplicial complex to trimesh object.

Parameters:
vertex_position_namestr, default=”position”

The simplex attribute name that contains the vertex positions.

Returns:
trimesh.Trimesh

The trimesh object corresponding to this simplicial complex.

up_laplacian_matrix(rank: int, signed: bool = True, weight: str | None = None, index: bool = False)[source]#

Compute the up Laplacian matrix of the simplicial complex.

Parameters:
rankint

Rank of the up Laplacian matrix.

signedbool

Whether to return the signed or unsigned up laplacian.

weightstr, optional

The name of the simplex attribute to use as weights for the hodge laplacian matrix. If None, the matrix is binary.

indexbool, default=False

Whether to return the indices of the rows and columns of the laplacian.

Returns:
row_indices, col_indiceslist

List identifying the rows and columns of the laplacian matrix with the simplices of this complex. Only returned if index is True.

up_laplacianscipy.sparse.csr.csr_matrix

The upper laplacian matrix of this simplicial complex.

Examples

>>> SC = tnx.SimplicialComplex([[1, 2, 3], [2, 3, 5], [0, 1]])
>>> SC.up_laplacian_matrix(1)