fat_graph#

Fat graph.

class surface_dynamics.topology.fat_graph.FatGraph(vp=None, fp=None, max_num_dart=None, mutable=False, check=True)#

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) and fp (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)#

Insert an isolated edge.

By default the created edge is isolated. If arguments i or j is provided its start or respectively end

INPUT:

  • i, j – optional darts. If specified, the new darts of the edge will be inserted after i and j 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)#

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)#

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)#

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)#

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()#

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)#

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)#

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()#

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(mutable, *args)#

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()#

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)#

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)#

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

TESTS:

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, 2])
Traceback (most recent call last):
...
ValueError: the length of b must be the number of faces
edges()#

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()#

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()#

EXAMPLES:

sage: from surface_dynamics import FatGraph
sage: FatGraph('()', '()').face_degrees()
[0]
face_max_degree()#
face_min_degree()#
face_permutation(copy=True)#

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()#

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)#

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)#

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()#

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)#

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()#

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()#

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)#

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)#

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()#

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)#

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()#

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)#

Move the dart i so that it comes right after j around the corresponding vertex.

If j is None, then i becomes the only dart adjacent to a new vertex. This requires i to be adjacent to two distinct faces.

INPUT:

  • i (integer) – a dart

  • j (integer or None) – (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()#

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()#

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()#

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()#

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)#

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)#

Remove an edge.

If the edge has the same face on both sides, then the genus drops by 1. Inverse operation of split_face() or trisect_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
remove_face_trisection(x)#

Remove the face trisection at x.

TESTS:

sage: from surface_dynamics.topology.fat_graph import FatGraph

sage: vp = '(0,2,1,3)'
sage: fp = '(0,2,1,3)'
sage: for i, j, k in [(0, 2, 1), (0, 0, 3), (0, 3, 3), (0, 3, 0), (0, 0, 0)]:
....:     cm = FatGraph(vp, fp, 8, mutable=True)
....:     cm.trisect_face(i, j, k)
....:     cm.remove_face_trisection(4)
....:     assert cm == FatGraph(vp, fp)
rotate(e)#

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)#

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)#

Insert an edge between the darts i and j to split the face.

This method only work when i and j 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 two swap(): swap(j, k) and swap(i, k + 1) where (k, k+1) is the added edge. Though, contrarily to swap() 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)#

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 two swap(): swap(j, k), and swap(i, k + 1) where (k, k+1) is the added edge. Though, contrarily to swap() 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)#

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()#

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)#

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()#

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]
vertex_max_degree()#
vertex_min_degree()#
vertex_permutation(copy=True)#

Return the vertex permutation as a list.

EXAMPLES:

sage: from surface_dynamics.topology.fat_graph import FatGraph
sage: FatGraph('(0,4,2)(1,5,3)').vertex_permutation()
[4, 5, 0, 1, 2, 3]
vertices()#

Return the vertices as a list of lists.

EXAMPLES:

sage: from surface_dynamics.topology.fat_graph import FatGraph
sage: FatGraph('(0,4,2,5,1,3)', '(0,5)(1,3,4,2)').vertices()
[[0, 4, 2, 5, 1, 3]]
surface_dynamics.topology.fat_graph.list_extrems(l, n)#

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)
surface_dynamics.topology.fat_graph.num_and_weighted_num(it)#