homology
#
Simplicial complex, homology of surfaces and translation surfaces
In this module are implemented simple homology computation for translation surfaces. There are three main classes:
RibbonGraph
: decomposition of a surface into polygons. The combinatorics is stored as a triple of permutations v (vertices), e (edges), f (faces) so that the product vef is the identity. The domain of the permutations correspond to the half edges or darts. The permutation e is an involution so that e(i) is the other half of the edge starting at i. The fixed points of e corresponds to edge glued to themselves. The permutation v is obtained by turning around a vertex, while f turning around a face.RibbonGraphWithAngles
: a ribbon graph with an additional angle structure.RibbonGraphWithHolonomies
: a ribbon graph with an additional holonomy structure on its edges.
EXAMPLES:
sage: from surface_dynamics import *
To create a ribbon graph you just need to fix two of the permutations v, e, f:
sage: R = RibbonGraph(vertices='(0,1,4,3)(5,2)',edges='(0,3)(1,2)(4,5)')
sage: R
Ribbon graph with 2 vertices, 3 edges and 3 faces
The vertices, edges and faces are by definition the cycles of the permutation.
Calling the method vertices()
, edges()
or
faces()
gives you access to these cycles:
sage: R.vertices()
[[0, 1, 4, 3], [2, 5]]
sage: R.edges()
[[0, 3], [1, 2], [4, 5]]
sage: R.faces()
[[0, 4, 2], [1, 5], [3]]
Given a half edge (i.e. a dart), you can get the index of the vertex, edge or
face it belongs with the methods dart_to_vertex()
,
dart_to_edge()
and dart_to_edge()
:
sage: R.dart_to_vertex(1)
0
sage: 1 in R.vertices()[0]
True
sage: R.dart_to_vertex(2)
1
sage: 2 in R.vertices()[1]
True
sage: R.dart_to_edge(3)
0
sage: R.dart_to_edge(4)
2
sage: R.dart_to_face(4)
0
sage: R.dart_to_face(3)
2
To initialize a ribbon graph with angles, you have to input the standard data to initialize a ribbon graph plus a list of positive rational numbers which corresponds to the angles between darts (more precisely, the number at position i is the angle between i and v(i)):
sage: e = '(0,1)(2,3)'
sage: f = '(0,2,1,3)'
sage: a = [1/2,1/2,1/2,1/2]
sage: r = RibbonGraphWithAngles(edges=e,faces=f,angles=a)
sage: r.spin_parity()
1
TODO:
Oriented ribbon graphs are exactly separatrix diagrams. An orientation on a Ribbon graphs is a choice of orientation for each vertex so that each face has all of its edges oriented in the same direction as we go along the boundaries.
- class surface_dynamics.flat_surfaces.homology.RibbonGraph(vertices=None, edges=None, faces=None, connected=True, check=True)#
Generic class for Ribbon graph.
A Ribbon graph (or fat graph or combinatorial map) is a graph embedded in a surface. This class uses representation as a triple
(v,e,f)
of permutations such that vef = 1 and the action of the group generated by v,e,f acts transitively in the domain. The cycles ofv
are considered as vertices, the ones ofe
are considered as edges and the ones off
as the faces. Each element of the domain is a half-edge which is called a dart. A dart is also associated to an oriented edge.The domain of the permutations must be a subset of [0, …, N-1] for some N. The edges are always considered to be (i, ~i) where i is the bit complement of the integer i (~0 = -1, ~-3 = 2). So that an edge always has a canonical representative number given by the non-negative version.
A dense ribbon graph has the following attributes
total_darts - non negative integer - the total number darts
num_darts - non negative integer - the number of active darts
active_darts - bitset - list of lengths _total_darts with True or False. The position i is True if i is an active dart.
vertices, vertices_inv - list - partial permutations of [0,N] which are inverse of each other
vertex_cycles - the cycles of the partial permutation vertices
dart_to_vertex_index
edges, edges_inv - list - partial permutations of [0,N] which are inverse of each other
edge_cycles - the cycles of the partial permutation edge
dart_to_edge_index
faces, faces_inv - list - partial permutations of [0,N] which are inverse of each other
face_cycles - the cycles of the partial permutation faces
dart_to_face_index
EXAMPLES:
sage: from surface_dynamics import * sage: RibbonGraph([],[],[]) Ribbon graph with 1 vertex, 0 edge and 1 face sage: RibbonGraph('()','(0,1)','(0,1)') Ribbon graph with 2 vertices, 1 edge and 1 face sage: G = RibbonGraph('(0,3)(1,2)','(0,1)(2,3)','(0,2)(1,3)') sage: G Ribbon graph with 2 vertices, 2 edges and 2 faces sage: G.darts() [0, 1, 2, 3] sage: G.genus() 0 sage: G = RibbonGraph(edges='(0,2)(1,3)(4,6)(5,7)',faces='(0,1,2,3,4,5,6,7)') sage: G Ribbon graph with 1 vertex, 4 edges and 1 face sage: G.darts() [0, 1, 2, 3, 4, 5, 6, 7] sage: G.genus() 2 sage: G = RibbonGraph(vertices='(0,2,3,6)(1,4,5,7)') sage: G Ribbon graph with 2 vertices, 4 edges and 4 faces sage: G.edges() [[0, 1], [2, 3], [4, 5], [6, 7]] sage: G.genus() 0
- add_extra_darts(n)#
Add extra darts to support a total of
2 n
darts.
- automorphism_group(fix_vertices=False, fix_edges=False, fix_faces=False)#
Return the automorphism group of this ribbon graph.
The automorphism is simply the intersection of the centralizers of the defining permutations.
EXAMPLES:
sage: from surface_dynamics import * sage: r = RibbonGraph('(1,6,5)(2,3,4)', '(1,2)(3,4)(5,6)', '(1,4,2,5)(3)(6)') sage: r.automorphism_group().cardinality() 2 sage: r.automorphism_group(fix_faces=True).cardinality() 1 sage: r = RibbonGraph('(1,4,5)(2,6,3)', '(1,2)(3,4)(5,6)', '(1,3)(2,5)(4,6)') sage: r.automorphism_group().cardinality() 6
Examples in genus 1:
sage: r = RibbonGraph('(1,5,4)(6,3,2)', '(1,2)(3,4)(5,6)', '(1,3,5,2,4,6)') sage: A = r.automorphism_group() sage: A Subgroup ... sage: A.cardinality() 6 sage: r.automorphism_group(fix_faces=True) == A True sage: r = RibbonGraph('(1,3,2,4)', '(1,2)(3,4)', '(1,3,2,4)') sage: r.automorphism_group().cardinality() 4
- boundaries()#
Return the list of cycles which are boundaries.
A cycle is a boundary if it bounds a face.
EXAMPLES:
sage: from surface_dynamics import * sage: r = RibbonGraph('(1,2,3)(4,5,6)','(1,2)(3,4)(5,6)') sage: r.boundaries() [[(1, 2)], [(2, 1), (3, 4), (6, 5), (4, 3)], [(5, 6)]] sage: r = RibbonGraph('(1,2,3)(4,5)(6,7,8)',edges='(1,2)(3,4)(5,6)(7,8)') sage: r.boundaries() [[(1, 2)], [(2, 1), (3, 4), (5, 6), (8, 7), (6, 5), (4, 3)], [(7, 8)]]
- collapse(spanning_tree=None)#
Return a ribbon graph callapsed along a spanning tree.
The resulting graph is on the same surface as the preceding but has only one vertex. It could be used twice to provide a polygonal representation with one vertex and one face.
EXAMPLES:
sage: from surface_dynamics import * sage: R = RibbonGraph(vertices='(0,1,2,5)(3,7)(4,10,9)(6,11,12)(8,13)') sage: R.genus() 1 sage: R.num_vertices() 5 sage: R.num_edges() 7 sage: R.num_faces() 2 sage: R2 = R.collapse() sage: R2 Ribbon graph with 1 vertex, 3 edges and 2 faces sage: R Ribbon graph with 5 vertices, 7 edges and 2 faces sage: R3 = R2.dual().collapse().dual() sage: R3 Ribbon graph with 1 vertex, 2 edges and 1 face
- cycle_basis(intersection=False, verbose=False)#
Returns a base of oriented cycles of the Ribbon graph modulo boundaries.
If
intersection
is set to True then the method also returns the intersection matrix of the cycles.EXAMPLES:
sage: from surface_dynamics import * sage: r = RibbonGraph('(1,2,3)(4,5,6)','(1,2)(3,4)(5,6)') sage: r.cycle_basis() [] sage: r = RibbonGraph('(1,2,3)(4,5)(6,7,8)',edges='(1,2)(3,4)(5,6)(7,8)') sage: r.cycle_basis() [] sage: r = RibbonGraph('(1,4,5)(2,3)(6,7,8)',edges='(1,2)(3,4)(5,6)(7,8)') sage: r.cycle_basis() [] sage: e = '(1,3)(2,4)(5,7)(6,8)' sage: f = '(1,2,3,4,5,6,7,8)' sage: r = RibbonGraph(edges=e,faces=f) sage: r.cycle_basis() [[[1, 3]], [[2, 4]], [[5, 7]], [[6, 8]]] sage: f = '(0,10,13)(6,17,11)(2,14,7)(15,12,3)(16,20,19)(18,1,9)(4,22,21)(23,8,5)' sage: e = tuple((i,i+1) for i in range(0,24,2)) sage: r = RibbonGraph(edges=e,faces=f); r Ribbon graph with 2 vertices, 12 edges and 8 faces sage: c,m = r.cycle_basis(intersection=True) sage: c [[(0, 1), [4, 5]], [[8, 9]], [[12, 13]], [[14, 15], (1, 0)]] sage: m [ 0 1 0 0] [-1 0 0 0] [ 0 0 0 1] [ 0 0 -1 0]
- dart_to_edge(i, orientation=False)#
Returns the edge the darts
i
belongs to.If orientation is set to
True
then the output is a 2-tuple(e,o)
wheree
is the index of the edge ando
is its orientation as+1
or-1
.
- dart_to_face(i)#
- dart_to_vertex(i)#
Return the vertex on which the dart
i
is attached.
- darts()#
Return the list of darts
- dual()#
Returns the dual Ribbon graph.
The dual ribbon graph of (v,e,f) is (f^{-1}, e, v^{-1}).
EXAMPLES:
sage: from surface_dynamics import * sage: r = RibbonGraph(edges='(0,1)',faces='(0)(1)'); r Ribbon graph with 1 vertex, 1 edge and 2 faces sage: r.dual() Ribbon graph with 2 vertices, 1 edge and 1 face
- edge_orbit(i)#
Return the orbit of the dart
i
under the permutation that defines the edges.
- edge_perm()#
Return the permutation that define the edges.
- edges()#
Return the set of edges.
- euler_characteristic()#
Returns the Euler characteristic of the embedded surface.
The Euler characteristic of a surface complex is V - E + F, where V is the number of vertices, E the number of edges and F the number of faces.
EXAMPLES:
sage: from surface_dynamics import * sage: r = RibbonGraph(edges='(0,1)(2,3)(4,5)',faces='(0,2,4)(1)(3,5)') sage: r.euler_characteristic() 2 sage: r = RibbonGraph(edges='(0,1)(2,3)',faces='(0,2,1,3)') sage: r.euler_characteristic() 0
- face_orbit(i)#
Return the orbit of
i
under the permutation associated to faces.
- face_perm()#
Return the permutation that defines the face.
- faces()#
Return the list of faces.
- genus()#
Return the genus of the surface associated to this Ribbon graph.
EXAMPLES:
sage: from surface_dynamics import * sage: R = RibbonGraph(vertices='(1)(2)',edges='(1,2)') sage: R.genus() 0 sage: e='(1,3)(2,4)' sage: f='(1,2,3,4)' sage: RibbonGraph(edges=e,faces=f).genus() 1 sage: e='(1,3)(2,4)(5,7)(6,8)' sage: f='(1,2,3,4,5,6,7,8)' sage: RibbonGraph(edges=e,faces=f).genus() 2 sage: e='(1,3)(2,4)(5,7)(6,8)(9,11)(10,12)' sage: f='(1,2,3,4,5,6,7,8,9,10,11,12)' sage: RibbonGraph(edges=e,faces=f).genus() 3
- is_connected(force_computation=False)#
- is_cycle(c)#
Test whether
c
is a cycle.A path is a sequence of oriented edges such that each edge starts where the preceding one ends. A cycle is a path which starts where it ends.
- is_plane()#
Returns true if and only if the ribbon graph belongs in a sphere. In other words if it has genus 0.
EXAMPLES:
sage: from surface_dynamics import * sage: r = RibbonGraph(vertices='(0)(1)',edges='(0,1)') sage: r.is_plane() True sage: r = RibbonGraph(vertices='(0,1)',edges='(0,1)') sage: r.is_plane() True sage: r = RibbonGraph(edges='(0,1)(2,3)',faces='(0,2)(1,3)') sage: r.is_plane() True sage: r = RibbonGraph(edges='(0,1)(2,3)',faces='(0,2,1,3)') sage: r.is_plane() False
- is_plane_tree()#
Returns True if and only if the ribbon graph is a planar tree. In other words, it has genus 0 and only one face.
EXAMPLES:
sage: from surface_dynamics import * sage: r = RibbonGraph(vertices='(0)(1)',edges='(0,1)') sage: r.is_plane_tree() True sage: r = RibbonGraph(vertices='(0)(1,2,4)(3)(5)',edges='(0,1)(2,3)(4,5)') sage: r.is_plane_tree() True sage: r = RibbonGraph(vertices='(0,1)',edges='(0,1)') sage: r.is_plane_tree() False sage: r.is_plane() True
- is_triangulated()#
Returns True if the surface is triangulated. In other words, faces consist only of the product of 3-cycles.
EXAMPLES:
sage: from surface_dynamics import * sage: r = RibbonGraph(edges='(0,1)(2,3)(4,5)',faces='(0,2,4)(1,5,3)') sage: r.is_triangulated() True sage: r = RibbonGraph(edges='(0,1)(2,3)',faces='(0,2,1,3)') sage: r.is_triangulated() False
- length_rational_fraction(var='b')#
Return the generating series for the number of lengths with the given boundaries
- monodromy_group()#
Return the group generated by the three defining permutations.
EXAMPLES:
sage: from surface_dynamics import RibbonGraph sage: r1 = RibbonGraph(vertices='(0)(2)(1,3,4,5)', edges='(0,1)(2,3)(4,5)') sage: G1 = r1.monodromy_group() sage: G1 Subgroup ... sage: G1.is_isomorphic(SymmetricGroup(5)) True sage: r2 = RibbonGraph(vertices='(0)(2)(1,3,4)(5,6,7)', edges='(0,1)(2,3)(4,5)(6,7)') sage: G2 = r2.monodromy_group() sage: G2 Subgroup ... sage: G2.is_isomorphic(PSL(2,7)) True sage: r3 = RibbonGraph(vertices='(1)(5)', edges='(1,5)', faces='(1,5)') sage: r3.monodromy_group() Subgroup ...
- num_darts()#
Returns the number of darts.
- num_edges()#
Returns the number of edges.
- num_faces()#
Return the number of faces.
- num_vertices()#
Returns the number of vertices.
- relabel(perm=None)#
perm is a of range(0,N)
If
perm
is None, relabel the darts on 0,2M keeping the relative order of the darts.
- spanning_tree()#
Return a spanning tree
OUTPUT:
spanning tree as a DiGraph
remaining edges as 2-tuples
(i,e[i])
EXAMPLES:
sage: from surface_dynamics import * sage: R = RibbonGraph('(1,2,3)','(1,2)(3,4)') sage: R Ribbon graph with 2 vertices, 2 edges and 2 faces sage: T,o = R.spanning_tree() sage: T Digraph on 2 vertices sage: T.edges(sort=True) [(0, 1, (3, 4))] sage: o [(1, 2)] sage: R = RibbonGraph('(1,2,3)(4,5,6)','(1,2)(3,4)(5,6)') sage: R Ribbon graph with 2 vertices, 3 edges and 3 faces sage: T,o = R.spanning_tree() sage: T Digraph on 2 vertices sage: T.edges(sort=True) [(0, 1, (3, 4))] sage: o [(1, 2), (5, 6)] sage: e = '(1,3)(5,7)(2,4)(6,8)' sage: f = '(1,2,3,4,5,6,7,8)' sage: R = RibbonGraph(edges=e, faces=f) sage: T,o = R.spanning_tree() sage: T Digraph on 1 vertex sage: o [[1, 3], [2, 4], [5, 7], [6, 8]]
- vertex_orbit(i)#
Return the orbit of
i
under the permutation that define the vertices.
- vertex_perm()#
Returns the permutation that define the vertices.
- vertices()#
Return the list of vertices as cycles decomposition of the vertex permutation.
- class surface_dynamics.flat_surfaces.homology.RibbonGraphWithAngles(vertices=None, edges=None, faces=None, angles=None)#
A Ribbon graph with angles between edges
Currently angles can only be rational multiples of pi.
TODO:
allows any kind of angles by providing a sum for the total and considering each angle as a (projective) portion of the total angle.
- angle_at_vertex(v)#
Angle at a vertex (coefficient of pi)
- angle_at_vertices()#
Return the list of angles at a vertex.
- angle_between_darts(d1, d2)#
Return the angle between the darts
d1
andd2
- has_trivial_holonomy()#
Test whether self has trivial holonomy representation
EXAMPLES:
sage: from surface_dynamics import * sage: e = '(0,1)(2,3)' sage: f = '(0,2,1,3)' sage: a = [1/2,1/2,1/2,1/2] sage: r = RibbonGraphWithAngles(edges=e,faces=f,angles=a) sage: r.has_trivial_holonomy() True sage: e = '(0,1)(2,3)(4,5)' sage: f = '(0,2,4)(1,5,3)' sage: a = [1/3,7/15,1/5,1/5,7/15,1/3] sage: r = RibbonGraphWithAngles(edges=e,faces=f,angles=a) sage: r.has_trivial_holonomy() False
- holonomy_representation()#
Return the holonomy representation in SO(2) as two lists.
The first list correspond to cycles around vertices, while the second correspond to a cycle basis that generate homology.
EXAMPLES:
sage: from surface_dynamics import * sage: e = '(0,1)(2,3)' sage: f = '(0,2,1,3)' sage: a = [1/2,1/2,1/2,1/2] sage: r = RibbonGraphWithAngles(edges=e,faces=f,angles=a) sage: r.holonomy_representation() ([0], [0, 0])
The standard cube:
sage: e = tuple((i,i+1) for i in range(0,24,2)) sage: f = '(0,20,7,10)(16,22,19,21)(2,9,5,23)(14,3,17,1)(12,8,15,11)(18,4,13,6)' sage: a = [1/2]*24 sage: r = RibbonGraphWithAngles(edges=e,faces=f,angles=a) sage: r.holonomy_representation() ([3/2, 3/2, 3/2, 3/2, 3/2, 3/2, 3/2, 3/2], [])
Two copies of a triangle:
sage: e = '(0,1)(2,3)(4,5)' sage: f = '(0,2,4)(1,5,3)' sage: a = [1/2,1/6,1/3,1/3,1/6,1/2] sage: r = RibbonGraphWithAngles(edges=e,faces=f,angles=a) sage: r.holonomy_representation() ([1, 1/2, 1/2], []) sage: a = [1/3,7/15,1/5,1/5,7/15,1/3] sage: r = RibbonGraphWithAngles(edges=e,faces=f,angles=a) sage: r.holonomy_representation() ([2/3, 2/3, 2/3], [])
- spin_parity(check=True, verbose=False)#
Return the spin parity of the Ribbon graph with angles.
The surface should be holonomy free and with odd multiple of 2 pi angles. Implements the formula of [Joh80].
EXAMPLES:
sage: from surface_dynamics import *
We first consider the case of the torus:
sage: e = '(0,1)(2,3)' sage: f = '(0,2,1,3)' sage: a = [1/2,1/2,1/2,1/2] sage: r = RibbonGraphWithAngles(edges=e,faces=f,angles=a) sage: r.spin_parity() 1
Then the case of genus 2 surface (with an angle of 6pi):
sage: e = '(0,1)(2,3)(4,5)(6,7)' sage: f = '(0,2,4,3,6,1,7,5)' sage: a = [1/2,1/2,1,1/2,1/2,1,3/2,1/2] sage: r = RibbonGraphWithAngles(edges=e,faces=f,angles=a) sage: r.spin_parity() 1 sage: e = '(0,1)(2,3)(4,5)(6,7)' sage: f = '(0,2,4,6,1,3,5,7)' sage: a = [1/2,1/2,1,1,1,1,1/2,1/2] sage: r = RibbonGraphWithAngles(edges=e,faces=f,angles=a) sage: r.spin_parity() 1 sage: e = '(0,1)(2,3)(4,5)(6,7)' sage: f = '(0,2,4,6,1,3,5,7)' sage: a = [3/4]*8 sage: r = RibbonGraphWithAngles(edges=e,faces=f,angles=a) sage: r.spin_parity() 1
In genus 3 two spin parities occur for one conical angle 10pi:
sage: e = '(0,1)(2,3)(4,5)(6,7)(8,9)(10,11)' sage: f1 = '(0,4,6,8,10,2,1,9,11,5,7,3)' sage: f2 = '(0,4,6,8,10,2,1,5,7,9,11,3)' sage: a = [1/2,1/2,1/2,1/2] + [1]*8 sage: r1 = RibbonGraphWithAngles(edges=e,faces=f1,angles=a) sage: r1.spin_parity() 1 sage: r2 = RibbonGraphWithAngles(edges=e,faces=f2,angles=a) sage: r2.spin_parity() 0
- winding(c)#
Return winding number along the cycle
c
.This is NOT well defined because it depends on the way we choose to pass on the left or on the right at singularity.
- class surface_dynamics.flat_surfaces.homology.RibbonGraphWithHolonomies(vertices=None, edges=None, faces=None, holonomies=None)#
A Ribbon graph with holonomies.
For now
- surface_dynamics.flat_surfaces.homology.angle(v)#
Return the argument of the vector
v
.