cohomology

Absolute and relative (simplicial) cohomology of surfaces.

EXAMPLES:

The absolute cohomology of the regular octagon:

sage: from flatsurf import translation_surfaces, SimplicialCohomology
sage: S = translation_surfaces.regular_octagon()
sage: H = SimplicialCohomology(S)

A basis of cohomology:

sage: H.gens()
[{B[(0, 1)]: 1}, {B[(0, 2)]: 1}, {B[(0, 3)]: 1}, {B[(0, 0)]: 1}]

The absolute cohomology of the unfolding of the (3, 4, 13) triangle:

sage: from flatsurf import Polygon, similarity_surfaces
sage: P = Polygon(angles=[3, 4, 13])
sage: S = similarity_surfaces.billiard(P).minimal_cover(cover_type="translation")
sage: H = SimplicialCohomology(S)
sage: len(H.gens())
16

The relative cohomology, relative to the vertices:

sage: S = S.erase_marked_points()  # optional: pyflatsurf  # random output due to deprecation warnings
sage: H = SimplicialCohomology(S, relative=S.vertices())  # optional: pyflatsurf
sage: len(H.gens())  # optional: pyflatsurf
17
flatsurf.geometry.cohomology.SimplicialCohomology(surface, k=1, coefficients=None, relative=None, implementation='dual', category=None)[source]

Return the k-th simplicial cohomology group of surface.

INPUT:

  • surface – a surface

  • k – an integer (default: 1)

  • coefficients – a ring (default: the reals); consider cohomology with coefficients in this ring

  • relative – a set (default: the empty set); if non-empty, then relative cohomology with respect to this set is constructed.

  • implementation – a string (default: "dual"); the algorithm used to compute the cohomology groups. Currently only "dual" is supported, i.e., the groups are computed as duals of the generic homology groups from SageMath.

  • category – a category; if not specified, a category for the cohomology group is chosen automatically depending on coefficients.

class flatsurf.geometry.cohomology.SimplicialCohomologyClass(parent, values)[source]

A cohomology class.

INPUT:

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialCohomology
sage: S = translation_surfaces.regular_octagon()
sage: H = SimplicialCohomology(S)

sage: f, _, _, _ = H.gens()

sage: from flatsurf.geometry.cohomology import SimplicialCohomologyClass
sage: isinstance(f, SimplicialCohomologyClass)
True
class flatsurf.geometry.cohomology.SimplicialCohomologyGroup(surface, k, coefficients, relative, implementation, category)[source]

The k-th simplicial cohomology group of the surface with coefficients.

INPUT:

  • surface – a finite type surface without boundary

  • k – an integer

  • coefficients – a ring

  • relative – a subset of points of the surface

  • implementation – a string; the algorithm used to compute the cohomology, only "dual is supported at the moment.

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialCohomology
sage: T = translation_surfaces.square_torus()

sage: SimplicialCohomology(T)
H¹(Translation Surface in H_1(0) built from a square)
Element

alias of SimplicialCohomologyClass

gens()[source]

Return generators of this cohomology.

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialCohomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialCohomology(T)
sage: H.gens()
[{B[(0, 1)]: 1}, {B[(0, 0)]: 1}]
homology()[source]

Return the homology of the underlying space (with integer coefficients).

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialCohomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialCohomology(T)
sage: H.homology()
H₁(Translation Surface in H_1(0) built from a square)
is_absolute()[source]

Return whether this is an absolute cohomology (and not a relative one).

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialCohomology
sage: T = translation_surfaces.square_torus()

sage: H = SimplicialCohomology(T)
sage: H.is_absolute()
True
surface()[source]

Return the surface over which this cohomology is defined.

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialCohomology
sage: T = translation_surfaces.square_torus()

sage: H = SimplicialCohomology(T)
sage: H.surface() is T
True