gl2r_orbit_closure

GL(2,R)-orbit closure of translation surfaces

EXAMPLES:

Let us first construct a Veech surface in the stratum H(2):

sage: from flatsurf import translation_surfaces
sage: from flatsurf import GL2ROrbitClosure

sage: x = polygen(QQ)
sage: K.<a> = NumberField(x^3 - 2, embedding=AA(2)**(1/3))
sage: S = translation_surfaces.mcmullen_L(1,1,1,a)
sage: O = GL2ROrbitClosure(S) # optional: pyflatsurf  # random output due to matplotlib warnings with some combinations of setuptools and matplotlib
sage: O.decomposition((1,2)).cylinders() # optional: pyflatsurf
[Cylinder with perimeter [...]]

The following is also a Veech surface. However the flow decomposition in directions with long cylinders might not discover them if a limit is set:

sage: S = translation_surfaces.mcmullen_genus2_prototype(4,2,1,1,1/4)
sage: l = S.base_ring().gen()
sage: O = GL2ROrbitClosure(S) # optional: pyflatsurf
sage: dec = O.decomposition((8*l - 25, 16), 9) # optional: pyflatsurf
sage: dec.undeterminedComponents() # optional: pyflatsurf
[Component with perimeter [...]]

Further refinement might change the status of undetermined components:

sage: import pyflatsurf # optional: pyflatsurf
sage: dec.decompose(10r) # optional: pyflatsurf
True
sage: dec.undeterminedComponents() # optional: pyflatsurf
[]

Veech surfaces have the property that any saddle connection direction is parabolic:

sage: S = translation_surfaces.veech_double_n_gon(5)
sage: O = GL2ROrbitClosure(S)  # optional: pyflatsurf
sage: all(d.parabolic() for d in O.decompositions_depth_first(3))  # optional: pyflatsurf
True

For surfaces in rank one loci, even though they are completely periodic, they are generally not parabolic:

sage: S = translation_surfaces.mcmullen_genus2_prototype(4,2,1,1,1/4)
sage: O = GL2ROrbitClosure(S)  # optional: pyflatsurf
sage: all((d.hasCylinder() == False) or d.parabolic() for d in O.decompositions(6))  # optional: pyflatsurf
False
sage: all((d.completelyPeriodic() == True) or (d.hasCylinder() == False) for d in O.decompositions(6))  # optional: pyflatsurf
True
class flatsurf.geometry.gl2r_orbit_closure.GL2ROrbitClosure(surface)[source]

Lower bound approximation to the tangent space of a GL(2,R)-orbit closure of a linear family of translation surfaces.

EXAMPLES:

sage: from flatsurf import polygons, similarity_surfaces
sage: from flatsurf import GL2ROrbitClosure  # optional: pyflatsurf

sage: T = polygons.triangle(3, 3, 5)
sage: S = similarity_surfaces.billiard(T)
sage: S = S.minimal_cover(cover_type="translation")
sage: GL2ROrbitClosure(S)  # optional: pyflatsurf
GL(2,R)-orbit closure of dimension at least 2 in H_5(4, 2^2) (ambient dimension 12)

Computing an orbit closure over an exact real ring with transcendental elements:

sage: from flatsurf import Polygon, EuclideanPolygonsWithAngles
sage: from pyexactreal import ExactReals  # optional: pyexactreal  # random output due to matplotlib warnings with some combinations of setuptools and matplotlib

sage: E = EuclideanPolygonsWithAngles((1, 5, 5, 5))
sage: R = ExactReals(E.base_ring())  # optional: pyexactreal
sage: slopes = E.slopes()
sage: T = Polygon(angles=(1, 5, 5, 5), edges=[slopes[0], R.random_element(1/4) * slopes[1]])  # optional: pyexactreal
sage: S = similarity_surfaces.billiard(T)  # optional: pyexactreal
sage: S = S.minimal_cover(cover_type="translation")  # optional: pyexactreal
sage: O = GL2ROrbitClosure(S); O  # optional: pyflatsurf, optional: pyexactreal
GL(2,R)-orbit closure of dimension at least 4 in H_7(4^3, 0) (ambient dimension 17)
sage: bound = E.billiard_unfolding_stratum('half-translation', marked_points=True).dimension()
sage: for decomposition in O.decompositions(1):  # long time, optional: pyflatsurf, optional: pyexactreal
....:     O.update_tangent_space_from_flow_decomposition(decomposition)
....:     if O.dimension() == bound: break
sage: O  # long time, optional: pyflatsurf, optional: pyexactreal
GL(2,R)-orbit closure of dimension at least 8 in H_7(4^3, 0) (ambient dimension 17)
absolute_dimension()[source]

EXAMPLES:

sage: from flatsurf import polygons, similarity_surfaces
sage: from flatsurf import GL2ROrbitClosure  # optional: pyflatsurf
sage: T = polygons.triangle(1,3,4)  # Veech octagon
sage: S = similarity_surfaces.billiard(T)
sage: S = S.minimal_cover("translation")
sage: O = GL2ROrbitClosure(S)  # optional: pyflatsurf
sage: O.absolute_dimension()  # optional: pyflatsurf
2

The triangular billiard (5,6,7) belongs to the canonical double cover of the stratum Q(5,3,0^3) in genus 3. The orbit is dense and we can check that the absolute dimension is indeed \(6 = 2 rank\):

sage: T = polygons.triangle(5,6,7)
sage: S = similarity_surfaces.billiard(T)
sage: S = S.minimal_cover("translation")
sage: O = GL2ROrbitClosure(S)  # optional: pyflatsurf
sage: for d in O.decompositions(5, 100):  # long time (3s)  # optional: pyflatsurf
....:     O.update_tangent_space_from_flow_decomposition(d)
....:     if O.dimension() == 9:
....:         break
sage: O.absolute_dimension()  # long time (above)  # optional: pyflatsurf
6
absolute_homology()[source]
ambient_stratum()[source]

Return the stratum of Abelian differentials this surface belongs to.

EXAMPLES:

sage: from flatsurf import Polygon, similarity_surfaces
sage: from flatsurf import GL2ROrbitClosure # optional: pyflatsurf
sage: T = Polygon(angles=(1, 3, 5))
sage: S = similarity_surfaces.billiard(T)
sage: S = S.minimal_cover(cover_type="translation")
sage: O = GL2ROrbitClosure(S) # optional: pyflatsurf
sage: O.ambient_stratum() # optional: pyflatsurf
H_3(4, 0^4)
base_ring()[source]

Return the underlying base ring

boundaries()[source]

Return the list of boundaries (ie sum of edges around a triangular face).

These are elements of H_1(S, Sigma; Z).

cylinder_circumference(component, A, sc_index, proj)[source]

Return the circumference of the cylinder component in the homology of the underlying surface.

INPUT:

  • component – a cylinder

  • A, sc_index, proj – the output of flow_decomposition_kontsevich_zorich_cocycle

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf import GL2ROrbitClosure  # optional: pyflatsurf

sage: S = translation_surfaces.veech_double_n_gon(5)
sage: O = GL2ROrbitClosure(S)  # optional: pyflatsurf
sage: for dec in O.decompositions_depth_first(1):  # optional: pyflatsurf
....:     break
sage: c0, c1 = dec.components() # optional: pyflatsurf
sage: kz = O.flow_decomposition_kontsevich_zorich_cocycle(dec) # optional: pyflatsurf
sage: O.cylinder_circumference(c0, *kz) # optional: pyflatsurf
(1, 1, 1, -1)
sage: O.cylinder_circumference(c1, *kz) # optional: pyflatsurf
(0, 0, 1, 0)
cylinder_deformation_subspace(decomposition)[source]

Return a subspace included in the tangent space to the GL(2,R)-orbit closure computed from the flow decomposition decomposition.

From A. Wright cylinder deformation Theorem.

decomposition(v, limit=-1)[source]
decompositions(bound, limit=-1, bfs=False)[source]
decompositions_breadth_first(bound, limit=-1)[source]
decompositions_depth_first(bound, limit=-1)[source]
dimension()[source]

Return the current complex dimension of the GL(2,R)-orbit closure.

Note that this is not the dimension of the orbit closure but only a lower bound. It is always at least 2 (coming from a GL(2,R)-orbit). The current tangent space could be refined via update_tangent_space_from_flow_decomposition().

EXAMPLES:

sage: from flatsurf import Polygon, similarity_surfaces
sage: from flatsurf import GL2ROrbitClosure # optional: pyflatsurf
sage: T = Polygon(angles=(1, 3, 5))
sage: S = similarity_surfaces.billiard(T)
sage: S = S.minimal_cover(cover_type="translation")
sage: O = GL2ROrbitClosure(S) # optional: pyflatsurf
sage: O.dimension() # optional: pyflatsurf
2
sage: bound = T.category().billiard_unfolding_stratum('half-translation', marked_points=True).dimension()
sage: for decomposition in O.decompositions(1):  # long time, optional: pyflatsurf
....:     if O.dimension() == bound: break
....:     O.update_tangent_space_from_flow_decomposition(decomposition)
....:     print(O.dimension())
3
5
5
7
9
10
field_of_definition()[source]

Return the field of definition of the current subspace.

Warning

This involves the computation of the echelon form of the matrix. It might be rather expensive if the computation of the tangent space is not terminated.

EXAMPLES:

sage: from flatsurf import Polygon, similarity_surfaces, EuclideanPolygonsWithAngles
sage: from flatsurf import GL2ROrbitClosure  # optional: pyflatsurf
sage: from pyexactreal import ExactReals  # optional: pyexactreal
sage: E = EuclideanPolygonsWithAngles((1, 5, 5, 5))
sage: R = ExactReals(E.base_ring())  # optional: pyexactreal
sage: slopes = E.slopes()
sage: T = Polygon(angles=(1, 5, 5, 5), edges=[slopes[0], R.random_element(1/4) * slopes[1]])  # optional: pyexactreal
sage: S = similarity_surfaces.billiard(T)  # optional: pyexactreal
sage: S = S.minimal_cover(cover_type="translation")  # optional: pyexactreal
sage: O = GL2ROrbitClosure(S); O  # optional: pyflatsurf, optional: pyexactreal
GL(2,R)-orbit closure of dimension at least 4 in H_7(4^3, 0) (ambient dimension 17)
sage: O.field_of_definition() # optional: pyflatsurf, optional: pyexactreal
Number Field in c0 with defining polynomial x^2 - 2 with c0 = 1.414213562373095?
sage: bound = E.billiard_unfolding_stratum('half-translation', marked_points=True).dimension()
sage: for decomposition in O.decompositions(1):  # long time, optional: pyflatsurf, optional: pyexactreal
....:     if O.dimension() == bound: break
....:     O.update_tangent_space_from_flow_decomposition(decomposition)
sage: O.field_of_definition()  # long time, optional: pyflatsurf, optional: pyexactreal
Rational Field

sage: T = Polygon(angles=(1, 3, 5))
sage: S = similarity_surfaces.billiard(T)
sage: S = S.minimal_cover(cover_type="translation")
sage: O = GL2ROrbitClosure(S) # optional: pyflatsurf
sage: O.field_of_definition() # optional: pyflatsurf
Number Field in c0 with defining polynomial x^3 - 3*x - 1 with c0 = 1.879385241571817?
sage: bound = T.category().billiard_unfolding_stratum('half-translation', marked_points=True).dimension()
sage: for decomposition in O.decompositions(1):  # long time, optional: pyflatsurf
....:     if O.dimension() == bound: break
....:     O.update_tangent_space_from_flow_decomposition(decomposition)
sage: O.field_of_definition()  # long time, optional: pyflatsurf
Rational Field
flow_decomposition_kontsevich_zorich_cocycle(decomposition)[source]

Base change from the homology of decomposition to the underlying surface.

Return a triple A, sc_index, proj where

  • A is a \(d \times d\) matrix where \(d\) is the dimension of the ambient stratum (ie the dimension of the homology relative to the singularities)

  • sc_index is a dictionary whose keys are saddle connections appearing in the various components of this decomposition and the corresponding value is an integral index

  • proj is a projection matrix from vectors expressed as a sum of saddle connections appearing in this component to the chosen representative of cohomology with respect to this component

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf import GL2ROrbitClosure  # optional: pyflatsurf

sage: S = translation_surfaces.square_torus()
sage: O = GL2ROrbitClosure(S)  # optional: pyflatsurf
sage: for dec in O.decompositions_depth_first(3):  # optional: pyflatsurf
....:     kz = O.flow_decomposition_kontsevich_zorich_cocycle(dec) # optional: pyflatsurf
....:     print(kz[0])
[ 0  1]
[-1 -1]
[ 0  1]
[-1 -2]
[-1 -1]
[ 1  0]
[-1 -1]
[ 2  1]
[-2 -1]
[ 3  1]
[-1 -1]
[ 3  2]
[1 0]
[0 1]
[ 1  0]
[-1  1]
holonomy(v)[source]

Return the holonomy of v (with respect to the chosen homology basis)

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf import GL2ROrbitClosure # optional: pyflatsurf
sage: x = polygen(QQ)
sage: K.<a> = NumberField(x^3 - 2, embedding=AA(2)**(1/3))
sage: S = translation_surfaces.mcmullen_L(1,1,1,a)
sage: O = GL2ROrbitClosure(S) # optional: pyflatsurf
sage: edges = O._surface.edges() # optional: pyflatsurf
sage: F = FreeModule(ZZ, len(edges)) # optional: pyflatsurf
sage: all(O.V2(O.holonomy(O.proj * F.gen(i))).vector == O.V2(O._surface.fromHalfEdge(e.positive())).vector for i, e in enumerate(edges)) # optional: pyflatsurf
True
holonomy_dual(v)[source]

Return the holonomy of the dual of v

is_teichmueller_curve(bound, limit=-1)[source]

Return False when the program can find a direction which is either completely periodic with incomensurable moduli or a direction with at least one cylinder and at least one minimal component.

EXAMPLES:

sage: from flatsurf import polygons, similarity_surfaces
sage: from flatsurf import GL2ROrbitClosure  # optional: pyflatsurf
sage: for a in range(1,6):  # long time (8s), optional: pyflatsurf
....:     for b in range(a,6):
....:         for c in range(b,6):
....:             if a + b + c > 7 or gcd([a,b,c]) != 1:
....:                 continue
....:             T = polygons.triangle(a, b, c)
....:             S = similarity_surfaces.billiard(T)
....:             S = S.minimal_cover(cover_type="translation")
....:             O = GL2ROrbitClosure(S)
....:             if O.is_teichmueller_curve(3, 50) != False:
....:                 print(a,b,c)
1 1 1
1 1 2
1 1 4
1 2 2
1 2 3
1 3 3
lift(v)[source]

Given a vector in the “spanning set basis” return a vector on the full basis of edges.

EXAMPLES:

sage: from flatsurf import polygons, translation_surfaces, similarity_surfaces
sage: from flatsurf import GL2ROrbitClosure  # optional: pyflatsurf

sage: S = translation_surfaces.mcmullen_genus2_prototype(4,2,1,1,0)
sage: O = GL2ROrbitClosure(S)  # optional: pyflatsurf
sage: u0,u1 = O.tangent_space_basis()  # optional: pyflatsurf
sage: v0 = O.lift(u0)  # optional: pyflatsurf
sage: v1 = O.lift(u1)  # optional: pyflatsurf
sage: span([v0, v1])  # optional: pyflatsurf
Vector space of degree 9 and dimension 2 over Real Embedded Number Field in l with defining polynomial x^2 - x - 8 with l = 3.372281323269015?
Basis matrix:
[ 1  0  -1   (1/8*l+7/8 ~ 1.2965352)  (-1/8*l+1/8 ~ -0.29653517)  -1   (5/8*l-5/8 ~ 1.4826758)  (-1/2*l+3/2 ~ -0.18614066)  (-5/8*l+13/8 ~ -0.48267583)]
[ 0  1  -1  (1/4*l-1/4 ~ 0.59307033)  (-1/4*l+1/4 ~ -0.59307033)   0  (1/4*l-1/4 ~ 0.59307033)                           0   (-1/4*l+1/4 ~ -0.59307033)]

This can be used to deform the surface:

sage: T = polygons.triangle(3,4,13)
sage: S = similarity_surfaces.billiard(T)
sage: S = S.minimal_cover("translation").erase_marked_points() # long time (3s, #122), optional: pyflatsurf
sage: O = GL2ROrbitClosure(S)  # long time (above), optional: pyflatsurf
sage: for d in O.decompositions(4, 20):  # long time (2s, #124), optional: pyflatsurf
....:     O.update_tangent_space_from_flow_decomposition(d)
....:     if O.dimension() == 4:
....:         break
sage: d1,d2,d3,d4 = [O.lift(b) for b in O.tangent_space_basis()]  # long time (above), optional: pyflatsurf
sage: dreal = d1/132 + d2/227 + d3/1280 - d4/13201  # long time (above), optional: pyflatsurf
sage: dimag = d1/141 - d2/233 + d4/1230 + d4/14250  # long time (above), optional: pyflatsurf
sage: d = [O.V2((x,y)).vector for x,y in zip(dreal,dimag)]  # long time (above), optional: pyflatsurf
sage: S2 = O._surface + d  # long time (6s), optional: pyflatsurf

sage: O2 = GL2ROrbitClosure(S2.surface())  # long time (above), optional: pyflatsurf
sage: for d in O2.decompositions(1, 20):  # long time (25s, #124), optional: pyflatsurf
....:     O2.update_tangent_space_from_flow_decomposition(d)
tangent_space_basis()[source]
update_tangent_space_from_flow_decomposition(decomposition)[source]

Update the current tangent space by using the cylinder deformation vectors from decomposition.

EXAMPLES:

sage: from flatsurf import polygons, similarity_surfaces
sage: from flatsurf import GL2ROrbitClosure  # optional: pyflatsurf

sage: T = polygons.triangle(1, 2, 5)
sage: S = similarity_surfaces.billiard(T)
sage: S = S.minimal_cover(cover_type="translation")
sage: O = GL2ROrbitClosure(S)  # optional: pyflatsurf
sage: for d in O.decompositions(1):  # long time (1s), optional: pyflatsurf
....:     O.update_tangent_space_from_flow_decomposition(d)
sage: assert O.dimension() == 2  # long time (above), optional: pyflatsurf
update_tangent_space_from_vector(v)[source]