fat_graph_exhaustive_generation#

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

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)

TESTS:

sage: from surface_dynamics.topology.fat_graph_exhaustive_generation import FatGraphs
sage: FatGraphs(g=0, nf=1, nv=1).list()
[FatGraph('()', '()')]
sage: FatGraphs(g=1, nf=1, nv=1).list()
[FatGraph('(0,2,1,3)', '(0,2,1,3)')]
sage: FatGraphs(g=0, nf=2, nv=1).list()
[FatGraph('(0,1)', '(0)(1)')]
sage: FatGraphs(g=0, nf=1, nv=2).list()
[FatGraph('(0)(1)', '(0,1)')]
sage: FatGraphs(g=0, ne=1).list()
[FatGraph('(0,1)', '(0)(1)'), FatGraph('(0)(1)', '(0,1)')]

sage: F3 = FatGraphs(g=0, nf_max=4, vertex_min_degree=3)
sage: for fg in F3.list(): assert all(d >= 3 for d in fg.vertex_degrees()), fg
sage: F3.cardinality_and_weighted_cardinality()
(3, 7/6)
sage: F4 = FatGraphs(g=0, nf_max=4, vertex_min_degree=4)
sage: for fg in F4.list(): assert all(d >= 4 for d in fg.vertex_degrees()), fg
sage: F4.cardinality_and_weighted_cardinality()
(1, 1/2)

sage: FatGraphs(g=1, ne=3).list()
[FatGraph('(0,5,4,2,1,3)', '(0,2,1,3,4)(5)'),
 FatGraph('(0,5,2,1,4,3)', '(0,2,4)(1,3,5)'),
 FatGraph('(0,5,2,1,3,4)', '(0,2,1,4)(3,5)'),
 FatGraph('(0,4,2,1,3)(5)', '(0,2,1,3,4,5)'),
 FatGraph('(0,4,1,3)(2,5)', '(0,4,2,1,3,5)'),
 FatGraph('(0,4,3)(1,5,2)', '(0,2,4,1,3,5)')]
cardinality_and_weighted_cardinality(filter=None)#
list()#

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

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)')
weighted_cardinality(filter=None)#
class surface_dynamics.topology.fat_graph_exhaustive_generation.FatGraphsTrace(filename=None, verbosity=0)#

A class to trace the execution of the fat graphs generation.

It is mostly used for debugging/profiling/illustration purposes.

add_edge(s0, s1)#
add_vertex(s, aut_grp, depth)#
canonical_edge(parent, cm, aut_grp, caller)#
graphviz_tree(filename=None)#
non_canonical_edge(parent, cm, caller)#
root(cm, aut_grp)#
summary(filename=None)#
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)#

TESTS:

sage: from surface_dynamics.topology.fat_graph_exhaustive_generation import FatGraphs_g_nf_nv
sage: F = FatGraphs_g_nf_nv(g=1, nf=1, nv=1)
doctest:warning
...
DeprecationWarning: FatGraphs_g_nf_nv is deprecated. Use FatGraphs
sage: F.list()
[FatGraph('(0,2,1,3)', '(0,2,1,3)')]
class surface_dynamics.topology.fat_graph_exhaustive_generation.ListCallback(mutable=False)#
list()#
class surface_dynamics.topology.fat_graph_exhaustive_generation.StackCallback(gmin, gmax, fmin, fmax, emin, emax, vmin, vmax, vertex_min_degree, callback, filter)#
run()#
surface_dynamics.topology.fat_graph_exhaustive_generation.augment1(cm, aut_grp, g, callback)#

Given a unicellular map cm with a single vertex and automorphism group aut_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)#

Given a map cm with a single vertex and automorphism group aut_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)#

Given a map cm, its automorphism group aut_grp and a minimum degree min_degree, iterate through all the canonical extensions of cm that are obtained by splitting depth 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.

surface_dynamics.topology.fat_graph_exhaustive_generation.minmax(l)#
surface_dynamics.topology.fat_graph_exhaustive_generation.num_and_weighted_num(it)#