simplicial: simplicial complexes, homology, topology and persistent homology

Date

2021-01-13 (updated 2023-12-13)

What for

Simplified package for understanding how to compute homology, persistent homology and how to visualize simplicial complexes and alike.

Warning

Work in progress (very much)

INSTALL:

First install (download, or exec)

$ wget https://www.dlfer.xyz/var/simplicial.py

from a shell, or

>>> !wget https://www.dlfer.xyz/var/simplicial.py

from a jupyter notebook cell). Check it first if it is safe.

This is a simple module, maybe to explain simplicial complexes and simplicial homology.

USAGE:

From jupyter/python console, then import it

>>> from simplicial import *

EXAMPLES:

>>> X=simplicial_sphere(['a','b'])
>>> print(X)
Simplicial Complex of dim 0:
 vertices: ['a', 'b']
 simplices:
  0: [['a'], ['b']]
>>> betti_numbers(X)
[2]

Join of simplices

>>> s0=simplicial_sphere(["A","B"])
>>> t0=simplicial_sphere([0,1,2])
>>> print(s0+t0)
They must be surfaces!
None
>>> print(s0*t0)
Simplicial Complex of dim 2:
 vertices: [0, 1, 2, 3, 4]
 simplices:
  0: [[0], [1], [2], [3], [4]]
  1: [[2, 3], [2, 4], [3, 4], [0, 2], [0, 3], [0, 4], [1, 2], [1, 3], [1, 4]]
  2: [[0, 2, 3], [0, 2, 4], [0, 3, 4], [1, 2, 3], [1, 2, 4], [1, 3, 4]]

Cartesian product

>>> print(s0 % t0)
Simplicial Complex of dim 1:
 vertices: [('A', 0), ('A', 1), ('A', 2), ('B', 0), ('B', 1), ('B', 2)]
 simplices:
  0: [[('A', 0)], [('A', 1)], [('A', 2)], [('B', 0)], [('B', 1)], [('B', 2)]]
  1: 6 simplices
>>> K= simplicial_sphere(1)
>>> print(K)
Simplicial Complex of dim 1:
 vertices: [0, 1, 2]
 simplices:
  0: [[0], [1], [2]]
  1: [[0, 1], [0, 2], [1, 2]]
>>> K=SimplicialComplex( ["A","B"], maximal_simplices={1: [["A","B"]] } , must_reindex=True)
>>> print(K)
Simplicial Complex of dim 1:
 vertices: ['A', 'B']
 simplices:
  0: [['A'], ['B']]
  1: [['A', 'B']]
>>> EPchar(K)
1
>>> s1=simplicial_sphere(1)
>>> T=s1 % s1
>>> K= ( T + T + T ) * s1
>>> print(K)
Simplicial Complex of dim 4:
 vertices: [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23]
 simplices:
  0: 24 simplices
  1: 141 simplices
  2: 338 simplices
  3: 375 simplices
  4: 150 simplices
>>> dims_c(K)
[24, 141, 338, 375, 150]
>>> betti_numbers(K)
****[1, 0, 0, 6, 1]

Connected sum.

>>> s1=simplicial_sphere(1)
>>> T=s1 % s1
>>> betti_numbers(T)
**[1, 2, 1]
>>> K=T+T+T
>>> betti_numbers(K)
**[1, 6, 1]
>>> K= (  ( s1 % s1 ) + simplicial_sphere(2) ) * s1
>>> betti_numbers(K)
****[1, 0, 0, 2, 1]
>>> K= disjoint_union(T,T)
>>> print(K)
Simplicial Complex of dim 2:
 vertices: [(0, 0), (0, 1), (0, 2), (1, 0), (1, 1), (1, 2), (2, 0), (2, 1), (2, 2), '(0, 0)+', '(0, 1)+', '(0, 2)+', '(1, 0)+', '(1, 1)+', '(1, 2)+', '(2, 0)+', '(2, 1)+', '(2, 2)+']
 simplices:
  0: 18 simplices
  1: 54 simplices
  2: 36 simplices
>>> betti_numbers(K)
**[2, 4, 2]
>>> EPchar(K)
0
>>> dims_c(K)
[18, 54, 36]
>>> K=simplicial_sphere([ "A","B","C" ] )
>>> L=simplicial_sphere([ "a","b","c" ])
>>> f=SimplicialMap({"A":"a","B":"c","C":"b"},K,L)
>>> g=SimplicialMap({"a":"A","b":"B","c":"B"},L,K)
>>> print(f)
SimplicialMap{'A': 'a', 'B': 'c', 'C': 'b'}
>>> print("is_simplicial f:", f.is_simplicial())
is_simplicial f: True
>>> print("is_simplicial_g:", g.is_simplicial())
is_simplicial_g: True
>>> print(f(["A","B"]) )
['a', 'c']
>>> print(f["C"] )
b
>>> print(f*g)
SimplicialMap{'a': 'a', 'b': 'c', 'c': 'c'}
>>> print(f*g*f)
SimplicialMap{'A': 'a', 'B': 'c', 'C': 'c'}
>>> print(homology_of_map(f))
[Matrix([[1]]), Matrix([[-1]])]
>>> print(homology_of_map(g))
[Matrix([[1]]), Matrix([[0]])]
>>> f=SimplicialMap({"A":"b","B":"b","C":"b"},K,L)
>>> print(f)
SimplicialMap{'A': 'b', 'B': 'b', 'C': 'b'}
>>> print(homology_of_map(f))
[Matrix([[1]]), Matrix([[0]])]
>>> K0=SimplicialComplex([0,1,2,3])
>>> K1=SimplicialComplex([0,1,2,3],maximal_simplices={1:[[0,1]],0:[[2],[3]]} )
>>> K2=SimplicialComplex([0,1,2,3],maximal_simplices={1:[[0,1],[0,2]],0: [ [3] ] } )
>>> K3=SimplicialComplex([0,1,2,3], maximal_simplices={1:[[0,1],[0,2],[1,2]], 0:[[3]] } )
>>> K4=SimplicialComplex([0,1,2,3],maximal_simplices={2:[[0,1,2]],0:[[3]] } )
>>> K5=K4
>>> filtration=[K0,K1,K2,K3,K4,K5]
>>> id_map={0:0,1:1,2:2,3:3}
>>> list_of_maps=[ SimplicialMap(id_map,filtration[j],filtration[j+1]) for j in range(len(filtration)-1)]
>>> for f in list_of_maps:
...     print("H:", homology_of_map(f))
H: [Matrix([
[1, 1, 0, 0],
[0, 0, 1, 0],
[0, 0, 0, 1]])]
H: [Matrix([
[1, 1, 0],
[0, 0, 1]]), Matrix(0, 0, [])]
H: [Matrix([
[1, 0],
[0, 1]]), Matrix(1, 0, [])]
H: [Matrix([
[1, 0],
[0, 1]]), Matrix(0, 1, [])]
H: [Matrix([
[1, 0],
[0, 1]]), Matrix(0, 0, []), Matrix(0, 0, [])]
>>> PH=persistent_homology(list_of_maps)
>>> print("PH:", PH)
PH: {0: [(0, 5), (0, 1), (0, 2), (0, 5)], 1: [(3, 4)], 2: []}

CONTENT:

class simplicial.SimplicialComplex(vertices, simplices=None, maximal_simplices=None, must_reindex=False, must_sort=True)

If just the set of verties: 0-dim simplicial complex. Otherwise, simplices, which is a dictionary of integer indices i.e. simplices. Otherwise, just maximal_simplices, and simplices are populated. if must_reindex=True, simplices in terms of symbols will be converted to simplices in terms of integer indices

operations:

  • connected_sum +

  • join *

  • cartesian_product %

  • disjoint_union

>>> K=SimplicialComplex(list(range(7)),maximal_simplices={1:[ sorted([x % 7, (x+1)% 7]) for x in range(7)] } )
>>> print(K)
Simplicial Complex of dim 1:
 vertices: [0, 1, 2, 3, 4, 5, 6]
 simplices:
  0: [[0], [1], [2], [3], [4], [5], [6]]
  1: [[0, 1], [1, 2], [2, 3], [3, 4], [4, 5], [5, 6], [0, 6]]
simplex_vertices(ind_simplex)

return the list of vertices of a index_simplex (i.e. a list of indices of vertices)

simplex_indices(list_of_vertices)

return the list of vertices corresponding to a list of indices

check()

check if it is really a simplicial complex

>>> K=SimplicialComplex(["a","b",2,3,4],maximal_simplices={1: [[0,1]]})
>>> print("K:", K)
K: Simplicial Complex of dim 1:
 vertices: ['a', 'b', 2, 3, 4]
 simplices:
  0: [['a'], ['b'], [2], [3], [4]]
  1: [['a', 'b']]
>>> print("check:", K.check() )
check: True
rename_vertices(other)

append the string other to the names of the vertices

class simplicial.SimplicialMap(f, K, L)

simple class to store the data of a simplicial map f : K -> L

is_simplicial()

check if it is really a simplicial map

in_chains()

chain map induced by f, represented as matrix, for each k = 0… dim K

simplicial.simplicial_sphere(vertici)

sfera simpliciale: il bordo di un simplesso standard di dimension n con vertici la lista vertici

simplicial.join_of_complexes(K, L, change_vertices=False)

make sure they do not have overlapping names…

simplicial.cartesian_product(K, L)

Cartesian product of two Simplicial Complexes

simplicial.connected_sum(K, L)

Connected sum of K and L : assume both have dim=2 and indices are ranges 0… n

simplicial.disjoint_union(K, L)

disjoint union of K and L

simplicial.display_2d(K, show_labels=False)

display_2d a 1-dimensional simplicial complex (in the plane). if vertices are points. Better in jupyter.

simplicial.display_3d(K, show_labels=False)

display a 3d simplicial complex using k3d, in jupyter

to install it: sudo pip install k3d sudo jupyter nbextension install –py k3d sudo jupyter nbextension enable –py k3d

simplicial.view(K, show_labels=False)

Display a euclidean simplicial complex

simplicial.boundary_operators(K)

gives the list list of boundary operators partial_k

simplicial.betti_numbers(K)

Betti numbers: b[k] dim of boundaries; c[k] dim of chains; z[k] dim of cycles

simplicial.sign_ordered(seq)

0 (zero) if it is not injective, otherwise sign of the permutation of the sequence seq of integers

simplicial.homology(K)

generators and projectors of each homology group it returns the 5-tuple L,shifted(D),R,r,H where

L,D,R are sequences of matrices, for k=0 .. dim(K), such that such that for each k the matrix D[k] is diagonal and

L[k] * boundary_operator[k] * R[k] = D[k]

r is the list of ranks of boundary_operator[:]

H is a list of dicts Each H[k] is a dict with two keys: gens and hcm.

H[k][‘gens’] is the list of generators of the k-th homology H[k][‘hcm’] is the homology class matrix. It is the projection from cycles to homology in base given by generators.

>>> S=simplicial_sphere(2)
>>> L,D,R,r,H = homology(S)
>>> for k in range(dim(S)+1):
...     print("k=", k, "\nL,D,R = ", L[k],D[k],R[k],"\nrank: ", r[k],"\n", H[k])
k= 0
L,D,R =  None None Matrix([[1, 0, 0, 0], [0, 1, 0, 0], [0, 0, 1, 0], [0, 0, 0, 1]])
rank:  0
 {'gens': [Matrix([
[0],
[0],
[0],
[1]])], 'hcm': Matrix([[1, 1, 1, 1]])}
k= 1
L,D,R =  Matrix([[-1, 0, 0, 0], [-1, -1, 0, 0], [-1, -1, -1, 0], [1, 1, 1, 1]]) Matrix([[0, 0, 0, 1, 0, 0], [0, 0, 0, 0, 1, 0], [0, 0, 0, 0, 0, 1], [0, 0, 0, 0, 0, 0]]) Matrix([[1, 0, 0, 1, -1, 0], [-1, 1, 0, 0, 1, -1], [0, -1, 0, 0, 0, 1], [1, -1, 1, 0, 0, 0], [0, 1, -1, 0, 0, 0], [0, 0, 1, 0, 0, 0]])
rank:  3
 {'gens': [], 'hcm': Matrix(0, 6, [])}
k= 2
L,D,R =  Matrix([[1, 0, 0, 0, 0, 0], [0, 1, 0, 0, 0, 0], [0, 0, 1, 0, 0, 0], [0, 0, 0, 1, 0, 0], [0, 0, 0, 0, 1, 0], [0, 0, 0, 0, 0, 1]]) Matrix([[0, 1, 0, 0], [0, 0, 1, 0], [0, 0, 0, 1], [0, 0, 0, 0], [0, 0, 0, 0], [0, 0, 0, 0]]) Matrix([[1, 1, -1, 1], [-1, 0, 1, -1], [1, 0, 0, 1], [-1, 0, 0, 0]])
rank:  3
 {'gens': [Matrix([
[ 1],
[-1],
[ 1],
[-1]])], 'hcm': Matrix([[0, 0, 0, -1]])}
simplicial.homology_of_map(f)

f is SimplicialMap: comput the matrix of its homology

>>> K=simplicial_sphere([ "A","B","C" ] )
>>> L=simplicial_sphere([ "a","b","c" ])
>>> f=SimplicialMap({"A":"a","B":"c","C":"b"},K,L)
>>> print(homology_of_map(f))
[Matrix([[1]]), Matrix([[-1]])]
simplicial.persistent_homology(list_of_maps)

very simple computation of persistent homology of a sequence of induced chain maps, using a naive row-echelon-form reduction of each homology morphism. It returns a list of pairs corresponding to time-of-birth and time-of-death, one for each generator of the homology.

>>> filtration=[
... SimplicialComplex([0,1,2]),
... SimplicialComplex([0,1,2],maximal_simplices={1:[[0,1]] } ),
... SimplicialComplex([0,1,2],maximal_simplices={1:[[0,1],[0,2]] } ) ,
... simplicial_sphere([0,1,2]) ,
... SimplicialComplex([0,1,2], maximal_simplices={2:[[0,1,2]] }) ,
... ]
>>> id_map={0:0,1:1,2:2}
>>> list_of_maps=[ SimplicialMap(id_map,filtration[j],filtration[j+1]) for j in
... range(len(filtration)-1)]
>>> PH=persistent_homology(list_of_maps)
>>> print("PH:", PH)
PH: {0: [(0, 4), (0, 1), (0, 2)], 1: [(3, 4)]}

It means that for the 0-th homology there are 3 generators, one survives all times, one dies soon, the other after two steps. In degree 1, there is just one generator, born at k=3 and dead in k=4.

class simplicial.Point(coords, label=None)
Just a point of a euclidean space of dimension 2 or 3 or more.

It may have a label.

>>> P=Point((1,0),label="A")
>>> print(P)
A=(1, 0)
>>> Q=Point((2,3),label="Q")
>>> print(Q)
Q=(2, 3)
>>> Z=Point((1,13,3))
>>> print(Z)
(1, 13, 3)
class simplicial.Embed(elastic_constant=1.0)

Class for finding a suitable euclidean embedding of simplicial complexes

optimize = <module 'scipy.optimize' from '/usr/lib/python3/dist-packages/scipy/optimize/__init__.py'>
Uij(P, Q, are_endpoints=False)

interaction of P and Q: edges elastic attract, all repulsive

grad_Uij(P, Q, are_endpoints=False)

grad of interaction of P and Q, wrt P

euclidean_embedding(K, dim=2)

return the Euclidean Embedding of the abstract complex K, using suitable interaction potentials (vertices and edges). it needs to import scipy.optimize.

OTHER:

simplicial.subsets(seq, k)

funzione ricorsiva/generatore per restituire i sottoinsiemi di k elementi della lista seq.

simplicial.EPchar(K)

Euler-Poincaré characteristic

>>> K=simplicial_sphere(3)
>>> EPchar(K)
0
>>> K=simplicial_sphere(2)
>>> EPchar(K)
2
simplicial.dim(K)

dimension of the simplicial complex K

simplicial.dims_c(K)

list of the numbers of simplices in each dimension

FOR MATRICES:

simplicial.row_echelon_form(M, get_pivots=False)
return the row echelon form E (almost reduced) of matrix A, together

with an invertible matrix L such that LM=E; if get_spivots is True, output the list of pivots as well.

>>> K=simplicial_sphere(2)
>>> print(row_echelon_form(boundary_operators(K)[1] ))
(Matrix([
[1, 1, 1, 0, 0, 0],
[0, 1, 1, 1, 1, 0],
[0, 0, 1, 0, 1, 1],
[0, 0, 0, 0, 0, 0]]), Matrix([
[-1,  0,  0, 0],
[-1, -1,  0, 0],
[-1, -1, -1, 0],
[ 1,  1,  1, 1]]))
simplicial.LDR(M)

return a diagonal matrix D with 1 on the diagonal, and 0 later, such that LMR = D, with L and R invertible

simplicial.shift_Z(c, r)

return the cxc matrix shifting the first r columns at the end, by right matrix multiplication

DOWNLOAD