conversion

Conversion of sage-flatsurf objects to libflatsurf/pyflatsurf and vice versa.

Ideally, there should be no need to call the functionality in this module directly. Interaction with libflatsurf/pyflatsurf should be handled transparently by the sage-flatsurf objects. Even for authors of sage-flatsurf there should essentially never be a need to use this module directly since objects should provide a pyflatsurf() method that returns a Conversion to libflatsurf/pyflatsurf.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion  # random output due to deprecation warnings in cppyy
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
sage: conversion  # optional: pyflatsurf
Conversion from Triangulation of Translation Surface in H_2(2) built from 2 regular pentagons to FlatTriangulationCombinatorial(...) with vectors ...
class flatsurf.geometry.pyflatsurf.conversion.Conversion(domain, codomain)[source]

Generic base class for a conversion from sage-flatsurf to pyflatsurf.

INPUT:

  • domain – the domain of this conversion, can be None, when the domain cannot be represented that easily with a sage-flatsurf object.

  • codomain – the codomain of this conversion, can be None, when the codomain cannot be represented that easily with libflatsurf/pyflatsurf object.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
__call__(x)[source]

Return the conversion at an element of domain() and return the corresponding pyflatsurf object.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: from flatsurf.geometry.surface_objects import SurfacePoint
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf

sage: p = SurfacePoint(S, (0, 2), (0, 1/2))
sage: conversion(p)  # optional: pyflatsurf
((-1/4*a^2 + 1/2*a + 1/2 ~ 0.54654802), (1/4*a^2 - 3/4 ~ 0.15450850), (1/4 ~ 0.25000000)) in (-6, 8, 9)
__dict__ = mappingproxy({'__module__': 'flatsurf.geometry.pyflatsurf.conversion', '__doc__': '\n    Generic base class for a conversion from sage-flatsurf to pyflatsurf.\n\n    INPUT:\n\n    - ``domain`` -- the domain of this conversion, can be ``None``, when the\n      domain cannot be represented that easily with a sage-flatsurf object.\n\n    - ``codomain`` -- the codomain of this conversion, can be ``None``, when\n      the codomain cannot be represented that easily with\n      libflatsurf/pyflatsurf object.\n\n    EXAMPLES::\n\n        sage: from flatsurf import translation_surfaces\n        sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion\n        sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()\n        sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf\n\n    TESTS::\n\n        sage: from flatsurf.geometry.pyflatsurf.conversion import Conversion\n        sage: isinstance(conversion, Conversion)  # optional: pyflatsurf\n        True\n\n    ', '__init__': <function Conversion.__init__>, 'to_pyflatsurf': <classmethod(<function Conversion.to_pyflatsurf>)>, 'from_pyflatsurf': <classmethod(<function Conversion.from_pyflatsurf>)>, 'to_pyflatsurf_from_elements': <classmethod(<function Conversion.to_pyflatsurf_from_elements>)>, 'domain': <function Conversion.domain>, 'codomain': <function Conversion.codomain>, '__call__': <function Conversion.__call__>, 'section': <function Conversion.section>, '__repr__': <function Conversion.__repr__>, '__eq__': <function Conversion.__eq__>, '__hash__': <function Conversion.__hash__>, '__dict__': <attribute '__dict__' of 'Conversion' objects>, '__weakref__': <attribute '__weakref__' of 'Conversion' objects>, '__annotations__': {}})
__eq__(other)[source]

Return whether this conversion is indistinguishable from other.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: from flatsurf.geometry.surface_objects import SurfacePoint
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion1 = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
sage: conversion2 = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf

sage: conversion1 == conversion2  # optional: pyflatsurf
True
__hash__()[source]

Return a hash value for this conversion that is compatible with __eq__().

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: from flatsurf.geometry.surface_objects import SurfacePoint
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion1 = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
sage: conversion2 = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf

sage: hash(conversion1) == hash(conversion2)  # optional: pyflatsurf
Traceback (most recent call last):
...
TypeError: unhashable type: 'FlatTriangulationConversion'
__init__(domain, codomain)[source]
__module__ = 'flatsurf.geometry.pyflatsurf.conversion'
__repr__()[source]

Return a printable representation of this conversion.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: from flatsurf.geometry.surface_objects import SurfacePoint
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
Conversion from Triangulation of Translation Surface in H_2(2) built from 2 regular pentagons to FlatTriangulationCombinatorial(...) with vectors ...
__weakref__

list of weak references to the object

codomain()[source]

Return the codomain of this conversion, a pyflatsurf object.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
sage: conversion.codomain()  # optional: pyflatsurf
FlatTriangulationCombinatorial(vertices = ...) with vectors {...}
domain()[source]

Return the domain of this conversion, a sage-flatsurf object.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
sage: conversion.domain()  # optional: pyflatsurf
Triangulation of Translation Surface in H_2(2) built from 2 regular pentagons
classmethod from_pyflatsurf(codomain, domain=None)[source]

Return a Conversion from domain to codomain.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
sage: FlatTriangulationConversion.from_pyflatsurf(conversion.codomain())  # optional: pyflatsurf
Conversion from ...
section(y)[source]

Return the conversion of an element of codomain() and return the corresponding sage-flatsurf object.

This is the inverse of __call__().

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: from flatsurf.geometry.surface_objects import SurfacePoint
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf

sage: p = SurfacePoint(S, (0, 2), (0, 1/2))
sage: q = conversion(p)  # optional: pyflatsurf
sage: conversion.section(q)  # optional: pyflatsurf
Point (0, 1/2) of polygon (0, 2)
classmethod to_pyflatsurf(domain, codomain=None)[source]

Return a Conversion from domain to the codomain.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
classmethod to_pyflatsurf_from_elements(elements, codomain=None)[source]

Return a Conversion that converts the sage-flatsurf elements to codomain.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: conversion = RingConversion.to_pyflatsurf_from_elements([1, 2, 3])  # optional: gmpxxyy  # random output due to cppyy deprecation warnings
sage: conversion  # optional: gmpxxyy
Conversion from Integer Ring to __gmp_expr<__mpz_struct[1],__mpz_struct[1]>
class flatsurf.geometry.pyflatsurf.conversion.FlatTriangulationConversion(domain, codomain, label_to_half_edge)[source]

Converts a sage-flatsurf surface to a FlatTriangulation object and vice-versa.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
__annotations__ = {}
__call__(x)[source]

Return the image of x under this conversion.

INPUT:

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: from flatsurf.geometry.surface_objects import SurfacePoint
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf

We map a point:

sage: p = SurfacePoint(S, (0, 2), (0, 1/2))
sage: conversion(p)  # optional: pyflatsurf
((-1/4*a^2 + 1/2*a + 1/2 ~ 0.54654802), (1/4*a^2 - 3/4 ~ 0.15450850), (1/4 ~ 0.25000000)) in (-6, 8, 9)

We map a half edge:

sage: conversion(((0, 0), 0))  # optional: pyflatsurf
1
__eq__(other)[source]

Return whether this conversion is indistinguishable from other.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: from flatsurf.geometry.surface_objects import SurfacePoint
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion1 = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
sage: conversion2 = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf

sage: conversion1 == conversion2  # optional: pyflatsurf
True
__hash__ = None
__init__(domain, codomain, label_to_half_edge)[source]

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: S = translation_surfaces.square_torus().triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf

sage: isinstance(conversion, FlatTriangulationConversion)  # optional: pyflatsurf
True
__module__ = 'flatsurf.geometry.pyflatsurf.conversion'
classmethod from_pyflatsurf(codomain, domain=None)[source]

Return a Conversion from domain to codomain.

INPUT:

  • codomain – a FlatTriangulation

  • domain – a sage-flatsurf surface or None (default: None); if None, the corresponding surface is constructed.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
sage: FlatTriangulationConversion.from_pyflatsurf(conversion.codomain())  # optional: pyflatsurf
Conversion from ...
ring_conversion()[source]

Return the conversion that maps the base ring of the domain of this conversion to the base ring of the codomain of this conversion.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: from flatsurf.geometry.surface_objects import SurfacePoint
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
sage: conversion.ring_conversion()  # optional: pyflatsurf
Conversion from Number Field in a with defining polynomial y^4 - 5*y^2 + 5 with a = 1.902113032590308? to NumberField(a^4 - 5*a^2 + 5, [...])
section(y)[source]

Return the preimage of y under this conversion.

INPUT:

  • y – an object defined in the codomain, e.g., a pyflatsurf Point

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: from flatsurf.geometry.surface_objects import SurfacePoint
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf

We roundtrip a point:

sage: p = SurfacePoint(S, (0, 2), (0, 1/2))
sage: q = conversion(p)  # optional: pyflatsurf
sage: conversion.section(q) == p  # optional: pyflatsurf
True

We roundtrip a half edge:

sage: half_edge = conversion(((0, 0), 0))  # optional: pyflatsurf
sage: conversion.section(half_edge)  # optional: pyflatsurf
((0, 0), 0)
classmethod to_pyflatsurf(domain, codomain=None)[source]

Return a Conversion from domain to the codomain.

INPUT:

  • domain – a sage-flatsurf surface

  • codomain – a FlatTriangulation or None (default: None); if None, the corresponding FlatTriangulation is constructed.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
vector_space_conversion()[source]

Return the conversion maps two-dimensional vectors over the base ring of the domain to Vector<T> for the codomain.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion
sage: from flatsurf.geometry.surface_objects import SurfacePoint
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
sage: conversion.vector_space_conversion()  # optional: pyflatsurf
Conversion from Vector space of dimension 2 over Number Field in a with defining polynomial y^4 - 5*y^2 + 5 with a = 1.902113032590308? to flatsurf::Vector<eantic::renf_elem_class>
class flatsurf.geometry.pyflatsurf.conversion.RingConversion(domain, codomain)[source]

A conversion between a SageMath ring and a C/C++ ring.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: conversion = RingConversion.to_pyflatsurf_from_elements([1, 2, 3])  # optional: gmpxxyy
__annotations__ = {}
__module__ = 'flatsurf.geometry.pyflatsurf.conversion'
classmethod from_pyflatsurf(codomain, domain=None)[source]

Return a RingConversion that maps domain to codomain.

INPUT:

  • codomain – a libflatsurf/pyflatsurf type or ring

  • domain – a SageMath ring, or None (default: None); if None, the ring is determined automatically.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: conversion = RingConversion.to_pyflatsurf(QQ)  # optional: gmpxxyy

sage: RingConversion.from_pyflatsurf(conversion.codomain())  # optional: gmpxxyy
Conversion from Rational Field to __gmp_expr<__mpq_struct[1],__mpq_struct[1]>
classmethod from_pyflatsurf_from_elements(elements, domain=None)[source]

Return a RingConversion that can map domain to the ring of elements.

INPUT:

  • elements – a sequence of ring elements that libflatsurf/pyflatsurf understands

  • domain – a SageMath ring, or None (default: None); if None, the ring is determined automatically.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion

sage: import gmpxxyy  # optional: gmpxxyy
sage: conversion = RingConversion.from_pyflatsurf_from_elements([gmpxxyy.mpz()])  # optional: gmpxxyy
classmethod from_pyflatsurf_from_flat_triangulation(flat_triangulation, domain=None)[source]

Return a RingConversion that can map domain to the ring over which flat_triangulation is defined.

INPUT:

  • flat_triangulation – a libflatsurf FlatTriangulation

  • domain – a SageMath ring, or None (default: None); if None, the ring is determined automatically.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion, RingConversion
sage: S = translation_surfaces.veech_double_n_gon(5).triangulate().codomain()
sage: flat_triangulation = FlatTriangulationConversion.to_pyflatsurf(S).codomain()  # optional: pyflatsurf
sage: conversion = RingConversion.from_pyflatsurf_from_flat_triangulation(flat_triangulation)  # optional: pyflatsurf

Note that this conversion does not roundtrip back to the same SageMath ring. An e-antic ring only has a single variable name but a SageMath number field has a variable name and a (potentially different) variable name in the defining polynomial:

sage: conversion.domain() is S.base_ring()  # optional: pyflatsurf
False
sage: conversion.domain()  # optional: pyflatsurf
Number Field in a with defining polynomial x^4 - 5*x^2 + 5 with a = 1.902113032590308?
sage: S.base_ring()
Number Field in a with defining polynomial y^4 - 5*y^2 + 5 with a = 1.902113032590308?

We can explicitly specify the domain to create the same conversion again:

sage: conversion = RingConversion.from_pyflatsurf_from_flat_triangulation(flat_triangulation, domain=S.base_ring())  # optional: pyflatsurf
sage: conversion.domain() is S.base_ring()  # optional: pyflatsurf
True
classmethod to_pyflatsurf(domain, codomain=None)[source]

Return a RingConversion that converts the SageMath ring domain to something that libflatsurf/pyflatsurf can understand.

INPUT:

  • domain – a ring

  • codomain – a C/C++ type or None (default: None); if None, the corresponding type is constructed.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: RingConversion.to_pyflatsurf(QQ)  # optional: gmpxxyy
Conversion from Rational Field to __gmp_expr<__mpq_struct[1],__mpq_struct[1]>
classmethod to_pyflatsurf_from_elements(elements, codomain=None)[source]

Return a RingConversion than can map from the elements to the codomain.

INPUT:

  • elements – a sequence of SageMath ring elements

  • codomain – a libflatsurf parent ring or None (default: None); if None, the ring is determined automatically.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion

sage: import cppyy  # optional: gmpxxyy
sage: conversion = RingConversion.to_pyflatsurf_from_elements([1])  # optional: gmpxxyy
class flatsurf.geometry.pyflatsurf.conversion.RingConversion_algebraic(domain, codomain)[source]

A conversion from the algebraic numbers in SageMath AA to an e-antic real embedded number field.

EXAMPLES:

There is no general notion of algebraic numbers in libflatsurf yet, so we can deduce such a conversion for a finite number of elements that span a finite number field:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: conversion = RingConversion.to_pyflatsurf(domain=AA)  # optional: pyeantic
Traceback (most recent call last):
...
NotImplementedError: ...

sage: conversion = RingConversion.to_pyflatsurf_from_elements(elements=[AA(1), AA(sqrt(2))])  # optional: pyeantic
__annotations__ = {}
__call__(x)[source]

Return the image of x under this conversion.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: conversion = RingConversion.to_pyflatsurf_from_elements(elements=[AA(sqrt(2))])  # optional: pyeantic

sage: conversion(AA(sqrt(2)))  # optional: pyeantic
(a ~ 1.4142136)
__init__(domain, codomain)[source]
__module__ = 'flatsurf.geometry.pyflatsurf.conversion'
class flatsurf.geometry.pyflatsurf.conversion.RingConversion_eantic(domain, codomain)[source]

A conversion from a SageMath number field to an e-antic real embedded number field.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: conversion = RingConversion.to_pyflatsurf(domain=QuadraticField(2))  # optional: pyeantic

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion_eantic
sage: isinstance(conversion, RingConversion_eantic)  # optional: pyeantic
True
__annotations__ = {}
__call__(x)[source]

Return the image of x under this conversion.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: domain = QuadraticField(2)
sage: conversion = RingConversion.to_pyflatsurf(domain)  # optional: pyeantic
sage: conversion(domain.gen())  # optional: pyeantic
(a ~ 1.4142136)
__module__ = 'flatsurf.geometry.pyflatsurf.conversion'
section(y)[source]

Return the preimage of y under this conversion.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: domain = QuadraticField(2)
sage: conversion = RingConversion.to_pyflatsurf(domain)  # optional: pyeantic
sage: gen = conversion(domain.gen())  # optional: pyeantic
sage: conversion.section(gen)  # optional: pyeantic
a
class flatsurf.geometry.pyflatsurf.conversion.RingConversion_exactreal(domain, codomain)[source]

A conversion from a pyexactreal SageMath ring to a exact-real ring.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: from pyexactreal import ExactReals  # optional: pyexactreal
sage: conversion = RingConversion.to_pyflatsurf(domain=ExactReals(QQ))  # optional: pyexactreal
Traceback (most recent call last):
...
NotImplementedError: cannot deduce the exact real module that corresponds to this generic ring of exact reals since there is no generic exact-real ring without a fixed set of generators in libexactreal yet

sage: from pyexactreal import QQModule, RealNumber  # optional: pyexactreal
sage: M = QQModule(RealNumber.rational(1))  # optional: pyexactreal
sage: conversion = RingConversion.to_pyflatsurf(domain=ExactReals(QQ), codomain=M)  # optional: pyexactreal

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion_exactreal
sage: isinstance(conversion, RingConversion_exactreal)  # optional: pyexactreal
True
__annotations__ = {}
__call__(x)[source]

Return the image of x under this conversion.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion_exactreal
sage: from pyexactreal import QQModule, RealNumber, ExactReals  # optional: pyexactreal
sage: domain = ExactReals(QQ)  # optional: pyexactreal
sage: M = QQModule(RealNumber.rational(1))  # optional: pyexactreal
sage: conversion = RingConversion_exactreal._create_conversion(domain=domain, codomain=M)  # optional: pyexactreal
sage: conversion(domain(1))  # optional: pyexactreal
1
__module__ = 'flatsurf.geometry.pyflatsurf.conversion'
section(y)[source]

Return the preimage of y under this conversion.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion_exactreal
sage: from pyexactreal import QQModule, RealNumber, ExactReals  # optional: pyexactreal
sage: domain = ExactReals(QQ)  # optional: pyexactreal
sage: M = QQModule(RealNumber.random())  # optional: pyexactreal
sage: conversion = RingConversion_exactreal._create_conversion(domain=domain, codomain=M)  # optional: pyexactreal
sage: conversion.section(M.gen(0R))  # optional: pyexactreal
(...)
class flatsurf.geometry.pyflatsurf.conversion.RingConversion_gmp(domain, codomain)[source]

Conversion between SageMath integers and rationals and GMP mpz and mpq.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: conversion = RingConversion.to_pyflatsurf(domain=ZZ)  # optional: gmpxxyy

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion_gmp
sage: isinstance(conversion, RingConversion_gmp)  # optional: gmpxxyy
True
__annotations__ = {}
__call__(x)[source]

Return the image of x under this conversion.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: domain = QQ
sage: conversion = RingConversion.to_pyflatsurf(QQ)  # optional: gmpxxyy
sage: x = 1/3
sage: y = conversion(x)  # optional: gmpxxyy
sage: conversion.section(y)  # optional: gmpxxyy
1/3
__module__ = 'flatsurf.geometry.pyflatsurf.conversion'
section(y)[source]

Return the preimage of y under this conversion.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: domain = QQ
sage: conversion = RingConversion.to_pyflatsurf(domain)  # optional: gmpxxyy
sage: y = conversion(1/3)  # optional: gmpxxyy
sage: conversion.section(y)  # optional: gmpxxyy
1/3
class flatsurf.geometry.pyflatsurf.conversion.RingConversion_int(domain, codomain)[source]

Conversion between SageMath integers and machine long long integers.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: conversion = RingConversion.to_pyflatsurf(domain=int)  # optional: gmpxxyy

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion_int
sage: isinstance(conversion, RingConversion_int)  # optional: gmpxxyy
True
__annotations__ = {}
__call__(x)[source]

Return the image of x under this conversion.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: domain = int
sage: conversion = RingConversion.to_pyflatsurf(int)  # optional: gmpxxyy
sage: conversion(1R)  # optional: gmpxxyy
1
__module__ = 'flatsurf.geometry.pyflatsurf.conversion'
section(y)[source]

Return the preimage of y under this conversion.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import RingConversion
sage: domain = int
sage: conversion = RingConversion.to_pyflatsurf(domain)  # optional: gmpxxyy
sage: y = conversion(3R)  # optional: gmpxxyy
sage: conversion.section(y)  # optional: gmpxxyy
3
class flatsurf.geometry.pyflatsurf.conversion.VectorSpaceConversion(domain, codomain, ring_conversion=None)[source]

Converts vectors in a SageMath vector space into libflatsurf Vector<T>s.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import VectorSpaceConversion
sage: conversion = VectorSpaceConversion.to_pyflatsurf(QQ^2)  # optional: pyflatsurf
__annotations__ = {}
__call__(vector)[source]

Return a Vector<T> corresponding to the SageMath vector.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import VectorSpaceConversion
sage: conversion = VectorSpaceConversion.to_pyflatsurf(QQ^2)  # optional: pyflatsurf
sage: conversion(vector(QQ, [1, 2]))  # optional: pyflatsurf
(1, 2)
__init__(domain, codomain, ring_conversion=None)[source]
__module__ = 'flatsurf.geometry.pyflatsurf.conversion'
classmethod from_pyflatsurf_from_elements(elements, domain=None, ring_conversion=None)[source]

Return a Conversion that converts the pyflatsurf elements into vectors in domain.

INPUT:

  • domain – a SageMath VectorSpace or None; if None, it is determined automatically.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import VectorSpaceConversion
sage: codomain = VectorSpaceConversion.to_pyflatsurf(QQ^2).codomain()  # optional: pyflatsurf

sage: VectorSpaceConversion.from_pyflatsurf_from_elements([codomain()])  # optional: pyflatsurf
Conversion from Vector space of dimension 2 over Rational Field to flatsurf::Vector<__gmp_expr<__mpq_struct[1],__mpq_struct[1]>...>
section(vector)[source]

Return the SageMath vector corresponding to the pyflatsurf vector.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import VectorSpaceConversion
sage: conversion = VectorSpaceConversion.to_pyflatsurf(QQ^2)  # optional: pyflatsurf
sage: v = vector(QQ, [1, 2])  # optional: pyflatsurf
sage: conversion.section(conversion(v)) == v  # optional: pyflatsurf
True
classmethod to_pyflatsurf(domain, codomain=None, ring_conversion=None)[source]

Return a Conversion from domain to codomain.

INPUT:

  • domain – a SageMath free module

  • codomain – a libflatsurf Vector<T> type or None (default: None); if None, the type is determined automatically.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import VectorSpaceConversion
sage: conversion = VectorSpaceConversion.to_pyflatsurf(QQ^2)  # optional: pyflatsurf
classmethod to_pyflatsurf_from_elements(elements, codomain=None)[source]

Return a conversion that can convert the SageMath elements to codomain.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import VectorSpaceConversion
sage: VectorSpaceConversion.to_pyflatsurf_from_elements([vector((1, 2))])  # optional: pyflatsurf
Conversion from Ambient free module of rank 2 over the principal ideal domain Integer Ring to ...
sage: VectorSpaceConversion.to_pyflatsurf_from_elements([vector((1, 2)), vector((1/2, 2/3))])  # optional: pyflatsurf
Conversion from Vector space of dimension 2 over Rational Field to ...
flatsurf.geometry.pyflatsurf.conversion.from_pyflatsurf(T)[source]

Given T a flatsurf::FlatTriangulation from libflatsurf/pyflatsurf, return a sage-flatsurf translation surface.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import to_pyflatsurf, from_pyflatsurf # optional: pyflatsurf
sage: S = translation_surfaces.veech_double_n_gon(5) # optional: pyflatsurf
sage: T = from_pyflatsurf(to_pyflatsurf(S)) # optional: pyflatsurf
sage: T  # optional: pyflatsurf
Translation Surface in H_2(2) built from 6 isosceles triangles
flatsurf.geometry.pyflatsurf.conversion.sage_ring(surface)[source]

Return the SageMath ring over which the pyflatsurf surface surface can be constructed in sage-flatsurf.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: from flatsurf.geometry.pyflatsurf.conversion import to_pyflatsurf, sage_ring # optional: pyflatsurf
sage: S = to_pyflatsurf(translation_surfaces.veech_double_n_gon(5)) # optional: pyflatsurf  # random output due to matplotlib warnings with some combinations of setuptools and matplotlib
sage: sage_ring(S) # optional: pyflatsurf
Number Field in a with defining polynomial x^4 - 5*x^2 + 5 with a = 1.902113032590308?
flatsurf.geometry.pyflatsurf.conversion.to_pyflatsurf(S)[source]

Given S a translation surface from sage-flatsurf return a flatsurf::FlatTriangulation from libflatsurf/pyflatsurf.

flatsurf.geometry.pyflatsurf.conversion.to_sage_ring(x)[source]

Given a coordinate of a flatsurf::Vector, return a SageMath element from which from_pyflatsurf() can eventually construct a translation surface.

EXAMPLES:

sage: from flatsurf.geometry.pyflatsurf.conversion import to_sage_ring  # optional: pyflatsurf
sage: to_sage_ring(1R).parent()  # optional: pyflatsurf
Integer Ring

GMP coordinate types:

sage: import cppyy  # optional: pyflatsurf
sage: import pyeantic  # optional: pyflatsurf
sage: to_sage_ring(cppyy.gbl.mpz_class(1)).parent()  # optional: pyflatsurf
Integer Ring
sage: to_sage_ring(cppyy.gbl.mpq_class(1, 2)).parent()  # optional: pyflatsurf
Rational Field

e-antic coordinate types:

sage: import pyeantic  # optional: pyflatsurf
sage: K = pyeantic.eantic.renf_class.make("a^3 - 3*a + 1", "a", "0.34 +/- 0.01", 64R)  # optional: pyflatsurf
sage: to_sage_ring(K.gen()).parent()  # optional: pyflatsurf
Number Field in a with defining polynomial x^3 - 3*x + 1 with a = 0.3472963553338607?

exact-real coordinate types:

sage: from pyexactreal import QQModule, RealNumber  # optional: pyflatsurf
sage: M = QQModule(RealNumber.random())   # optional: pyflatsurf
sage: to_sage_ring(M.gen(0R)).parent()  # optional: pyflatsurf
Real Numbers as (Rational Field)-Module