Smooth projective curves over a field

Let \(k\) be a field and \(F/k\) a finitely generated field extension of transcendence degree one (i.e. a ‘function field over \(k\)’). Then there exists a smooth projective curve \(X\) over \(Spec(k)\) with function field \(F\), unique up to unique isomorphism. The set of closed points on \(X\) are in natural bijection with the set of discrete valuations on \(F\) which are trivial on \(k\). See

    1. Hartshorne, Algebraic Geometry, Theorem I.6.9.

The classes in this module provide the definition and some basic functionality for such curves.

A curve \(X\) is defined via its function field \(F_X\). Points are represented by the corresponding valuations on \(F_X\), and no smooth projective model of \(X\) is actually computed. However, we do compute a list of ‘coordinate functions’ \(x_1,..,x_n\) which separate all points, meaning that the closure of the rational map from \(X\) to projective space of dimension \(n\) is injective. Then a (closed) point \(x\) on \(X\) can also be represented by the tupel \((f_1(x),..,f_n(x))\). This is useful to test for equality of points.

A function field in Sage is always realized as a simple separable extension of a rational function field. Geometrically, this means that the curve \(X\) is implicitly equipped with a finite separable morphism \(\phi:X\to\mathbb{P}^1_k\) to the projective line over the base field \(k\).

The base field \(k\) is called the constant base field of the curve, and it is part of the data. We do not assume that the extension \(F_X/k\) is regular (i.e. that \(k\) is algebraically closed in \(F_X\)). Geometrically this means that the curve \(X\) may not be absolutely irreducibel as a \(k\)-scheme. The field of constants of \(X\) is defined as the algebraic closure of \(k\) inside \(F_X\). It is a finite extension \(k_c/k\). If we regard \(X\) as a curve over its fields of constants then it becomes absolute irreducible.

It would be interesting to have an efficient algorithm for computing the field of constants, but it seems that this has not been implemented in Sage yet. To compute the genus of \(X\) it is necessary to know at least the degree \([k_c:k]\). If \(k\) is a finite field, it is actually easy to compute \(k_c\). If \(k\) is a number field we use a probabilistic algorithm for computing the degree \([k_c:k]\), by reducing the curve modulo several small primes.

Currently, the function field \(F\) defining a smooth projective curve must be a simple finite extension of a rational function field, i.e. of the form

\[F = k(x)[y]/(G)\]

where \(G\) is an irreducible polynomial over \(k(x)\). If not explicitly stated otherwise, it is assumed that \(k\) is the constant base field of the curve \(X\). If \(k\) is a finite field, then one may also declare any subfield \(k_0\) of \(k\) to be the constant base field. Geometrically, this means that we consider \(X\) as a curve over \({\rm Spec}(k_0)\). In any case, the field of constants of \(X\) is a finite extension of \(k\).

AUTHORS:

  • Stefan Wewers (2016-11-11): initial version

EXAMPLES:

sage: from mclf import *
sage: K = GF(2)
sage: FX.<x> = FunctionField(K)
sage: R.<T> = FX[]
sage: FY.<y> = FX.extension(T^2+x^2*T+x^5+x^3)
sage: Y = SmoothProjectiveCurve(FY)
sage: Y
the smooth projective curve with Function field in y defined by y^2 + x^2*y + x^5 + x^3
sage: Y.genus()
1
sage: Y.zeta_function()
(2*T^2 + T + 1)/(2*T^2 - 3*T + 1)

Over finite fields, we are allowed to specify the constant base field:

sage: K = GF(4)
sage: F.<x> = FunctionField(K)
sage: X = SmoothProjectiveCurve(F, k=GF(2))
sage: X
the smooth projective curve with Rational function field in x over Finite Field in z2 of size 2^2 and constant base field Finite Field of size 2
sage: X.field_of_constants()
Finite Field in z2 of size 2^2

A curve may also be defined by an irreducible bivariate polynomial:

sage: A.<x,y> = QQ[]
sage: X = SmoothProjectiveCurve(y^2 - x^3 - 1)
sage: X
the smooth projective curve with Function field in y defined by y^2 - x^3 - 1

If the curve is defined over a number field, we can find a prime of good reduction, and compute the reduction:

sage: v_p = X.prime_of_good_reduction()
sage: v_p
5-adic valuation

The result is a discrete (p-adic) valuation on the constant base field. The reduction is a smooth projective curve of the same genus:

sage: Xb = X.good_reduction(v_p)
sage: Xb
the smooth projective curve with Function field in y defined by y^2 + 4*x^3 + 4
sage: Xb.genus()
1

Todo

  • allow to specify the constant base field in a more flexible way
  • write more doctests !!
  • implement a general and deterministic algorithm for computing the field of constants (and not just the degree)
  • the residue field of a point should be explicitly an extension of the constant base field.
  • treat the base curve \(X\) as a curve, not just as a function field
  • realize morphisms between curves, in particular the canonical map to \(X\)
class mclf.curves.smooth_projective_curves.PointOnSmoothProjectiveCurve(X, v)

Bases: SageObject

A closed point on a smooth projective curve.

A point on a curve \(X\) is identified with the corresponding valuation \(v_x\) on the function field \(F\) of \(X\).

Alternatively, a point \(x\) on \(X\) can be represented by the vector

\[[v_x(f_1),\ldots, v_x(f_n)]\]

where \(f_1,\ldots,f_n\) is a list of coordinate functions, i.e. rational functions which define an injective map from \(X\) into \(\mathbb{P}^1\times\ldots\times\mathbb{P}^1\).

We use the latter representation to check for equality of points.

absolute_degree()

Return the absolute degree of self.

The absolute degree of a point \(x\) on a curve \(X\) over \(k\) is the degree of the extension \(k(x)/k\).

Here \(k\) is the constant base field of the curve, which may not be equal to the field of constants.

coordinates()

Return the coordinate tupel of the point.

NOTE:

for a curve over a number field and for a point whose residue field is of high degree, this can be very slow. It would be better to implement this function in a lazy way, for instance as an iterator.

curve()

Return the underlying curve of the point.

degree()

Return the degree of self.

The degree of a point \(x\) on a curve \(X\) over \(k\) is the degree of the residue field \(k(x)\) as an extension of the field of constants of \(X\). The latter may be a proper extension of the base field \(k\)!

is_equal(P)

Check whether this point is equal to P.

INPUT:

  • P – a point on the curve underlying this point

OUTPUT:

True is \(P\) is equal to self, False otherwise.

Currently, the check for equality is done using the ‘coordinates’ of the points. This may be very slow. It would probably be better to test the equality of the underlying valuations. But here we can’t rely on Sage, so this would require a hack.

EXAMPLES:

sage: import mclf.curves.smooth_projective_curves
sage: from mclf.curves.smooth_projective_curves import SmoothProjectiveCurve
sage: F0.<x> = FunctionField(GF(3))
sage: R.<y> = F0[]
sage: F.<y> = F0.extension(y^2 - (x+1)*x^2)
sage: Y = SmoothProjectiveCurve(F)
sage: v0 = F0.valuation(x)
sage: fiber = Y.fiber(v0)
sage: fiber[0].is_equal(fiber[1])
False
order(f)

Return the order of the function in the point.

This is the same as self.valuation()(f).

residue_field()

Return the residue field of the point.

valuation()

Return the valuation corresponding to the point.

value(f)

Return the value of the function f in the point.

If f has a pole then Infinity is returned.

class mclf.curves.smooth_projective_curves.SmoothProjectiveCurve(F, k=None, assume_regular=False)

Bases: SageObject

Return the smooth projective curve with function field \(F\).

INPUT:

  • F – a function field, or an irreducible bivariate polynomial over a field
  • k – a field which has a natural embedding into the constant base field of \(F\), such that the constant base field is a finite extension of k (or None).
  • assume_regular – a boolean (default: False)

OUTPUT:

the smooth projective curve \(X\) with function field \(F\). If \(F\) is an irreducible bivariate polynomial, we use the function field with two generators and relation \(F\).

If \(k\) is given, then \(X\) is considered as a \(k\)-scheme. If \(k\) is not given then we use the field of constants of \(F\) instead.

NOTE:

At the moment, \(k\) should only be different from the constant base field of \(F\) if \(k\) is finite (because it is then easy to compute the degree of the degree of the constant base field of \(F\) over \(k\)).

compute_separable_model()

Compute a separable model of the curve (if necessary).

OUTPUT: None

This function only has to be called only once. It then decides whether or not the function field of the curve is given as a separable extension of the base field or not. If it is not separable then we compute a separable model, which is a tripel \((Y_1,\phi, q)\) where

  • \(Y_1\) is a smooth projective curve over the same constant base field \(k\) as the curve \(Y\) itself, and which is given by a separable extension,
  • \(\phi\) is a field homomorphism from the function field of \(Y_1\) into the function field of \(Y\) corresponding to a purely inseparable extension,
  • \(q\) is the degree of the extension given by \(\phi\), i.e. the degree of inseparability of the map \(Y\to Y_1\) given by \(\phi\). Note that \(q\) is a power of the characteristic of \(k\).
constant_base_field()

Return the constant base field.

coordinate_functions()

Return a list of coordinate functions.

By ‘list of coordinate functions’ we mean elements \(f_i\) in the function field, such that the map

\[x \mapsto (f_1(x),\ldots, f_n(x))\]

from \(X\) to \((\mathbb{P}^1)^n\) is injective.

Note that this map may not be an embedding, i.e. image of this map may not be a smooth model of the curve.

count_points(d)

Return number of points of degree less or equal to d.

INPUT:

  • d – an interger \(\geq 1\)

OUTPUT:

a list N, where N[i] is the number of points on self of absolute degree \(i\), for \(i=1,..,d\).

Recall that the absolute degree of a point is the degree of the residue field of the point over the constant base field (not over the field of constants).

This is a very slow realization and should be improved at some point.

covering_degree()

Return the degree of this curve as a covering of the projective line.

degree(D)

Return the degree of the divisor D.

Note that the degree of \(D\) is defined relative to the field of constants of the curve.

degree_of_inseparability()

Return the degree of inseparability of this curve.

OUTPUT: positive integer, which is a power of the characteristic of the function field of this curve.

divisor_of_poles(f)

Return the divisor of poles of f.

divisor_of_zeroes(f)

Return the divisor of zeroes of f.

fiber(v0)

Return the set of points lying above a point on the projective line.

INPUT:

  • v0 – a function field valuation on the rational function field

OUTPUT:

a list containing the points on this curve \(Y\) corresponding to extensions of v0 to the function field of \(Y\).

field_of_constants()

Return the field of constants of this curve.

If \(F\) is the function field of the curve and \(k\) the constant base field, then the field of constants is the algebraic closure of \(k\) in \(F\).

For the moment, this is implemented only if the constant base field is a finite field.

EXAMPLES:

sage: from mclf import *
sage: F.<x> = FunctionField(GF(2))
sage: R.<y> = F[]
sage: G = (y+x)^4 + (y + x) + 1
sage: F1.<y> = F.extension(G)
sage: Y1 = SmoothProjectiveCurve(F1)
sage: Y1.field_of_constants()
Finite Field in z4 of size 2^4
field_of_constants_degree()

Return the degree of the field of constants over the constant base field.

If \(F\) is the function field of the curve and \(k\) the constant base field, then the field of constants \(k_c\) is the algebraic closure of \(k\) in \(F\).

If \(k\) is a finite field then we actually compute the field of constants, and the result is provably correct. If \(k\) is a number field, then we use a heuristic method: we find at least \(10\) different primes of \(k\) for which the reduction of the defining equation remains irreducible, and then we apply the method for finite fields to the reduced equation. The result is very likely the true degree of the field of constants, and if the result is equal to \(1\) then it is provably correct.

EXAMPLES:

sage: from mclf import *
sage: k = GF(2^3)
sage: F.<x> = FunctionField(k)
sage: R.<y> = F[]
sage: G = (y+x)^4 + (y+x) + 1
sage: F1.<y> = F.extension(G)
sage: Y1 = SmoothProjectiveCurve(F1, GF(2))
sage: Y1.field_of_constants_degree()
12
sage: F.<x> = FunctionField(QQ)
sage: R.<y> = F[]
sage: G = y^4 + x*y + 1
sage: F2.<y> = F.extension(G)
sage: Y2 = SmoothProjectiveCurve(F2)
sage: Y2.field_of_constants_degree()
1
sage: R.<y> = F[]
sage: G = (y+x)^3 + (y+x) + 1
sage: F3.<y> = F.extension(G)
sage: Y3 = SmoothProjectiveCurve(F3)
sage: Y3.field_of_constants_degree()
3
sage: Y3.genus()
0

Todo

  • implement a deterministic algorithm for number fields
function_field()

Return the function field of the curve self.

genus(use_reduction=True)

Return the genus of the curve.

INPUT:

  • use_reduction – a boolean (default: True)

OUTPUT: the genus of this curve.

The genus of the curve is defined as the dimension of the cohomology group \(H^1(X,\mathcal{O}_X)\), as a vector space over the field of constants `k_c`.

The genus \(g\) of the curve \(X\) is computed using the Riemann-Hurwitz formula, applied to the cover \(X\to\mathbb{P}^1\) corresponding to the underlying realization of the function field of \(X\) as a finite separable extension of a rational function field. See:

  • Hartshorne, Algebraic Geometry, Corollary IV.2.4

If the constant base field is finite, we compute the degree of the ‘ramification divisor’. If it is not, we assume that the characteristic is zero, and we use the ‘tame’ Riemann Hurwitz Formula.

If the curve is defined over a number field, and use_reduction is True (the default) then the genus of a reduction of this curve to some prime of good reduction is computed. This may be consderably faster.

EXAMPLES:

sage: from mclf import *
sage: F0.<x> = FunctionField(GF(2))
sage: R.<T> = F0[]
sage: G = T^2 + T + x^3 + x + 1
sage: F.<y> = F0.extension(G)
sage: Y = SmoothProjectiveCurve(F)
sage: Y.genus()
1
sage: G = T^2 + x^3 + x + 1
sage: F.<y> = F0.extension(G)
sage: Y = SmoothProjectiveCurve(F)
sage: Y.genus()
0
good_reduction(v_p)

Return the reduction of this curve at a prime of good reduction.

INPUT:

  • v_p – a discrete valuation on the constant base field \(K\)

OUTPUT: the reduction of this curve with respect to \(v_p\), assuming that this is again a smooth projective curve of the same genus. Otherwise, an error is raised.

Note that we just reduce the plane equation for this curve with respect to \(v_p\). This is a very naive notion of good reduction. If it works, then the curve does indeed have good reduction at \(v_p\), and the result is correct.

is_separable()

Check whether this curve is represented as a separable cover of the projective line.

phi()

Return the natural embedding of the function field of the separable model into the function field of this curve.

OUTPUT: a field homomorphism

plane_equation()

Return the plane equation of this curve.

OUTPUT:

A polynomial in two variables over the constant base field which defines the plane model of this curve, where the first variable corresponds to the base field F0.

point(v)

Returns the point on the curve corresponding to v.

INPUT:

  • v – a discrete valuation on the function field of the curve

OUTPUT:

The point on the curve corresponding to v. The valuation \(v\) must be trivial on the constant base field.

points_with_coordinates(a)

Return all points with given coordinates.

INPUT:

  • a – a tupel of coordinates, of lenght \(n\), at most the number of coordinate functions of the curve

OUTPUT:

a list containing all points on the curve whose first \(n\) coordinate values agree with a.

potential_branch_divisor()

Return list of valuations containing the branch locus.

OUTPUT:

A list of pairs \((v,d)', where \) runs over a list of valuations of the base field \(F_0 = K(x)\) containing all the valuations corresponding to a branch point of the cover of curves, and \(d\) is the degree of \(v\).

prime_of_good_reduction()

Return a prime ideal where this curve has good reduction.

OUTPUT:

We assume that the constant base field \(K\) is a number field. We return a discrete valuation \(v\) on \(K\) such that the following holds:

  • all the coefficients of the plane equation \(G(x,y)=0\) of this curve are \(v\)-integral
  • the reduction of \(G\) to the residue field of \(v\) is irreducible and defines a plane curve with the same genus as the original curve.

Note that this implies that \(v\) is inert in the field of constants of the curve.

EXAMPLES:

sage: from mclf import *
sage: R.<T> = QQ[]
sage: K.<zeta> = NumberField(T^2+T+1)
sage: A.<x,y> = K[]
sage: X = SmoothProjectiveCurve(y^3 - y - x^2 + 1 + zeta)
sage: X.prime_of_good_reduction()
5-adic valuation
sage: X.good_reduction(_)
the smooth projective curve with Function field in y defined by y^3 + 4*y + 4*x^2 + u1 + 1
principal_divisor(f)

Return the principal divisor of f.

INPUT:

  • f – a nonzero element of the function field of self

OUTPUT: the principal divisor \(D =(f)\). This is a list of pairs \((P, m)\), where \(P\) is a point and \(m\) is an integer.

ramification_divisor()

Return the ramification divisor of self.

OUTPUT:

The ramification divisor of the curve, if the curve is given by a separable model. Otherwise, an error is raised.

So the function field of of the curve is a finite separable extension of a rational function field. Geometrically, this means that the curve \(X\) is represented as a separable cover of the projective line. The ramification divisor of this cover is supported in the set of ramification points of this cover. Sheaf theoretically, the divisor represents the sheaf of relative differentials. See:

  • Hartshorne, Algebraic Geometry, Definition IV.2.
random_point()

Return a random closed point on the curve.

rational_function_field()

Return the rational function field underlying the function field of \(X\).

By definition, the function field \(F_X\) of the curve \(X\) is a finite separable extension of a rational function field \(k(x)\), where \(k\) is the base field of \(X\).

separable_model()

Return the separable model of this curve.

OUTPUT: a smooth projective curve over the same constant base field

The separable model of this curve \(Y\) is a curve \(Y_s\) defined over the same constant base field and whose defining equation realizes \(Y_s\) as a finite separable cover of the projective line. It comes equipped with a finite, purely inseparable morphism \(Y\to Y_s\). In particular, \(Y_s\) has the same genus as \(Y\).

The inclusion of function fields \(\phi:F(Y_s)\to F(Y)\) can be accessed via the method phi(), the degree of the extension \(Y/Y_s\) via the method degree_of_inseparability.

singular_locus()

Return the singular locus of the affine model of the curve.

OUTPUT:

a list of discrete valuations on the base field \(k(x)\) which represent the \(x\)-coordinates of the points where the affine model of the curve given by the defining equation of the function field may be singular.

structure_map()

Return the canonical map from this curve to the projective line.

zeta_function(var_name='T')

Return the Zeta function of the curve.

For any scheme \(X\) of finite type over \(\mathbb{Z}\), the arithmetic zeta funtion of \(X\) is defined as the product

\[\zeta(X,s) := \prod_x \frac{1}{1-N(x)^{-s}},\]

where \(x\) runs over over all closed points of \(X\) and \(N(x)\) denotes the cardinality of the residue field of \(x\).

If \(X\) is a smooth projective curve over a field with \(q\) elements, then \(\zeta(X,s) = Z(X,q^{-s})\), where \(Z(X,T)\) is a rational function in \(T\) of the form

\[Z(X,T) = \frac{P(T)}{(1-T)(1-qT)},\]

for a polynomial \(P\) of degree \(2g\), with some extra properties reflecting the Weil conjectures. See:

  • Hartshorn, Algebraic Geometry, Appendix C, Section 1.

Note that that this makes only sense if the constant base field of self is finite, and that \(Z(X,T)\) depends on the choice of the constant base field (unlike the function \(\zeta(X,s)\)!).

mclf.curves.smooth_projective_curves.absolute_degree(K)

Return the absolute degree of a (finite) field.

mclf.curves.smooth_projective_curves.compute_value(v, f)

Return the value of f at v.

INPUT:

  • v – a function field valuation on \(F\)
  • f – an element of \(F\)

OUTPUT: The value of f at the point corresponding to v.

This is either an element of the residue field of the valuation v (which is a finite extension of the base field of \(F\)), or \(\infty\).

mclf.curves.smooth_projective_curves.e_f_of_valuation(v)

Return the ramification index of this valuation.

INPUT:

  • v – a function field valuation on an extension of a rational function field

OUTPUT: the ramification index of \(v\) over the base field

mclf.curves.smooth_projective_curves.extension_degree(K, L, check=False)

Return the degree of the field extension.

INPUT:

  • K, L – two fields, where K has a natural embedding into L
  • check (default: False) – boolean

OUTPUT:

the degree \([L:K]\)

At the moment, this works correctly only if \(K\) and \(L\) are finite extensions of a common base field. It is not checked whether \(K\) really maps to \(L\).

mclf.curves.smooth_projective_curves.extension_of_finite_field(K, n)

Return a field extension of this finite field of degree n.

INPUT:

  • K – a finite field
  • n – a positive integer

OUTPUT: a field extension of \(K\) of degree \(n\).

This function is useful if \(K\) is constructed as an explicit extension \(K = K_0[x]/(f)\); then K.extension(n) is not implemented.

Note

This function should be removed once trac.sagemath.org/ticket/26103 has been merged.

mclf.curves.smooth_projective_curves.field_of_constant_degree_of_polynomial(G, return_field=False)

Return the degree of the field of constants of a polynomial.

INPUT:

  • G – an irreducible monic polynomial over a rational function field
  • return_field – a boolean (default:\(False\))

OUTPUT: the degree of the field of constants of the function field defined by G. If return_field is True then the actual field of constants is returned. This is currently implemented for finite fields only.

This is a helper function for SmoothProjectiveCurve.field_of_constants_degree.

mclf.curves.smooth_projective_curves.make_finite_field(k)

Return the finite field isomorphic to this field.

INPUT:

  • k – a finite field

OUTPUT: a triple \((k_1,\phi,\psi)\) where \(k_1\) is a ‘true’ finite field, \(\phi\) is an isomorphism from \(k\) to \(k_1\) and \(\psi\) is the inverse of \(\phi\).

This function is useful when \(k\) is constructed as a tower of extensions with a finite field as a base field.

Note

This function should be removed once trac.sagemath.org/ticket/26103 has been merged.

mclf.curves.smooth_projective_curves.separate_points(coordinate_functions, valuations)

Add new coordinate functions to separate a given number of points.

INPUT:

  • coordinate_functions – a list of elements of a function field \(F\)
  • valuations – a list of function field valuations on \(F\)

OUTPUT: enlarges the list coordinate_functions in such a way that the lists [value(v,x) for x in coordinate_functions], where v runs through valuations, are pairwise distinct.

mclf.curves.smooth_projective_curves.separate_two_points(v1, v2)

Return a rational function which separates two points

INPUT:

  • v1, v2 – discrete, nonequivalent valuations on a common function field \(F\)

OUTPUT:

An element \(f\) of \(F\) which takes distinct values at the two points corresponding to v1 and v2.

mclf.curves.smooth_projective_curves.sum_of_divisors(D1, D2)

Return the sum of the divisors D1 and D2.

INPUT:

  • D1, D2: divisors on the same curve \(X\)

OUTPUT:

\(D_1\) is replaced by the sum \(D_1+D_2\) (note that this changes \(D_1\)!).

Here a divisor \(D\) is given by a dictionary with entries (a:(P,m)), where a is a coordinate tupel, P is a point on \(X\) with coordinates a and m is the multiplicity of P in \(D\).