Origamis

Individual origamis and Teichmüller curves

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1)

Dense origamis

An origami is a pair (r, u) of permutations up to conjugacy. The canonical representative of an origami is computed by an independent C program in normal_form.c.

A pillowcase cover is a quadruple (g_0, g_1, g_2, g_3) of permutations such that the product g_0 g_1 g_2 g_3 is the identity. It is sometimes called a 4-constellation. The canonical representative of a pillowcase cover is computed by an independent C program in normal_form.c.

class surface_dynamics.flat_surfaces.origamis.origami_dense.Origami_dense_pyx

Bases: object

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 405)

Origami or square-tiled surface.

An origami is a flat surface which is a covering of a once-punctured torus. It can be described either by a pair of permutations, up to simultaneous conjugacy in the symmetric group, or by a subgroup of finite index of the free group on two generators, which is the fundamental group of the once-punctured torus.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: Origami([2, 1, 3], [3, 2, 1])
(1,2)(3)
(1,3)(2)
absolute_period_generators()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3333)

Return a generating set of the absolute periods of this origami.

To each curve on an origami, we can associate its holonomy (that is an element of ZZ times ZZ). This function returns a generating set of the module generated by holonomies of closed curves.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2,3,4)(5,6)', '(1,5)(2,6)')
sage: o.absolute_period_generators()
[(2, 0), (2, 0), (0, 1), (0, 1)]
as_graph()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1979)

Return the graph associated to self

The graph associated to an origami is the graph on [1, …, N] for which the edges correspond to the action of the permutations r and u.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: G = o.as_graph(); G
Looped multi-digraph on 3 vertices
sage: G.vertices(sort=True)
[0, 1, 2]
sage: G.edges(sort=True)
[(0, 1, 'r'), (0, 2, 'u'), (1, 0, 'r'), (1, 1, 'u'), (2, 0, 'u'), (2, 2, 'r')]
automorphism_group()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2824)

Return the automorphism group of the origami as a permutation group.

The automorpism group of a translation surface is the set of affine diffeomorphisms which have a trivial linear part. For an origami, it corresponds combinatorially to the centralizer of the group generated by the permutations r and u that define this origami.

EXAMPLES:

sage: from surface_dynamics.all import Origami

The L with 3 squares has no automorphisms:

sage: o = Origami('(1,2)', '(1,3)')
sage: G = o.automorphism_group()
sage: G.order()
1

This 4-square square-tiled surface in H(1, 1) has non trivial automorphism for which the quotient is a torus with two squares:

sage: o = Origami('(1,2)(3,4)', '(2,3)')
sage: G = o.automorphism_group()
sage: G.order()
2
sage: oo = o.quotient(G)
sage: oo
(1,2)
(1)(2)
connected_components()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2029)

Return the list of connected origami that composes this origami.

cover(sr, su, check=True, as_tuple=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2694)

Build the (ramified) cover of this origami by sr and su.

INPUT:

  • sr, su - two list of N permutations where N is the number of squares of this origami

  • check - whether or not check the input

  • as_tuple - assume that sr and su are list of tuples of the same length and corresponds to permutations of [0, …, d-1] (much more efficient in time)

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)'); o
(1,2)(3)
(1,3)(2)
sage: o.cover(['(1,2)', '', ''], ['', '', ''])
(1,5,4,2)(3)(6)
(1,3)(2)(4,6)(5)
cylinder_decomposition()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3121)

Return the cylinder decomposition of the standard form of this origami.

OUTPUT:

A list of cylinders where each cylinder is a 6-tuple (bot, top, w, h, bot_twist, top_twist) where

  • bot and top are list of right squares adjacent to singularities (the order is in the direction of the permutation r of the origami)

  • w and h are width and height of the cylinder

  • bot_twist and top_twist are twist between the minimum square number and the minimum square number adjacent to a singularity.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: r  = '(1,2,3,4)(5,6)'
sage: u0 = '(1,5)(2,6)(3)(4)'
sage: u1 = '(1,5,2,6,3,4)'

sage: Origami(r, u0).cylinder_decomposition()
[([(1, 3)], [(3, 1)], 2, 1, 0, 0),
 ([(3, 1), (5, 5)], [(5, 5), (1, 3)], 4, 1, 0, 0)]
sage: Origami(r, u1).cylinder_decomposition()
[([(1, 6)], [(3, 1)], 2, 1, 0, 0),
 ([(3, 1), (5, 4)], [(1, 6), (5, 4)], 4, 1, 0, 1)]
cylinder_diagram(data=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3218)

Return the cylinder diagram corresponding to the horizontal direction.

If data is True, also return the list of lengths of separatrices, the heights of cylinders and the twists.

The cylinder diagram of a completely periodic surface encodes the combinatorics of cylinders and saddle connections. From a cylinder diagram and its metric data, it is possible to build the surface back.

INPUT:

  • data - boolean (default: False) - if True, return the cylinder diagrams, the lengths, the heights and twists.

EXAMPLES:

sage: from surface_dynamics.all import Origami

The two examples in the stratum H(2):

sage: o1 = Origami('(1,2,3)', '(2,3)')
sage: o1.stratum()
H_2(2)
sage: c1 = o1.cylinder_diagram()
sage: c1.ncyls()
1
sage: c1.nseps()
3
sage: o2 = Origami('(1,2)', '(1,3)')
sage: o2.stratum()
H_2(2)
sage: c2 = o2.cylinder_diagram()
sage: c2.ncyls()
2
sage: c2.nseps()
3

sage: r = (1, 2, 0, 4, 5, 6, 3)
sage: u = (1, 2, 3, 5, 4, 0, 6)
sage: o = Origami(r, u, as_tuple=True)
sage: c, lengths, heights, twists = o.cylinder_diagram(True)
sage: c.cylcoord_to_origami(lengths, heights, twists) == o
True
genus()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3390)

Return the genus of the origami

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)', '(1,3)')
sage: o.genus()
2
sage: o = Origami('(1,2)(3,4)', '(1,3)')
sage: o.genus()
2
gl2z_edges()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1642)

Return the GL(2,ZZ) action on the GL(2,ZZ)-orbit of this origami.

The action is returned as a pair of dictionaries (l_edges, i_edges) representing the action of the following two generators for GL(2, ZZ):

\[\begin{split}L = \begin{pmatrix}1&1\\0&1\end{pmatrix} I = \begin{pmatrix}0&1\\1&0\end{pmatrix} \quad \text{and} \quad\end{split}\]

BEWARE: Do not modify these dictionaries!

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)', '(1,3)')
sage: l, i = o.gl2z_edges()
sage: for oo in sorted(l): print("(%s,%s) -> (%s,%s)" % (oo.r(), oo.u(), l[oo].r(), l[oo].u()))
((2,3),(1,2)) -> ((2,3),(1,2,3))
((2,3),(1,2,3)) -> ((2,3),(1,2))
((1,2,3),(2,3)) -> ((1,2,3),(2,3))
horizontal_cycle_representatives()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1306)

EXAMPLES:

sage: from surface_dynamics import Origami
sage: Origami('(1,2)', '(1,3)').horizontal_cycle_representatives()
[0, 2]
horizontal_symmetry()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1251)

Return the origami (r, u^{-1}).

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2,3)(5,6)', '(1,4)(3,5,7)')
sage: oh = o.horizontal_symmetry(); oh
(1,2,3)(4)(5,6)(7)
(1,4)(2)(3,7,5)(6)

We check that it commutes with vertical symmetry of the cylinder diagram:

sage: ch1 = o.cylinder_diagram().horizontal_symmetry()
sage: ch2 = oh.cylinder_diagram()
sage: ch1.is_isomorphic(ch2)
True
horizontal_twist(width=1, cylinder=None)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1352)

Return the origami (r, r^{-width}u) which is obtained by the action of an horizontal twist of width k on this origami.

INPUT:

  • width - integer (default: 1) - the width of the twist

  • cylinder - integer or None (default None) - if not None performs a twist only in the band of squares that contain i. In that case, it is not a SL(2,R)-deformation.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2,3,4,5,6)', '(1,7)')
sage: o
(1,2,3,4,5,6)(7)
(1,7)(2)(3)(4)(5)(6)
sage: o.horizontal_twist()
(1,2,3,4,5,6)(7)
(1,6,5,4,3,2,7)
sage: o.horizontal_twist(-1)
(1,2,3,4,5,6)(7)
(1,2,3,4,5,6,7)
sage: o.horizontal_twist(6) == o
True

sage: S = SymmetricGroup(5)
sage: r = S('(1,2,3)(4,5)')
sage: u = S('(1,4)(2,5)')
sage: o = Origami(r, u)
sage: o.horizontal_twist(-1, 1) == o.horizontal_twist(2, 1) == Origami(r, S('(1,2,3)') * u)
True
sage: o.horizontal_twist(-2, 1) == o.horizontal_twist(1, 1) == Origami(r, S('(1,2,3)')^2 * u)
True
sage: o.horizontal_twist(-3, 1) == o.horizontal_twist(-3, 1) == o
True
intermediate_covers(degree=None)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2517)

Return the list of intermediate covers of this origami.

If degree is specified, only intermediate covers of given degree are returned.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(3,4,5)', '(1,2,3)(4,6,7)(5,8,9)')
sage: for oo in o.intermediate_covers():
....:    print(oo.nb_squares())
....:    print(oo)
....:    print("- - - - - -")
1
(1)
(1)
- - - - - -
3
(1)(2)(3)
(1,2,3)
- - - - - -
9
(1)(2)(3,4,5)(6)(7)(8)(9)
(1,2,3)(4,6,7)(5,8,9)
- - - - - -
sage: o.intermediate_covers(degree=3)
[(1)(2)(3)
(1,2,3)]
inverse()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1180)

Return the origami (r^{-1}, u^{-1}) which corresponds to the action of -Id on the origami.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2,3)', '(1,2)')
sage: o
(1,2,3)
(1,2)(3)
sage: o.inverse()
(1,3,2)
(1,2)(3)

We check that it commutes with the inverse operation on cylinder diagrams:

sage: o = Origami('(1,2,3)(5,6)', '(1,4)(3,5,7)')
sage: os = o.inverse(); os
(1,3,2)(4)(5,6)(7)
(1,4)(2)(3,7,5)(6)

sage: c = o.cylinder_diagram()
sage: cs1 = o.cylinder_diagram().inverse()
sage: cs2 = os.cylinder_diagram()
sage: cs1.is_isomorphic(cs2)
True
is_connected()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2009)

Check whether the origami is connected or not

It is equivalent to ask whether the group generated by r and u acts transitively on the {1, dots, n}.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.is_connected()
True
sage: o = Origami('(1,2)(3,4)', '(1,2)', check=False)
sage: o.is_connected()
False
is_hyperelliptic(stratum=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2347)

Return True if this origami is hyperelliptic

If stratum is set to True, then return also the corresponding stratum of quadratic differentials this origami is a cover from.

EXAMPLES:

sage: from surface_dynamics.all import Origami, origamis

sage: o = Origami('(1,2)', '(1,3)')
sage: o.is_hyperelliptic()
True
sage: o.is_hyperelliptic(stratum=True)
(True, Q_0(1, -1^5))

sage: o = origamis.Podium([3, 3, 2, 1])
sage: o.is_hyperelliptic()
False
sage: o.is_hyperelliptic(stratum=True)
(False, None)
is_isomorphic(other, certificate=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2049)

Isomorphism test

EXAMPLES:

sage: from surface_dynamics import Origami

sage: o1 = Origami('(1,2)', '(1,3)')
sage: o2 = Origami('(1,2)', '(2,3)')
sage: o3 = Origami('(1,3)', '(1,2)')
sage: o1.is_isomorphic(o2) and o2.is_isomorphic(o1)
True
sage: o1.is_isomorphic(o3) and o3.is_isomorphic(o1)
True
sage: o2.is_isomorphic(o3) and o3.is_isomorphic(o2)
True


sage: from surface_dynamics.misc.permutation import perm_conjugate
sage: for a,b in [(o1,o2),(o2,o1),(o1,o3),(o3,o1),(o2,o3)]:
....:     m = a.is_isomorphic(b, certificate=True)[1]
....:     assert perm_conjugate(a.r_tuple(), m) == list(b.r_tuple())
....:     assert perm_conjugate(a.u_tuple(), m) == list(b.u_tuple())
is_normal()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2667)

Tests if this origami is a normal cover of the torus

An origami is normal if the subgroup of F_2 that defines the cover is normal. It is equivalent to say that the order of the automorphism group equals the number of squares or that the automorphism group acts transitively on the squares.

EXAMPLES:

sage: from surface_dynamics.all import Origami, origamis

sage: o = Origami('(1,2)', '(1,3)')
sage: o.is_normal()
False
sage: o.is_regular()
False

sage: o = origamis.Escalator(4)
sage: o.is_normal()
True
sage: o.is_normal() == o.is_regular()
True
is_orientation_cover()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2168)

Return true if the origami is an orientation cover of a quadratic differential.

It is equivalent to say that -1 is in the Veech group.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.is_orientation_cover()
True

sage: r = '(1,2,3,4,5,6,7,8,9,10)'
sage: u = '(3,5,7,9,4,8,10)'
sage: o = Origami(r, u)
sage: o.is_orientation_cover()
False
is_primitive(return_base=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2460)

An origami is primitive if it does not cover an other origami.

An origami is primitive if the action of the monodromy group has no non trivial block.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.is_primitive()
True
sage: o = Origami('(1,2)(3,4)', '(1,3,5,6)(2,4)')
sage: o.is_primitive()
False
is_quadratic_cover()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2168)

Return true if the origami is an orientation cover of a quadratic differential.

It is equivalent to say that -1 is in the Veech group.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.is_orientation_cover()
True

sage: r = '(1,2,3,4,5,6,7,8,9,10)'
sage: u = '(3,5,7,9,4,8,10)'
sage: o = Origami(r, u)
sage: o.is_orientation_cover()
False
is_quasi_primitive()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2484)

An origami is quasi primitive if it is reduced and all intermediate covers are of genus 1.

SEE ALSO:

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)(3,4)', '(1,3)')
sage: o.is_primitive()
False
sage: o.is_quasi_primitive()
True

sage: o = Origami('(1,2,3,4)(5)(6)', '(1,5)(3,6)')
sage: o.is_primitive()
False
sage: o.is_quasi_primitive()
False
is_quasi_regular()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2650)

An origami (r, u) is quasi-regular if the normal closure of the commutator c = rur^{-1}u^{-1} is contained in the automorphism group.

Equivalently, a quasi-regular origami is a translation surface which is a normal cover of a torus ramified over several rational points.

is_reduced()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1012)

Test of reducibility

An origami is reduced, if it is not a ramified cover of a bigger torus with only one ramification point. In other terms, it is equivalent to say that the period of the origami generates ZZ^2.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)', '(1,3)')
sage: o.is_reduced()
True
sage: o = Origami('(1,2,3,4)(5,6)', '(1,5)(2,6)')
sage: o.is_reduced()
False
sage: o = Origami('(1,2)(3,4)', '(1,3,5,6)(2,4)')
sage: o.is_reduced()
False
sage: o = Origami('(1,2,3,4)(5,6)', '(1,5)(2,6)')
sage: o.is_reduced()
False
sage: o = Origami('(1,2,3,4)(5,6)', '(1,5)(2,6)')
sage: o.is_reduced()
False

sage: o = Origami('(1)(2)(3)', '(1,2,3)')
sage: o.is_reduced()
False
is_regular()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2622)

An origami is regular if its automorphism group acts transitively on the squares.

EXAMPLES:

sage: from surface_dynamics.all import Origami, origamis

sage: o = origamis.EierlegendeWollmilchsau()
sage: o.is_regular()
True
sage: o.is_normal()
True
sage: o = Origami('(1,3,2,4,5,6)', '(1,5)')
sage: o.is_regular()
False
sage: o.is_normal()
False
lattice_of_absolute_periods()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 957)

Return (a, t, u) for a standard basis of the lattice of absolute periods.

The corresponding standard basis is ((a, 0), (t, u)).

The lattice of absolute periods of an origami is the sublattice of ZZ^2 generated by the holonomy vectors of loops drawn on it. Any sublattice of ZZ^2 has a standard basis consisting of a horizontal vector (a, 0) and a nonhorizontal vector (t, u), where a, t, u are integers satisfying 0 <= t < a and 0 < t.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)(3,4)', '(2,3)')
sage: o.lattice_of_absolute_periods()
(2, 0, 1)

sage: r = '(1,2)(3,4)(5,6)(7,8,9,10)(11,12)'
sage: u = '(1,3,5,7)(2,4,6,8)(9,11,10,12)'
sage: o = Origami(r, u)
sage: o.lattice_of_absolute_periods()
(2, 1, 2)
lattice_of_periods()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 930)

Return (a, t, u) for a standard basis of the lattice of periods.

The corresponding standard basis is ((a, 0), (t, u)).

The lattice of periods of an origami is the sublattice of ZZ^2 generated by the holonomy vectors of its saddle connections. Any sublattice of ZZ^2 has a standard basis consisting of a horizontal vector (a, 0) and a nonhorizontal vector (t, u), where a, t, u are integers satisfying 0 <= t < a and 0 < t.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.lattice_of_periods()
(1, 0, 1)
sage: r = '(1,2,3,4,5,6,7,8,9)(10,11,12)(13,14,15,16,17,18,19,20,21)'
sage: u = '(1,14,21,19,8,3,10,5,12,4,11,6)(2,15,16,17,18,7)(9,13,20)'
sage: oy = Origami(r, u)
sage: oy.lattice_of_periods()
(3, 2, 1)
lattice_of_quotients(verbose=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2573)

Return the lattice of quotients of this origami.

The set of quotients of an origami contain a maximal element (itself) and a minimal element (the 1-torus). More generally, it is organised as a lattice.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: P = o.normal_cover().lattice_of_quotients(); P
Finite lattice containing 6 elements
sage: for p in P:
....:     print("%d %s" % (p.nb_squares(), p.stratum_component()))
6 H_3(2^2)^odd
3 H_2(2)^hyp
3 H_2(2)^hyp
3 H_2(2)^hyp
2 H_1(0)^hyp
1 H_1(0)^hyp

sage: o = Origami('(1,2)(3,4)', '(1,3)')
sage: o.lattice_of_quotients()
Finite lattice containing 3 elements
lyapunov_exponents_approx(nb_iterations=20480, nb_experiments=4, only_mean=True, seed=None, nb_vectors=None, nb_vectors_p=None, nb_vectors_m=None, involution=None)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1778)

Approximation of Lyapunov exponents of the Kontsevich-Zorich cocycle.

An origami defines a Teichmuller curve in the moduli space of translation surfaces. The Kontsevich-Zorich cocycle above this Teichmuller curve (for the Haar measure) has the following form

\[1 = \lambda_1, \lambda_2, ... \lambda_g, -\lambda_g, ..., -\lambda_2, -\lambda_1 = -1\]

This function return the approximations of lambda_2, lambda_3, …, lambda_g) as a list.

INPUT:

  • nb_iterations - integer (default: 2**17) - the number of iterations performed in the algorithm

  • nb_experiments - integer (default: 4) - the number of experiments to perform

  • only_mean - boolean (default: True) - if True, return the list of mean exponents, otherwise return a list of lists.

  • nb_vectors - integer (default: genus-1) - the number of vectors to consider

  • involution - permutation or boolean - if True or an inan involution for the origami with derivative either 1.

  • nb_vectors_p, nb_vectors_m - if involution is not None, then it will be interpreted as the number of + and - vectors to consider.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)(3,4)(5,6)', '(2,3)(4,5)')
sage: lexp = o.lyapunov_exponents_approx(nb_iterations=2**21)
sage: lexp  # random
[0.6666, 0.3333]
sage: len(lexp)
2

sage: o = Origami('(1,2)(3,4)(5,6)(7,8)(9,10)', '(2,3)(4,5)(6,7)(8,9)')
sage: s = SymmetricGroup(10)('(1,10)(2,9)(3,8)(4,7)(5,6)')
sage: lexp = o.lyapunov_exponents_approx(involution=s, nb_iterations=2**21)
sage: lexp  # random
([0.6000, 0.2000],
 [0.8000, 0.4000])
sage: print(len(lexp[0]), len(lexp[1]))
2 2
mirror()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1281)

Return the origami (u, r)

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2,3,4)', '(1,5)'); o
(1,2,3,4)(5)
(1,5)(2)(3)(4)
sage: o.mirror()
(1,5)(2)(3)(4)
(1,2,3,4)(5)
monodromy(relative=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2774)

Return the monodromy group of the origami.

The monodromy group of an origami is the group generated by the permutations r and u from which it is defined.

INPUT:

  • relative – if True return the monodromy relative to the largest torus over which this origami is a cover (possibly ramified over several points)

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: G = o.monodromy()
sage: G
Permutation Group with generators [(1,2), (1,3)]
sage: G.order()
6
nb_squares()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 627)

Return the number of squares.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.nb_squares()
3
nb_vertices()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3621)

Return the number of singularities of this origami

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.nb_vertices()
1
sage: o = Origami('(1,2,3)(4,5,6)', '(3,4)')
sage: o.nb_vertices()
2
normal_cover()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2758)

Return the normal cover of this origami.

num_cylinders()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3201)

Return the number of cylinders of this origami.

EXAMPLES:

sage: from surface_dynamics.all import Origami, origamis

sage: Origami('(1,2)', '(1,3)').num_cylinders()
2
sage: origamis.Stair(5).num_cylinders()
3
sage: Origami('(1,2)', '(1)(2)').num_cylinders()
1
optimal_degree()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 984)

Return the optimal degree of this origami.

The optimal degree of self is the degree of the map to the largest torus.

Any origami X to T factors as i circ pi_{opt} where i is an isogeny. The optimal degree is the degree of pi_{opt}.

EXAMPLES:

sage: from surface_dynamics.all import Origami, origamis
sage: E = origamis.EierlegendeWollmilchsau()
sage: E.optimal_degree()
2

sage: o = Origami('(1,2)(3,4)', '(2,3)')
sage: o.optimal_degree()
2

sage: o = Origami('(1,2)', '(1,3)')
sage: o.optimal_degree()
3
orientation_data(points=False, verbose=0)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2193)

Return the list of quadratic stratum and ramification data associated to the orientation quotients of this origami. If the origami is primitive, then there is at most one orientation quotient.

Each element of the list is a 3-tuple containing:

  • a quadratic stratum

  • the list of degrees of zeros which are ramified in the covering (consider only integer points)

  • the partition of half-integers points which are mapped to poles (middle of squares, horizontal edges, vertical edges)

INPUT:

  • points - boolean (default: False) - return singularitiy tuples and not only degrees

EXAMPLES:

sage: from surface_dynamics.all import Origami

The stratum H(2) contains two families of primitive origamis for an odd number of squares. Every surface in H(2) is a covering of a quadratic differential in Q(1,-1^5). The ramification data gives an invariant for those families:

sage: o = Origami('(1,2,3,4,5)', '(2,1)')
sage: o.stratum_component()
H_2(2)^hyp
sage: o.orientation_data()
[(Q_0(1, -1^5), (2,), (1, 3, 1))]
sage: o = Origami('(1,2,3)', '(1,4,5)')
sage: o.stratum_component()
H_2(2)^hyp
sage: o.orientation_data()
[(Q_0(1, -1^5), (2, 0, 0), (1, 1, 1))]

sage: o = Origami('(1)(2)(3,4)(5,6,7)', '(1,2,3)(4,5)(6)(7)')
sage: o.stratum_component()
H_3(4)^hyp
sage: o.orientation_data()
[(Q_0(3, -1^7), (0, 4, 0), (3, 1, 1))]
sage: o = Origami('(1)(2)(3)(4,5)(6,7)', '(1,2,3,4)(5,6,7)')
sage: o.stratum_component()
H_3(4)^hyp
sage: o.orientation_data()
[(Q_0(3, -1^7), (4,), (3, 1, 3))]
sage: o = Origami('(1)(2)(3)(4,5)(6,7)', '(1,2,3,4)(5,6)(7)')
sage: o.stratum_component()
H_3(4)^hyp
sage: o.orientation_data()
[(Q_0(3, -1^7), (4,), (5, 1, 1))]

sage: o = Origami('(1,2,4)(3,6,5)', '(1,3)(2,5)(4,6)')
sage: for q, _, _ in sorted(o.orientation_data()): print(q)
Q_1(4, -1^4)
Q_1(4, -1^4)
Q_1(4, -1^4)
Q_0(1^2, -1^6)
period_generators()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 849)

Return a list of periods that generate the lattice of periods.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,3,6)(2,5,7)(4)', '(1,2,4,3,5,6,7)')
sage: sorted(o.period_generators())
[(-1, 2), (0, 1), (0, 2), (1, 0), (1, 0), (2, 0)]

sage: r = '(1,17,6,18,10,8)(11,2,12,5,7,9)(15,16)(3,13)(14,4)'
sage: u = '(1,11,15,3,14,10,7,6,12)(17,2,16,13,4,8,9,18,5)'
sage: o = Origami(r, u)
sage: o.stratum()
H_2(2)
sage: sorted(o.period_generators())
[(-2, 2), (0, 2), (0, 3), (2, 0), (2, 0), (4, 0)]

sage: o = Origami('(1)(2)(3)', '(1,2,3)')
sage: o.period_generators()
[(1, 0), (0, 3)]

sage: o = Origami('(1,2)', '(1,2)')
sage: o.period_generators()
[(2, 0), (-1, 1)]
pgl2z_edges()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1715)

Return the PGL(2,ZZ)-action on the PGL(2,ZZ)-orbit of this origami.

The action is returned as a pair of dictionaries (l, i) representing the projective action of the matrices:

L = [1 1] [0 1]

and

I = [0 1] [1 0]

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)', '(1,3)')
sage: l, i = o.pgl2z_edges()
sage: len(l), len(i)
(3, 3)
plot(**args)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3000)

Plot self.

The positions of each square follow a naive algorithm. If you believe that a better picture exists look at the method .set_positions()

EXAMPLES:

sage: from surface_dynamics.all import origamis
sage: o = origamis.Escalator(3)
sage: o.plot()  # not tested (problem with matplotlib font caches)
Graphics object consisting of 71 graphics primitives
psl2z_edges()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1742)

Return the PSL(2,ZZ)-action on the PSL(2,ZZ)-orbit of this origami.

The action is returned as a triple of dictionaries (l, r, s) representing the projective action of the matrices:

L = [1 1] [0 1]

R = [0 1] [1 1]

and

S = [0 -1] [1 0]

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)', '(1,3)')
sage: l, r, s = o.psl2z_edges()
sage: len(l), len(r), len(s)
(3, 3, 3)
quotient(H=None)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2864)

Return a quotient of self by the group H.

The group H must be a subgroup of the automorphism group of this origami. If H is None, it is set by default to the full automorphism group.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)(3,4)', '(2,3)')
sage: G = o.automorphism_group()
sage: G.order()
2
sage: oo = o.quotient(G)
sage: print(oo)
(1,2)
(1)(2)
sage: oo.genus()
1
r()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 750)

Return the right permutation of the origami.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2,3)', '(1,2)')
sage: o.r()
(1,2,3)
r_inv()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 763)

Return the inverse of the right permutation of the origami.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2,3)', '(1,2)')
sage: o.r_inv()
(1,3,2)
r_inv_tuple()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 681)

Return the inverse of the right permutation as a tuple on {0, …, n-1}

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2,3)', '(1,2)')
sage: o.r_inv_tuple()
(2, 0, 1)
r_orbit(i)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 656)

Return the orbit of i under the r permutation as a list.

EXAMPLES:

sage: from surface_dynamics import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.r_orbit(0)
[0, 1]
sage: o.r_orbit(1)
[1, 0]
sage: o.r_orbit(2)
[2]
r_tuple()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 641)

Return the right permutation of the origami as a tuple on {0, …, n-1}

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.r_tuple()
(1, 0, 2)
reduce()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2384)

Return a reduced origami isomorphic (up to SL(2, QQ) action) to that origami.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)(3,4)', '(1,3,5,6)(2,4)')
sage: o.lattice_of_periods()
(1, 0, 2)
sage: o.reduce()
(1,2)(3)
(1,3)(2)

sage: o = Origami('(1,2)(3,4,5,6)', '(1,3,5)(2,4,6)')
sage: o.lattice_of_periods()
(2, 0, 1)
sage: o.reduce()
(1)(2,3)
(1,2,3)

sage: o = Origami('(1,2)(3,4,5,6)', '(1,3,4,5)(2,6)')
sage: o.lattice_of_periods()
(2, 1, 1)
sage: o.reduce()
(1)(2,3)
(1,2,3)
relabel(return_map=False, inplace=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1139)

Relabel self

INPUT:

  • return_map – return the labelization

  • inplace – modify self, the default is False. It might be dangerous to set it True as an origami is hashable and the hash is modified by this operation.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2,3,4)(5,6)', '(2,5)(3,6)')
sage: o2, p = o.relabel(return_map=True)
sage: old_r = o.r_tuple()
sage: old_u = o.u_tuple()
sage: new_r = o2.r_tuple()
sage: new_u = o2.u_tuple()
sage: all(p[old_r[i]] == new_r[p[i]] for i in range(6))
True
sage: all(p[old_u[i]] == new_u[p[i]] for i in range(6))
True
rename(name)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2767)

set_positions(pos=None)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2949)

Set positions of the squares for plotting.

INPUT:

  • pos - list of pairs of coordinates, one for each square

show()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3111)

Show a picture of this origami.

sl2z_edges()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1682)

Return the action of SL(2, Z) on the SL(2, Z)-orbit of this origami.

The action is returned as a triple of dictionaries (l, r, s) representing the action of the matrices L, R, S, where:

L = [1 1] [0 1]

R = [1 0] [1 1]

S = R * ~L * R = ~L * R * ~L = [0 -1] [1 0]

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1)(2,3)(4,5)', '(1,2,3,4)(5)')
sage: l, r, s = o.sl2z_edges()
sage: len(l), len(r), len(s)
(12, 12, 12)
sage: all(l[s[l[s[l[s[o]]]]]] == s[s[o]] for o in l)
True
stratum(fake_zeros=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3358)

Stratum of this origami.

INPUT:

  • fake_zeros – boolean (default False) whether the stratum should include zeros of order 0 corresponding to the regular corners of the origami.

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.stratum()
H_2(2)

sage: r = '(1,2,3,4,5,6,7,8,9,10,11,12)'
sage: u = '(1,2,3,4,12)(5,7)(6)(8,10)(9)(11)'
sage: o = Origami(r, u)
sage: o.stratum()
H_4(3, 2, 1)
sage: o.stratum(True)
H_4(3, 2, 1, 0^3)
stratum_component(fake_zeros=False, verbose=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2094)

Return the component of stratum this origami belongs to.

INPUT:

  • fake_zeros – (boolean, default False) whether the stratum should include 0 degree order

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.stratum_component()
H_2(2)^hyp

sage: r = '(1)(2)(3,4)(5,6,7,8)(9,10)'
sage: u = '(1,2,3,5,6,10)(4,9)(7,8)'
sage: Origami(r, u).stratum_component()
H_4(2^3)^even
sage: Origami(u, r).stratum_component()
H_4(2^3)^even

sage: r = '(1,2,3,4,5)(6,7,8,9,10)'
sage: u = '(1,6)(2,10)(3,9)(4,8)(5,7)'
sage: o = Origami(r, u)
sage: o.stratum_component()
H_5(4^2)^odd

sage: r = '(1,2,3,4,5,6,7,8,9,10,11,12,13)'
sage: u1 = '(1,2,3,4,13)(5,12)(6,11)(7,10)(8,9)'
sage: u2 = '(1,2,3,4,13)(5,6)(7)(8,9)(10,11)(12)'
sage: u3 = '(1,2,3,4,13)(5,12,10,9,6,11,8,7)'
sage: o1 = Origami(r, u1)
sage: print(o1.stratum_component(), o1.stratum_component(True))
H_5(4^2)^hyp H_5(4^2, 0^3)^hyp
sage: o2 = Origami(r, u2)
sage: print(o2.stratum_component(), o2.stratum_component(True))
H_5(4^2)^odd H_5(4^2, 0^3)^odd
sage: o3 = Origami(r, u3)
sage: print(o3.stratum_component(), o3.stratum_component(True))
H_5(4^2)^even H_5(4^2, 0^3)^even
sum_of_lyapunov_exponents()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3517)

Return the sum of Lyapunov exponents for this origami

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.sum_of_lyapunov_exponents()
4/3
sage: o = Origami('(1,2)(3,4)', '(2,3)')
sage: o.sum_of_lyapunov_exponents()
3/2
teichmueller_curve

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3480)

Return the teichmueller curve of this origami

The result is cached for future usage.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)', '(1,3)')
sage: t = o.teichmueller_curve()
sage: t
Teichmueller curve of the origami
(1)(2,3)
(1,2)(3)
sage: t.sum_of_lyapunov_exponents()
4/3
sage: for o in t.cusp_representatives():
....:     print(o[0])
....:     print(o[1])
(1)(2,3)
(1,2)(3)
2
(1,2,3)
(1)(2,3)
1
to_standard_form(return_map=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1102)

Return an isomorphic origami in standard form.

INPUT:

  • return_map - boolean (default: False) - if True return the associated mapping

EXAMPLES:

sage: from surface_dynamics.all import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: oo, m = o.to_standard_form(return_map=True)
sage: ~m * o.r() * m == oo.r()
True
sage: ~m * o.u() * m == oo.u()
True


sage: o = Origami('(1,2,3)(4,5)(6)(7,8,9)', '(1,6,3,8,2,7,4,9)')
sage: oo, m = o.to_standard_form(return_map=True)
sage: ~m * o.r() * m == oo.r()
True
sage: ~m * o.u() * m == oo.u()
True
translation_group()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 2824)

Return the automorphism group of the origami as a permutation group.

The automorpism group of a translation surface is the set of affine diffeomorphisms which have a trivial linear part. For an origami, it corresponds combinatorially to the centralizer of the group generated by the permutations r and u that define this origami.

EXAMPLES:

sage: from surface_dynamics.all import Origami

The L with 3 squares has no automorphisms:

sage: o = Origami('(1,2)', '(1,3)')
sage: G = o.automorphism_group()
sage: G.order()
1

This 4-square square-tiled surface in H(1, 1) has non trivial automorphism for which the quotient is a torus with two squares:

sage: o = Origami('(1,2)(3,4)', '(2,3)')
sage: G = o.automorphism_group()
sage: G.order()
2
sage: oo = o.quotient(G)
sage: oo
(1,2)
(1)(2)
u()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 776)

Return the up permutation of the origami.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)', '(1,2,3,4)')
sage: o.u()
(1,2,3,4)
u_inv()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 789)

Return the inverse of the up permutation of the origami.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)', '(1,2,3,4)')
sage: o.u_inv()
(1,4,3,2)
u_inv_tuple()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 733)

Return the inverse of the up permutation as a tuple on {0, …, n-1}

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)', '(1,2,3)')
sage: o.u_inv_tuple()
(2, 0, 1)
u_orbit(i)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 708)

Return the orbit of i under the u permutation as a list.

EXAMPLES:

sage: from surface_dynamics import Origami

sage: o = Origami('(1,2)', '(1,3)')
sage: o.u_orbit(0)
[0, 2]
sage: o.u_orbit(1)
[1]
sage: o.u_orbit(2)
[2, 0]
u_tuple()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 694)

Return the up permutation of the origami as a tuple on {0, …, n-1}

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)', '(1,3)')
sage: o.u_tuple()
(2, 1, 0)
veech_group()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3406)

Return the Veech group of this origami.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2)', '(1,3)')
sage: G = o.veech_group()
sage: G
Arithmetic subgroup with permutations of right cosets
 S2=(2,3)
 S3=(1,2,3)
 L=(1,2)
 R=(1,3)

Most geometric information on the quotient of the upper half plane by the Veech group can be recovered from G:

sage: G.index()
3
sage: G.genus()
0
sage: G.ncusps()
2
sage: G.nu2()
1
sage: G.nu3()
0

As well as some arithmetic information:

sage: G.is_congruence()
True
sage: G.generalised_level()
2

Note that the fact of being congruent is rather exceptional:

sage: o = Origami('(1,2,3,4)', '(4,5)')
sage: o.veech_group().is_congruence()
False

sage: o = Origami('(1,2,3,4,5)', '(5,6)')
sage: o.veech_group().is_congruence()
False

sage: o = Origami('(1,2,3,4,5,6)', '(6,7)')
sage: o.veech_group().is_congruence()
False
vertex_degrees(fake_zeros=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3573)

Return the list of degree of the vertices

INPUT:

  • fake_zeros – boolean (default False) whether we return also the 2 pi angle degrees

EXAMPLES:

sage: from surface_dynamics.all import Origami, origamis

sage: o = Origami('(1,2,3,4,5)', '(1,5,3,2,4)')
sage: o.vertex_degrees()
[4]

sage: o = origamis.ProjectiveLine(5)
sage: o.vertex_degrees()
[2, 2]

sage: r = '(1,2,3,4,5,6,7,8,9,10,11,12)'
sage: u = '(1,2,3,4,12)(5,7)(6)(8,10)(9)(11)'
sage: o = Origami(r, u)
sage: o.vertex_degrees()
[3, 2, 1]
sage: o.vertex_degrees(fake_zeros=True)
[3, 2, 1, 0, 0, 0]
vertical_cycle_representatives()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1329)

EXAMPLES:

sage: from surface_dynamics import Origami
sage: Origami('(1,2)', '(1,3)').vertical_cycle_representatives()
[0, 1]
vertical_symmetry()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1221)

Return the origami (r^{-1}, u).

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2,3)(5,6)', '(1,4)(3,5,7)')
sage: ov = o.vertical_symmetry(); ov
(1,3,2)(4)(5,6)(7)
(1,4)(2)(3,5,7)(6)

We check that it commutes with vertical symmetry of the cylinder diagram:

sage: cv1 = o.cylinder_diagram().vertical_symmetry()
sage: cv2 = ov.cylinder_diagram()
sage: cv1.is_isomorphic(cv2)
True
vertical_twist(width=1, cylinder=None)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 1437)

Return the origami (ru^{-k}, u) obtained by the action of a vertical twist of width k on this origami.

INPUT:

  • width - integer (default: 1) - the width of the twist

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: o = Origami('(1,2,3,4)', '(4,5,6,7)')
sage: o.vertical_twist()
(1,2,3,4,7,6,5)
(1)(2)(3)(4,5,6,7)
sage: o.vertical_twist(-1)
(1,2,3,4,5,6,7)
(1)(2)(3)(4,5,6,7)
sage: o.vertical_twist(4) == o
True

sage: S = SymmetricGroup(5)
sage: r = S('(1,4)(2,5)')
sage: u = S('(1,2,3)(4,5)')
sage: o = Origami(r, u)
sage: o.vertical_twist(-1, 1) == o.vertical_twist(2, 1) == Origami(S('(1,2,3)') * r, u)
True
sage: o.vertical_twist(-2, 1) == o.vertical_twist(1, 1) == Origami(S('(1,2,3)')^2 * r, u)
True
sage: o.vertical_twist(-3, 1) == o.vertical_twist(-3, 1) == o
True
vertices(register_automorphism_action=False)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3556)

INPUT:

  • register_automorphism_action - (default is False) whether the action of the automorphism group of the origami is registered on the vertices

EXAMPLES:

sage: from surface_dynamics.all import origamis
sage: o = origamis.Escalator(4)
sage: o.vertices()
[vertex (1, 5), vertex (2, 6), vertex (3, 7), vertex (4, 8)]
widths_and_heights()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 802)

Return the list of widths and heights of cylinder.

EXAMPLES:

sage: from surface_dynamics.all import Origami
sage: Origami('(1,2)', '(1,3)').widths_and_heights()
[(1, 1), (2, 1)]
sage: Origami('(1,2)(3,4)', '(1,3,5)(2,4)').widths_and_heights()
[(1, 1), (2, 2)]
sage: Origami('(1,2)', '(1,3,4)').widths_and_heights()
[(1, 2), (2, 1)]
sage: Origami('(1,2)(3,4)', '(1,3,5,6)(2,4)').widths_and_heights()
[(1, 2), (2, 2)]
class surface_dynamics.flat_surfaces.origamis.origami_dense.PillowcaseCover_dense_pyx

Bases: object

degree()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 4061)

The degree of the covering.

EXAMPLES:

sage: from surface_dynamics.all import PillowcaseCover

sage: p = PillowcaseCover([2, 1, 3, 4], [3, 2, 1, 4], [4, 2, 3, 1])
sage: p.degree()
4
g_tuple(i)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 4034)

Return a tuple.

EXAMPLES:

sage: from surface_dynamics.all import PillowcaseCover

sage: p = PillowcaseCover([2, 1, 3, 4], [3, 2, 1, 4], [4, 2, 3, 1])
sage: p.g_tuple(0)
(1, 0, 2, 3)
sage: p.g_tuple(1)
(2, 1, 0, 3)
sage: p.g_tuple(2)
(3, 1, 2, 0)
sage: p.g_tuple(3)
(3, 0, 1, 2)

sage: p.g_tuple(4)
Traceback (most recent call last):
...
IndexError: the index i (= 4) must be in {0, 1, 2, 3}
nb_pillows()

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3880)

Return the number of pillows.

EXAMPLES:

sage: from surface_dynamics.all import PillowcaseCover
sage: PillowcaseCover([0], [0], [0], [0], as_tuple=True).nb_pillows()
1
sage: PillowcaseCover([2, 1, 3, 4], [1, 3, 2, 4], [4, 2, 3, 1]).nb_pillows()
4
surface_dynamics.flat_surfaces.origamis.origami_dense.lattice(vectors)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 111)

Return a triple (a, t, u) where ((a, 0), (t, u)) is a basis for the integer lattice generated by vectors.

EXAMPLES:

sage: from surface_dynamics.flat_surfaces.origamis.origami_dense import lattice

sage: lattice([(2, 0), (3, 1)])
(2, 1, 1)
sage: lattice([(2, 0), (3, 1), (2, 1)])
(1, 0, 1)
sage: lattice([(1, 2), (3, 6), (5, 10)])
Traceback (most recent call last):
...
ValueError: all input vectors are collinear
surface_dynamics.flat_surfaces.origamis.origami_dense.origami_from_gap_permutations(r, u)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 82)

surface_dynamics.flat_surfaces.origamis.origami_dense.origami_unpickle(r, u, pos, name)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 396)

surface_dynamics.flat_surfaces.origamis.origami_dense.sl2z_orbits(origamis, n, limit)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 3799)

Action of the matrices l, r, s

EXAMPLES:

sage: from surface_dynamics.all import Stratum
sage: from surface_dynamics.flat_surfaces.origamis.origami_dense import sl2z_orbits
sage: C = Stratum([2,2], k=1).odd_component()
sage: origamis = C.origamis(10)
sage: len(origamis)
8955
sage: slorbits = sl2z_orbits(origamis, 10, 0)
sage: len(slorbits)
19
sage: sum(len(o[0]) for o in slorbits)
8955
sage: for n in range(6, 10):
....:     print(len(C.origamis(n)))
69
270
1260
3384
surface_dynamics.flat_surfaces.origamis.origami_dense.sl_orbit_from_gl_orbit(o, L, I)

File: surface_dynamics/flat_surfaces/origamis/origami_dense.pyx (starting at line 319)

Compute the SL(2, Z) orbit of the origami o knowing the action of GL(2, Z).

TODO: this has nothing to do with origamis… but rather with the action of SL(2, Z), GL(2, Z) and PSL(2, Z).

INPUT:

  • o - an origami

  • L - the action of the matrix l

  • I - the action of the matrix i

EXAMPLES:

On the following example, the SL(2, Z) action has two orbits whereas the GL(2, Z) action has only one:

sage: l_edges = {0:1,1:9,2:8,3:0,4:13,5:10,6:5,7:4,8:11,9:3,10:6,11:2,12:7,13:12}
sage: i_edges = {0:10,1:5,2:12,3:4,4:3,5:1,6:11,7:9,8:13,9:7,10:0,11:6,12:2,13:8}
sage: from surface_dynamics.flat_surfaces.origamis.origami_dense import sl_orbit_from_gl_orbit
sage: s0 = sl_orbit_from_gl_orbit(0, l_edges, i_edges)
sage: print(sorted(s0[0].items()))
[(0, 1), (1, 9), (2, 8), (3, 0), (8, 11), (9, 3), (11, 2)]
sage: print(sorted(s0[1].items()))
[(0, 11), (1, 0), (2, 9), (3, 8), (8, 2), (9, 3), (11, 1)]
sage: print(sorted(s0[2].items()))
[(0, 2), (1, 8), (2, 0), (3, 9), (8, 1), (9, 3), (11, 11)]
sage: s4 = sl_orbit_from_gl_orbit(4, l_edges, i_edges)
sage: print(sorted(s4[0].items()))
[(4, 13), (5, 10), (6, 5), (7, 4), (10, 6), (12, 7), (13, 12)]
sage: print(sorted(s4[1].items()))
[(4, 10), (5, 7), (6, 12), (7, 4), (10, 5), (12, 13), (13, 6)]
sage: print(sorted(s4[2].items()))
[(4, 7), (5, 13), (6, 6), (7, 4), (10, 12), (12, 10), (13, 5)]
sage: s1 = sl_orbit_from_gl_orbit(1, l_edges, i_edges)
sage: s0 == s1, s0 == s4
(True, False)

Teichmueller curves of Origamis.

class surface_dynamics.flat_surfaces.origamis.teichmueller_curve.Cusp(parent, origami, slope)[source]

Bases: SageObject

A cusp in a Teichmueller curve.

  • width of the cusp

  • cylinder decomposition

  • lengths

  • heights

  • a representative

cylinder_diagram(data=False)[source]

Return the cylinder diagram associated to this cusp.

slope()[source]

Return one slope

width()[source]

Return the width of the cusp

class surface_dynamics.flat_surfaces.origamis.teichmueller_curve.TeichmuellerCurve[source]

Bases: SageObject

surface_dynamics.flat_surfaces.origamis.teichmueller_curve.TeichmuellerCurveOfOrigami(origami)[source]

Return the teichmueller curve for an origami

class surface_dynamics.flat_surfaces.origamis.teichmueller_curve.TeichmuellerCurveOfOrigami_class(mapping, inv_mapping, l_edges, r_edges, s2_edges, s3_edges)[source]

Bases: TeichmuellerCurve

an_element()

Returns an origami in this Teichmueller curve.

Beware that the labels of the initial origami may have changed.

EXAMPLES:

sage: from surface_dynamics import *

sage: o = Origami('(1,2)','(1,3)')
sage: t = o.teichmueller_curve()
sage: t.origami()
(1)(2,3)
(1,2)(3)
cusp_representative_iterator()[source]

Iterator over the cusp of self.

Each term is a couple (o,w) where o is a representative of the cusp (an origami) and w is the width of the cusp (an integer).

cusp_representatives()[source]
orbit_graph(s2_edges=True, s3_edges=True, l_edges=False, r_edges=False, vertex_labels=True)[source]

Return the graph of action of PSL on this origami

INPUT:

  • return_map - return the list of origamis in the orbit

origami()[source]

Returns an origami in this Teichmueller curve.

Beware that the labels of the initial origami may have changed.

EXAMPLES:

sage: from surface_dynamics import *

sage: o = Origami('(1,2)','(1,3)')
sage: t = o.teichmueller_curve()
sage: t.origami()
(1)(2,3)
(1,2)(3)
plot_graph()[source]

Plot the graph of the veech group.

The graph corresponds to the action of the generators on the cosets determined by the Veech group in PSL(2,ZZ).

simplicial_action_generators()[source]

Return action of generators of the Veech group on homology

stratum()[source]

Returns the stratum of the Teichmueller curve

EXAMPLES:

sage: from surface_dynamics import *

sage: o = Origami('(1,2)','(1,3)')
sage: t = o.teichmueller_curve()
sage: t.stratum()
H_2(2)
sum_of_lyapunov_exponents()[source]

Returns the sum of Lyapunov exponents for this origami

EXAMPLES:

sage: from surface_dynamics import *

Let us consider few examples in H(2) for which the sum is independent of the origami:

sage: o = Origami('(1,2)','(1,3)')
sage: o.stratum()
H_2(2)
sage: o.sum_of_lyapunov_exponents()
4/3
sage: o = Origami('(1,2,3)(4,5,6)','(1,5,7)(2,6)(3,4)')
sage: o.stratum()
H_2(2)
sage: o.sum_of_lyapunov_exponents()
4/3

This is true as well for the stratum H(1,1):

sage: o = Origami('(1,2)','(1,3)(2,4)')
sage: o.stratum()
H_2(1^2)
sage: o.sum_of_lyapunov_exponents()
3/2
sage: o = Origami('(1,2,3,4,5,6,7)','(2,6)(3,7)')
sage: o.stratum()
H_2(1^2)
sage: o.sum_of_lyapunov_exponents()
3/2

ALGORITHM:

Kontsevich-Zorich formula

veech_group()[source]

Return the veech group of the Teichmueller curve

surface_dynamics.flat_surfaces.origamis.teichmueller_curve.TeichmuellerCurvesOfOrigamis(origamis, assume_normal_form=False, limit=0, verbose_level=0)[source]

Return a set of Teichmueller curve from a set of origamis

INPUT:

  • origamis - iterable of origamis

  • assume_normal_form - whether or not assume that each origami is in normal form

  • limit - integer (default: 0) - if positive, then stop the computation if the size of an orbit is bigger than limit.

Database

Database of reduced origamis.

Origamis are cover of the torus over one point. The group SL(2,ZZ) act on them (as their mapping class group). This file contains the implementation of a database of these SL(2,Z)-orbits. To do concrete computations with origamis see ?.

The database of origamis contain information about SL(2,ZZ) orbits of the origamis that are arithmetic Teichmueller curves.

EXAMPLES:

sage: from surface_dynamics import *

The most useful way to explore the database through queries. We develop several examples below and refer to the documentation in the method OrigamiDatabase.query() for the complete specifications.

Let us start by finding the two exceptional surfaces whose Teichmueller curves have a complete degenerate Lyapunov spectrum: the Eierlegende Wollmilchsau and the Ornythorinque:

sage: D = OrigamiDatabase()

sage: q = D.query(sum_of_L_exp=1)
sage: q.number_of()
2
sage: o0, o1 = q.list()
sage: o0
(1,2,3,4)(5,6,7,8)
(1,5,3,7)(2,8,4,6)
sage: o1
(1,2,3,4,5,6)(7,8,9,10,11,12)
(1,7,5,9,3,11)(2,8,4,12,6,10)
sage: o0.is_isomorphic(origamis.CyclicCover([1,1,1,1]))
True
sage: o1.is_isomorphic(origamis.CyclicCover([1,1,1,3]))
True

We show two methods to get information on a given query: number_of() and list(). To simply display the result on the screen in a table, one can use show() of OrigamiQuery:

sage: q.show()
Origami
--------------------------------
r=1234 5678  u=1537 2846
r=123456 789abc  u=17593b 284c6a

We can display much more information by modifying the displayed columns:

sage: q.cols('nb_squares', 'stratum', 'component', 'regular', 'quasi_regular')
sage: q.show()
Nb squares  Stratum   Comp.  Regular   Quasi regular
----------------------------------------------------
8           H_3(1^4)  c      True      True
12          H_4(2^3)  even   False     True

And the complete list of columns in the database is obtained through:

sage: D.cols()
['representative',
 'stratum',
 'component',
 'primitive',
 'quasi_primitive',
 'orientation_cover',
 'hyperelliptic',
 'regular',
 'quasi_regular',
...
 'relative_monodromy_nilpotent',
 'relative_monodromy_gap_primitive_id',
 'orientation_stratum',
 'orientation_genus',
 'pole_partition',
 'automorphism_group_order',
 'automorphism_group_name']

Now, we do some simple counting of the number of primitive arithmetic Teichmueller curves. Let us first check the classification of arithmetic Teichmueller curves in H(2) from Hubert-Lelievre and McMullen:

sage: A = Stratum([2], k=1)
sage: for n in range(3, 17):
....:     q = D.query(stratum=A, nb_squares=n)
....:     print("%2d %d"%(n, q.number_of()))
 3 1
 4 1
 5 2
 6 1
 7 2
 8 1
 9 2
10 1
11 2
12 1
13 2
14 1
15 2
16 1

And look at the conjecture of Delecroix-Lelievre in the stratum H(1,1):

sage: A = Stratum([1,1], k=1)
sage: for n in range(4,20):
....:     q = D.query(stratum=A, nb_squares=n, primitive=True)
....:     print("%2d %d"%(n, q.number_of()))
 4 1
 5 1
 6 2
 7 2
 8 2
 9 2
10 2
11 2
12 2
13 2
14 2
15 2
16 2
17 2
18 2
19 2

You can get an overview of the content of the database:

sage: D.info()
genus 2
=======
 H_2(2)^hyp        : 205 T. curves (up to 139 squares)
 H_2(1^2)^hyp      : 452 T. curves (up to 80 squares)
<BLANKLINE>
genus 3
=======
 H_3(4)^hyp        : 163 T. curves (up to 51 squares)
 H_3(4)^odd        : 118 T. curves (up to 41 squares)
 ...
genus 6
=======
 H_6(10)^hyp       :  46 T. curves (up to 15 squares)
 H_6(10)^odd       :   4 T. curves (up to 11 squares)
 H_6(10)^even      :  33 T. curves (up to 12 squares)
<BLANKLINE>
<BLANKLINE>
Total: 4688 Teichmueller curves

Here is a last example of the list of regular origamis (i.e. such that their group of translation acts transitively on the set of squares):

sage: q = D.query(regular=True)
sage: q.cols('nb_squares', 'stratum', 'automorphism_group_name')
sage: q.show()
Nb squares  Stratum   Automorphism
----------------------------------
6           H_3(2^2)  S3
8           H_3(1^4)  D8
8           H_3(1^4)  Q8
10          H_5(4^2)  D10
12          H_4(1^6)  A4

AUTHOR:

  • Vincent Delecroix (2011-2014): initial version

class surface_dynamics.flat_surfaces.origamis.origami_database.EnhancedSQLQuery(query_string, database=None)[source]

Bases: SQLQuery

surface_dynamics.flat_surfaces.origamis.origami_database.L_exp_approx_to_data(t)

Convert a tuple of real numbers with same precision into a string.

The output string is a list of numbers written in base 36 (0, 1, …, 9, a, b, …, z) separated by space ‘ ‘. The first number is the precision of the real field. Then each real number consists of three numbers as sign, mantissa, exponent.

EXAMPLES:

sage: from surface_dynamics import *

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: odb.real_tuple_to_data((1.23,-4.0))
'1h 1 1ijk9vqiyzy -1g -1 18ce53un18g -1e'

We may check that the first part consists of the precision:

sage: Integer('1h', 36)
53
sage: RR.precision()
53

And then of the two real numbers we input:

sage: sign = Integer('1', 36)
sage: mantissa = Integer('1ijk9vqiyzy', 36)
sage: exponent = Integer('-1g', 36)
sage: RR(sign * mantissa * 2 ** exponent)
1.23000000000000

sage: sign = Integer('-1', 36)
sage: mantissa = Integer('18ce53un18g', 36)
sage: exponent = Integer('-1e', 36)
sage: RR(sign * mantissa * 2 ** exponent)
-4.00000000000000
class surface_dynamics.flat_surfaces.origamis.origami_database.OrigamiDatabase(dblocation=None, read_only=True, force_creation=False, old_version=False)[source]

Bases: SQLDatabase

Database of arithmetic Teichmueller curves.

EXAMPLES:

sage: from surface_dynamics import *

To query the database the main method is meth:query:

sage: D = OrigamiDatabase()
sage: q = D.query(genus=3, nb_squares=12)
sage: q.number_of()
146
sage: l = q.list()
sage: o = l[0]
sage: o.genus()
3
sage: o.nb_squares()
12

For the precise usage of the method query, see the documentation of that function. Here is an example to show how to found the Eierlegende Wollmilchsau in the database from its property of complete degenerate spectrum:

sage: q = D.query(sum_of_L_exp=1, stratum=Stratum([1,1,1,1],k=1))
sage: len(q)
1
sage: o = q.list()[0]
sage: o.is_isomorphic(origamis.EierlegendeWollmilchsau())
True

We check the classification of arithmetic Teichmueller curves in H(2) from Hubert-Lelievre and McMullen:

sage: A = Stratum([2], k=1)
sage: for n in range(3, 15):
....:     q = D.query(stratum=A, nb_squares=n)
....:     print("%2d %d"%(n, q.number_of()))
 3 1
 4 1
 5 2
 6 1
 7 2
 8 1
 9 2
10 1
11 2
12 1
13 2
14 1

To know all entries of the database look at meth:cols, meth:info or meth:help. The latter also provides a summary description of the columns.

build(comp, N, force_computation=False, verbose=False)[source]

Update the database for the given component comp up to N squares.

Note that in order to work the database should not be in read only mode.

INPUT:

  • comp - stratum or component of stratum

  • N - integer

  • force_computation - force computation of data which may be yet in the database.

  • verbose - boolean - if True, print useful interactive information during the process.

EXAMPLES:

sage: from surface_dynamics import *

sage: import os
sage: import tempfile
sage: with tempfile.TemporaryDirectory() as tmpdir:
....:     db_name = os.path.join(tmpdir, 'my_db.db')
....:     D = OrigamiDatabase(db_name, read_only=False)
....:     D.build(Stratum([4],k=1).odd_component(), 8)  # optional: gap_packages
....:     D.info()                                       # optional: gap_packages
genus 2
=======
<BLANKLINE>
genus 3
=======
 H_3(4)^odd   :  8 T. curves (up to  7 squares)
<BLANKLINE>
<BLANKLINE>
 Total: 8 Teichmueller curves
cols()[source]

Returns the skeleton of self (which is the list of possible entries).

EXAMPLES:

sage: from surface_dynamics import *

sage: O = OrigamiDatabase()
sage: cols = O.cols()
sage: "representative" in cols
True
sage: len(cols)
45
help(cols=None)[source]

Print some help relative to the columns of the database.

If cols is provided then gives help only for these columns.

info(genus=None, dimension=None, print_all=False)[source]

Print the list of connected components and the number of squares up to which the database is filled.

INPUT:

  • genus - integer (default: None) - if not None, print only info for that genus.

  • dimension - Integer (default: None) - if not None, print only info for that dimension.

  • print_all - boolean (default: False) - print also the components for which nothing has been computed yet.

EXAMPLES:

sage: from surface_dynamics import *

sage: O = OrigamiDatabase()
sage: O.info()
genus 2
=======
...
max_nb_squares(comp=None)[source]

Returns the maximum number of squares for which Teichmueller curves have been computed.

If comp is None (default), then returns the biggest integer for which the database contains all data up to that integer. If comp is a stratum or a component of stratum, then returns the maximum number of squares computed for that stratum.

EXAMPLES:

sage: from surface_dynamics import *

sage: O = OrigamiDatabase()
sage: O.max_nb_squares(Stratum([2],k=1))
139
sage: O.max_nb_squares()
11
query(*query_list, **kwds)[source]

From a list of restriction, returns a list of possible entry in the database. By default, returns only the found origamis.

Where to find the possible entries self.cols() then a sign among ‘=’ (equality), ‘<>’ (difference), ‘<’, ‘>’, ‘<=’, ‘>=’ (comparisons).

EXAMPLES:

sage: from surface_dynamics import *

sage: D = OrigamiDatabase()
sage: for o in D.query(stratum=Stratum([1,1], k=1), nb_squares=6):
....:     print("%s\n---------------" % o)
(1)(2)(3,4,5,6)
(1,2,3)(4,5,6)
---------------
(1)(2)(3)(4,5,6)
(1,2,3,4)(5,6)
---------------
(1)(2)(3,4)(5)(6)
(1,2,3)(4,5,6)
---------------
(1)(2)(3)(4,5)(6)
(1,2,3,4)(5,6)
---------------
(1)(2,3)(4,5)(6)
(1,2,4)(3,5,6)
---------------
rebuild(q=None, local_data=False, lyapunov_exponents=False, global_data=False, nb_iterations=4096, nb_experiments=5, verbose=False)[source]

Rebuild some of the data for the origami in the query q.

INPUT:

  • q - a query

  • local_data - boolean - whether or not we rebuild local data.

  • lyapunov_exponents - boolean - whether or not rebuild lyapunov exponents approximation (time consuming).

  • global_data - boolean - whether or not rebuild global data (time and memory consuming).

  • nb_experiments, nb_iterations - integers - option for the computation of Lyapunov exponents.

  • verbose - boolean - if True displays nice information in real time.

update(other, replace=False, verbose=False)[source]

Update the content of this database with the content of another one.

INPUT:

  • other - a query, an origami database or a path to an origami database

  • replace - boolean - whether or not replace entries which are yet in the database.

  • verbose - boolean - if True, displays information during the transfer.

EXAMPLES:

sage: from surface_dynamics import *

sage: import os
sage: import tempfile
sage: with tempfile.TemporaryDirectory() as tmpdir:
....:     db1_name = os.path.join(tmpdir, 'the_first_one.db')
....:     db2_name = os.path.join(tmpdir, 'the_second_one.db')
....:     D1 = OrigamiDatabase(db1_name, read_only=False)
....:     D2 = OrigamiDatabase(db2_name, read_only=False)
....:     D1.build(Stratum([1,1],k=1).unique_component(), 7)  # optional: gap_packages
....:     D2.build(Stratum([2],k=1).unique_component(), 5)    # optional: gap_packages
....:     D2.update(D1)                                        # optional: gap_packages
....:     D2.info()                                            # optional: gap_packages
genus 2
=======
 H_2(2)^hyp  :   2 T. curves (up to  4 squares)
 H_2(1^2)^hyp:   8 T. curves (up to  6 squares)
<BLANKLINE>
<BLANKLINE>
  Total: 10 Teichmueller curves
class surface_dynamics.flat_surfaces.origamis.origami_database.OrigamiQuery(db, query_string, cols, **kwds)[source]

Bases: object

Origami database query.

A query for an instance of OrigamiDatabase. This class nicely wraps the SQLQuery class located in surface_dynamics.databases.database.py to make the query constraints intuitive and with as many pre-definitions as possible. (i.e.: since it has to be a OrigamiDatabase, we already know the table structure and types; and since it is immutable, we can treat these as a guarantee).

cols(*cols)[source]

Get or modify columns.

In a sql query, it corresponds to the clause ‘SELECT’.

EXAMPLES:

sage: from surface_dynamics import *

sage: O = OrigamiDatabase()
sage: q = O.query()
sage: q.cols()
['representative']
sage: q.cols("nb_squares")
sage: q.cols()
['nb_squares']
sage: q.cols("veech_group_index","primitive")
sage: q.cols()
['veech_group_index', 'primitive']
database()[source]

Returns the database of that query.

EXAMPLES:

sage: from surface_dynamics import *

sage: D = OrigamiDatabase()
sage: q = D.query(stratum=Stratum([6], k=1))
sage: q.database()
Database of origamis
sage: q.database() is D
True
dict()[source]

Returns a list of dictionaries: col -> value.

EXAMPLES:

sage: from surface_dynamics import *

sage: D = OrigamiDatabase()
sage: q = D.query(stratum=Stratum([1,1],k=1), nb_squares=8)
sage: q.cols('teich_curve_genus')
sage: q.dict()
[{'teich_curve_genus': 1},
 {'teich_curve_genus': 1},
 {'teich_curve_genus': 0},
 {'teich_curve_genus': 0}]

sage: q.cols('pole_partition', 'primitive')
sage: q.dict()
[{'pole_partition': (0, 2, 2, 2), 'primitive': True},
 {'pole_partition': (2, 0, 2, 2), 'primitive': True},
 {'pole_partition': (0, 0, 2, 4), 'primitive': False},
 {'pole_partition': (0, 0, 2, 4), 'primitive': False}]
get_query_string()[source]

Output the query string in sql format.

This is the method where the attribute of this object are translated into a sql query.

EXAMPLES:

sage: from surface_dynamics import *

sage: D = OrigamiDatabase()
sage: q = D.query(('stratum','=',Stratum([1,1], k=1)), ('nb_squares','<', 13))
sage: q.get_query_string()
"SELECT representative FROM origamis WHERE stratum='1 1' AND nb_squares<13 ORDER BY nb_squares ASC"
sage: q.order(("nb_squares",1),("pole_partition",-1))
sage: q.get_query_string()
"SELECT representative FROM origamis WHERE stratum='1 1' AND nb_squares<13 ORDER BY nb_squares ASC,pole_partition DESC"
list()[source]

Returns the list of entries of the query.

The output is either a list of objects if there is only one column for that query. Otherwise, it is a list of lists where each item is the entries of columns.

See also meth:dict to get a dictionary output.

EXAMPLES:

sage: from surface_dynamics import *

sage: S = OrigamiDatabase(read_only=False)
sage: q = S.query(('stratum','=',Stratum([1,1],k=1)),('nb_squares','=',6))
sage: q.list()
[(1)(2)(3,4,5,6)
(1,2,3)(4,5,6), (1)(2)(3)(4,5,6)
(1,2,3,4)(5,6), (1)(2)(3,4)(5)(6)
(1,2,3)(4,5,6), (1)(2)(3)(4,5)(6)
(1,2,3,4)(5,6), (1)(2,3)(4,5)(6)
(1,2,4)(3,5,6)]
number_of()[source]

Returns the number of entries in the database that satisfy the query.

order(*args)[source]

Get or modify order.

Order should be a list (col_name, +1) or (col_name, -1). First one means ascending and the second one descending. In a sql query, it corresponds to the clause ‘ORDER BY’.

EXAMPLES:

sage: from surface_dynamics import *

sage: O = OrigamiDatabase()
sage: q = O.query()
sage: q.order(("nb_squares",1),("pole_partition",1))
sage: q.get_query_string()
'SELECT representative FROM origamis ORDER BY nb_squares ASC,pole_partition ASC'
sage: q.order(("nb_squares",-1))
sage: q.get_query_string()
'SELECT representative FROM origamis ORDER BY nb_squares DESC'
show(**opts)[source]

Output a text array with the results of that query.

EXAMPLES:

sage: from surface_dynamics import *

There is a problem with show method the SQLQuery:

sage: O = OrigamiDatabase()
sage: q = O.query(("nb_squares","=",6))
sage: q.cols(("stratum","sum_of_L_exp","L_exp_approx"))
sage: q.show()
Stratum    Sum of L exp  L exp approx
---------------------------------------------------------------
...
sql_query()[source]

Returns the SQLQuery.

surface_dynamics.flat_surfaces.origamis.origami_database.are_skeleton_equal(sk1, sk2)[source]

Test whether the two skeleton sk1 and sk2 are equals.

EXAMPLES:

sage: from surface_dynamics.flat_surfaces.origamis.origami_database import are_skeleton_equal
sage: sk1 = {'table1': {'col1': {'sql': 'TEXT', 'unique': True}, 'col2': {'sql': 'INTEGER'}}}
sage: sk2 = {'table1': {'col1': {'sql': 'TEXT', 'unique': True}, 'col2': {'sql': 'INTEGER', 'unique': False}}}
sage: are_skeleton_equal(sk1,sk2)
True
surface_dynamics.flat_surfaces.origamis.origami_database.build_global_data(o, c=None)[source]

Compute the global data that are obtained from the Teichmueller curve of the origami o. If the Teichmueller curve c is not provided, then it is recomputed from scratch (and may be long).

EXAMPLES:

sage: from surface_dynamics import *

sage: o = Origami('(1,2)', '(1,3)')
sage: from surface_dynamics.flat_surfaces.origamis.origami_database import build_global_data
sage: data = build_global_data(o)
sage: for item in sorted(data): print("%-25s %s" % (item, data[item]))
max_hom_dim               2
max_nb_of_cyls            2
min_hom_dim               1
min_nb_of_cyls            1
minus_identity_invariant  True
sum_of_L_exp              4/3
teich_curve_genus         0
teich_curve_ncusps        2
teich_curve_nu2           1
teich_curve_nu3           0
veech_group_congruence    True
veech_group_index         3
veech_group_level         2
surface_dynamics.flat_surfaces.origamis.origami_database.build_local_data(o)[source]

Build local data for the origami o.

The output is a dictionary that is intended to be used to feed the database of origamis. The local data are geometrical aspects (stratum, genus, …) the monodromy group (primitivity, orientation cover, …) and the automorphisms of o.

See also build_lyapunov_exponents() to construct Lyapunov exponents and build_global_data().

EXAMPLES:

sage: from surface_dynamics import *

sage: o = Origami('(1,2)','(1,3)')
sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: data = odb.build_local_data(o)  # optional: gap_packages
sage: data['stratum']                 # optional: gap_packages
H_2(2)
sage: data['nb_squares']              # optional: gap_packages
3
surface_dynamics.flat_surfaces.origamis.origami_database.build_lyapunov_exponents(o, nb_iterations=65536, nb_experiments=10)[source]

Compute the lyapunov exponents for the origami o and update the database.

EXAMPLES:

sage: from surface_dynamics import *

sage: o = Origami('(1,2)', '(1,3)')
sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: data = odb.build_lyapunov_exponents(o)
sage: data       # abs tol .1
{'L_exp_approx': [0.333348091]}
surface_dynamics.flat_surfaces.origamis.origami_database.data_to_L_exp_approx(s)

Convert a string into a tuple of real numbers.

For the encoding convention, see meth:real_tuple_to_data.

EXAMPLES:

sage: from surface_dynamics import *

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: R = RealField(22)
sage: t = (R(1), R(pi))
sage: s = odb.real_tuple_to_data(t)
sage: tt = odb.data_to_real_tuple(s)
sage: tt == t
True
sage: tt[0].parent()
Real Field with 22 bits of precision
sage: tt[1].parent()
Real Field with 22 bits of precision
surface_dynamics.flat_surfaces.origamis.origami_database.data_to_integer_tuple(s)[source]

Convert a string into a tuple of integers.

For the encoding convention, see meth:integer_tuple_to_data.

EXAMPLES:

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: t = (-12, 351234123, 45)
sage: s = odb.integer_tuple_to_data(t)
sage: odb.data_to_integer_tuple(s) == t
True
surface_dynamics.flat_surfaces.origamis.origami_database.data_to_nb_cyls_spectrum(s)

Convert a string into a tuple of integers.

For the encoding convention, see meth:integer_tuple_to_data.

EXAMPLES:

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: t = (-12, 351234123, 45)
sage: s = odb.integer_tuple_to_data(t)
sage: odb.data_to_integer_tuple(s) == t
True
surface_dynamics.flat_surfaces.origamis.origami_database.data_to_orientation_stratum(s)[source]
surface_dynamics.flat_surfaces.origamis.origami_database.data_to_pole_partition(s)[source]
surface_dynamics.flat_surfaces.origamis.origami_database.data_to_real_tuple(s)[source]

Convert a string into a tuple of real numbers.

For the encoding convention, see meth:real_tuple_to_data.

EXAMPLES:

sage: from surface_dynamics import *

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: R = RealField(22)
sage: t = (R(1), R(pi))
sage: s = odb.real_tuple_to_data(t)
sage: tt = odb.data_to_real_tuple(s)
sage: tt == t
True
sage: tt[0].parent()
Real Field with 22 bits of precision
sage: tt[1].parent()
Real Field with 22 bits of precision
surface_dynamics.flat_surfaces.origamis.origami_database.data_to_representative(s)[source]

Convert data to representative.

For encoding convention, see meth:representative_to_data.

surface_dynamics.flat_surfaces.origamis.origami_database.data_to_small_positive_integer_tuple(s)[source]

Convert a string into a tuple of Integer.

For encoding convention, see meth:small_positive_integer_tuple_to_data.

EXAMPLES:

sage: from surface_dynamics import *

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: odb.data_to_small_positive_integer_tuple('14m')
(1, 4, 22)
sage: odb.data_to_small_positive_integer_tuple('')
()
surface_dynamics.flat_surfaces.origamis.origami_database.data_to_stratum(s)[source]

Convert a string into a stratum.

For encoding convention, see meth:stratum_to_data.

EXAMPLES:

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: odb.data_to_stratum('f f 2')
H_17(15^2, 2)
surface_dynamics.flat_surfaces.origamis.origami_database.format_pole_partition(p)[source]

Format into a nice readable string the pole partition.

EXAMPLES:

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: s = odb.pole_partition_to_data((0,3,5,1))
sage: odb.format_pole_partition(s)
'0 (3,5,1)'
surface_dynamics.flat_surfaces.origamis.origami_database.format_representative(s)[source]

Convert a string that encodes an origami into a human readable string.

The output form consists of the concatenation of the cycles without the parenthesis. In other words the permutation (1,5)(2)(3,7,4) will be convert into ‘15 374’.

EXAMPLES:

sage: from surface_dynamics import *

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: o = Origami('(1,2)','(1,3)')
sage: s = odb.representative_to_data(o)
sage: odb.format_representative(s)
'r=12  u=13'
surface_dynamics.flat_surfaces.origamis.origami_database.integer_tuple_to_data(t)[source]

Convert a tuple of arbitrary integers into a string.

The encoding consists in a string that are representation of the integers in base 36 separated by space.

EXAMPLES:

sage: from surface_dynamics import *

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: t = (0,-12435123,3)
sage: s = odb.integer_tuple_to_data(t)
sage: s
'0 -7ej03 3'
sage: list(map(lambda x: Integer(x,36), s.split(' ')))
[0, -12435123, 3]
surface_dynamics.flat_surfaces.origamis.origami_database.nb_cyls_spectrum_to_data(t)

Convert a tuple of arbitrary integers into a string.

The encoding consists in a string that are representation of the integers in base 36 separated by space.

EXAMPLES:

sage: from surface_dynamics import *

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: t = (0,-12435123,3)
sage: s = odb.integer_tuple_to_data(t)
sage: s
'0 -7ej03 3'
sage: list(map(lambda x: Integer(x,36), s.split(' ')))
[0, -12435123, 3]
surface_dynamics.flat_surfaces.origamis.origami_database.orientation_stratum_to_data(q)[source]
surface_dynamics.flat_surfaces.origamis.origami_database.pole_partition_to_data(t)

Convert a tuple of integers between 0 and 35 to a string.

The encoding consists of a string of the same length as t where each character is the representation of the number in base 36 (0, 1, …, 9, a, b, …, z).

EXAMPLES:

sage: from surface_dynamics import *

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: t = (0, 35, 25, 2, 12)
sage: s = odb.small_positive_integer_tuple_to_data(t)
sage: s
'0zp2c'
sage: list(map(lambda x: Integer(x, 36), s))
[0, 35, 25, 2, 12]
surface_dynamics.flat_surfaces.origamis.origami_database.real_tuple_to_data(t)[source]

Convert a tuple of real numbers with same precision into a string.

The output string is a list of numbers written in base 36 (0, 1, …, 9, a, b, …, z) separated by space ‘ ‘. The first number is the precision of the real field. Then each real number consists of three numbers as sign, mantissa, exponent.

EXAMPLES:

sage: from surface_dynamics import *

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: odb.real_tuple_to_data((1.23,-4.0))
'1h 1 1ijk9vqiyzy -1g -1 18ce53un18g -1e'

We may check that the first part consists of the precision:

sage: Integer('1h', 36)
53
sage: RR.precision()
53

And then of the two real numbers we input:

sage: sign = Integer('1', 36)
sage: mantissa = Integer('1ijk9vqiyzy', 36)
sage: exponent = Integer('-1g', 36)
sage: RR(sign * mantissa * 2 ** exponent)
1.23000000000000

sage: sign = Integer('-1', 36)
sage: mantissa = Integer('18ce53un18g', 36)
sage: exponent = Integer('-1e', 36)
sage: RR(sign * mantissa * 2 ** exponent)
-4.00000000000000
surface_dynamics.flat_surfaces.origamis.origami_database.representative_to_data(o)[source]

The maximum number of squares is 95. The encoding consists of the concatenation of the two permutations that define the origami.

surface_dynamics.flat_surfaces.origamis.origami_database.small_positive_integer_tuple_to_data(t)[source]

Convert a tuple of integers between 0 and 35 to a string.

The encoding consists of a string of the same length as t where each character is the representation of the number in base 36 (0, 1, …, 9, a, b, …, z).

EXAMPLES:

sage: from surface_dynamics import *

sage: import surface_dynamics.flat_surfaces.origamis.origami_database as odb
sage: t = (0, 35, 25, 2, 12)
sage: s = odb.small_positive_integer_tuple_to_data(t)
sage: s
'0zp2c'
sage: list(map(lambda x: Integer(x, 36), s))
[0, 35, 25, 2, 12]
surface_dynamics.flat_surfaces.origamis.origami_database.square_num_to_str(i)[source]
surface_dynamics.flat_surfaces.origamis.origami_database.stratum_to_data(h)[source]

Encode the stratum entry.