homology

Absolute and relative (simplicial) homology of surfaces.

EXAMPLES:

The absolute homology of the regular octagon:

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

A basis of homology, with generators written as (sums of) oriented edges:

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

The absolute homology 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: S.genus()
8
sage: H = SimplicialHomology(S)
sage: len(H.gens())
16

Relative homology, relative to the singularities of the surface:

sage: S = S.erase_marked_points()  # optional: pyflatsurf  # random output due to deprecation warnings
sage: H1 = SimplicialHomology(S, relative=S.vertices())  # optional: pyflatsurf
sage: len(H1.gens())  # optional: pyflatsurf
17

We can also form relative \(H_0\) and \(H_2\), though they are not overly interesting of course:

sage: H0 = SimplicialHomology(S, relative=S.vertices(), k=0)
sage: len(H0.gens())
0

sage: H2 = SimplicialHomology(S, relative=S.vertices(), k=2)
sage: len(H2.gens())
1

We create the homology class corresponding to the core curve of a cylinder (the interface here is terrible at the moment, see https://github.com/flatsurf/sage-flatsurf/issues/166):

sage: from flatsurf import Polygon, similarity_surfaces, SimplicialHomology, GL2ROrbitClosure

sage: P = Polygon(angles=[3, 4, 13])
sage: S = similarity_surfaces.billiard(P).minimal_cover(cover_type="translation").triangulate().codomain()

sage: from flatsurf.geometry.pyflatsurf.conversion import FlatTriangulationConversion  # optional: pyflatsurf
sage: conversion = FlatTriangulationConversion.to_pyflatsurf(S)  # optional: pyflatsurf
sage: T = conversion.codomain()  # optional: pyflatsurf
sage: O = GL2ROrbitClosure(T)  # optional: pyflatsurf

sage: D = O.decomposition((13, 37))  # optional: pyflatsurf
sage: cylinder = D.cylinders()[0]  # optional: pyflatsurf

sage: H = SimplicialHomology(S)
sage: core = sum(int(str(chain[edge])) * H(conversion.section(edge.positive())) for segment in cylinder.right() for chain in [segment.saddleConnection().chain()] for edge in T.edges())  # optional: pyflatsurf
sage: core  # optional: pyflatsurf  # random output, the chosen generators vary between operating systems
972725347814111665129717*B[((0, -1/2*c0, -1/2*c0^2 + 3/2), 2)] + 587352809047576581321682*B[((0, -1/2*c0^2 + 1, -1/2*c0^3 + 3/2*c0), 2)] + 60771110563809382932401*B[((0, -1/2*c0^2 + 1, 1/2*c0^3 - 3/2*c0), 2)] ...
flatsurf.geometry.homology.SimplicialHomology(surface, k=1, coefficients=None, relative=None, implementation='generic', category=None)[source]

Return the k-th simplicial homology group of surface.

INPUT:

  • surface – a surface

  • k – an integer (default: 1)

  • coefficients – a ring (default: the integer ring); consider the homology with coefficients in this ring

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

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

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

class flatsurf.geometry.homology.SimplicialHomologyClass(parent, chain)[source]

An element of a homology group.

INPUT:

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: S = translation_surfaces.regular_octagon()
sage: H0 = SimplicialHomology(S, k=0)
sage: g0 = H0.gens()[0]
sage: g0
B[Vertex 0 of polygon 0]

sage: H1 = SimplicialHomology(S, k=1)
sage: g1 = H1.gens()[0]
sage: g1
B[(0, 1)]

sage: H2 = SimplicialHomology(S, k=2)
sage: g2 = H2.gens()[0]
sage: g2
B[0]
algebraic_intersection(other)[source]

Return the algebraic intersection of this class of a closed curve with other.

INPUT:

EXAMPLES:

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

sage: H((0, 0)).algebraic_intersection(H((0, 1)))
1

sage: a = H((0, 1))
sage: b = 3 * H((0, 0)) + 5 * H((0, 2)) - 2 * H((0, 4))
sage: a.algebraic_intersection(b)
0

sage: a = 2 * H((0, 0)) + H((0, 1)) + 3 * H((0, 2)) + H((0, 3)) + H((0, 4)) + H((0, 5)) + H((0, 7))
sage: b = H((0, 0)) + 2 * H((0, 1)) + H((0, 2)) + H((0, 3)) + 2 * H((0, 4)) + 3 * H((0, 5)) + 4 * H((0, 6)) + 3 * H((0, 7))
sage: a.algebraic_intersection(b)
-6

sage: S = translation_surfaces.cathedral(1, 4)
sage: H = SimplicialHomology(S)
sage: a = H((0, 3))
sage: b = H((2, 1))
sage: a.algebraic_intersection(b)
0

sage: a = H((0, 3))
sage: b = H((3, 4)) + 3 * H((0, 3)) + 2 * H((0, 0)) - H((1, 7)) + 7 * H((2, 1)) - 2 * H((2, 2))
sage: a.algebraic_intersection(b)
2
coefficient(gen)[source]

Return the multiplicity of this class at a generator of homology.

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T)
sage: a,b = H.gens()
sage: a.coefficient(a)
1
sage: a.coefficient(b)
0
coefficients()[source]

Return the coefficients of this element in terms of the generators of homology.

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T)
sage: H.gens()[0].coefficients()
(1, 0)
sage: H.gens()[1].coefficients()
(0, 1)
surface()[source]

Return the surface on which this homology class in defined.

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T)
sage: h = H.gens()[0]
sage: h.surface()
Translation Surface in H_1(0) built from a square
class flatsurf.geometry.homology.SimplicialHomologyGroup(surface, k, coefficients, relative, implementation, category)[source]

The k-th simplicial homology 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 homology, only "generic" is supported at the moment which uses the generic homology machinery of SageMath.

EXAMPLES:

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

Surfaces must be immutable to compute their homology:

sage: T = MutableOrientedSimilaritySurface.from_surface(T)
sage: SimplicialHomology(T)
Traceback (most recent call last):
...
ValueError: surface must be immutable to compute homology
sage: T.set_immutable()
sage: SimplicialHomology(T)
H₁(Translation Surface in H_1(0) built from a square)
Element

alias of SimplicialHomologyClass

boundary(chain)[source]

Return the boundary of chain as an element of the chain_module() in lower dimension.

INPUT:

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T, k=0)
sage: c = H.chain_module().an_element(); c
2*B[Vertex 0 of polygon 0]
sage: H.boundary(c)
0
sage: H = SimplicialHomology(T, k=1)
sage: c = H.chain_module().an_element(); c
2*B[(0, 0)] + 2*B[(0, 1)]
sage: H.boundary(c)
0
sage: H = SimplicialHomology(T, k=2)
sage: c = H.chain_module().an_element(); c
2*B[0]
sage: H.boundary(c)
0
chain_module()[source]

Return the free module of simplicial chains of the surface i.e., formal sums of simplicies, e.g., formal sums of edges of the surface.

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T)
sage: H.chain_module()
Free module generated by {(0, 1), (0, 0)} over Integer Ring
change(k=None)[source]

Return this homology but in dimension k.

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T)
sage: H
H₁(Translation Surface in H_1(0) built from a square)

sage: H.change(k=0)
H₀(Translation Surface in H_1(0) built from a square)
degree()[source]

Return the degree \(k\) for this homology \(H_k\).

EXAMPLES:

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

sage: H = SimplicialHomology(T)
sage: H.degree()
1
gens()[source]

Return generators of homology.

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T)
sage: H.gens()
(B[(0, 1)], B[(0, 0)])
sage: H = SimplicialHomology(T, 0)
sage: H.gens()
(B[Vertex 0 of polygon 0],)
sage: H = SimplicialHomology(T, 2)
sage: H.gens()
(B[0],)
hom(f, codomain=None)[source]

Return the homomorphism of homology induced by f.

INPUT:

  • f – a morphism of surfaces or a matrix

  • codomain – the simplicial homology this morphism maps into

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: S = translation_surfaces.square_torus()
sage: A = S.affine_automorphism_group()
sage: M = matrix([[1, 2], [0, 1]])
sage: f = A.derivative().section()(M, check=False)

sage: from flatsurf import SimplicialHomology
sage: H = SimplicialHomology(S)
sage: g = H.hom(f)

sage: g.matrix()  # optional: pyflatsurf
[1 0]
[2 1]

sage: H.gens()
(B[(0, 1)], B[(0, 0)])
sage: [g(h) for h in H.gens()]  # optional: pyflatsurf
[2*B[(0, 0)] + B[(0, 1)], B[(0, 0)]]
is_absolute()[source]

Return whether this is absolute homology (and not relative to some set of points.)

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T)
sage: H.is_absolute()
True
ngens()[source]

Return the Betti number of this homology.

EXAMPLES:

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

sage: H = T.homology()
sage: H.ngens()
2
simplices()[source]

Return the simplices that form the generators of chain_module().

EXAMPLES:

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

In dimension 1, this is the set of edges:

sage: H = SimplicialHomology(T)
sage: H.simplices()
((0, 1), (0, 0))

In dimension 0, this is the set of vertices:

sage: H = SimplicialHomology(T, k=0)
sage: H.simplices()
(Vertex 0 of polygon 0,)

In dimension 2, this is the set of polygons:

sage: H = SimplicialHomology(T, k=2)
sage: H.simplices()
(0,)

In all other dimensions, there are no simplices:

sage: H = SimplicialHomology(T, k=12)
sage: H.simplices()
()
some_elements()[source]

Return some typical homology classes (for testing.)

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T)
sage: H.some_elements()
[0, B[(0, 1)], B[(0, 0)]]
surface()[source]

Return the surface of which this is the homology.

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T)
sage: H.surface() == T
True
symplectic_basis()[source]

Return a symplectic basis of generators of this homology group.

A symplectic basis is a basis of the form \((e_1, \ldots, e_n, f_1, \ldots, f_n)\) such that for the SimplicialHomologyClass.algebraic_intersection(), \(e_i \cdot e_j = f_i \cdot f_j = 0\) and \(e_i \cdot f_j = \delta_{ij}\) for all \(i\) and \(j\).

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T)
sage: H.symplectic_basis()
[B[(0, 0)], B[(0, 1)]]
zero()[source]

Return the zero homology class.

EXAMPLES:

sage: from flatsurf import translation_surfaces, SimplicialHomology
sage: T = translation_surfaces.square_torus()
sage: H = SimplicialHomology(T)
sage: H.zero()
0
class flatsurf.geometry.homology.SimplicialHomologyMorphismSpace(domain, codomain, category=None)[source]

The space of homomorphisms from the homology domain to codomain.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: S = translation_surfaces.square_torus()
sage: A = S.affine_automorphism_group()
sage: M = matrix([[1, 0], [0, 1]])
sage: f = A.derivative().section()(M, check=False)

sage: from flatsurf import SimplicialHomology
sage: H = SimplicialHomology(S)
sage: g = H.hom(f)
sage: G = g.parent()

Since these are homomorphisms in homology, they preserve the linear structure of the homology:

sage: g.category()
Category of endsets of modules over Integer Ring
an_element()[source]

Return some homomorphism in homology.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: S = translation_surfaces.square_torus()
sage: T = translation_surfaces.mcmullen_L(1, 1, 1, 1)

sage: End(S.homology()).an_element()
Generic endomorphism of H₁(Translation Surface in H_1(0) built from a square)
  Defn: [1 0]
        [0 1]

sage: Hom(S.homology(), T.homology()).an_element()
Generic morphism:
  From: H₁(Translation Surface in H_1(0) built from a square)
  To:   H₁(Translation Surface in H_2(2) built from 3 squares)
  Defn: [0 0]
        [0 0]
        [0 0]
        [0 0]
base_ring()[source]

Return the ring over which these homomorphisms are defined.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: S = translation_surfaces.square_torus()

sage: End(S.homology()).base_ring()
Integer Ring
identity()[source]

Return the identity homomorphism in this space (if it exists.)

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: S = translation_surfaces.square_torus()
sage: T = translation_surfaces.mcmullen_L(1, 1, 1, 1)

sage: End(S.homology()).identity()
Generic endomorphism of H₁(Translation Surface in H_1(0) built from a square)
  Defn: [1 0]
        [0 1]
zero()[source]

Return the zero homomorphism.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: S = translation_surfaces.square_torus()

sage: End(S.homology()).zero()
Generic endomorphism of H₁(Translation Surface in H_1(0) built from a square)
  Defn: [0 0]
        [0 0]
class flatsurf.geometry.homology.SimplicialHomologyMorphism_base[source]

Base class for all homomorphisms in homology.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: S = translation_surfaces.square_torus()
sage: A = S.affine_automorphism_group()
sage: M = matrix([[1, 0], [0, 1]])
sage: f = A.derivative().section()(M, check=False)

sage: from flatsurf import SimplicialHomology
sage: H = SimplicialHomology(S)
sage: g = H.hom(f)
matrix()[source]

Return the matrix describing this homomorphism on the generators of homology (as a multiplication from the left.)

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: S = translation_surfaces.mcmullen_L(1, 1, 1, 1)
sage: A = S.affine_automorphism_group()
sage: M = matrix([[1, 2], [0, 1]])
sage: f = A.derivative().section()(M, check=False)

sage: from flatsurf import SimplicialHomology
sage: H = SimplicialHomology(S)
sage: H.gens()
(B[(0, 1)], B[(0, 0)], B[(1, 1)], B[(2, 0)])
sage: g = H.hom(f)

sage: g.matrix()  # optional: pyflatsurf
[1 0 0 0]
[1 1 2 0]
[0 0 1 0]
[1 0 0 1]
class flatsurf.geometry.homology.SimplicialHomologyMorphism_induced(parent, morphism)[source]

A homomorphism of homology induced by a morphism of surfaces.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: S = translation_surfaces.square_torus()
sage: A = S.affine_automorphism_group()
sage: M = matrix([[1, 2], [0, 1]])
sage: f = A.derivative().section()(M, check=False)

sage: from flatsurf import SimplicialHomology
sage: H = SimplicialHomology(S)
sage: g = H.hom(f)
class flatsurf.geometry.homology.SimplicialHomologyMorphism_matrix(parent, matrix)[source]

A homomorphism of homology that is given by a matrix that describes the homomorphism on the generators.

EXAMPLES:

sage: from flatsurf import translation_surfaces
sage: S = translation_surfaces.mcmullen_L(1, 1, 1, 1)
sage: T = translation_surfaces.square_torus()

sage: f = S.homology().hom(matrix([[1, 2, 3, 4], [5, 6, 7, 8]]), codomain=T.homology())
sage: f
Generic morphism:
  From: H₁(Translation Surface in H_2(2) built from 3 squares)
  To:   H₁(Translation Surface in H_1(0) built from a square)
  Defn: [1 2 3 4]
        [5 6 7 8]