Riemann matrices and endomorphism rings of algebraic Riemann surfaces

This module provides a class, RiemannSurface, to model the Riemann surface determined by a plane algebraic curve over a subfield of the complex numbers.

A homology basis is derived from the edges of a Voronoi cell decomposition based on the branch locus. The pull-back of these edges to the Riemann surface provides a graph on it that contains a homology basis.

The class provides methods for computing the Riemann period matrix of the surface numerically, using a certified homotopy continuation method due to [Kr2016].

The class also provides facilities for computing the endomorphism ring of the period lattice numerically, by determining integer (near) solutions to the relevant approximate linear equations.

AUTHORS:

  • Alexandre Zotine, Nils Bruin (2017-06-10): initial version

  • Nils Bruin, Jeroen Sijsling (2018-01-05): algebraization, isomorphisms

EXAMPLES:

We compute the Riemann matrix of a genus 3 curve:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<x,y> = QQ[]
sage: f = x^4-x^3*y+2*x^3+2*x^2*y+2*x^2-2*x*y^2+4*x*y-y^3+3*y^2+2*y+1
sage: S = RiemannSurface(f,prec=100)
sage: M = S.riemann_matrix()

We test the usual properties, i.e., that the period matrix is symmetric and that the imaginary part is positive definite:

sage: all(abs(a) < 1e-20 for a in (M-M.T).list())
True
sage: iM = Matrix(RDF,3,3,[a.imag_part() for a in M.list()])
sage: iM.is_positive_definite()
True

We compute the endomorphism ring and check it has \(\ZZ\)-rank 6:

sage: A = S.endomorphism_basis(80,8)
sage: len(A) == 6
True

In fact it is an order in a number field:

sage: T.<t> = QQ[]
sage: K.<a> = NumberField(t^6 - t^5 + 2*t^4 + 8*t^3 - t^2 - 5*t + 7)
sage: all(len(a.minpoly().roots(K)) == a.minpoly().degree() for a in A)
True
exception sage.schemes.riemann_surfaces.riemann_surface.ConvergenceError

Bases: ValueError

Error object suitable for raising and catching when Newton iteration fails.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import ConvergenceError
sage: raise ConvergenceError("test")
Traceback (most recent call last):
...
ConvergenceError: test
sage: isinstance(ConvergenceError(),ValueError)
True
class sage.schemes.riemann_surfaces.riemann_surface.RiemannSurface(f, prec=53, certification=True, differentials=None)

Bases: object

Construct a Riemann Surface. This is specified by the zeroes of a bivariate polynomial with rational coefficients \(f(z,w) = 0\).

INPUT:

  • f – a bivariate polynomial with rational coefficients. The surface is interpreted as the covering space of the coordinate plane in the first variable.

  • prec – the desired precision of computations on the surface in bits (default: 53)

  • certification – a boolean (default: True) value indicating whether homotopy continuation is certified or not. Uncertified homotopy continuation can be faster.

  • differentials – (default: None). If specified, provides a list of polynomials \(h\) such that \(h/(df/dw) dz\) is a regular differential on the Riemann surface. This is taken as a basis of the regular differentials, so the genus is assumed to be equal to the length of this list. The results from the homology basis computation are checked against this value. Providing this parameter makes the computation independent from Singular. For a nonsingular plane curve of degree \(d\), an appropriate set is given by the monomials of degree up to \(d-3\).

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = w^2 - z^3 + 1
sage: RiemannSurface(f)
Riemann surface defined by polynomial f = -z^3 + w^2 + 1 = 0, with 53 bits of precision

Another Riemann surface with 100 bits of precision:

sage: S = RiemannSurface(f, prec=100); S
Riemann surface defined by polynomial f = -z^3 + w^2 + 1 = 0, with 100 bits of precision
sage: S.riemann_matrix()^6 #abs tol 0.00000001
[1.0000000000000000000000000000 - 1.1832913578315177081175928479e-30*I]

We can also work with Riemann surfaces that are defined over fields with a complex embedding, but since the current interface for computing genus and regular differentials in Singular presently does not support extensions of QQ, we need to specify a description of the differentials ourselves. We give an example of a CM elliptic curve:

sage: Qt.<t> = QQ[]
sage: K.<a> = NumberField(t^2-t+3,embedding=CC(0.5+1.6*I))
sage: R.<x,y> = K[]
sage: f = y^2+y-(x^3+(1-a)*x^2-(2+a)*x-2)
sage: S = RiemannSurface(f,prec=100,differentials=[1])
sage: A = S.endomorphism_basis()
sage: len(A)
2
sage: all( len(T.minpoly().roots(K)) > 0 for T in A)
True
cohomology_basis(option=1)

Compute the cohomology basis of this surface.

INPUT:

  • option – Presently, this routine uses Singular’s adjointIdeal and passes the option parameter on. Legal values are 1, 2, 3 ,4, where 1 is the default. See the Singular documentation for the meaning. The backend for this function may change, and support for this parameter may disappear.

OUTPUT:

This returns a list of polynomials \(g\) representing the holomorphic differentials \(g/(df/dw) dz\), where \(f(z,w)=0\) is the equation specifying the Riemann surface.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = z^3*w + w^3 + z
sage: S = RiemannSurface(f)
sage: S.cohomology_basis()
[1, w, z]
downstairs_edges()

Compute the edgeset of the Voronoi diagram.

OUTPUT:

A list of integer tuples corresponding to edges between vertices in the Voronoi diagram.

EXAMPLES:

Form a Riemann surface, one with a particularly simple branch locus:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = w^2 + z^3 - z^2
sage: S = RiemannSurface(f)

Compute the edges:

sage: S.downstairs_edges()
[(0, 1), (0, 5), (1, 4), (2, 3), (2, 4), (3, 5), (4, 5)]

This now gives an edgeset which one could use to form a graph.

Note

The numbering of the vertices is given by the Voronoi package.

downstairs_graph()

Return the Voronoi decomposition as a planar graph.

The result of this routine can be useful to interpret the labelling of the vertices.

OUTPUT:

The Voronoi decomposition as a graph, with appropriate planar embedding.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = w^2 - z^4 + 1
sage: S = RiemannSurface(f)
sage: S.downstairs_graph()
Graph on 11 vertices

Similarly one can form the graph of the upstairs edges, which is visually rather less attractive but can be instructive to verify that a homology basis is likely correctly computed.:

sage: G = Graph(S.upstairs_edges()); G
Graph on 22 vertices
sage: G.is_planar()
False
sage: G.genus()
1
sage: G.is_connected()
True
edge_permutations()

Compute the permutations of branches associated to each edge.

Over the vertices of the Voronoi decomposition around the branch locus, we label the fibres. By following along an edge, the lifts of the edge induce a permutation of that labelling.

OUTPUT:

A dictionary with as keys the edges of the Voronoi decomposition and as values the corresponding permutations.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = w^2 + z^2+1
sage: S = RiemannSurface(f)
sage: S.edge_permutations()
{(0, 2): (),
 (0, 4): (),
 (1, 2): (),
 (1, 3): (0,1),
 (1, 6): (),
 (2, 0): (),
 (2, 1): (),
 (2, 5): (0,1),
 (3, 1): (0,1),
 (3, 4): (),
 (4, 0): (),
 (4, 3): (),
 (5, 2): (0,1),
 (5, 7): (),
 (6, 1): (),
 (6, 7): (),
 (7, 5): (),
 (7, 6): ()}
endomorphism_basis(b=None, r=None)

Numerically compute a \(\ZZ\)-basis for the endomorphism ring.

Let \(\left(I | M \right)\) be the normalized period matrix (\(M\) is the \(g\times g\) riemann_matrix()). We consider the system of matrix equations \(MA + C = (MB + D)M\) where \(A, B, C, D\) are \(g\times g\) integer matrices. We determine small integer (near) solutions using LLL reductions. These solutions are returned as \(2g \times 2g\) integer matrices obtained by stacking \(\left(D | B\right)\) on top of \(\left(C | A\right)\).

INPUT:

  • b – integer (default provided). The equation coefficients are scaled by \(2^b\) before rounding to integers.

  • r – integer (default: b/4). Solutions that have all coefficients smaller than \(2^r\) in absolute value are reported as actual solutions.

OUTPUT:

A list of \(2g \times 2g\) integer matrices that, for large enough r and b-r, generate the endomorphism ring.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<x,y> = QQ[]
sage: S = RiemannSurface(x^3 + y^3 + 1)
sage: B = S.endomorphism_basis(); B #random
[
[1 0]  [ 0 -1]
[0 1], [ 1  1]
]
sage: sorted([b.minpoly().disc() for b in B])
[-3, 1]
homology_basis()

Compute the homology basis of the Riemann surface.

OUTPUT:

A list of paths \(L = [P_1, \dots, P_n]\). Each path \(P_i\) is of the form \((k, [p_1 ... p_m, p_1])\), where \(k\) is the number of times to traverse the path (if negative, to traverse it backwards), and the \(p_i\) are vertices of the upstairs graph.

EXAMPLES:

In this example, there are two paths that form the homology basis:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: g = w^2 - z^4 + 1
sage: S = RiemannSurface(g)
sage: S.homology_basis() #random
[[(1, [(3, 1), (5, 0), (9, 0), (10, 0), (2, 0), (4, 0),
    (7, 1), (10, 1), (3, 1)])],
 [(1, [(8, 0), (6, 0), (7, 0), (10, 0), (2, 0), (4, 0),
    (7, 1), (10, 1), (9, 1), (8, 0)])]]

In order to check that the answer returned above is reasonable, we test some basic properties. We express the faces of the downstairs graph as ZZ-linear combinations of the edges and check that the projection of the homology basis upstairs projects down to independent linear combinations of an even number of faces:

sage: dg = S.downstairs_graph()
sage: edges = dg.edges()
sage: E = ZZ^len(edges)
sage: edge_to_E = { e[:2]: E.gen(i) for i,e in enumerate(edges)}
sage: edge_to_E.update({ (e[1],e[0]): -E.gen(i) for i,e in enumerate(edges)})
sage: face_span = E.submodule([sum(edge_to_E[e] for e in f) for f in dg.faces()])
sage: def path_to_E(path):
....:     k,P = path
....:     return k*sum(edge_to_E[(P[i][0],P[i+1][0])] for i in range(len(P)-1))
sage: hom_basis = [sum(path_to_E(p) for p in loop) for loop in S.homology_basis()]
sage: face_span.submodule(hom_basis).rank()
2
sage: [sum(face_span.coordinate_vector(b))%2 for b in hom_basis]
[0, 0]
homomorphism_basis(other, b=None, r=None)

Numerically compute a \(\ZZ\)-basis for module of homomorphisms to a given complex torus.

Given another complex torus (given as the analytic Jacobian of a Riemann surface), numerically compute a basis for the homomorphism module. The answer is returned as a list of 2g x 2g integer matrices T=(D, B; C, A) such that if the columns of (I|M1) generate the lattice defining the Jacobian of the Riemann surface and the columns of (I|M2) do this for the codomain, then approximately we have (I|M2)T=(D+M2C)(I|M1), i.e., up to a choice of basis for \(\CC^g\) as a complex vector space, we we realize (I|M1) as a sublattice of (I|M2).

INPUT:

  • b – integer (default provided). The equation coefficients are scaled by \(2^b\) before rounding to integers.

  • r – integer (default: b/4). Solutions that have all coefficients smaller than \(2^r\) in absolute value are reported as actual solutions.

OUTPUT:

A list of \(2g \times 2g\) integer matrices that, for large enough r and b-r, generate the homomorphism module.

EXAMPLES:

sage: S1 = EllipticCurve("11a1").riemann_surface()
sage: S2 = EllipticCurve("11a3").riemann_surface()
sage: [m.det() for m in S1.homomorphism_basis(S2)]
[5]
homotopy_continuation(edge)

Perform homotopy continuation along an edge of the Voronoi diagram using Newton iteration.

INPUT:

  • edge – a tuple of integers indicating an edge of the Voronoi diagram

OUTPUT:

A list of complex numbers corresponding to the points which are reached when traversing along the direction of the edge. The ordering of these points indicates how they have been permuted due to the weaving of the curve.

EXAMPLES:

We check that continued values along an edge correspond (up to the appropriate permutation) to what is stored. Note that the permutation was originally computed from this data:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = z^3*w + w^3 + z
sage: S = RiemannSurface(f)
sage: edge1 = sorted(S.edge_permutations())[0]
sage: sigma = S.edge_permutations()[edge1]
sage: continued_values = S.homotopy_continuation(edge1)
sage: stored_values = S.w_values(S._vertices[edge1[1]])
sage: all( abs(continued_values[i]-stored_values[sigma(i)]) < 1e-8 for i in range(3))
True
make_zw_interpolator(upstairs_edge)

Given an upstairs edge for which continuation data has been stored, return a function that computes \(z(t),w(t)\) , where \(t\) in \([0,1]\) is a parametrization of the edge.

INPUT:

  • upstairs_edge – a pair of integer tuples indicating an edge on the upstairs graph of the surface

OUTPUT:

A tuple (g, d), where g is the function that computes the interpolation along the edge and d is the difference of the z-values of the end and start point.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = w^2 - z^4 + 1
sage: S = RiemannSurface(f)
sage: _ = S.homology_basis()
sage: g,d = S.make_zw_interpolator([(0,0),(1,0)]);
sage: all(f(*g(i*0.1)).abs() < 1e-13for i in range(10))
True
sage: abs((g(1)[0]-g(0)[0]) - d) < 1e-13
True
matrix_of_integral_values(differentials)

Compute the path integrals of the given differentials along the homology basis.

The returned answer has a row for each differential. If the Riemann surface is given by the equation \(f(z,w)=0\), then the differentials are encoded by polynomials g, signifying the differential \(g(z,w)/(df/dw) dz\).

INPUT:

  • differentials – a list of polynomials.

OUTPUT:

A matrix, one row per differential, containing the values of the path integrals along the homology basis of the Riemann surface.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<x,y> = QQ[]
sage: S = RiemannSurface(x^3 + y^3 + 1)
sage: B = S.cohomology_basis()
sage: m = S.matrix_of_integral_values(B)
sage: parent(m)
Full MatrixSpace of 1 by 2 dense matrices over Complex Field with 53 bits of precision
sage: (m[0,0]/m[0,1]).algdep(3).degree() #curve is CM, so the period is quadratic
2
monodromy_group()

Compute local monodromy generators of the Riemann surface.

For each branch point, the local monodromy is encoded by a permutation. The permutations returned correspond to positively oriented loops around each branch point, with a fixed base point. This means the generators are properly conjugated to ensure that together they generate the global monodromy. The list has an entry for every finite point stored in self.branch_locus, plus an entry for the ramification above infinity.

OUTPUT:

A list of permutations, encoding the local monodromy at each branch point.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z, w> = QQ[]
sage: f = z^3*w + w^3 + z
sage: S = RiemannSurface(f)
sage: G = S.monodromy_group(); G
[(0,1,2), (0,1), (0,2), (1,2), (1,2), (1,2), (0,1), (0,2), (0,2)]

The permutations give the local monodromy generators for the branch points:

sage: list(zip(S.branch_locus + [unsigned_infinity], G)) #abs tol 0.0000001
[(0.000000000000000, (0,1,2)),
 (-1.31362670141929, (0,1)),
 (-0.819032851784253 - 1.02703471138023*I, (0,2)),
 (-0.819032851784253 + 1.02703471138023*I, (1,2)),
 (0.292309440469772 - 1.28069133740100*I, (1,2)),
 (0.292309440469772 + 1.28069133740100*I, (1,2)),
 (1.18353676202412 - 0.569961265016465*I, (0,1)),
 (1.18353676202412 + 0.569961265016465*I, (0,2)),
 (Infinity, (0,2))]

We can check the ramification by looking at the cycle lengths and verify it agrees with the Riemann-Hurwitz formula:

sage: 2*S.genus-2 == -2*S.degree + sum(e-1 for g in G for e in g.cycle_type())
True
period_matrix()

Compute the period matrix of the surface.

OUTPUT:

A matrix of complex values.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = z^3*w + w^3 + z
sage: S = RiemannSurface(f, prec=30)
sage: M = S.period_matrix()

The results are highly arbitrary, so it is hard to check if the result produced is correct. The closely related riemann_matrix is somewhat easier to test.:

sage: parent(M)
Full MatrixSpace of 3 by 6 dense matrices over Complex Field with 30 bits of precision
sage: M.rank()
3
plot_paths()

Make a graphical representation of the integration paths.

This returns a two dimensional plot containing the branch points (in red) and the integration paths (obtained from the Voronoi cells of the branch points). The integration paths are plotted by plotting the points that have been computed for homotopy continuation, so the density gives an indication of where numerically sensitive features occur.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<x,y> = QQ[]
sage: S = RiemannSurface(y^2 - x^3 - x)
sage: S.plot_paths()
Graphics object consisting of 2 graphics primitives
plot_paths3d(thickness=0.01)

Return the homology basis as a graph in 3-space.

The homology basis of the surface is constructed by taking the Voronoi cells around the branch points and taking the inverse image of the edges on the Riemann surface. If the surface is given by the equation \(f(z,w)\), the returned object gives the image of this graph in 3-space with coordinates \(\left(\operatorname{Re}(z), \operatorname{Im}(z), \operatorname{Im}(w)\right)\).

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<x,y> = QQ[]
sage: S = RiemannSurface(y^2-x^3-x)
sage: S.plot_paths3d()
Graphics3d Object
riemann_matrix()

Compute the Riemann matrix.

OUTPUT:

A matrix of complex values.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = z^3*w + w^3 + z
sage: S = RiemannSurface(f, prec=60)
sage: M = S.riemann_matrix()

The Klein quartic has a Riemann matrix with values is a quadratic field:

sage: x = polygen(QQ)
sage: K.<a> = NumberField(x^2-x+2)
sage: all(len(m.algdep(6).roots(K)) > 0 for m in M.list())
True
rosati_involution(R)

Compute the Rosati involution of an endomorphism.

The endomorphism in question should be given by its homology representation with respect to the symplectic basis of the Jacobian.

INPUT:

  • R – integral matrix.

OUTPUT:

The result of applying the Rosati involution to R.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: A.<x,y> = QQ[]
sage: S = RiemannSurface(y^2 - (x^6 + 2*x^4 + 4*x^2 + 8), prec = 100)
sage: Rs = S.endomorphism_basis()
sage: S.rosati_involution(S.rosati_involution(Rs[1])) == Rs[1]
True
simple_vector_line_integral(upstairs_edge, differentials)

Perform vectorized integration along a straight path.

INPUT:

  • upstairs_edge – a pair of integer tuples corresponding to an edge of the upstairs graph.

  • differentials – a list of polynomials; a polynomial \(g\) represents the differential \(g(z,w)/(df/dw) dz\) where \(f(z,w)=0\) is the equation defining the Riemann surface.

OUTPUT:

A complex number, the value of the line integral.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = w^2 - z^4 + 1
sage: S = RiemannSurface(f); S
Riemann surface defined by polynomial f = -z^4 + w^2 + 1 = 0, with 53 bits of precision

Since we make use of data from homotopy continuation, we need to compute the necessary data:

sage: M = S.riemann_matrix()
sage: differentials = S.cohomology_basis()
sage: S.simple_vector_line_integral([(0,0),(1,0)], differentials) #abs tol 0.00000001
(1.14590610929717e-16 - 0.352971844594760*I)

Note

Uses data that “homology_basis” initializes.

symplectic_automorphism_group(endo_basis=None, b=None, r=None)

Numerically compute the symplectic automorphism group as a permutation group.

INPUT:

  • endo_basis (default: None) – a \(\ZZ\)-basis of the endomorphisms of self, as obtained from endomorphism_basis(). If you have already calculated this basis, it saves time to pass it via this keyword argument. Otherwise the method will calculate it.

  • b – integer (default provided): as for homomorphism_basis(), and used in its invocation if (re)calculating said basis.

  • r – integer (default: b/4). as for homomorphism_basis(), and used in its invocation if (re)calculating said basis.

OUTPUT:

The symplectic automorphism group of the Jacobian of the Riemann surface. The automorphism group of the Riemann surface itself can be recovered from this; if the curve is hyperelliptic, then it is identical, and if not, then one divides out by the central element corresponding to multiplication by -1.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: A.<x,y> = QQ[]
sage: S = RiemannSurface(y^2 - (x^6 + 2*x^4 + 4*x^2 + 8), prec = 100)
sage: G = S.symplectic_automorphism_group()
sage: G.as_permutation_group().is_isomorphic(DihedralGroup(4))
True
symplectic_isomorphisms(other=None, hom_basis=None, b=None, r=None)

Numerically compute symplectic isomorphisms.

INPUT:

  • other (default: self) – the codomain, another Riemann surface.

  • hom_basis (default: None) – a \(\ZZ\)-basis of the homomorphisms from self to other, as obtained from homomorphism_basis(). If you have already calculated this basis, it saves time to pass it via this keyword argument. Otherwise the method will calculate it.

  • b – integer (default provided): as for homomorphism_basis(), and used in its invocation if (re)calculating said basis.

  • r – integer (default: b/4). as for homomorphism_basis(), and used in its invocation if (re)calculating said basis.

OUTPUT:

This returns the combinations of the elements of homomorphism_basis() that correspond to symplectic isomorphisms between the Jacobians of self and other.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<x,y> = QQ[]
sage: f = y^2 - (x^6 + 2*x^4 + 4*x^2 + 8)
sage: X = RiemannSurface(f, prec=100)
sage: P = X.period_matrix()
sage: g = y^2 - (x^6 + x^4 + x^2 + 1)
sage: Y = RiemannSurface(g, prec=100)
sage: Q = Y.period_matrix()
sage: Rs = X.symplectic_isomorphisms(Y)
sage: Ts = X.tangent_representation_numerical(Rs, other = Y)
sage: test1 = all(((T*P - Q*R).norm() < 2^(-80)) for [T, R] in zip(Ts, Rs))
sage: test2 = all(det(R) == 1 for R in Rs)
sage: test1 and test2
True
tangent_representation_algebraic(Rs, other=None, epscomp=None)

Compute the algebraic tangent representations corresponding to the homology representations in Rs.

The representations on homology Rs have to be given with respect to the symplectic homology basis of the Jacobian of self and other. Such matrices can for example be obtained via endomorphism_basis().

Let \(P\) and \(Q\) be the period matrices of self and other. Then for a homology representation \(R\), the corresponding tangential representation \(T\) satisfies \(T P = Q R\).

INPUT:

  • Rs – a set of matrices on homology to be converted to their tangent representations.

  • other (default: self) – the codomain, another Riemann surface.

  • epscomp – real number (default: 2^(-prec + 30)). Used to determine whether a complex number is close enough to a root of a polynomial.

OUTPUT:

The algebraic tangent representations of the matrices in Rs.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: A.<x,y> = QQ[]
sage: S = RiemannSurface(y^2 - (x^6 + 2*x^4 + 4*x^2 + 8), prec = 100)
sage: Rs = S.endomorphism_basis()
sage: Ts = S.tangent_representation_algebraic(Rs)
sage: Ts[0].base_ring().maximal_order().discriminant() == 8
True
tangent_representation_numerical(Rs, other=None)

Compute the numerical tangent representations corresponding to the homology representations in Rs.

The representations on homology Rs have to be given with respect to the symplectic homology basis of the Jacobian of self and other. Such matrices can for example be obtained via endomorphism_basis().

Let \(P\) and \(Q\) be the period matrices of self and other. Then for a homology representation \(R\), the corresponding tangential representation \(T\) satisfies \(T P = Q R\).

INPUT:

  • Rs – a set of matrices on homology to be converted to their tangent representations.

  • other (default: self) – the codomain, another Riemann surface.

OUTPUT:

The numerical tangent representations of the matrices in Rs.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: A.<x,y> = QQ[]
sage: S = RiemannSurface(y^2 - (x^6 + 2*x^4 + 4*x^2 + 8), prec = 100)
sage: P = S.period_matrix()
sage: Rs = S.endomorphism_basis()
sage: Ts = S.tangent_representation_numerical(Rs)
sage: all(((T*P - P*R).norm() < 2^(-80)) for [T, R] in zip(Ts, Rs))
True
upstairs_edges()

Compute the edgeset of the lift of the downstairs graph onto the Riemann surface.

OUTPUT:

An edgeset between vertices (i, j), where i corresponds to the i-th point in the Voronoi diagram vertices, and j is the j-th w-value associated with that point.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = w^2 + z^3 - z^2
sage: S = RiemannSurface(f)
sage: edgeset = S.upstairs_edges()
sage: len(edgeset) == S.degree*len(S.downstairs_edges())
True
sage: {(v[0],w[0]) for v,w in edgeset} == set(S.downstairs_edges())
True
w_values(z0)

Return the points lying on the surface above z0.

INPUT:

  • z0 – (complex) a point in the complex z-plane.

OUTPUT:

A set of complex numbers corresponding to solutions of \(f(z_0,w) = 0\).

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface
sage: R.<z,w> = QQ[]
sage: f = w^2 - z^4 + 1
sage: S = RiemannSurface(f)

Find the w-values above the origin, i.e. the solutions of \(w^2 + 1 = 0\):

sage: S.w_values(0) # abs tol 1e-14
[-1.00000000000000*I, 1.00000000000000*I]
class sage.schemes.riemann_surfaces.riemann_surface.RiemannSurfaceSum(L)

Bases: sage.schemes.riemann_surfaces.riemann_surface.RiemannSurface

Represent the disjoint union of finitely many Riemann surfaces.

Rudimentary class to represent disjoint unions of Riemann surfaces. Exists mainly (and this is the only functionality actually implemented) to represents direct products of the complex tori that arise as analytic Jacobians of Riemann surfaces.

INPUT:

  • L – list of RiemannSurface objects

EXAMPLES:

sage: _.<x> = QQ[]
sage: SC = HyperellipticCurve(x^6-2*x^4+3*x^2-7).riemann_surface(prec=60)
sage: S1 = HyperellipticCurve(x^3-2*x^2+3*x-7).riemann_surface(prec=60)
sage: S2 = HyperellipticCurve(1-2*x+3*x^2-7*x^3).riemann_surface(prec=60)
sage: len(SC.homomorphism_basis(S1+S2))
2
period_matrix()

Return the period matrix of the surface.

This is just the diagonal block matrix constructed from the period matrices of the constituents.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface, RiemannSurfaceSum
sage: R.<x,y> = QQ[]
sage: S1 = RiemannSurface(y^2-x^3-x-1)
sage: S2 = RiemannSurface(y^2-x^3-x-5)
sage: S = RiemannSurfaceSum([S1,S2])
sage: S1S2 = S1.period_matrix().block_sum(S2.period_matrix())
sage: S.period_matrix() == S1S2[[0,1],[0,2,1,3]]
True
riemann_matrix()

Return the normalized period matrix of the surface.

This is just the diagonal block matrix constructed from the Riemann matrices of the constituents.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface, RiemannSurfaceSum
sage: R.<x,y> = QQ[]
sage: S1 = RiemannSurface(y^2-x^3-x-1)
sage: S2 = RiemannSurface(y^2-x^3-x-5)
sage: S = RiemannSurfaceSum([S1,S2])
sage: S.riemann_matrix() == S1.riemann_matrix().block_sum(S2.riemann_matrix())
True
sage.schemes.riemann_surfaces.riemann_surface.bisect(L, t)

Find position in a sorted list using bisection.

Given a list \(L = [(t_0,...),(t_1,...),...(t_n,...)]\) with increasing \(t_i\), find the index i such that \(t_i <= t < t_{i+1}\) using bisection. The rest of the tuple is available for whatever use required.

INPUT:

  • L – A list of tuples such that the first term of each tuple is a real number between 0 and 1. These real numbers must be increasing.

  • t – A real number between \(t_0\) and \(t_n\).

OUTPUT:

An integer i, giving the position in L where t would be in

EXAMPLES:

Form a list of the desired form, and pick a real number between 0 and 1:

sage: from sage.schemes.riemann_surfaces.riemann_surface import bisect
sage: L = [(0.0, 'a'), (0.3, 'b'), (0.7, 'c'), (0.8, 'd'), (0.9, 'e'), (1.0, 'f')]
sage: t = 0.5
sage: bisect(L,t)
1

Another example which demonstrates that if t is equal to one of the t_i, it returns that index:

sage: L = [(0.0, 'a'), (0.1, 'b'), (0.45, 'c'), (0.5, 'd'), (0.65, 'e'), (1.0, 'f')]
sage: t = 0.5
sage: bisect(L,t)
3
sage.schemes.riemann_surfaces.riemann_surface.differential_basis_baker(f)

Compute a differential bases for a curve that is nonsingular outside (1:0:0),(0:1:0),(0:0:1)

Baker’s theorem tells us that if a curve has its singularities at the coordinate vertices and meets some further easily tested genericity criteria, then we can read off a basis for the regular differentials from the interior of the Newton polygon spanned by the monomials. While this theorem only applies to special plane curves it is worth implementing because the analysis is relatively cheap and it applies to a lot of commonly encountered curves (e.g., curves given by a hyperelliptic model). Other advantages include that we can do the computation over any exact base ring (the alternative Singular based method for computing the adjoint ideal requires the rationals), and that we can avoid being affected by subtle bugs in the Singular code.

None is returned when f does not describe a curve of the relevant type. If f is of the relevant type, but is of genus \(0\) then [] is returned (which are both False values, but they are not equal).

INPUT:

  • \(f\) – a bivariate polynomial

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import differential_basis_baker
sage: R.<x,y> = QQ[]
sage: f = x^3+y^3+x^5*y^5
sage: differential_basis_baker(f)
[y^2, x*y, x*y^2, x^2, x^2*y, x^2*y^2, x^2*y^3, x^3*y^2, x^3*y^3]
sage: f = y^2-(x-3)^2*x
sage: differential_basis_baker(f) is None
True
sage: differential_basis_baker(x^2+y^2-1)
[]
sage.schemes.riemann_surfaces.riemann_surface.integer_matrix_relations(M1, M2, b=None, r=None)

Determine integer relations between complex matrices.

Given two square matrices with complex entries of size g, h respectively, numerically determine an (approximate) ZZ-basis for the 2g x 2h matrices with integer entries of the shape (D, B; C, A) such that B+M1*A=(D+M1*C)*M2. By considering real and imaginary parts separately we obtain \(2gh\) equations with real coefficients in \(4gh\) variables. We scale the coefficients by a constant \(2^b\) and round them to integers, in order to obtain an integer system of equations. Standard application of LLL allows us to determine near solutions.

The user can specify the parameter \(b\), but by default the system will choose a \(b\) based on the size of the coefficients and the precision with which they are given.

INPUT:

  • M1 – square complex valued matrix

  • M2 – square complex valued matrix of same size as M1

  • b – integer (default provided). The equation coefficients are scaled by \(2^b\) before rounding to integers.

  • r – integer (default: b/4). The vectors found by LLL that satisfy the scaled equations to within \(2^r\) are reported as solutions.

OUTPUT:

A list of 2g x 2h integer matrices that, for large enough \(r\), \(b-r\), generate the ZZ-module of relevant transformations.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import integer_matrix_relations
sage: M1=M2=matrix(CC,2,2,[sqrt(d) for d in [2,-3,-3,-6]])
sage: T=integer_matrix_relations(M1,M2)
sage: id=parent(M1)(1)
sage: M1t=[id.augment(M1) * t for t in T]
sage: [((m[:,:2]^(-1)*m)[:,2:]-M2).norm() < 1e-13 for m in M1t]
[True, True]
sage.schemes.riemann_surfaces.riemann_surface.numerical_inverse(C)

Compute numerical inverse of a matrix via LU decomposition

INPUT:

  • C – A real or complex invertible square matrix

EXAMPLES:

sage: C = matrix(CC,3,3,[-4.5606e-31 + 1.2326e-31*I,
....: -0.21313 + 0.24166*I,
....: -3.4513e-31 + 0.16111*I,
....: -1.0175 + 9.8608e-32*I,
....: 0.30912 + 0.19962*I,
....: -4.9304e-32 + 0.39923*I,
....: 0.96793 - 3.4513e-31*I,
....: -0.091587 + 0.19276*I,
....: 3.9443e-31 + 0.38552*I])
sage: from sage.schemes.riemann_surfaces.riemann_surface import numerical_inverse
sage: max(abs(c) for c in (C^(-1)*C-C^0).list()) < 1e-10
False
sage: max(abs(c) for c in (numerical_inverse(C)*C-C^0).list()) < 1e-10
True
sage.schemes.riemann_surfaces.riemann_surface.voronoi_ghost(cpoints, n=6, CC=Complex Double Field)

Convert a set of complex points to a list of real tuples \((x,y)\), and appends n points in a big circle around them.

The effect is that, with n >= 3, a Voronoi decomposition will have only finite cells around the original points. Furthermore, because the extra points are placed on a circle centered on the average of the given points, with a radius 3/2 times the largest distance between the center and the given points, these finite cells form a simply connected region.

INPUT:

  • cpoints – a list of complex numbers

OUTPUT:

A list of real tuples \((x,y)\) consisting of the original points and a set of points which surround them.

EXAMPLES:

sage: from sage.schemes.riemann_surfaces.riemann_surface import voronoi_ghost
sage: L = [1 + 1*I, 1 - 1*I, -1 + 1*I, -1 - 1*I]
sage: voronoi_ghost(L) # abs tol 1e-6
[(1.0, 1.0),
 (1.0, -1.0),
 (-1.0, 1.0),
 (-1.0, -1.0),
 (2.121320343559643, 0.0),
 (1.0606601717798216, 1.8371173070873836),
 (-1.060660171779821, 1.8371173070873839),
 (-2.121320343559643, 2.59786816870648e-16),
 (-1.0606601717798223, -1.8371173070873832),
 (1.06066017177982, -1.8371173070873845)]