Surface topology and geometry¶
Separatrix diagrams¶
Separatrix diagrams and cylinder diagrams
A separatrix diagram is a couple of permutation (bot,top)
that have the same
number of cycles in their cycle decompositions. A cylinder diagram is a
separatrix diagram together with a bijection between the cycles of bot
and
top
.
A cylinder diagram encodes the combinatorics of cylinder decomposition of a completely periodic direction in a translation surface. If we adjoin coordinates to this combinatorial datum, we have a complete description of the underlying surface. In the case of arithmetic curves, the coordinates can be taken to be rational numbers.
This representation of a surface is used in various constructions:
square tiled surfaces
Thurston-Veech construction of pseudo-Anosov diffeomorphism
description of the cusp of Teichmueller curves
EXAMPLES:
sage: from surface_dynamics import *
Separatrix diagrams:
sage: s = SeparatrixDiagram('(0,1,2)(3,4)(5,6,7)','(0,4,1,2)(3,7)(5,6)')
sage: s
(0,1,2)(3,4)(5,6,7)-(0,4,1,2)(3,7)(5,6)
sage: s.bot_cycle_tuples()
[(0, 1, 2), (3, 4), (5, 6, 7)]
sage: s.top_cycle_tuples()
[(0, 4, 1, 2), (3, 7), (5, 6)]
Cylinder diagrams:
sage: c = CylinderDiagram([((0,),(4,)),((1,2),(0,1,3)),((3,4),(2,))])
sage: print(c)
(0)-(4) (1,2)-(0,1,3) (3,4)-(2)
sage: print(c.separatrix_diagram())
(0)(1,2)(3,4)-(0,1,3)(2)(4)
They can also be built from separatrix diagram:
sage: s = SeparatrixDiagram('(0,1,2)(3,4)(5,6,7)','(0,4,1,2)(3,7)(5,6)')
sage: s
(0,1,2)(3,4)(5,6,7)-(0,4,1,2)(3,7)(5,6)
sage: s.to_cylinder_diagram([(0,1),(1,0),(2,2)])
(0,1,2)-(3,7) (3,4)-(0,4,1,2) (5,6,7)-(5,6)
- class surface_dynamics.flat_surfaces.separatrix_diagram.CylinderDiagram(data, check=True)[source]¶
Bases:
SeparatrixDiagram
Separatrix diagram with pairing.
Each cylinder is stored as a couple (bot,top) for which the orientation is as follows:
+--------------------+ | <-- top -- | | | | | | -- bot --> | +--------------------+
INPUT:
data
- list of 2-tuples - matching of bottom-top pairs
EXAMPLES:
sage: from surface_dynamics import *
We first build the simplest cylinder diagram which corresponds to a torus:
sage: CylinderDiagram([((0,),(0,))]) (0)-(0)
The same initialized from a string:
sage: CylinderDiagram('(0)-(0)') (0)-(0)
The following initialize a cylinder diagram with two cylinder which gives a surface of genus 2 with one singularity of degree 2:
sage: CylinderDiagram([((0,1),(0,2)),((2,),(1,))]) (0,1)-(0,2) (2)-(1)
ALGORITHM:
A cylinder is represented by a couple (i,j) where i is the min in bot and j is the min in top. The data _top_to_cyl and _bot_to_cyl corresponds to the association of a separatrix to the corresponding 2-tuple. The twist coordinate correspond to the shift between those two indices.
- an_origami()[source]¶
Return one origami with this diagram cylinder if any.
EXAMPLES:
sage: from surface_dynamics import * sage: cyl = CylinderDiagram('(0,1)-(0,2) (2,3)-(1,3)') sage: cyl.an_origami() (1,2)(3,4) (1,3,4,2)
- automorphism_group(order=False)[source]¶
Return the automorphism group
INPUT:
order
- boolean (default: False) - whether or not return the order of the group
EXAMPLES:
sage: from surface_dynamics import * sage: cyl = CylinderDiagram('(0,1)-(0,2) (2,3)-(1,3)') sage: cyl.automorphism_group() Permutation Group with generators [(0,3)(1,2)]
- bot_to_cyl(j)[source]¶
Return the cylinder above the
j
-th separatrix.EXAMPLES:
sage: from surface_dynamics import * sage: c = CylinderDiagram('(0,2,4)-(1,3,5) (1,5)-(0) (3)-(2,4)') sage: c (0,2,4)-(1,3,5) (1,5)-(0) (3)-(2,4) sage: c.bot_to_cyl(0) ((0, 2, 4), (1, 3, 5)) sage: c.bot_to_cyl(1) ((1, 5), (0,)) sage: c.bot_to_cyl(3) ((3,), (2, 4))
- canonical_label(inplace=False, return_map=False)[source]¶
Return a cylinder diagram with canonical labels.
Note
The canonical label might change depending on your Sage version and the optional packages available.
EXAMPLES:
sage: from surface_dynamics import CylinderDiagram sage: c1 = CylinderDiagram('(0,5,4)-(0,3,2,1) (1,3,2)-(4,5)') sage: c1.canonical_label() # random (0,4,5)-(1,3,2,5) (1,3,2)-(0,4) sage: c2 = c1.relabel([2,4,5,1,3,0]) sage: c3 = c1.relabel([5,3,0,1,4,2]) sage: c1.canonical_label() == c2.canonical_label() == c3.canonical_label() True
- cylcoord_to_origami(lengths, heights, twists=None)[source]¶
Convert coordinates of the cylinders into an origami.
INPUT:
lengths
- lengths of the separatricesheights
- heights of the cylinderstwists
- twists for cylinders
EXAMPLES:
sage: from surface_dynamics import * sage: c = CylinderDiagram([((0,1),(1,2)),((2,),(0,))]) sage: c.stratum() H_2(2) sage: c.cylcoord_to_origami([1,1,1],[1,1]).stratum() H_2(2) sage: o1 = c.cylcoord_to_origami([2,1,2],[1,1],[1,0]) sage: o1 = o1.relabel() sage: o2 = c.cylcoord_to_origami([2,1,2],[1,1],[0,1]) sage: o2 = o2.relabel() sage: o3 = c.cylcoord_to_origami([2,1,2],[1,1],[1,1]) sage: o3 = o3.relabel() sage: all(o.stratum() == Stratum([2], k=1) for o in [o1,o2,o3]) True sage: o1 == o2 or o1 == o3 or o3 == o1 False
If the lengths are not compatible with the cylinder diagram a ValueError is raised:
sage: c.cylcoord_to_origami([1,2,3],[1,1]) Traceback (most recent call last): ... ValueError: lengths are not compatible with cylinder equations
- cylcoord_to_origami_iterator(lengths, heights)[source]¶
Convert coordinates of the cylinders into an origami.
INPUT:
lengths
- lengths of the separatricesheights
- heights of the cylinders
OUTPUT:
iterator over all possible origamis with those lengths and heights…
EXAMPLES:
sage: from surface_dynamics import * sage: cyl = CylinderDiagram('(0,1,2)-(3,1,2) (3)-(0)') sage: for o in cyl.cylcoord_to_origami_iterator((1,1,1,1),(1,1)): ....: print(o) (1,2,3)(4) (1,4)(2,3) (1,2,3)(4) (1,2,4)(3) (1,2,3)(4) (1,3,4)(2)
The number of origamis generated is just the product of the widths:
sage: sum(1 for _ in cyl.cylcoord_to_origami_iterator((2,1,1,2),(3,2))) 8
- cylinder_graph()[source]¶
Return the cylinder graph.
The cylinder graph is the graph whose vertex set are the cylinders and for each saddle connection there is a directed edge from the adjacent cylinders. The multiplicities are encoded in the labels.
EXAMPLES:
sage: from surface_dynamics import CylinderDiagram sage: c = CylinderDiagram('(0,1,5)-(2,5) (2)-(0,1,3) (3,4)-(4)') sage: c.cylinder_graph() Looped digraph on 3 vertices sage: c.cylinder_graph().edges(sort=True) [(0, 0, 1), (0, 1, 1), (1, 0, 2), (1, 2, 1), (2, 2, 1)] sage: c = CylinderDiagram('(0,1,3,5)-(2,5,3) (2,4)-(0,4,1)') sage: c.cylinder_graph().edges(sort=True) [(0, 0, 2), (0, 1, 1), (1, 0, 2), (1, 1, 1)]
- cylinders()[source]¶
Return the cylinders as a list of pairs
(bot, top)
.EXAMPLES:
sage: from surface_dynamics import * sage: c = CylinderDiagram('(0,2,4)-(1,3,5) (1,5)-(0) (3)-(2,4)') sage: c (0,2,4)-(1,3,5) (1,5)-(0) (3)-(2,4) sage: c.cylinders() [((0, 2, 4), (1, 3, 5)), ((1, 5), (0,)), ((3,), (2, 4))]
- dual_graph()[source]¶
The dual graph of the stable curve at infinity in the horizontal direction.
This graph is defines as follows. Cut each horizontal cylinder along a circumference, then the vertices are the equivalence class of half cylinder modulo the relation “linked by a saddle connection” and the edges are the circumferences.
EXAMPLES:
sage: from surface_dynamics import *
We consider the three diagrams of the stratum H(1,1):
sage: c1 = CylinderDiagram('(0,1,2,3)-(0,1,2,3)') sage: c1.stratum() H_2(1^2) sage: c1.dual_graph() Looped multi-graph on 1 vertex sage: c2 = CylinderDiagram('(0,1)-(1,2) (2,3)-(0,3)') sage: c2.stratum() H_2(1^2) sage: c2.dual_graph() Looped multi-graph on 1 vertex sage: c3 = CylinderDiagram('(0,1)-(2,3) (2)-(0) (3)-(1)') sage: c3.stratum() H_2(1^2) sage: c3.dual_graph() Looped multi-graph on 2 vertices
- horizontal_symmetry()[source]¶
Return the cylinder diagram obtained by reflecting the cylinder configuration along the horizontal axis.
EXAMPLES:
sage: from surface_dynamics import * sage: c = CylinderDiagram('(0,3,4)-(0,3,5) (1,2,5)-(1,2,4)') sage: c.horizontal_symmetry() (0,5,3)-(0,4,3) (1,4,2)-(1,5,2) sage: c.separatrix_diagram().horizontal_symmetry() == c.horizontal_symmetry().separatrix_diagram() True sage: A = Stratum([2,2], k=1) sage: all(c.horizontal_symmetry().stratum() == A for c in A.cylinder_diagrams()) True
- inverse()[source]¶
Return the inverse cylinder diagram.
The inverse of a cylinder diagram is the cylinder diagram in which all cylinders have been reversed. It corresponds to the multiplication by -1 on the underlying Abelian differential.
Combinatorially the operation is b0-t0 … bk-tk becomes t0-b0 … tk-bk
EXAMPLES:
sage: from surface_dynamics import * sage: c = CylinderDiagram('(0,1)-(0,2) (3,5,4)-(1,4,6) (2,6)-(3,5)') sage: c (0,1)-(0,2) (2,6)-(3,5) (3,5,4)-(1,4,6) sage: c.inverse() (0,2)-(0,1) (1,4,6)-(3,5,4) (3,5)-(2,6)
The operation can also be defined at the level of the separatrix diagrams and the two operation commutes:
sage: c.separatrix_diagram().inverse() == c.inverse().separatrix_diagram() True
The inversion can also be seen as the composition of the horizontal and vertical symmetries:
sage: c.horizontal_symmetry().vertical_symmetry() == c.inverse() True sage: c.vertical_symmetry().horizontal_symmetry() == c.inverse() True
The inversion is an involution on cylinder diagrams:
sage: all(cc.inverse().inverse() == cc for cc in Stratum([4], k=1).cylinder_diagrams()) # long time True
- is_hyperelliptic(verbose=False)[source]¶
Test of hyperellipticity
Each stratum of Abelian differentials as up to three connected components. For the strata H(2g-2) and H(g-1,g-1) there is a special component called hyperelliptic in which all translation surfaces (X,omega) in that component are such that X is hyperelliptic.
This function returns True if and only if the cylinder diagrams correspond to a decomposition of a surface associated to the hyperelliptic components in H(2g-2) or H(g-1,g-1).
EXAMPLES:
sage: from surface_dynamics import *
In genus 2, strata H(2) and H(1,1), all surfaces are hyperelliptic:
sage: for c in Stratum([2], k=1).cylinder_diagrams(): ....: print("%d %s" % (c.ncyls(), c.is_hyperelliptic())) 1 True 2 True sage: for c in Stratum([1,1], k=1).cylinder_diagrams(): ....: print("%d %s" % (c.ncyls(), c.is_hyperelliptic())) 1 True 2 True 2 True 3 True
In higher genera, some of them are, some of them are not:
sage: C = Stratum([4], k=1).cylinder_diagrams() sage: len(C) 15 sage: sum(c.is_hyperelliptic() for c in C) 5 sage: C = Stratum([2,2], k=1).cylinder_diagrams() sage: len(C) 41 sage: sum(c.is_hyperelliptic() for c in C) 11
- is_isomorphic(other, return_map=False)[source]¶
Test whether this cylinder diagram is isomorphic to
other
.EXAMPLES:
sage: from surface_dynamics import * sage: c1 = CylinderDiagram('(0,2,1)-(3,4,5) (3)-(1) (4)-(2) (5)-(0)') sage: c2 = CylinderDiagram('(0,2,1)-(3,5,4) (3)-(1) (4)-(2) (5)-(0)') sage: c1.is_isomorphic(c2) False sage: c3 = CylinderDiagram('(1,3,5)-(2,4,0) (2)-(5) (4)-(3) (0)-(1)') sage: c1.is_isomorphic(c3) True sage: ans, m = c1.is_isomorphic(c3, return_map=True) sage: assert ans is True and c1.relabel(m) == c3 sage: c4 = CylinderDiagram('(4,2,3)-(1,5,0) (1)-(3) (0)-(2) (5)-(4)') sage: c2.is_isomorphic(c4) True sage: ans, m = c2.is_isomorphic(c4, return_map=True) sage: assert ans is True and c2.relabel(m) == c4 sage: c1 = CylinderDiagram('(0,5)-(3,4) (1,4)-(2,5) (2)-(0) (3)-(1)') sage: c2 = CylinderDiagram('(0,5)-(3,4) (1,4)-(2,5) (2)-(1) (3)-(0)') sage: c1.is_isomorphic(c2) False sage: c1 = CylinderDiagram('(0,1)-(4,5) (2,4)-(0,3) (3,5)-(1,2)') sage: c2 = CylinderDiagram('(0,1)-(2,3) (2,4)-(1,4) (3,5)-(0,5)') sage: c3 = CylinderDiagram('(0,5)-(2,5) (1,4)-(0,4) (2,3)-(1,3)') sage: c1.is_isomorphic(c2) False sage: c1.is_isomorphic(c3) False sage: c2.is_isomorphic(c3) False
- lengths_cone()[source]¶
Return the rational polyhedron corresponding to the set of length with the given fixed heights.
-> one can obtain ehrhart series for each of them! It tells us that we have a nice asymptotics… and the asymptotics is simply given by the volume of this polytope (up to the ignored twists parameters)!
- lengths_minimal_solution()[source]¶
Return an integral solution for the lengths equation with minimal sum
- matrix_relation()[source]¶
Return the matrix of relation on the lengths of the separatrices.
The output matrix has size ncyls times nseps.
EXAMPLES:
sage: from surface_dynamics import *
For a one cylinder diagram, there is no relations:
sage: cyl = CylinderDiagram('(0,1,2,3)-(0,1,2,3)') sage: cyl.matrix_relation() [0 0 0 0]
Here is an example in the stratum H(2):
sage: cyl = CylinderDiagram('(0,1)-(0,2) (2)-(1)') sage: cyl.stratum() H_2(2) sage: cyl.matrix_relation() [ 0 1 -1] [ 0 -1 1]
- origami_iterator(n)[source]¶
Iteration over all origamis with n squares.
INPUT:
n
- positive integer - the number of squares
EXAMPLES:
sage: from surface_dynamics import * sage: cyl = CylinderDiagram('(0,1,2)-(3,1,2) (3)-(0)') sage: for o in cyl.origami_iterator(4): ....: print(o) ....: print(o.stratum()) ....: print(o.nb_squares()) (1,2,3)(4) (1,4)(2,3) H_2(1^2) 4 (1,2,3)(4) (1,2,4)(3) H_2(1^2) 4 (1,2,3)(4) (1,3,4)(2) H_2(1^2) 4
- origamis(n=None)[source]¶
Return the set of origamis having
n
squares.If
n
is None then return the origamis with less number of squares.EXAMPLES:
sage: from surface_dynamics import * sage: cyl = CylinderDiagram('(0,1,2)-(0,1,3) (3)-(2)') sage: o5 = cyl.origamis(5) sage: o5[0] (1,2,3,4)(5) (1,5,4,2,3) sage: o5[1].nb_squares() 5 sage: o5[2].stratum_component() H_2(1^2)^hyp
- relabel(perm, inplace=False)[source]¶
EXAMPLES:
sage: from surface_dynamics import CylinderDiagram sage: c = CylinderDiagram('(0,3)-(3,7) (1,6,4)-(0,2,4,6,8) (2,8,7,5)-(1,5)') sage: c.relabel([3,1,2,6,4,5,0,7,8]) (0,4,1)-(0,8,3,2,4) (2,8,7,5)-(1,5) (3,6)-(6,7) sage: c (0,3)-(3,7) (1,6,4)-(0,2,4,6,8) (2,8,7,5)-(1,5) sage: c.relabel([3,1,2,6,4,5,0,7,8], inplace=True) (0,4,1)-(0,8,3,2,4) (2,8,7,5)-(1,5) (3,6)-(6,7) sage: c (0,4,1)-(0,8,3,2,4) (2,8,7,5)-(1,5) (3,6)-(6,7)
- separatrix_diagram()[source]¶
Return the underlying separatrix diagram
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1)(2,3,4)','(0,3)(1,4,2)'); s (0,1)(2,3,4)-(0,3)(1,4,2) sage: c = s.to_cylinder_diagram([(0,1),(1,0)]); c (0,1)-(1,4,2) (2,3,4)-(0,3) sage: c.separatrix_diagram() == s True
- smallest_integer_lengths()[source]¶
Check if there is a integer solution that satisfy the cylinder conditions.
If there is a solution, the function returns a list a pair
(total_length, list_of_lengths)
that consists of the sum of the length of the separatrices together with the list of lengths. Otherwise, returns False.EXAMPLES:
sage: from surface_dynamics import * sage: c = CylinderDiagram('(0,1)-(0,2) (2,3)-(1,3)') sage: c.smallest_integer_lengths() (4, [1, 1, 1, 1]) sage: c = CylinderDiagram('(0,1,2)-(3) (3)-(0) (4)-(1,2,4)') sage: c.smallest_integer_lengths() False sage: c = CylinderDiagram('(0,1)-(0,5) (2)-(3) (3,6)-(8) (4,8)-(6,7) (5)-(2,4) (7)-(1)') sage: c.smallest_integer_lengths() (13, [1, 2, 1, 1, 1, 2, 1, 2, 2])
- spin_parity()[source]¶
Return the spin parity of any surface that is built from this cylinder diagram.
EXAMPLES:
sage: from surface_dynamics import * sage: c = CylinderDiagram('(0,1,2,3,4)-(0,1,2,3,4)') sage: c.spin_parity() 0 sage: c = CylinderDiagram('(0,1,2,3,4)-(0,1,4,2,3)') sage: c.spin_parity() 1 sage: c = CylinderDiagram('(0,2,6,1)-(0,8,1,9,2,5,7,4) (3,7,4,8,9,5)-(3,6)') sage: c.spin_parity() 0 sage: c = CylinderDiagram('(3,7,4,8,9,5)-(0,8,1,9,2,5,7,4) (0,2,6,1)-(3,6)') sage: c.spin_parity() 1
- stratum_component()[source]¶
Return the connected component of stratum of
self
.EXAMPLES:
sage: from surface_dynamics import * sage: CylinderDiagram('(0,1)-(0,2) (2)-(1)').stratum_component() H_2(2)^hyp sage: c = CylinderDiagram('(0,3,2,1)-(1,4,3,2) (4,7,6,5)-(0,7,6,5)') sage: c.stratum_component() H_4(3^2)^hyp sage: c = CylinderDiagram('(0,1,4)-(1,6,7) (2,5,3)-(0,2,4) (6)-(5) (7)-(3)') sage: c.stratum_component() H_4(3^2)^nonhyp sage: c = CylinderDiagram('(0,6)-(1,7) (1,5,4,3,2)-(2,6,5,4,3) (7,9,8)-(0,9,8)') sage: c.stratum_component() H_5(4^2)^hyp sage: c = CylinderDiagram('(0,2,6,1)-(0,8,1,9,2,5,7,4) (3,7,4,8,9,5)-(3,6)') sage: c.stratum_component() H_5(4^2)^even sage: c = CylinderDiagram('(3,7,4,8,9,5)-(0,8,1,9,2,5,7,4) (0,2,6,1)-(3,6)') sage: c.stratum_component() H_5(4^2)^odd
- symmetries()[source]¶
Return a triple
(horiz_sym, vert_sym, inv_sym)
EXAMPLES:
sage: from surface_dynamics import * sage: c = CylinderDiagram('(0,1)-(2,3,5) (2,3,4)-(1) (5)-(0,4)') sage: c.symmetries() (False, True, False) sage: c.horizontal_symmetry().is_isomorphic(c) False sage: c.vertical_symmetry().is_isomorphic(c) True sage: c.inverse().is_isomorphic(c) False sage: CylinderDiagram('(0,2,1,4,3)-(0,4,2,1,3)').symmetries() (True, True, True) sage: CylinderDiagram('(0,3)-(0,4,5) (1,4,2)-(1,3) (5)-(2)').symmetries() (False, False, True) sage: CylinderDiagram('(0,2,3,1)-(0,2,1,4) (4)-(3)').symmetries() (True, False, False) sage: CylinderDiagram('(0,1,2)-(0,3,1,4) (3,4)-(2)').symmetries() (False, True, False) sage: CylinderDiagram('(0,1)-(0,3,4) (2,3)-(1) (4)-(2)').symmetries() (False, False, False)
- tikz_string(cyl_height=0.8, cyl_sep=1.3, cyl_pos=None, sep_labels=None)[source]¶
INPUT:
cyl_height
- (optional) height of cylinderscyl_sep
- (optional) vertical space between cylinderscyl_pos
- (optional) position of left corners of cylinders. If provided,cyl_sep
is ignored.sep_labels
- labels of the separatrices
- to_directed_graph()[source]¶
Return a labeled directed graph that encodes the cylinder diagram.
EXAMPLES:
sage: from surface_dynamics import * sage: c = CylinderDiagram('(0,1,5)-(2,5) (2)-(0,1,3) (3,4)-(4)'); c (0,1,5)-(2,5) (2)-(0,1,3) (3,4)-(4) sage: G = c.to_directed_graph(); G Looped multi-digraph on 6 vertices sage: G.edges(sort=True) [(0, 1, 'b'), (0, 1, 't'), (0, 2, 'c'), (0, 5, 'c'), (1, 2, 'c'), (1, 3, 't'), (1, 5, 'b'), (1, 5, 'c'), (2, 0, 'c'), (2, 1, 'c'), (2, 2, 'b'), (2, 3, 'c'), (2, 5, 't'), (3, 0, 't'), (3, 4, 'b'), (3, 4, 'c'), (4, 3, 'b'), (4, 4, 'c'), (4, 4, 't'), (5, 0, 'b'), (5, 2, 'c'), (5, 2, 't'), (5, 5, 'c')]
- to_ribbon_graph()[source]¶
Return a ribbon graph
A ribbon graph is a graph embedded in an oriented surface such that its complement is a union of topological discs. To a cylinder diagram we associate the graph which consists of separatrices together with a choice of one vertical edge in each cylinder.
The edges of the ribbon graph are labeled by
(i,nseps+i)
for separatrices and by(2(nseps+j),2(nseps+j)+1)
for vertical in cylinders.EXAMPLES:
sage: from surface_dynamics import * sage: C = CylinderDiagram([((0,1),(0,2)),((2,),(1,))]) sage: C.stratum() H_2(2) sage: R = C.to_ribbon_graph(); R Ribbon graph with 1 vertex, 5 edges and 2 faces sage: l,m = R.cycle_basis(intersection=True) sage: m.rank() == 2 * C.genus() True
- top_to_cyl(j)[source]¶
Return the cylinder below the
j
-th separatrix.EXAMPLES:
sage: from surface_dynamics import * sage: c = CylinderDiagram('(0,2,4)-(1,3,5) (1,5)-(0) (3)-(2,4)') sage: c.top_to_cyl(0) ((1, 5), (0,)) sage: c.top_to_cyl(2) ((3,), (2, 4))
- vertical_symmetry()[source]¶
Return the cylinder diagram obtained by reflecting the cylinder configuration along the vertical axis.
EXAMPLES:
sage: from surface_dynamics import * sage: c = CylinderDiagram('(0,3,4)-(0,3,5) (1,2,5)-(1,2,4)') sage: c.vertical_symmetry() (0,4,3)-(0,5,3) (1,5,2)-(1,4,2) sage: c.separatrix_diagram().vertical_symmetry() == c.vertical_symmetry().separatrix_diagram() True sage: A = Stratum([2,2], k=1) sage: all(c.vertical_symmetry().stratum() == A for c in A.cylinder_diagrams()) True
- volume_contribution()[source]¶
Return the volume contribution of this cylinder diagram as a generalized multiple zeta values.
EXAMPLES:
sage: from surface_dynamics import * sage: c0, c1 = Stratum([2], k=1).cylinder_diagrams() sage: v0 = c0.volume_contribution() # optional: latte_int sage: v0 # optional: latte_int (1/3)/((w)^4) sage: v0.integral_sum_as_mzv() # optional: latte_int 1/3*ζ(4) sage: v1 = c1.volume_contribution() # optional: latte_int sage: v1 # optional: latte_int (2/3)/((w1)*(w0 + w1)^3) + (1/3)/((w1)^2*(w0 + w1)^2) sage: v1.integral_sum_as_mzv() # optional: latte_int 2/3*ζ(1,3) + 1/3*ζ(2,2) sage: for c in Stratum([1,1], k=1).cylinder_diagrams(): # optional: latte_int ....: print(c, c.volume_contribution().integral_sum_as_mzv()) (0,3,1,2)-(0,3,1,2) 1/6*ζ(5) (0)-(1) (1,2,3)-(0,2,3) 1/3*ζ(2,3) + 1/3*ζ(3,2) (0,3)-(1,3) (1,2)-(0,2) ζ(1,4) + 1/3*ζ(2,3) (0,1)-(2,3) (2)-(1) (3)-(0) 1/3*ζ(1,3) + 1/3*ζ(2,2) - 1/3*ζ(2,3) - 1/3*ζ(3,2) + 1/3*ζ(4) - 1/3*ζ(5) sage: sum(c.volume_contribution() for c in Stratum([2,1,1], k=1).cylinder_diagrams(1)).integral_sum_as_mzv() # optional: latte_int 7/180*ζ(8)
Detailed contribution of 2 cylinder diagrams:
sage: cyls = Stratum([2,1,1], k=1).cylinder_diagrams(2) sage: sum(cyls[k].volume_contribution() for k in [2,7,8,21,22]).integral_sum_as_mzv() # optional: latte_int 13/630*ζ(5,3) + 13/252*ζ(6,2) sage: sum(cyls[k].volume_contribution() for k in [0,11,19,20]).integral_sum_as_mzv() # optional: latte_int 1/21*ζ(4,4) + 4/63*ζ(5,3) sage: sum(cyls[k].volume_contribution() for k in [3,10]).integral_sum_as_mzv() # optional: latte_int 2/35*ζ(3,5) + 3/70*ζ(4,4) sage: sum(cyls[k].volume_contribution() for k in [13,23,24]).integral_sum_as_mzv() # optional: latte_int 2/21*ζ(2,6) + 4/105*ζ(3,5) sage: sum(cyls[k].volume_contribution() for k in [9]).integral_sum_as_mzv() # optional: latte_int 1/21*ζ(1,7) + 1/126*ζ(2,6) sage: sum(cyls[k].volume_contribution() for k in [5]).integral_sum_as_mzv() # optional: latte_int 2/105*ζ(3,5) + 1/70*ζ(4,4) sage: sum(cyls[k].volume_contribution() for k in [6,26]).integral_sum_as_mzv() # optional: latte_int 1/7*ζ(1,7) + 1/14*ζ(2,6) + 1/21*ζ(3,5) + 1/28*ζ(4,4) + 2/105*ζ(5,3) sage: sum(cyls[k].volume_contribution() for k in [1,14]).integral_sum_as_mzv() # optional: latte_int 2/7*ζ(1,7) + 1/7*ζ(2,6) + 2/21*ζ(3,5) + 3/70*ζ(4,4) sage: sum(cyls[k].volume_contribution() for k in [17,27]).integral_sum_as_mzv() # optional: latte_int 2/7*ζ(1,7) + 1/7*ζ(2,6) + 4/105*ζ(3,5) sage: sum(cyls[k].volume_contribution() for k in [4,25]).integral_sum_as_mzv() # optional: latte_int 1/7*ζ(1,7) + 1/42*ζ(2,6) sage: sum(cyls[k].volume_contribution() for k in [12,28]).integral_sum_as_mzv() # optional: latte_int 4/21*ζ(1,7) + 2/21*ζ(2,6) + 8/315*ζ(3,5) sage: sum(cyls[k].volume_contribution() for k in [15,16]).integral_sum_as_mzv() # optional: latte_int 4/21*ζ(1,7) + 2/21*ζ(2,6) + 8/315*ζ(3,5) sage: sum(cyls[k].volume_contribution() for k in [18]).integral_sum_as_mzv() # optional: latte_int 1/36*ζ(7) - 1/36*ζ(8)
- weighted_adjacency_matrix(canonicalize=True)[source]¶
Return the weighted adjacency matrix of this cylinder diagram.
There is a vertex per cylinder and the weight from cyl1 to cyl2 is the number of saddle connections which are on top of cyl1 and in the bottom of cyl2.
EXAMPLES:
sage: from surface_dynamics import Stratum sage: for cd in Stratum([1,1], k=1).cylinder_diagrams(): ....: print(cd.weighted_adjacency_matrix()) [4] [0 1] [1 2] [1 1] [1 1] [0 0 1] [0 0 1] [1 1 0]
- widths_and_heights_iterator(n, height_one=False)[source]¶
Iterate over the possible integer widths and heights of the cylinders for which the corresponding translation surface has area
n
.At each iteration, the output is a pair of
(lengths, heights)
. You can then usecylcoord_to_origami()
to build the corresponding origami.INPUT:
height_one
– (boolean defaultFalse
) whether to return only coordinates with cylinder height one.
EXAMPLES:
sage: from surface_dynamics import * sage: cyl = CylinderDiagram([((0,1),(0,2)),((2,),(1,))]) sage: cyl (0,1)-(0,2) (2)-(1) sage: it = cyl.widths_and_heights_iterator(10) sage: l,h = next(it) sage: print(l) (2, 1, 1) sage: print(h) [3, 1] sage: cyl.cylcoord_to_origami(l,h) (1,2,3)(4,5,6)(7,8,9)(10) (1,4,7)(2,5,8)(3,6,9,10) sage: it = cyl.widths_and_heights_iterator(10, height_one=True) sage: l,h = next(it) sage: print(l) (8, 1, 1) sage: print(h) [1, 1] sage: cyl.cylcoord_to_origami(l,h) (1,2,3,4,5,6,7,8,9)(10) (1)(2)(3)(4)(5)(6)(7)(8)(9,10)
- widths_generating_series(var='w')[source]¶
Generating series of the number of lengths solutions given the widths.
WARNING: when a triangulation is involved, the generating series ignore some lower dimensional polytopes!!
INPUT:
var
- (optional) an optional name for the variable
EXAMPLES:
sage: from surface_dynamics import CylinderDiagram sage: c = CylinderDiagram('(0,1)-(0,2) (2)-(1)') sage: c.widths_generating_series() # optional: latte_int (1)/((1 - w0)*(1 - w0*w1)) sage: c = CylinderDiagram('(0)-(2) (1,2,3)-(4,5) (4)-(3) (5)-(0,1)') sage: c.widths_generating_series() # optional: latte_int (1)/((1 - w1*w3)*(1 - w1*w2)*(1 - w0*w1*w3)) sage: c = CylinderDiagram('(0,1,3)-(0,2,5) (2,4)-(1,3) (5)-(4)') sage: c.widths_generating_series() # optional: latte_int (1)/((1 - w0)*(1 - w0*w1)*(1 - w0*w1*w2)^2) + (1)/((1 - w0)*(1 - w0*w1)^2*(1 - w0*w1*w2))
- class surface_dynamics.flat_surfaces.separatrix_diagram.QuadraticCylinderDiagram(arg1, arg2=None)[source]¶
Bases:
SageObject
Cylinder diagram for quadratic differentials
Cylinder diagrams are encoded as a Ribbon graph together with a pairing of faces (in particular the number of faces must be even).
EXAMPLES:
sage: from surface_dynamics import *
If you start with strings, the cylinders are preserved but the names of saddle connections are changed:
sage: QuadraticCylinderDiagram('(4,4,5)-(6,6,1) (2,3,2,0)-(1,0,5,3)') (0,0,1)-(2,2,3) (4,5,4,6)-(3,6,1,5)
- cylcoord_to_pillowcase_cover(lengths, heights, twists=None, verbose=False)[source]¶
Convert coordinates of the cylinders into a pillowcase cover.
The pillow is considered as made of two 1 x 1 squares in order to avoid denominators in the input
lengths
,heights
andtwists
.INPUT:
lengths
- positive integers - lengths of the separatricesheights
- positive integers - heights of the cylinderstwists
- (optional) non-negative integers - twists. The twist is measured as the difference in horizontal coordinates between the smallest element in bottom and the smallest in element in top.
OUTPUT: a pillowcase cover
EXAMPLES:
sage: from surface_dynamics import *
Some pillows in Q(-1^4):
sage: q = QuadraticCylinderDiagram('(0,0)-(1,1)') sage: q.cylcoord_to_pillowcase_cover([1,1],[1],[0]) g0 = (1) g1 = (1) g2 = (1) g3 = (1) sage: q.cylcoord_to_pillowcase_cover([1,1],[1],[1]) g0 = (1) g1 = (1) g2 = (1) g3 = (1) sage: q.cylcoord_to_pillowcase_cover([2,2],[1],[0]) g0 = (1)(2) g1 = (1,2) g2 = (1,2) g3 = (1)(2) sage: q.cylcoord_to_pillowcase_cover([2,2],[1],[1]) g0 = (1)(2) g1 = (1,2) g2 = (1)(2) g3 = (1,2) sage: q.cylcoord_to_pillowcase_cover([2,2],[1],[2]) == q.cylcoord_to_pillowcase_cover([2,2],[1],[2]) True sage: q.cylcoord_to_pillowcase_cover([3,3],[1],[0]) g0 = (1)(2,3) g1 = (1,3)(2) g2 = (1,3)(2) g3 = (1)(2,3) sage: q.cylcoord_to_pillowcase_cover([3,3],[1],[1]) g0 = (1)(2,3) g1 = (1,3)(2) g2 = (1)(2,3) g3 = (1,2)(3) sage: q.cylcoord_to_pillowcase_cover([1,1],[2],[0]) g0 = (1)(2) g1 = (1)(2) g2 = (1,2) g3 = (1,2) sage: q.cylcoord_to_pillowcase_cover([1,1],[2],[1]) == q.cylcoord_to_pillowcase_cover([1,1],[2],[0]) True sage: q.cylcoord_to_pillowcase_cover([2,2],[2],[1]) g0 = (1)(2)(3,4) g1 = (1,2)(3)(4) g2 = (1,3)(2,4) g3 = (1,4)(2,3)
Two one cylinder examples in Q(2^2):
sage: q1 = QuadraticCylinderDiagram('(0,1,0,1)-(2,3,2,3)') sage: q2 = QuadraticCylinderDiagram('(0,1,0,2)-(3,1,3,2)') sage: q1.cylcoord_to_pillowcase_cover([2,2,2,2],[1],[0]) g0 = (1,2,3,4) g1 = (1,3)(2,4) g2 = (1,3)(2,4) g3 = (1,4,3,2) sage: q1.cylcoord_to_pillowcase_cover([2,2,2,2],[1],[1]) g0 = (1,2,3,4) g1 = (1,3)(2,4) g2 = (1,4,3,2) g3 = (1,3)(2,4) sage: p = q1.cylcoord_to_pillowcase_cover([2,6,4,4],[3],[1]) sage: p.stratum() Q_2(2^2) sage: p.nb_pillows() 24
One two cylinders example in Q(2^2):
sage: q = QuadraticCylinderDiagram('(0,1)-(2,3) (0,3)-(1,2)') sage: p = q.cylcoord_to_pillowcase_cover([1,1,1,1], [2,2], [0,1]) sage: p g0 = (1,4,2,3) g1 = (1,3,2,4) g2 = (1,2)(3,4) g3 = (1,2)(3,4) sage: p.stratum() Q_2(2^2) sage: q.cylcoord_to_pillowcase_cover([1,3,1,3], [2,2], [0,1]).stratum() Q_2(2^2)
- cylinders(dart=False)[source]¶
Cylinders of self
Return a list of pairs
(bot, top)
where, by convention the bottom corresponds to the face with smaller index in the list of faces.EXAMPLES:
sage: from surface_dynamics import * sage: from surface_dynamics.flat_surfaces.separatrix_diagram import QuadraticCylinderDiagram sage: rg = RibbonGraph(edges='(0,1)(2,3)(4,5)(6,7)', faces='(0,1)(2,4,5)(3)(6,7)', connected=False) sage: q = QuadraticCylinderDiagram(rg, '(0,1)(2,3)') sage: q.cylinders() [((0, 0), (1, 2, 2)), ((1,), (3, 3))] sage: q.cylinders(True) [((0, 1), (2, 4, 5)), ((3,), (6, 7))]
- lengths_cone()[source]¶
Return the polytope of admissible lengths.
EXAMPLES:
sage: from surface_dynamics import * sage: rg = RibbonGraph(edges='(0,1)(2,3)(4,5)(6,7)', faces='(0,1)(2,4,5)(3)(6,7)', connected=False) sage: q = QuadraticCylinderDiagram(rg, '(0,1)(2,3)') sage: q.stratum() Q_0(1, -1^5) sage: q.lengths_cone().rays_list() [[1, 0, 1, 0], [1, 2, 0, 1]] sage: rg = RibbonGraph(edges='(0,1)(2,3)(4,5)(6,7)', faces='(0,2,3)(1)(4,6,7)(5)', connected=False) sage: q = QuadraticCylinderDiagram(rg, '(0,2)(1,3)') sage: rg = RibbonGraph(edges='(1,2)(3,4)(5,6)(7,8)(9,10)(11,12)', faces='(1,2)(3,4,6)(5)(8)(7,9,10)(11,12)', connected=False) sage: q = QuadraticCylinderDiagram(rg, '(0,1)(2,4)(3,5)') sage: q.stratum() Q_0(1^2, -1^6) sage: L = q.lengths_cone() sage: L A 3-dimensional polyhedron in QQ^6 defined as the convex hull of 1 vertex and 3 rays
- ncyls()¶
Number of cylinders.
- nseps()¶
Number of edges
- stratum()[source]¶
Return the stratum of quadratic differentials associated to this cylinder diagram
EXAMPLES:
sage: from surface_dynamics import * sage: QuadraticCylinderDiagram('(0,0)-(1,1)').stratum() Q_0(-1^4) sage: QuadraticCylinderDiagram('(0,0)-(1,1,2,2,3,3)').stratum() Q_0(1, -1^5) sage: QuadraticCylinderDiagram('(0,2,3,2)-(1,0,1,3)').stratum() Q_2(2^2) sage: QuadraticCylinderDiagram('(0,1)-(2,3) (0)-(4,4) (1)-(5,5) (2)-(6,6) (3)-(7,7)').stratum() Q_0(2^2, -1^8)
- widths_generating_series(var='w')[source]¶
Generating series of the number of saddle connection lengths of this quadratic differential cylinder diagram.
Warning
When a triangulation is involved, the generating series ignore some lower dimensional polytopes that are counted twice!
EXAMPLES:
sage: from surface_dynamics import * sage: q = QuadraticCylinderDiagram('(0,1,2,3,3)-(0,4,4,2,1)') sage: q.widths_generating_series() # optional -- latte_int (1)/((1 - w)^3*(1 - w^2)) sage: q = QuadraticCylinderDiagram('(0,0,1,1,2,2)-(3,3,4,4)') sage: q.widths_generating_series() # optional -- latte_int (3)/((1 - w^2)^4) sage: q = QuadraticCylinderDiagram('(0,0,1)-(2,2,3) (1,4)-(3,4)') sage: q.widths_generating_series() # optional -- latte_int (1)/((1 - w1)*(1 - w0*w1)*(1 - w0^2)) sage: q = QuadraticCylinderDiagram('(0,0,1,2,3)-(1,4,4,5,6) (2,5,7,7,8)-(3,6,8,9,9)') sage: F = q.widths_generating_series() # optional -- latte_int sage: q = QuadraticCylinderDiagram('(0,0,1,2,3)-(1,4,4,5,6) (2,5,7,8,8)-(3,6,7,9,9)') sage: F = q.widths_generating_series() # optional -- latte_int
- class surface_dynamics.flat_surfaces.separatrix_diagram.SeparatrixDiagram(data, top=None, check=True, copy=True)[source]¶
Bases:
SageObject
Separatrix diagram of oriented foliation.
A separatrix diagram is a 2-tuple of permutations
(bot,top)
such thatbot
andtop
share the same number of cycles.bot (resp. top) has to be thought a bottom (resp. top) of a potential face as in the following:
-- bot --> ------------------- <-- top --
The order for bot and top is chosen in such a way that it corresponds to the orientation of a face.
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,2)(1,3,4)','(0,4)(2,1,3)') sage: print(s) (0,2)(1,3,4)-(0,4)(1,3,2) sage: print(s.stratum()) H_3(4)
- automorphism_group(implementation='graph')[source]¶
Return the automorphism group of self.
That is the centralizer of the permutations top and bottom.
INPUT:
implementation
- either graph or gap
EXAMPLES:
sage: from surface_dynamics import * sage: S = SeparatrixDiagram('(0,3,1,4,2)','(0,1,2,3,4)') sage: G1 = S.automorphism_group(implementation='graph'); G1 Permutation Group with generators [(0,1,2,3,4)] sage: G2 = S.automorphism_group(implementation='gap'); G2 Subgroup ... sage: G1.is_isomorphic(G2) True
- bot()[source]¶
The bot permutation as a list from 0 to nseps-1
Warning: the output list should not be modified
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0)(1,2)','(0,1)(2)') sage: s.bot() [0, 2, 1]
- bot_cycle_string()[source]¶
Return the cycles of the top permutation as a string.
EXAMPLES:
sage: from surface_dynamics import * sage: S = SeparatrixDiagram('(0,2)(3,4)','(0)(1,2,3)') sage: S.bot_cycle_string() '(0,2)(1)(3,4)'
- bot_cycle_tuples()[source]¶
Return the cycles of the bottom permutation as a list of tuples.
EXAMPLES:
sage: from surface_dynamics import * sage: S = SeparatrixDiagram('(0,2)(3,4)','(0)(1,2,3)') sage: S.bot_cycle_tuples() [(0, 2), (1,), (3, 4)]
- bot_orbit(i)[source]¶
Return the orbit of i under the bot permutation
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1)(2,5)(3,4,6)','(0,1,5)(2,3,6)(4)') sage: s.bot_orbit(0) (0, 1) sage: s.bot_orbit(4) (3, 4, 6)
- bot_perm()[source]¶
Return the bot as a permutation (element of a group)
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0)(1,2)','(0,1)(2)') sage: s.bot_perm() (2,3)
- canonical_label(inplace=False)[source]¶
Relabel self according to some canonical labels.
The result is cached.
INPUT:
inplace
- boolean (default:True
) - if True modify self if not return a new separatrix diagram.
EXAMPLES:
sage: from surface_dynamics import * sage: bot = '(0,1,3,6,7,5)(2,4)(8)(9)' sage: top = '(0)(1,2)(3,4,5)(6,7,8,9)' sage: s = SeparatrixDiagram(bot,top) sage: s.canonical_label() (0)(1)(2,3,4,5,6,7)(8,9)-(0,1,2,3)(4,7,9)(5)(6,8)
- cylinder_diagram_iterator(connected=True, up_to_symmetry=True)[source]¶
Construct all cylinder diagrams from given separatrix diagram (i.e. a pair of permutations).
INPUT:
connected
- boolean (default: True) - if true, returns only connected cylinder diagrams.up_to_symmetry
- boolean (default: True) take care of the horizontal and vertical symmetries.
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1)(2,3)(4,5)','(1,2)(3,4)(5,0)') sage: for c in s.cylinder_diagram_iterator(): print(c) # random (0,5)-(0,4) (1,4)-(1,3) (2,3)-(2,5) (0,3)-(0,5) (1,2)-(1,4) (4,5)-(2,3) (0,5)-(3,4) (1,4)-(0,2) (2,3)-(1,5) sage: sum(1 for _ in s.cylinder_diagram_iterator(up_to_symmetry=False)) 6 sage: sum(1 for _ in s.cylinder_diagram_iterator(up_to_symmetry=True)) 3
Here is an example with some symmetry:
sage: s = SeparatrixDiagram('(0)(1)(2,3)(4,5,6)-(0,1)(2,4)(3,5)(6)') sage: s.vertical_symmetry().canonical_label() == s True sage: C1 = [CylinderDiagram('(0,1)-(4) (2,4,3)-(5,6) (5)-(0,2) (6)-(1,3)'), ....: CylinderDiagram('(0,3,1)-(0,6) (2,6)-(4,5) (4)-(1) (5)-(2,3)'), ....: CylinderDiagram('(0,1)-(0,4) (2,3,4)-(5,6) (5)-(2) (6)-(1,3)')] sage: C2 = s.cylinder_diagrams() sage: assert len(C1) == len(C2) sage: for (c1, c2) in zip(C1, C2): assert c1.is_isomorphic(c2)
- cylinder_diagrams(connected=True, up_to_symmetry=True)[source]¶
Return the list of cylinder diagrams associated to this separatrix diagram.
We warn that the cylinder diagram may be renumeroted in the output list (in order to prevent repetitions). If you care about numerotation the option
up_to_symmetry
should be set to False.INPUT:
connected
- boolean (default: True)up_to_symmetry
- boolean (default: True)
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0)(1)(2)','(0)(1)(2)') sage: for c in s.cylinder_diagrams(connected=True): print(c) (0)-(2) (1)-(0) (2)-(1) sage: for c in s.cylinder_diagrams(connected=False): print(c) (0)-(0) (1)-(1) (2)-(2) (0)-(1) (1)-(0) (2)-(2) (0)-(2) (1)-(0) (2)-(1) sage: s = SeparatrixDiagram('(0,1)(2)','(0)(1,2)') sage: C1 = [CylinderDiagram('(0,1)-(0,2) (2)-(1)')] sage: C2 = s.cylinder_diagrams() sage: assert len(C1) == len(C2) sage: for (c1, c2) in zip(C1, C2): assert c1.is_isomorphic(c2)
In the example below, there is no isomorphism problem for the cylinder diagram generation as the separatrix diagram admit no automorphism:
sage: s = SeparatrixDiagram('(0,3)(1,4,5)(2)','(0)(1,2)(3,4,5)') sage: s.automorphism_group() Permutation Group with generators [()] sage: C1 = [CylinderDiagram('(0,3,1)-(5) (2,5)-(3,4) (4)-(0,2,1)'), ....: CylinderDiagram('(0,1,2)-(0,1,5) (3,5)-(2,4) (4)-(3)'), ....: CylinderDiagram('(0,2,3)-(2,5) (1,4)-(0,1,3) (5)-(4)')] sage: C2 = s.cylinder_diagrams() sage: C3 = s.cylinder_diagrams(up_to_symmetry=False) sage: assert len(C1) == len(C2) == len(C3) sage: for (c1, c2, c3) in zip(C1, C2, C3): assert c1.is_isomorphic(c2) and c1.is_isomorphic(c3)
- degree()[source]¶
Return the degree (number of separatrices) of this separatrix diagram.
EXAMPLES:
sage: from surface_dynamics import * sage: S = SeparatrixDiagram('(0,1)(2,3)','(1,3,2)(0)') sage: S.degree() 4
- euler_characteristic()[source]¶
Return the Euler characteristic
EXAMPLES:
sage: from surface_dynamics import * sage: SeparatrixDiagram('(0)','(0)').euler_characteristic() 0 sage: CylinderDiagram([((0,),(0,))]).euler_characteristic() 0 sage: CylinderDiagram([((0,1),(0,2)), ((2,),(1,))]).euler_characteristic() -2
- genus()[source]¶
Return the genus
EXAMPLES:
sage: from surface_dynamics import * sage: CylinderDiagram([((0,),(0,))]).genus() 1 sage: CylinderDiagram([((0,1),(0,1))]).genus() 1 sage: CylinderDiagram([((0,1,2),(0,1,2))]).genus() 2 sage: CylinderDiagram([((0,1,2,3),(0,1,2,3))]).genus() 2 sage: CylinderDiagram([((0,1,2,3,4),(0,1,2,3,4))]).genus() 3
- homological_dimension_of_cylinders()[source]¶
Returns the dimension in the first homology group of the span of waist curves of horizontal cylinders.
EXAMPLES:
sage: from surface_dynamics import *
Homological dimension in the stratum H(2):
sage: c = CylinderDiagram('(0,1,2)-(0,1,2)') sage: c.stratum() H_2(2) sage: c.homological_dimension_of_cylinders() 1 sage: c = CylinderDiagram('(0,1)-(1,2) (2)-(0)') sage: c.stratum() H_2(2) sage: c.homological_dimension_of_cylinders() 2
Homological dimensions for cylinder diagrams in H(1,1):
sage: c = CylinderDiagram('(0,1,2,3)-(0,1,2,3)') sage: c.stratum() H_2(1^2) sage: c.homological_dimension_of_cylinders() 1 sage: c = CylinderDiagram('(0,1)-(0,2) (2,3)-(1,3)') sage: c.stratum() H_2(1^2) sage: c.homological_dimension_of_cylinders() 2 sage: c = CylinderDiagram('(0,1,2)-(1,2,3) (3)-(0)') sage: c.stratum() H_2(1^2) sage: c.homological_dimension_of_cylinders() 2 sage: c = CylinderDiagram('(0,1)-(2,3) (2)-(0) (3)-(1)') sage: c.stratum() H_2(1^2) sage: c.homological_dimension_of_cylinders() 2
- homologous_cylinders()[source]¶
Return the list of homologous cylinders.
OUTPUT: a list of lists. Each sublist is an equivalence class of > 1 homologous cylinders.
EXAMPLES:
sage: from surface_dynamics import CylinderDiagram sage: c = CylinderDiagram('(0,7,1,2)-(3,6,4,5) (3,6,4,5)-(0,7,1,2)') sage: c.homologous_cylinders() [[0, 1]] sage: c = CylinderDiagram('(0,1)-(2,3) (2)-(0) (3)-(1)') sage: c.homologous_cylinders() [] sage: c = CylinderDiagram('(0,2,1)-(9) (3,6,4,5)-(7,10,8) (7,9,8)-(3,6,4,5) (10)-(0,2,1)') sage: c.homologous_cylinders() [[0, 3], [1, 2]]
- horizontal_symmetry()[source]¶
Return the horizontal symmetric of this separatrix diagram.
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1,2,3)(4,5)','(1,2,3)(4,5,0)') sage: sh = s.horizontal_symmetry() sage: print(sh) (0,5,4)(1,3,2)-(0,3,2,1)(4,5) sage: c = s.cylinder_diagrams(up_to_symmetry=False)[0].horizontal_symmetry() sage: ch = sh.cylinder_diagrams(up_to_symmetry=False)[0] sage: c.is_isomorphic(ch) True
- incoming_edges_perm()[source]¶
Permutation associated to turning around vertices in trigonometric order.
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1)','(2,3)') sage: s.incoming_edges_perm() [1, 0, 3, 2] sage: s = SeparatrixDiagram('(0,5,2)(1,3,4)(6,7,8)','(0,3,7,8)(1,5)(2,4,6)') sage: s.incoming_edges_perm() [7, 0, 8, 2, 5, 4, 3, 1, 6] sage: c = CylinderDiagram('(0,5,2,1,4,3)-(0,4,2,1,5,3)') sage: c.incoming_edges_perm() [4, 5, 1, 0, 3, 2]
- inverse()[source]¶
Return the inverse of this separatrix diagram, that is the one we obtain after application of -Id.
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1,2)(3,4,5,6,7,8)-(0,1,3,5,7)(2,4,6,8)') sage: s.inverse() (0,1,3,5,7)(2,4,6,8)-(0,1,2)(3,4,5,6,7,8) sage: s.horizontal_symmetry().vertical_symmetry() == s.inverse() True sage: s.vertical_symmetry().horizontal_symmetry() == s.inverse() True
- is_in_normal_form()[source]¶
Test normal form
Return True if self is in normal form and False otherwise.
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1,2)(3,4,5)(6,7,8)','(0,3,7,8)(1,5)(2,4,6)') sage: s.is_in_normal_form() False sage: s.canonical_label().is_in_normal_form() True
- is_isomorphic(other, return_map=False)[source]¶
Test whether this separatrix diagram is isomorphic to
other
.EXAMPLES:
sage: from surface_dynamics import * sage: bot = [1,2,0,3] sage: top = [1,0,3,2] sage: s = SeparatrixDiagram(bot,top); s (0,1,2)(3)-(0,1)(2,3) sage: m = [3,0,1,2] sage: bot2 = [0]*4 sage: top2 = [0]*4 sage: for i in range(4): ....: bot2[m[i]] = m[bot[i]] ....: top2[m[i]] = m[top[i]] sage: ss = SeparatrixDiagram(bot2,top2) sage: s.is_isomorphic(ss) True sage: m = [1,2,0,3] sage: for i in range(4): ....: bot2[m[i]] = m[bot[i]] ....: top2[m[i]] = m[top[i]] sage: ss = SeparatrixDiagram(bot2,top2) sage: s.is_isomorphic(ss) True
- ncyls()[source]¶
Return the number of cylinders of this separatrix diagram.
EXAMPLES:
sage: from surface_dynamics import * sage: S = SeparatrixDiagram('(0,1)(2,3)','(1,3,2)(0)') sage: S.ncyls() 2
- nseps()¶
Return the degree (number of separatrices) of this separatrix diagram.
EXAMPLES:
sage: from surface_dynamics import * sage: S = SeparatrixDiagram('(0,1)(2,3)','(1,3,2)(0)') sage: S.degree() 4
- num_edges()¶
Return the degree (number of separatrices) of this separatrix diagram.
EXAMPLES:
sage: from surface_dynamics import * sage: S = SeparatrixDiagram('(0,1)(2,3)','(1,3,2)(0)') sage: S.degree() 4
- outgoing_edges_perm()[source]¶
Permutation associated to turning around vertices in trigonometric order.
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1)','(2,3)') sage: s.outgoing_edges_perm() [1, 0, 3, 2] sage: s = SeparatrixDiagram('(0,5,2)(1,3,4)(6,7,8)','(0,3,7,8)(1,5)(2,4,6)') sage: s.outgoing_edges_perm() [6, 2, 1, 5, 0, 8, 7, 4, 3] sage: c = CylinderDiagram('(0,5,2,1,4,3)-(0,4,2,1,5,3)') sage: c.outgoing_edges_perm() [5, 4, 1, 0, 2, 3]
- profile()[source]¶
Return the angles around each vertex
EXAMPLES:
sage: from surface_dynamics import * sage: a = Stratum([1,1,0], k=1) sage: s = a.separatrix_diagrams()[0] sage: s.profile() [2, 2, 1]
- relabel(perm, inplace=False)[source]¶
Relabel self according to the permutation
perm
.EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0)(2,3,4)','(0,3,2)(1)') sage: s (0)(1)(2,3,4)-(0,3,2)(1)(4) sage: s.relabel(perm=[1,0,2,3,4]) (0)(1)(2,3,4)-(0)(1,3,2)(4) sage: s.relabel(perm=[1,2,0,3,4]) (0,3,4)(1)(2)-(0,1,3)(2)(4)
- saddle_connections_graph(mutable=False)[source]¶
Return the fat graph (or ribbon graph) made by the saddle connections.
The return graph is a
FatGraph
. The saddle connection labelled i on this diagram gets labels 2i and 2i+1 in the graph (there one label per half-edge in the fat graph). The even labels correspond to half-edges in the bottom of cylinders while the odd ones correspond to the top.EXAMPLES:
sage: from surface_dynamics import Stratum sage: H11 = Stratum([1,1], k=1).unique_component() sage: for cd in H11.cylinder_diagrams(): ....: fg = cd.saddle_connections_graph() ....: print(cd.ncyls(), [comp.genus() for comp in fg.connected_components()]) 1 [1] 2 [0] 2 [0] 3 [0, 0]
- stratum()[source]¶
Return the Abelian stratum this separatrix diagram belongs to.
EXAMPLES:
sage: from surface_dynamics import * sage: SeparatrixDiagram('(0)(1)(2)','(0)(1)(2)').stratum() H_1(0^3) sage: SeparatrixDiagram('(0,1)(2)','(0,2)(1)').stratum() H_2(2)
- symmetries()[source]¶
Return a triple of boolean
(horiz_sym, vert_sym, inverse_sym)
which correspond to the symmetry ofself
.EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1,2)(3,4,5)-(0,1)(2,3,4,5)') sage: s.symmetries() (False, True, False) sage: s.horizontal_symmetry().is_isomorphic(s) False sage: s.vertical_symmetry().is_isomorphic(s) True sage: s.inverse().is_isomorphic(s) False sage: s = SeparatrixDiagram('(0,1,3,5)(2,4)-(0,4,1,5)(2,3)') sage: s.symmetries() (True, False, False) sage: s.horizontal_symmetry().is_isomorphic(s) True sage: s.vertical_symmetry().is_isomorphic(s) False sage: s.inverse().is_isomorphic(s) False sage: s = SeparatrixDiagram('(0,1,3,5)(2,4)-(0,3,2,1)(5,4)') sage: s.symmetries() (False, False, True) sage: s.horizontal_symmetry().is_isomorphic(s) False sage: s.vertical_symmetry().is_isomorphic(s) False sage: s.inverse().is_isomorphic(s) True sage: s = SeparatrixDiagram('(0)(1,2,3,4,5)-(0,1,2,5,3)(4)') sage: s.symmetries() (False, False, False) sage: s.horizontal_symmetry().is_isomorphic(s) False sage: s.vertical_symmetry().is_isomorphic(s) False sage: s.inverse().is_isomorphic(s) False
- to_cylinder_diagram(pairing)[source]¶
Return a cylinder diagram with the given pairing
The pairing should be a list of 2-tuples of integer.
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1,3)(2,4)','(0,2)(1,4,3)'); s (0,1,3)(2,4)-(0,2)(1,4,3) sage: s.to_cylinder_diagram([(0,0),(1,1)]) (0,1,3)-(0,2) (2,4)-(1,4,3) sage: s.to_cylinder_diagram([(1,1),(0,0)]) (0,1,3)-(0,2) (2,4)-(1,4,3) sage: s.to_cylinder_diagram([(0,1),(1,0)]) (0,1,3)-(1,4,3) (2,4)-(0,2) sage: s.to_cylinder_diagram([(1,0),(0,1)]) (0,1,3)-(1,4,3) (2,4)-(0,2)
- to_directed_graph()[source]¶
Return a graph that encodes this separatrix diagram.
The vertices correspond to separatrix and the edges are of two types
‘b’ neighbor corresponds to the right neighbors on the bottom permutation
‘t’ edges correspond to the neighbor of the top permutation
EXAMPLES:
sage: from surface_dynamics import * sage: S = SeparatrixDiagram('(0,1)(2,3,4)','(0,3,2)(1,4)') sage: G = S.to_directed_graph(); G Looped multi-digraph on 5 vertices sage: G.vertices(sort=True) [0, 1, 2, 3, 4] sage: G.edges(sort=True) [(0, 1, 'b'), (0, 3, 't'), (1, 0, 'b'), (1, 4, 't'), (2, 0, 't'), (2, 3, 'b'), (3, 2, 't'), (3, 4, 'b'), (4, 1, 't'), (4, 2, 'b')]
- top()[source]¶
Return the top permutation of self as a list.
Warning: the output should not be modified
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1,3)(2,4)','(0,4)(1,2,3)') sage: s.top() [4, 2, 3, 1, 0]
- top_cycle_string()[source]¶
Return the cycle of the top permutation as a string.
EXAMPLES:
sage: from surface_dynamics import * sage: S = SeparatrixDiagram('(0,2)(3,4)','(0)(1,2,3)') sage: S.top_cycle_string() '(0)(1,2,3)(4)'
- top_cycle_tuples()[source]¶
Return the cycle of the top permutation as a list of tuples.
EXAMPLES:
sage: from surface_dynamics import * sage: S = SeparatrixDiagram('(0,2)(3,4)','(0)(1,2,3)') sage: S.top_cycle_tuples() [(0,), (1, 2, 3), (4,)]
- top_orbit(i)[source]¶
Return the orbit of
i
under the top permutation.EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1)(2,5)(3,4,6)','(0,1,5)(2,3,6)(4)') sage: s.top_orbit(0) (0, 1, 5) sage: s.top_orbit(6) (2, 3, 6)
- top_perm()[source]¶
Return the top as a permutation
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0)(1,2)','(1)(0,2)') sage: s.top_perm() (1,3)
- vertical_symmetry()[source]¶
Return the vertical symmetric of this separatrix diagram.
EXAMPLES:
sage: from surface_dynamics import * sage: s = SeparatrixDiagram('(0,1,2,3)(4,5)','(1,2,3)(4,5,0)') sage: sv = s.vertical_symmetry() sage: print(sv) (0,3,2,1)(4,5)-(0,5,4)(1,3,2) sage: c = sv.cylinder_diagrams(up_to_symmetry=False)[0].vertical_symmetry() sage: cv = sv.cylinder_diagrams(up_to_symmetry=False)[0] sage: c.is_isomorphic(cv) True
- vertices()[source]¶
Return a list of pairs
(vout, vin)
where each pair represent a vertex andvout
,vin
are respectively the labels of outgoing/incoming separatrices at this vertex.EXAMPLES:
sage: from surface_dynamics import CylinderDiagram sage: from surface_dynamics.misc.permutation import perm_cycles sage: c = CylinderDiagram('(0,1)-(0,2,6,7,4,5,3,8) (2,7,3,4,8,6,5)-(1)') sage: c.vertices() [([0, 1, 8, 7], [0, 4, 2, 1]), ([2, 4, 5], [3, 6, 5]), ([3, 6], [7, 8])] sage: perm_cycles(c.outgoing_edges_perm()) [[0, 1, 8, 7], [2, 4, 5], [3, 6]] sage: perm_cycles(c.incoming_edges_perm()) [[0, 4, 2, 1], [3, 6, 5], [7, 8]] sage: c = CylinderDiagram('(0,1,3,4,2)-(0,1,5,7,6) (5,6)-(4) (7)-(2,3)') sage: perm_cycles(c.outgoing_edges_perm()) [[0, 3], [1, 6], [2, 4], [5, 7]] sage: perm_cycles(c.incoming_edges_perm()) [[0, 5], [1, 2], [3, 4], [6, 7]]
- surface_dynamics.flat_surfaces.separatrix_diagram.cyclic_direction(x, y, z)[source]¶
Returns 1 or -1 depending on the cyclic ordering of (x,y,z)
- surface_dynamics.flat_surfaces.separatrix_diagram.hyperelliptic_cylinder_diagram_iterator(a, verbose=False)[source]¶
Return an iterator over cylinder diagrams of Abelian differentials that double covers Q((a-2), -1^(a+2)).
The generator is up to isomorphism.
TODO:
An optimization could be obtained by considering the generation of k-subsets of {1,…,n} up to the cyclic symmetry of the tree.
INPUT:
a
- integer - angle of the conical singularity of the quadratic differential.verbose
- integer (default: 0) - output various information during the iteration (mainly for debug).
EXAMPLES:
sage: from surface_dynamics import * sage: from surface_dynamics.flat_surfaces.separatrix_diagram import hyperelliptic_cylinder_diagram_iterator sage: it = hyperelliptic_cylinder_diagram_iterator(3) sage: c = next(it); c.is_isomorphic(CylinderDiagram('(0,1)-(0,2) (2)-(1)')) True sage: c.stratum_component() H_2(2)^hyp sage: hyp = Stratum([2,2], k=1).hyperelliptic_component() sage: all(c.stratum_component() == hyp for c in hyperelliptic_cylinder_diagram_iterator(6)) True
- surface_dynamics.flat_surfaces.separatrix_diagram.move_backward(i, v, g01, g23)[source]¶
Helper function to build pillowcase covers
EXAMPLES:
sage: from surface_dynamics.flat_surfaces.separatrix_diagram import move_backward sage: g23 = [1,2,0] sage: g01 = [2,0,1] sage: i,v = 0,0 sage: for _ in range(6): ....: i,v = move_backward(i, v, g01, g23) ....: print("%d %d" % (i,v)) 2 1 2 0 1 1 1 0 0 1 0 0 sage: i,v = 0,2 sage: for _ in range(6): ....: i,v = move_backward(i, v, g01, g23) ....: print("%d %d" % (i,v)) 1 3 1 2 2 3 2 2 0 3 0 2
- surface_dynamics.flat_surfaces.separatrix_diagram.move_forward(i, v, g01, g23)[source]¶
Helper function to build pillowcase covers
EXAMPLES:
sage: from surface_dynamics.flat_surfaces.separatrix_diagram import move_forward sage: g23 = [1,2,0] sage: g01 = [2,0,1] sage: i,v = 0,0 sage: for _ in range(6): ....: i,v = move_forward(i, v, g01, g23) ....: print("%d %d" % (i,v)) 0 1 1 0 1 1 2 0 2 1 0 0 sage: i,v = 0,2 sage: for _ in range(6): ....: i,v = move_forward(i, v, g01, g23) ....: print("%d %d" % (i,v)) 0 3 2 2 2 3 1 2 1 3 0 2
- surface_dynamics.flat_surfaces.separatrix_diagram.orientation_cover(alpha, phi, a, verbose=0)[source]¶
Build the cylinder diagram of Abelian differentials that double covers it.
A quadratic differrential separatrix diagram is given by three permutations
sigma: the permutation of 1/2-separatrices around vertices
- alpha: the permutation of 1/2-separatrices that describe the separatrices
(it is a fixed point free involution)
phi: the permutation of 1/2-separatrices that describe the cycles.
INPUT:
alpha
– permutationphi
– permutationa
– number of half separatrices
EXAMPLES:
sage: from surface_dynamics import * sage: from surface_dynamics.flat_surfaces.separatrix_diagram import orientation_cover sage: alpha = [3, 2, 1, 0, 5, 4, 7, 6] sage: phi = [3, 1, 0, 2, 5, 4, 7, 6] sage: orientation_cover(alpha,phi,3) (0,2)-(0,1) (1)-(2)
- surface_dynamics.flat_surfaces.separatrix_diagram.separatrix_diagram_fast_iterator(profile, ncyls=None)[source]¶
Iterator over separatrix diagram with given
profile
Return a list of 3-tuples
[bot, top, s]
wherebot
andtop
are list on 0, …, nseps-1 that corresponds to a separatrix diagram with profileprofile
whiles
is the element conjugacy class corresponding to the profile which equalsbot * top
.If ncyls is not None, it should be a list of integers from which the number of cylinders is considered.
Warning: each isomorphism class of separatrix diagram is output more than once in general. If you want a unique representative in each isomorphism class you may consider the method separatrix_diagram_iterator instead.
Uses bounds from [Nav08].
EXAMPLES:
sage: from surface_dynamics import * sage: from surface_dynamics.flat_surfaces.separatrix_diagram import separatrix_diagram_fast_iterator sage: for s in sorted(separatrix_diagram_fast_iterator([3])): print(s) ([0, 2, 1], [1, 0, 2], [(0, 1, 2)]) ([1, 2, 0], [1, 2, 0], [(0, 2, 1)]) ([2, 1, 0], [1, 0, 2], [(0, 2, 1)]) sage: for s in sorted(separatrix_diagram_fast_iterator([2,2])): print(s) ([0, 1, 3, 2], [1, 0, 2, 3], [(0, 1), (2, 3)]) ([0, 2, 3, 1], [1, 2, 0, 3], [(0, 1), (2, 3)]) ([1, 2, 3, 0], [1, 2, 3, 0], [(0, 2), (1, 3)]) ([1, 3, 2, 0], [1, 2, 0, 3], [(0, 2), (1, 3)]) ([2, 3, 0, 1], [1, 0, 3, 2], [(0, 3), (1, 2)]) ([3, 1, 0, 2], [1, 2, 0, 3], [(0, 3), (1, 2)]) ([3, 2, 1, 0], [1, 0, 3, 2], [(0, 2), (1, 3)])
- surface_dynamics.flat_surfaces.separatrix_diagram.separatrix_diagram_iterator(profile, ncyls=None)[source]¶
Iterator over separatrix diagram with given
profile
and number of cylinders.Warning: to prevent isomorphism class to be output twice the function implement a cache mechanism. If you intend to iterate through a huge class of separatrix_diagram and do not care about isomorphism problem use separatrix_diagram_fast_iterator instead.
EXAMPLES:
sage: from surface_dynamics import * sage: from surface_dynamics.flat_surfaces.separatrix_diagram import separatrix_diagram_iterator sage: for s in sorted(separatrix_diagram_iterator([1,1])): print(s) (0,1)-(0,1) (0)(1)-(0)(1) sage: for s in sorted(separatrix_diagram_iterator([3])): print(s) (0,1,2)-(0,1,2) (0)(1,2)-(0,1)(2) sage: for s in sorted(separatrix_diagram_iterator([2,2])): print(s) (0,1,2,3)-(0,1,2,3) (0)(1,2,3)-(0,1,2)(3) (0,1)(2,3)-(0,2)(1,3) (0)(1)(2,3)-(0,1)(2)(3) sage: sum(1 for s in separatrix_diagram_iterator([3,2,2])) 64
- surface_dynamics.flat_surfaces.separatrix_diagram.simplex_count(rays)[source]¶
EXAMPLES:
sage: from surface_dynamics.flat_surfaces.separatrix_diagram import simplex_count sage: rays = [(0,1,1), (1,0,1), (1,1,0)] sage: simplex_count(rays) 1 sage: rays = [(0,1,1), (1,0,1), (1,1,0)] sage: simplex_count(rays) 1 sage: rays = [(0,1,1,1),(1,0,1,1),(1,1,0,1),(1,1,1,0)] sage: simplex_count(rays) 1
- surface_dynamics.flat_surfaces.separatrix_diagram.two_non_connected_perms_canonical_labels(bot, top)[source]¶
EXAMPLES:
sage: from surface_dynamics.flat_surfaces.separatrix_diagram import two_non_connected_perms_canonical_labels sage: two_non_connected_perms_canonical_labels([3,2,1,0],[0,1,2,3]) ([1, 0, 3, 2], [0, 1, 2, 3])
Homology of cylinder decomposition.
- class surface_dynamics.flat_surfaces.twist_space.TwistSpace(cd)[source]¶
Bases:
object
Subspace of relative homology generated by horizontal saddle connections.
This subspace contains the core curves of cylinders. The subspace spanned by the latter is an isotropic subspace of absolute homology.
EXAMPLES:
sage: from surface_dynamics import CylinderDiagram sage: from surface_dynamics.flat_surfaces.twist_space import TwistSpace sage: cd = CylinderDiagram("(0,1)-(0,5) (2)-(4) (3,4)-(1) (5)-(2,3)") sage: tw = TwistSpace(cd) sage: tw TwistSpace('(0,1)-(0,5) (2)-(4) (3,4)-(1) (5)-(2,3)') sage: tw.homologous_cylinders() [[2, 3]] sage: tw.cylinder_dimension() 3
In this example, we get that the last two cylinders (3,4)-(1) and (5)-(2,3) are homologous.
- cylinder_core_curves()[source]¶
Iterate through the core curves of cylinders as homology elements.
EXAMPLES:
sage: from surface_dynamics.flat_surfaces.separatrix_diagram import CylinderDiagram sage: from surface_dynamics.flat_surfaces.twist_space import TwistSpace sage: cd = CylinderDiagram('(0,7,1,2)-(3,6,4,5) (3,6,4,5)-(0,7,1,2)') sage: tw = TwistSpace(cd) sage: list(tw.cylinder_core_curves()) [(0, 0, 1, 1, 1, 1, 0), (0, 0, 1, 1, 1, 1, 0)] sage: cd = CylinderDiagram('(0,7,3,4)-(0,5,3,6) (1,5,2,6)-(1,7,2,4)') sage: tw = TwistSpace(cd) sage: list(tw.cylinder_core_curves()) [(1, 0, 0, 1, 1, 1, 0), (0, 1, 1, 0, 1, 1, 0)]
- cylinder_dimension()[source]¶
Return the dimension of the span of core curves.
See also the method
homological_dimension_of_cylinders
forCylinderDiagram
.EXAMPLES:
sage: from surface_dynamics import Stratum sage: from surface_dynamics.flat_surfaces.twist_space import TwistSpace sage: for cd in Stratum([1,1,1,1], k=1).cylinder_diagrams(): ....: assert TwistSpace(cd).cylinder_dimension() == cd.homological_dimension_of_cylinders()
- homologous_cylinders()[source]¶
Return the list of homologous cylinders.
OUTPUT: a list of lists. Each sublist is an equivalence class of > 1 homologous cylinders.
EXAMPLES:
sage: from surface_dynamics import Stratum sage: from surface_dynamics.flat_surfaces.twist_space import TwistSpace sage: for cd in Stratum([1,1,1,1], k=1).cylinder_diagrams(2): ....: hom_cyls = TwistSpace(cd).homologous_cylinders() ....: if hom_cyls: ....: print(cd) ....: print(hom_cyls) (0,7,1,2)-(3,6,4,5) (3,6,4,5)-(0,7,1,2) [[0, 1]]
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)[source]¶
Bases:
SageObject
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
- automorphism_group(fix_vertices=False, fix_edges=False, fix_faces=False)[source]¶
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()[source]¶
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)[source]¶
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)[source]¶
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)[source]¶
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
.
- dual()[source]¶
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
- euler_characteristic()[source]¶
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
- genus()[source]¶
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_cycle(c)[source]¶
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()[source]¶
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()[source]¶
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()[source]¶
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')[source]¶
Return the generating series for the number of lengths with the given boundaries
- monodromy_group()[source]¶
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 ...
- relabel(perm=None)[source]¶
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()[source]¶
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]]
- class surface_dynamics.flat_surfaces.homology.RibbonGraphWithAngles(vertices=None, edges=None, faces=None, angles=None)[source]¶
Bases:
RibbonGraph
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.
- has_trivial_holonomy()[source]¶
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()[source]¶
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)[source]¶
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
- class surface_dynamics.flat_surfaces.homology.RibbonGraphWithHolonomies(vertices=None, edges=None, faces=None, holonomies=None)[source]¶
Bases:
RibbonGraph
A Ribbon graph with holonomies.
For now
Fat graphs (or combinatorial maps)¶
Fat graph.
- class surface_dynamics.topology.fat_graph.FatGraph(vp=None, fp=None, max_num_dart=None, mutable=False, check=True)[source]¶
Bases:
object
Fat graph or graph embedded in an orientable surface.
A fat graph is a graph embedded in an orientable surface. It is encoded as a pair of permutations
vp
(for vertex permutation) andfp
(for face permutation). The encoding is done by first associating an integer to each oriented edge such that the orientation reversing maps is 0 leftrightarrow 1, 2 leftrightarrow 3, etc. Then the vertex permutation is the permutation of the oriented edge which corresponds to move counterclockwise to the next outgoing edge at each vertex. The face permutation is the permutation of the oriented edge which corresponds to follow the arrow inside the face on its left side (so that faces also run counterclockwise).EXAMPLES:
The once punctured torus:
sage: from surface_dynamics import FatGraph sage: vp = '(0,2,1,3)' sage: fp = '(0,2,1,3)' sage: FatGraph(vp, fp) FatGraph('(0,2,1,3)', '(0,2,1,3)')
Actually it is enough to specify one of the vertex permutation or face permutation as one determines the other:
sage: vp = '(0,3,1)(4)(2,5,6,7)' sage: fp = '(0,3,7,5,4,2)(1)(6)' sage: F0 = FatGraph(vp=vp, fp=fp) sage: F1 = FatGraph(fp=fp) sage: F2 = FatGraph(vp=vp, fp=fp) sage: F3 = FatGraph(vp=vp) sage: F0 == F1 and F0 == F2 and F0 == F3 True
- add_edge(i=None, j=None, check=True)[source]¶
Insert an isolated edge.
By default the created edge is isolated. If arguments
i
orj
is provided its start or respectively endINPUT:
i
,j
– optional darts. If specified, the new darts of the edge will be inserted afteri
andj
at the corresponding vertices.
EXAMPLES:
sage: from surface_dynamics import FatGraph sage: cm = FatGraph('', '', mutable=True, max_num_dart=4) sage: cm.add_edge() 0 sage: cm._check() sage: cm FatGraph('(0)(1)', '(0,1)') sage: cm.add_edge() 2 sage: cm._check() sage: cm FatGraph('(0)(1)(2)(3)', '(0,1)(2,3)') sage: cm = FatGraph('(0,1)', '(0)(1)', mutable=True) sage: cm._realloc(8) sage: cm.add_edge(0) 2 sage: cm._check() sage: cm FatGraph('(0,2,1)(3)', '(0,2,3)(1)') sage: cm.add_edge(0) 4 sage: cm._check() sage: cm FatGraph('(0,4,2,1)(3)(5)', '(0,2,3,4,5)(1)') sage: cm.add_edge(3) 6 sage: cm._check() sage: cm FatGraph('(0,4,2,1)(3,6)(5)(7)', '(0,2,6,7,3,4,5)(1)') sage: cm = FatGraph('(0,1)', '(0)(1)', mutable=True) sage: cm._realloc(8) sage: cm.add_edge(0, 1) 2 sage: cm FatGraph('(0,2,1,3)', '(0,2,1,3)') sage: cm.add_edge(0, 3) 4 sage: cm FatGraph('(0,4,2,1,3,5)', '(0,2,1,5)(3,4)') sage: cm.add_edge(1, 2) 6
- add_loop(check=True)[source]¶
Insert an isolated loop.
EXAMPLES:
sage: from surface_dynamics import FatGraph sage: cm = FatGraph('', '', mutable=True, max_num_dart=4) sage: cm.add_loop() 0 sage: cm.add_loop() 2 sage: cm FatGraph('(0,1)(2,3)', '(0)(1)(2)(3)')
- automorphism_face_action(g)[source]¶
Return the action of
g
on faces.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: r = FatGraph('(0,5,4)(1,2,3)', '(0,3,1,4)(2)(5)') sage: A = r.automorphism_group() sage: r.automorphism_face_action(A.gens()[0]) [0, 2, 1]
- automorphism_group(fix_vertices=False, fix_edges=False, fix_faces=False)[source]¶
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: from surface_dynamics.misc.permutation import perm_conjugate
The four unicellular map with 4 edges in genus 2:
sage: cm0 = FatGraph.from_unicellular_word([0,1,0,1,2,3,2,3]) sage: cm1 = FatGraph.from_unicellular_word([0,1,0,2,1,3,2,3]) sage: cm2 = FatGraph.from_unicellular_word([0,1,0,2,3,1,2,3]) sage: cm3 = FatGraph.from_unicellular_word([0,1,2,3,0,1,2,3]) sage: for cm in [cm0, cm1, cm2, cm3]: ....: P = cm.automorphism_group() ....: print(P.group_cardinality()) ....: vp = cm.vertex_permutation() ....: fp = cm.face_permutation() ....: for a in P.gens(): ....: pp = perm_conjugate(vp, a) ....: assert pp == vp, (vp, pp) 2 1 1 8 sage: cm = FatGraph.from_unicellular_word([0,1,2,3,0,4,1,2,3,4]) sage: cm.automorphism_group().group_cardinality() 2
An example with two faces:
sage: vp = '(0,9,6,5,7,4,8,2,1,3)' sage: fp = '(0,2,1,3,8)(4,6,5,7,9)' sage: cm = FatGraph(vp, fp) sage: cm.automorphism_group() PermutationGroupOrbit(10, [(0,4)(1,5)(2,6)(3,7)(8,9)])
One can compute face stabilizer:
sage: r = FatGraph('(0,5,4)(1,2,3)', '(0,3,1,4)(2)(5)') sage: r.automorphism_group().group_cardinality() 2 sage: r.automorphism_group(fix_faces=True).group_cardinality() 1 sage: r = FatGraph('(0,4,3)(5,2,1)', '(0,2,4,1,3,5)') sage: A = r.automorphism_group() sage: A.group_cardinality() 6 sage: r.automorphism_group(fix_faces=True) == A True sage: r = FatGraph('(0,2,1,3)', '(0,2,1,3)') sage: r.automorphism_group().group_cardinality() 4
- automorphism_vertex_action(g)[source]¶
Return the action of
g
on vertices.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: r = FatGraph('(0,5,4)(1,2,3)', '(0,3,1,4)(2)(5)') sage: A = r.automorphism_group() sage: r.automorphism_vertex_action(A.gens()[0]) [1, 0]
- connected_components()[source]¶
Return the list of connected components.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph(vp='(0,2)(1,3)').connected_components() [FatGraph('(0,2)(1,3)', '(0,3)(1,2)')] sage: FatGraph(vp='(0,1)(2,3)').connected_components() [FatGraph('(0,1)', '(0)(1)'), FatGraph('(0,1)', '(0)(1)')]
- contract_edge(i)[source]¶
Contract an edge between two distinct zeros.
Inverse operation of
split_vertex()
except that here we allow vertices of degree one.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: vp = '(0,2,1,3)' sage: fp = '(0,2,1,3)' sage: vp02 = '(0,4,1,3)(2,5)' sage: fp02 = '(0,4,2,1,3,5)' sage: cm = FatGraph(vp02, fp02, mutable=True) sage: cm.contract_edge(4) sage: cm == FatGraph(vp, fp) True sage: cm = FatGraph(vp02, fp02, mutable=True) sage: cm.contract_edge(5) sage: cm == FatGraph(vp, fp) True sage: vp01 = '(0,4,3)(1,5,2)' sage: fp01 = '(0,2,4,1,3,5)' sage: cm = FatGraph(vp01, fp01, mutable=True) sage: cm.contract_edge(4) sage: cm == FatGraph(vp, fp) True sage: cm = FatGraph(vp01, fp01, mutable=True) sage: cm.contract_edge(5) sage: cm == FatGraph(vp, fp) True sage: vp03 = '(0,4)(1,3,5,2)' sage: fp03 = '(0,2,1,4,3,5)' sage: cm = FatGraph(vp03, fp03, mutable=True) sage: cm.contract_edge(4) sage: cm == FatGraph(vp, fp) True sage: cm = FatGraph(vp03, fp03, mutable=True) sage: cm.contract_edge(5) sage: cm == FatGraph(vp, fp) True
Degree 1 vertices:
sage: cm = FatGraph('(0,2)(1)(3)', '(0,1,2,3)', mutable=True) sage: cm.contract_edge(2) sage: cm FatGraph('(0)(1)', '(0,1)') sage: cm2 = FatGraph('(0,2)(1)(3)', '(0,1,2,3)', mutable=True) sage: cm2.contract_edge(3) sage: cm == cm2 True
- copy(mutable=None)[source]¶
Return a copy of this fat graph.
INPUT:
mutable
– (optional) whether the copy must be mutable
EXAMPLES:
sage: from surface_dynamics import FatGraph sage: F = FatGraph.from_unicellular_word([0,1,0,2,3,4,1,4,3,2]) sage: G = F.copy() sage: G._check()
- darts()[source]¶
Iterator through the list of darts.
EXAMPLES:
sage: from surface_dynamics import FatGraph sage: list(FatGraph('(0,3,5)(1,4,2)').darts()) [0, 1, 2, 3, 4, 5]
- disjoint_union(*args, mutable=False)[source]¶
Return the union of
self
with the graphs provided as arguments.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: fg0 = FatGraph('(0,4,2,5)(1,6)(3,7)') sage: fg1 = FatGraph('(0,3,5)(1,2,4)') sage: fg2 = FatGraph('(0,1)') sage: fg = fg0.disjoint_union(fg1, fg2) sage: fg FatGraph('(0,4,2,5)(1,6)(3,7)(8,11,13)(9,10,12)(14,15)', '(0,6,3,4,2,7,1,5)(8,12,11,9,13,10)(14)(15)') sage: fg.connected_components() [FatGraph('(0,4,2,5)(1,6)(3,7)', '(0,6,3,4,2,7,1,5)'), FatGraph('(0,3,5)(1,2,4)', '(0,4,3,1,5,2)'), FatGraph('(0,1)', '(0)(1)')]
- dual()[source]¶
Return the dual fat graph.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: F = FatGraph(fp='(0)(1)') sage: F.dual() sage: F FatGraph('(0)(1)', '(0,1)') sage: F._check() sage: s = '20_i31027546b98jchedfag_23146758ab9igdhfejc0' sage: F = FatGraph.from_string(s) sage: F.dual() sage: F._check() sage: F.dual() sage: F._check() sage: F == FatGraph.from_string(s) True
- edge_flip(i)[source]¶
Return the half-edge that together with
i
makes an edge.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: F = FatGraph('(0,2,1,3)') sage: F.edge_flip(0) 1 sage: F.edge_flip(1) 0 sage: F.edge_flip(4) Traceback (most recent call last): ... ValueError: dart i=4 out of range
- edge_lengths_polytope(b, min_length=0)[source]¶
Return the polytope of edge lengths where the input
b
specifies the length of faces.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: vp = '(0,4,2,5,1,3)' sage: fp = '(0,5)(1,3,4,2)' sage: cm = FatGraph(vp, fp) sage: cm.edge_lengths_polytope([3, 5], min_length=1).vertices_list() [[1, 1, 1, 1, 2, 2], [2, 2, 1, 1, 1, 1]] sage: cm.edge_lengths_polytope([3, 5], min_length=0).vertices_list() [[0, 0, 1, 1, 3, 3], [3, 3, 1, 1, 0, 0]] sage: vp = '(0,3,4)(1,2,6)(5)(7)' sage: fp = '(0,6,7,2)(1,4,5,3)' sage: cm = FatGraph(vp, fp) sage: cm.edge_lengths_polytope([5, 7]) A 2-dimensional polyhedron in QQ^8 defined as the convex hull of 3 vertices
- edges()[source]¶
Return the edges.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph('(0,4,2,5,1,3)', '(0,5)(1,3,4,2)').edges() [[0, 1], [2, 3], [4, 5]]
- euler_characteristic()[source]¶
Return the Euler characteristic of the associated 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.topology.fat_graph import FatGraph sage: vp = '(0,4,2,5,1,3)' sage: fp = '(0,5)(1,3,4,2)' sage: cm = FatGraph(vp, fp) sage: cm.euler_characteristic() 0
A non-connected example (Euler characteristic is additive):
sage: cm.disjoint_union(cm, cm, cm).euler_characteristic() 0
- face_degrees()[source]¶
EXAMPLES:
sage: from surface_dynamics import FatGraph sage: FatGraph('()', '()').face_degrees() [0]
- face_permutation(copy=True)[source]¶
Return the face permutation as a list.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph('(0,4,2)(1,5,3)').face_permutation() [3, 2, 5, 4, 1, 0]
- faces()[source]¶
Return the faces.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph('(0,4,2,5,1,3)', '(0,5)(1,3,4,2)').faces() [[0, 5], [1, 3, 4, 2]]
- static from_string(s)[source]¶
Build a fat graph from a serialized string.
See also
to_string()
.EXAMPLES:
sage: from surface_dynamics import FatGraph sage: s = '20_i31027546b98jchedfag_23146758ab9igdhfejc0' sage: F = FatGraph.from_string(s) sage: F.to_string() == s True sage: FatGraph.from_string('0__') FatGraph('()', '()')
- static from_unicellular_word(w)[source]¶
Build a fat graph from a word on the letters {0, …, n-1} where each letter appears exactly twice.
EXAMPLES:
sage: from surface_dynamics import FatGraph sage: FatGraph.from_unicellular_word([0,1,0,2,3,4,1,4,3,2]) FatGraph('(0,4)(1,3,9,2)(5,6)(7,8)', '(0,2,1,4,6,8,3,9,7,5)') sage: FatGraph.from_unicellular_word([0,1,2,0,3,2,4,1,3,4]) FatGraph('(0,8,4,3,9,6)(1,5,7,2)', '(0,2,4,1,6,5,8,3,7,9)')
- genus()[source]¶
Return the genus of this graph.
If the graph is not connected, then the answer is the sum of the genera of the components.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: vp0 = '(0,4,2,5,1,3)' sage: fp0 = '(0,5)(1,3,4,2)' sage: cm0 = FatGraph(vp0, fp0) sage: cm0.genus() 1 sage: vp1 = '(0,6,5,7,2,1,3,4)' sage: fp1 = '(0,2,1,4,6,5,3,7)' sage: cm1 = FatGraph(vp1, fp1) sage: cm1.genus() 2
Non-connected examples:
sage: cm0.disjoint_union(cm0, cm1).genus() 4 sage: cm1.disjoint_union(cm0).genus() 3
- integral_points(b, min_length=1)[source]¶
Return the edge lengths solution to the face lengths constraint
b
.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: vp = '(0,4,2,5,1,3)' sage: fp = '(0,5)(1,3,4,2)' sage: cm = FatGraph(vp, fp) sage: cm.integral_points((2, 4)) ((1, 1, 1, 1, 1, 1),) sage: cm.integral_points((3, 5)) ((1, 1, 1, 1, 2, 2), (2, 2, 1, 1, 1, 1)) sage: cm.integral_points((5, 11)) ((1, 1, 3, 3, 4, 4), (2, 2, 3, 3, 3, 3), (3, 3, 3, 3, 2, 2), (4, 4, 3, 3, 1, 1))
- is_connected()[source]¶
Return whether the graph is connected.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph(vp='(0,2)(1,3)').is_connected() True sage: FatGraph(vp='(0,1)(2,3)').is_connected() False
- is_cubic()[source]¶
Return whether all vertices have degree 3.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph('(0,4,1)(2,3,5)').is_cubic() True sage: FatGraph('(0,4,1)(2)(3)(5)').is_cubic() False
- is_face_bipartite(certificate=False)[source]¶
Return whether the faces admit a proper 2-coloring.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import * sage: vp = '(0,2,1,3)' sage: fp = '(0,2,1,3)' sage: F = FatGraph(vp, fp, 6, mutable=True) sage: F.is_face_bipartite() False sage: F.split_face(0, 1) sage: F.is_face_bipartite() True sage: vp = '(0,5,2,1,3,4)' sage: fp = '(0,2,1,4)(3,5)' sage: FatGraph(vp, fp).is_face_bipartite() False sage: from surface_dynamics import FatGraphs sage: F = FatGraphs(g=1, nf=3, nv=3, vertex_min_degree=3) sage: F.cardinality_and_weighted_cardinality(filter=lambda x,a: x.is_face_bipartite()) (3, 5/3)
- is_face_regular(d=None)[source]¶
Return whether all faces have the same degree.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph(fp='(0,5,2,3)(1,6,7,4)').is_face_regular() True
- is_triangulation()[source]¶
Return whether the underlying graph is a triangulation.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph(fp='(0,4,1)(2,3,5)').is_triangulation() True sage: FatGraph(fp='(0,4,1)(2)(3)(5)').is_triangulation() False
- is_vertex_regular(d=None)[source]¶
Return whether all vertices have the same degree.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph('(0,5,2,3)(1,6,7,4)').is_vertex_regular() True
- monodromy_group()[source]¶
Return the group generated by the vertex and face permutations.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: r1 = FatGraph('(0)(2)(1,3,4,5)') sage: G1 = r1.monodromy_group() sage: G1 Subgroup ... sage: G1.is_isomorphic(SymmetricGroup(5)) True sage: r2 = FatGraph('(0)(2)(1,3,4)(5,6,7)') sage: G2 = r2.monodromy_group() sage: G2 Subgroup ... sage: G2.is_isomorphic(PSL(2,7)) True
- move_dart(i, j=None, check=True)[source]¶
Move the dart
i
so that it comes right afterj
around the corresponding vertex.If
j
isNone
, theni
becomes the only dart adjacent to a new vertex. This requiresi
to be adjacent to two distinct faces.INPUT:
i
(integer) – a dartj
(integer orNone
) – (optional) dart
EXAMPLES:
sage: from surface_dynamics import FatGraph sage: cm = FatGraph('(0,3,5)(1,4,2)', mutable=True) sage: cm.move_dart(0) sage: cm._check() sage: cm FatGraph('(0)(1,4,2)(3,5)', '(0,2,5,1)(3,4)') sage: cm.move_dart(1) Traceback (most recent call last): ... NotImplementedError: cylinder case: single face adjacent to i sage: cm = FatGraph('(0,3,5)(1,4,2)', mutable=True) sage: cm.move_dart(0, 2) sage: cm._check() sage: cm FatGraph('(0,1,4,2)(3,5)', '(0)(1,2,5)(3,4)') sage: cm = FatGraph('(0,3,5)(1,4,2)', mutable=True) sage: cm.move_dart(0, 1) sage: cm._check() sage: cm FatGraph('(0,4,2,1)(3,5)', '(0,2,5)(1)(3,4)')
- num_edges()[source]¶
Return the number of edges.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph('(0,4,2,5,1,3)', '(0,5)(1,3,4,2)').num_edges() 3 sage: FatGraph('()', '()').num_edges() 0
- num_faces()[source]¶
Return the number of faces.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph('(0,4,2,5,1,3)', '(0,5)(1,3,4,2)').num_faces() 2 sage: FatGraph('()', '()').num_faces() 1
- num_vertices()[source]¶
Return the number of vertices.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph('(0,4,2,5,1,3)', '(0,5)(1,3,4,2)').num_vertices() 1 sage: FatGraph('()', '()').num_vertices() 1
- profile()[source]¶
Return the pair of vertex and face profiles.
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: FatGraph('(0,4,2,5,1,3)', '(0,5)(1,3,4,2)').profile() ([6], [2, 4]) sage: FatGraph('()', '()').profile() ([0], [0])
- relabel(r)[source]¶
Relabel according to the permutation
r
.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: from surface_dynamics.misc.permutation import perm_on_list_inplace, perm_invert sage: cm = FatGraph('(0,11,9,14,7,10,12,4,2,1,3)(5,15,8,6)(13)', ....: '(0,2,1,3,4,6,14,5,12,13,10)(7,8,11)(9,15)') sage: r = [4,5,8,9,1,0,11,10,15,14,2,3,7,6,12,13] sage: cm.relabel(r) sage: cm._check()
- remove_edge(i)[source]¶
Remove an edge.
If the edge has the same face on both sides, then the genus drops by 1. Inverse operation of
split_face()
ortrisect_face()
.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: vp = '(0,2,1,3)' sage: fp = '(0,2,1,3)' sage: cm = FatGraph(vp, fp) sage: vp20 = '(0,5,4,2,1,3)' sage: fp20 = '(0,2,1,3,4)(5)' sage: cm2 = FatGraph(vp20, fp20, 6, mutable=True) sage: cm2.remove_edge(4) sage: cm2 == cm True sage: cm2 = FatGraph(vp20, fp20, 6, mutable=True) sage: cm2.remove_edge(5) sage: cm2 == cm True sage: vp10 = '(0,4,2,5,1,3)' sage: fp10 = '(0,5)(1,3,4,2)' sage: cm2 = FatGraph(vp10, fp10, mutable=True) sage: cm2.remove_edge(4) sage: cm2 == cm True sage: cm2 = FatGraph(vp10, fp10, mutable=True) sage: cm2.remove_edge(5) sage: cm2 == cm True sage: vp30 = '(0,4,2,1,5,3)' sage: fp30 = '(0,2,5)(1,3,4)' sage: cm2 = FatGraph(vp30, fp30, mutable=True) sage: cm2.remove_edge(4) sage: cm2 == cm True sage: cm2 = FatGraph(vp30, fp30, mutable=True) sage: cm2.remove_edge(5) sage: cm2 == cm True sage: vp00 = '(0,5,4,2,1,3)' sage: fp00 = '(0,2,1,3,4)(5)' sage: cm2 = FatGraph(vp00, fp00, mutable=True) sage: cm2.remove_edge(4) sage: cm2 == cm True sage: vp22 = '(0,2,5,4,1,3)' sage: fp22 = '(0,4,2,1,3)(5)' sage: cm2 = FatGraph(vp00, fp00, mutable=True) sage: cm2.remove_edge(4) sage: cm2 == cm True
- rotate(e)[source]¶
Rotate the edge
e
counterclockwise.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: fg = FatGraph("(0,5,3)(1,4,2)", "(0,2,5,1,3,4)", mutable=True) sage: fg.rotate(0) sage: fg FatGraph('(0,5,1,3)(2,4)', '(0,5,2,1,3,4)') sage: fg.rotate(5) sage: fg FatGraph('(0,5,1,4,3)(2)', '(0,5,1,3,2,4)') sage: vp = "(0,8,15)(1,16,12)(2,13,3)(4,11,5)(6,9,7)(10,17,14)" sage: fp = "(0,12,2,13,16,10,4,11,14,8,6,9)(1,15,17)(3)(5)(7)" sage: fg = FatGraph(vp, fp, mutable=True) sage: fg.rotate(0) sage: fg.rotate(1) sage: fg FatGraph('(0,14,10,17)(1,13,3,2)(4,11,5)(6,9,7)(8,15)(12,16)', '(0,2,13,16,10,4,11,14,8,6,9,15)(1,17,12)(3)(5)(7)')
- rotate_dual(e)[source]¶
Rotate the edge
e
counterclockwise in the dual graph.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: fg = FatGraph("(0,5,3)(1,4,2)", "(0,2,5,1,3,4)") sage: vp = "(0,8,15)(1,16,12)(2,13,3)(4,11,5)(6,9,7)(10,17,14)" sage: fp = "(0,12,2,13,16,10,4,11,14,8,6,9)(1,15,17)(3)(5)(7)" sage: fg = FatGraph(vp, fp, mutable=True) sage: fg.rotate_dual(0) sage: fg.rotate_dual(1) sage: fg FatGraph('(0,15,16)(1,12,8)(2,13,3)(4,11,5)(6,9,7)(10,17,14)', '(0,8,6,9,12,2,13,1,16,10,4,11,14)(3)(5)(7)(15,17)')
- split_face(i, j)[source]¶
Insert an edge between the darts
i
andj
to split the face.This method only work when
i
andj
belongs to the same face.One of the face will contains i, fp[i], …, (the x-face) and the other one will contain j, fp[j], … In the special case i=j, a monogon (= face with only one edge) is created.
The converse operation is implemented in
remove_edge()
.This operation is equivalent to first adding an isolated edge with
add_edge()
and perform twoswap()
:swap(j, k)
andswap(i, k + 1)
where(k, k+1)
is the added edge. Though, contrarily toswap()
this method does not perform relabelling and runs in constant time for labelled fat graphs.EXAMPLES:
The once punctured torus:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: vp = '(0,2,1,3)' sage: fp = '(0,2,1,3)' sage: vp20 = '(0,4,2,5,1,3)' sage: fp20 = '(0,5)(1,3,4,2)' sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: cm.split_face(2, 0) sage: cm == FatGraph(vp20, fp20) True sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: k = cm.add_edge() sage: cm.swap(0, k) sage: cm.swap(2, k + 1) sage: cm == FatGraph(vp20, fp20) True sage: vp10 = '(0,4,2,1,5,3)' sage: fp10 = '(0,2,5)(1,3,4)' sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: cm.split_face(1, 0) sage: cm == FatGraph(vp10, fp10) True sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: k = cm.add_edge() sage: cm.swap(0, k) sage: cm.swap(1, k + 1) sage: cm == FatGraph(vp10, fp10) True sage: vp30 = '(0,4,2,1,3,5)' sage: fp30 = '(0,2,1,5)(3,4)' sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: cm.split_face(3, 0) sage: cm == FatGraph(vp30, fp30) True sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: k = cm.add_edge() sage: cm.swap(0, k) sage: cm.swap(3, k + 1) sage: cm == FatGraph(vp30, fp30) True sage: vp00 = '(0,5,4,2,1,3)' sage: fp00 = '(0,2,1,3,4)(5)' sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: cm.split_face(0, 0) sage: cm == FatGraph(vp00, fp00) True sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: k = cm.add_edge() sage: cm.swap(0, k) sage: cm.swap(0, k + 1) sage: cm == FatGraph(vp00, fp00) True sage: vp22 = '(0,2,5,4,1,3)' sage: fp22 = '(0,4,2,1,3)(5)' sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: cm.split_face(2, 2) sage: cm == FatGraph(vp22, fp22) True sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: k = cm.add_edge() sage: cm.swap(2, k) sage: cm.swap(2, k + 1) sage: cm == FatGraph(vp22, fp22) True
A genus 2 surface:
sage: vp = '(0,3,6,8)(1,10,9,12,5)(2,7,4,11)(13)' sage: fp = '(0,5,7,3,11,1,8,10,4,12,13,9,6,2)' sage: cm = FatGraph(vp, fp, 21, mutable=True); cm._check() sage: cm.split_face(0, 1); cm._check() sage: cm.split_face(4, 13); cm._check() sage: cm.split_face(5, 14); cm._check() sage: cm FatGraph('(0,15,3,6,8)(1,14,18,10,9,12,5,19)(2,7,4,17,11)(13,16)', '(0,19,14)(1,8,10,17,13,9,6,2,15)(3,11,18,5,7)(4,12,16)') sage: cm.remove_edge(18); cm._check() sage: cm.remove_edge(16); cm._check() sage: cm.remove_edge(14); cm._check() sage: cm == FatGraph(vp, fp, 21) True
- split_vertex(i, j)[source]¶
Insert a new edge to split the vertex located at the darts i and j.
This operation keeps the genus constant. The inverse operation is implemented in
contract_edge()
.This operation is equivalent to first adding an isolated loop with
add_loop()
and perform twoswap()
:swap(j, k)
, andswap(i, k + 1)
where(k, k+1)
is the added edge. Though, contrarily toswap()
this method does not perform relabelling and runs in constant time for labelled fat graphs.EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: vp = '(0,2,1,3)' sage: fp = '(0,2,1,3)' sage: vp02 = '(0,4,1,3)(2,5)' sage: fp02 = '(0,4,2,1,3,5)' sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: cm.split_vertex(0, 2) sage: cm == FatGraph(vp02, fp02) True sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: k = cm.add_loop() sage: cm.swap(2, k) sage: cm.swap(0, k + 1) sage: cm == FatGraph(vp02, fp02) True sage: vp01 = '(0,4,3)(1,5,2)' sage: fp01 = '(0,2,4,1,3,5)' sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: cm.split_vertex(0, 1) sage: cm == FatGraph(vp01, fp01) True sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: k = cm.add_loop() sage: cm.swap(1, k) sage: cm.swap(0, k + 1) sage: cm == FatGraph(vp01, fp01) True sage: vp03 = '(0,4)(1,3,5,2)' sage: fp03 = '(0,2,1,4,3,5)' sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: cm.split_vertex(0, 3) sage: cm == FatGraph(vp03, fp03) True sage: cm = FatGraph(vp, fp, 6, mutable=True) sage: k = cm.add_loop() sage: cm.swap(3, k) sage: cm.swap(0, k + 1) sage: cm == FatGraph(vp03, fp03) True
- swap(i, j, check=True)[source]¶
Modify the fat graph by multiplying the vertex and face permutations by the transposition (i,j) respectively on left and right.
EXAMPLES:
sage: from surface_dynamics import FatGraph sage: g = FatGraph('(0,1,2,3)', '(0)(2)(3,1)', mutable=True) sage: g.swap(0, 1) sage: g FatGraph('(0,2,3)(1)', '(0,1,3)(2)') sage: g.swap(0, 2) sage: g FatGraph('(0,3)(1)(2)', '(0,1,3,2)') sage: g.swap(0, 3) sage: g FatGraph('(0)(1)(2)(3)', '(0,1)(2,3)')
- to_string()[source]¶
Serialization to string.
EXAMPLES:
sage: from surface_dynamics import FatGraph sage: FatGraph.from_unicellular_word([0,1,0,2,3,4,1,4,3,2]).to_string() '10_4319065872_2419608537' sage: FatGraph('', '').to_string() '0__'
- trisect_face(i, j, k)[source]¶
Insert a bridge
INPUT:
i
,j
,k
- dart in the same face in counter-clockwise order
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import FatGraph sage: vp = '(0,2,1,3)' sage: fp = '(0,2,1,3)' sage: cm = FatGraph(vp, fp, 8, mutable=True) sage: vp021 = '(0,7,2,6,5,1,4,3)' sage: fp021 = '(0,5,1,3,7,2,4,6)' sage: cm021 = FatGraph(vp021, fp021) sage: cm.trisect_face(0, 2, 1) sage: cm == cm021 True sage: cm = FatGraph(vp, fp, 10, mutable=True) sage: cm.trisect_face(0, 0, 3) sage: cm = FatGraph(vp, fp, 10, mutable=True) sage: cm.trisect_face(0, 3, 3) sage: cm = FatGraph(vp, fp, 10, mutable=True) sage: cm.trisect_face(0, 3, 0) sage: cm = FatGraph(vp, fp, 10, mutable=True) sage: cm.trisect_face(0, 0, 0)
- vertex_degrees()[source]¶
Return the degree of vertices.
EXAMPLES:
sage: from surface_dynamics import FatGraph sage: FatGraph('(0,3)(1,2,5)(4)').vertex_degrees() [2, 3, 1] sage: FatGraph('()', '()').vertex_degrees() [0]
- surface_dynamics.topology.fat_graph.list_extrems(l, n)[source]¶
EXAMPLES:
sage: from surface_dynamics.topology.fat_graph import list_extrems sage: list_extrems([3,5,3,17,7,2,1,19],4) (3, 17)
Exhaustive generation of fat graphs.
This is done following the McKay canonical augmentation. This module is experimental.
- class surface_dynamics.topology.fat_graph_exhaustive_generation.CountAndWeightedCount[source]¶
Bases:
object
- class surface_dynamics.topology.fat_graph_exhaustive_generation.FatGraphs(g=None, nf=None, ne=None, nv=None, vertex_min_degree=0, g_min=None, g_max=None, nf_min=None, nf_max=None, ne_min=None, ne_max=None, nv_min=None, nv_max=None)[source]¶
Bases:
object
Isomorphism classes of fat graphs with topological constraints.
EXAMPLES:
sage: from surface_dynamics import FatGraphs
Trees and their dual (maps with single vertex) in genus zero are counted by Catalan numbers:
sage: for n in range(2, 10): ....: ntrees1 = 2 * (n-1) * FatGraphs(g=0, nf=n, nv=1).weighted_cardinality() ....: ntrees2 = 2 * (n-1) * FatGraphs(g=0, nf=1, nv=n).weighted_cardinality() ....: assert catalan_number(n-1) == ntrees1 == ntrees2, (n, ntrees1, ntrees2)
Tutte formulas in genus 0 (combinatorial maps counted by number of edges):
sage: R.<t> = QQ[] sage: F = FatGraphs(g=0, ne_max=8) sage: poly = R.zero() sage: def update(cm, aut): ....: global poly ....: p = [cm.face_degrees()] ....: aut_card = 1 if aut is None else aut.group_cardinality() ....: poly += 2*cm.num_edges() // aut_card * t**cm.num_edges() sage: F.map_reduce(update) sage: poly 208494*t^7 + 24057*t^6 + 2916*t^5 + 378*t^4 + 54*t^3 + 9*t^2 + 2*t sage: sum(2 * 3**n / (n+2) / (n+1) * binomial(2*n,n) *t**n for n in range(1,8)) 208494*t^7 + 24057*t^6 + 2916*t^5 + 378*t^4 + 54*t^3 + 9*t^2 + 2*t
Genus zero with same number of vertices and faces:
sage: FatGraphs(g=0, nf=2, nv=2).cardinality_and_weighted_cardinality() (2, 5/4) sage: FatGraphs(g=0, nf=3, nv=3).cardinality_and_weighted_cardinality() (23, 41/2) sage: FatGraphs(g=0, nf=4, nv=4).cardinality_and_weighted_cardinality() (761, 8885/12)
Duality checks
sage: for g,nf,nv in [(0,2,3), (0,2,4), (0,3,4), ….: (1,1,2), (1,1,3), (1,1,4), (1,1,5), (1,2,3), (1,2,4), (1,3,4), ….: (2,1,2), (2,1,3)]: ….: n1 = FatGraphs(g=g, nf=nf, nv=nv).cardinality_and_weighted_cardinality() ….: n2 = FatGraphs(g=g, nf=nv, nv=nf).cardinality_and_weighted_cardinality() ….: assert n1 == n2 ….: print(g, nf, nv, n1) 0 2 3 (5, 11/3) 0 2 4 (14, 93/8) 0 3 4 (108, 103) 1 1 2 (3, 5/3) 1 1 3 (11, 35/4) 1 1 4 (46, 42) 1 1 5 (204, 385/2) 1 2 3 (180, 172) 1 2 4 (1198, 14065/12) 1 3 4 (18396, 18294) 2 1 2 (53, 483/10) 2 1 3 (553, 539)
Unicellular map with one vertex in genus 3:
sage: FatGraphs(g=3, nf=1, nv=1).cardinality_and_weighted_cardinality() (131, 495/4)
Minimum vertex degree bounds:
sage: for k in range(2,5): ....: F = FatGraphs(g=1, nf=2, nv=2, vertex_min_degree=1) ....: c1 = F.cardinality_and_weighted_cardinality(lambda cm,_: cm.vertex_min_degree() >= k) ....: G = FatGraphs(g=1, nf=2, nv=2, vertex_min_degree=k) ....: c2 = G.cardinality_and_weighted_cardinality() ....: assert c1 == c2 ....: print(c1) (14, 87/8) (8, 47/8) (4, 15/8) sage: for k in range(2,6): ....: F = FatGraphs(g=0, nf=5, nv=2, vertex_min_degree=1) ....: c1 = F.cardinality_and_weighted_cardinality(lambda cm,_: cm.vertex_min_degree() >= k) ....: G = FatGraphs(g=0, nf=5, nv=2, vertex_min_degree=k) ....: c2 = G.cardinality_and_weighted_cardinality() ....: assert c1 == c2 ....: print(c1) (28, 123/5) (21, 88/5) (13, 48/5) (7, 18/5)
Using ranges for vertex and face numbers:
sage: F = FatGraphs(g=1, nv_min=2, nv_max=4, nf_min=2, nf_max=4) sage: def check(cm,aut): ....: if cm.num_vertices() < 2 or \ ....: cm.num_vertices() >= 4 or \ ....: cm.num_faces() < 2 or \ ....: cm.num_faces() >= 4: ....: raise ValueError(str(cm)) sage: F.map_reduce(check) sage: for nf in [2,3]: ....: for nv in [2,3]: ....: c1 = F.cardinality_and_weighted_cardinality(lambda cm,_: cm.num_vertices() == nv and cm.num_faces() == nf) ....: c2 = FatGraphs(g=1, nf=nf, nv=nv).cardinality_and_weighted_cardinality() ....: assert c1 == c2 ....: print(nf, nv, c1) 2 2 (24, 167/8) 2 3 (180, 172) 3 2 (180, 172) 3 3 (2048, 6041/3)
- list()[source]¶
EXAMPLES:
sage: from surface_dynamics import FatGraphs sage: L21 = FatGraphs(g=0, nf=2, nv=1).list() sage: L21[0].num_faces() 2 sage: L21[0].num_vertices() 1 sage: L12 = FatGraphs(g=0, nf=1, nv=2).list() sage: L12[0].num_faces() 1 sage: L12[0].num_vertices() 2
- map_reduce(callback, filter=None)[source]¶
EXAMPLES:
sage: from surface_dynamics import FatGraphs sage: FatGraphs(g=1, nf=2, nv=2).map_reduce(lambda x,y: print(x)) FatGraph('(0,6,5,4,2,1,3)(7)', '(0,2,1,3,4,6,7)(5)') FatGraph('(0,6,2,1,3)(4,7,5)', '(0,2,1,3,6,4,7)(5)') FatGraph('(0,5,4,2,1,6,3)(7)', '(0,2,6,7,1,3,4)(5)') FatGraph('(0,7,3)(1,6,5,4,2)', '(0,2,7,1,3,4,6)(5)') FatGraph('(0,5,4,2,6,1,3)(7)', '(0,6,7,2,1,3,4)(5)') FatGraph('(0,7,1,3)(2,6,5,4)', '(0,7,2,1,3,4,6)(5)') FatGraph('(0,5,4,2,1,3,6)(7)', '(0,2,1,6,7,3,4)(5)') FatGraph('(0,7)(1,3,6,5,4,2)', '(0,2,1,7,3,4,6)(5)') FatGraph('(0,5,4,6,2,1,3)(7)', '(0,2,1,3,6,7,4)(5)') FatGraph('(0,5,4,6,1,3)(2,7)', '(0,6,2,1,3,7,4)(5)') FatGraph('(0,5,6,4,2,1,3)(7)', '(0,2,1,3,4)(5,6,7)') FatGraph('(0,5,6,2,1,3)(4,7)', '(0,2,1,3,6,4)(5,7)') FatGraph('(0,5,6,1,3)(2,7,4)', '(0,6,2,1,3,4)(5,7)') FatGraph('(0,5,6,3)(1,7,4,2)', '(0,2,6,1,3,4)(5,7)') FatGraph('(0,6,5,2,1,4,3)(7)', '(0,2,4,6,7)(1,3,5)') FatGraph('(0,6,2,1,4,3)(5,7)', '(0,2,4,7)(1,3,6,5)') FatGraph('(0,6,4,3)(1,7,5,2)', '(0,2,4,7)(1,3,5,6)') FatGraph('(0,6,5,2,1,3,4)(7)', '(0,2,1,4,6,7)(3,5)') FatGraph('(0,5,2,6,1,3,4)(7)', '(0,6,7,2,1,4)(3,5)') FatGraph('(0,5,2,6,3,4)(1,7)', '(0,7,2,6,1,4)(3,5)') FatGraph('(0,5,2,1,3,6,4)(7)', '(0,2,1,4)(3,5,6,7)') FatGraph('(0,5,2,1,3,6)(4,7)', '(0,2,1,6,4)(3,5,7)') FatGraph('(0,7,4)(1,3,6,5,2)', '(0,2,1,4,6)(3,5,7)') FatGraph('(0,5,7,4)(1,3,6,2)', '(0,2,1,4)(3,6,5,7)')
- class surface_dynamics.topology.fat_graph_exhaustive_generation.FatGraphsTrace(filename=None, verbosity=0)[source]¶
Bases:
object
A class to trace the execution of the fat graphs generation.
It is mostly used for debugging/profiling/illustration purposes.
- surface_dynamics.topology.fat_graph_exhaustive_generation.FatGraphs_g_nf_nv(g=None, nf=None, nv=None, vertex_min_degree=1, g_min=None, g_max=None, nf_min=None, nf_max=None, nv_min=None, nv_max=None)[source]¶
- class surface_dynamics.topology.fat_graph_exhaustive_generation.ListCallback(mutable=False)[source]¶
Bases:
object
- class surface_dynamics.topology.fat_graph_exhaustive_generation.StackCallback(gmin, gmax, fmin, fmax, emin, emax, vmin, vmax, vertex_min_degree, callback, filter)[source]¶
Bases:
object
- surface_dynamics.topology.fat_graph_exhaustive_generation.augment1(cm, aut_grp, g, callback)[source]¶
Given a unicellular map
cm
with a single vertex and automorphism groupaut_grp
, iterate through all its canonical extensions that are uniface-univertex maps of greater genus.This operation inserts two edges.
This augmentation function is sufficient to iterate through unicellular map.
- surface_dynamics.topology.fat_graph_exhaustive_generation.augment2(cm, aut_grp, depth, callback)[source]¶
Given a map
cm
with a single vertex and automorphism groupaut_grp
iterate through all its canonical extensions that are obtained by splitting one of its faces (by adding a single edge).Because of the chosen canonical labellings, we only need to consider the faces with maximal degree and split in such way that the secondly created face is still at least as big as the second biggest.
- surface_dynamics.topology.fat_graph_exhaustive_generation.augment3(cm, aut_grp, depth, min_degree, callback)[source]¶
Given a map
cm
, its automorphism groupaut_grp
and a minimum degreemin_degree
, iterate through all the canonical extensions ofcm
that are obtained by splittingdepth
times a vertices into two vertices.This augmentation add
depth
edges to the fat graph.In principle, because of the chosen canonical labellings, we only need to consider the vertices with maximal degree.