Points on the Jacobian

Points on Jac(C) are represented as divisors on C. They can be specified simply by giving points on C, or divisors on C, or in the Mumford representation, which is the way Magma returns them (and which it uses to store and manipulate them). Points can be added and subtracted when the base field is an exact field.

Representation of points on Jac(C): Let C be a hyperelliptic curve of genus g. A triple < a(x), b(x), d > specifies the divisor D of degree d on C defined by A(x, z) = 0, y = B(x, z) where A(x, z) is the degree d homogenisation of a(x), and B(x, z) is the degree (g + 1) homogenisation of b(x). Note that the equation y = B(x, z) makes sense projectively because y has weight g + 1 (as always for hyperelliptic curves in Magma).

The point on Jac(C) corresponding to < (a(x), b(x), d) > is then D minus a multiple of the divisor at infinity on C. So, when there is a single point P at infinity, we have D - d P. Otherwise, there is a Q-rational divisor P + ∞ + P - ∞ consisting of the two points at infinity; in this case d is required to be even and we have D - (d/2) (P + ∞ + P - ∞).

All points on Jac(C) can be expressed in this way, except when g is odd and there are no rational points at infinity (in which case the extra points can not be created in Magma, and arithmetic on points is not implemented). There is a uniquely determined "reduced" triple representing each point, which Magma uses to represent any point it encounters.

See the examples in the following section ("Creation of Points").

Technical details: In order to make sense, a triple < (a(x), b(x), d) > is required to satisfy:

(a)
a(x) is monic of degree at most g;
(b)
b(x) has degree at most g + 1, and a(x) divides b(x)2 + h(x)b(x) - f(x), where h(x) and f(x) are the defining polynomials of C;
(c)
d is a positive integer with deg(a(x)) ≤d ≤g + 1, such that the degree of b(x)2 + h(x)b(x) - f(x) is less than or equal to 2g + 2 - d + deg(a(x)).

For uniqueness of representation, in the case of one point at infinity, (the odd degree case, in particular), we require that d = deg(a(x)).

A triple < (a(x), b(x), d) > is reduced to a canonical representative as follows:

(a)
If d = deg(a(x)), we reduce b mod a, so deg(b(x)) < deg(a(x)).
(b)
If a = 1, we assume that deg(B(1, z)) = d.
(c)
A unique representative for b(x) is found such that the coefficient of xk in b is zero for deg(a(x)) ≤k ≤deg(a(x)) + g + 1 - d.

For certain models of C, not all rational points can be represented in the above form (with the restrictions on d) or the representation is not unique. The bad cases are g odd and (i) 0 or (ii) 2 rational points at infinity.

In case (i) we would have to allow divisors with deg(a(x)) = g + 1, which give linear pencils of equivalent divisors, and in case (ii) a Jacobian point with d = g + 1 has precisely two canonical representatives.

In case (i), there is no obvious canonical representative for these additional Jacobian points and currently Magma does not deal with them. Consequently, arithmetic is not implemented in this case.

However, in case (ii), the two representatives correspond to the two distinguished elements of a linear pencil of divisors which include contributions from one or the other point at infinity. By arbitrarily selecting one of these two infinite points as the default, a unique representative can be chosen. This is now performed internally by Magma, the default point being chosen at the time of creation of the Jacobian. Arithmetic is also implemented. Note that in this case, extra work is generally required in point addition to keep track of points at infinity and final reduction to the unique representation. If very fast addition is crucial when C is of this type, it is generally wise to use an isomorphic model with exactly one point at infinity, if possible, by moving a rational Weierstrass point to infinity.

Contents

Creation of Points

J ! 0 : JacHyp, RngIntElt -> JacHypPt
Id(J) : JacHyp -> JacHypPt
Identity(J) : JacHyp -> JacHypPt
The identity element on the Jacobian J.
J ! [a, b] : JacHyp, [ RngUPolElt ] -> JacHypPt
elt< J | a, b > : JacHyp, RngUPolElt, RngUPolElt -> JacHypPt
elt< J | [a, b] > : JacHyp, [ RngUPolElt ] -> JacHypPt
elt< J | a, b, d > : JacHyp, RngUPolElt, RngUPolElt, RngIntElt -> JacHypPt
elt< J | [a, b], d > : JacHyp, [ RngUPolElt ], RngIntElt -> JacHypPt
The point on the Jacobian J defined by the polynomials a and b and the positive integer d; if not specified then d is taken to be deg(a).
P - Q : PtHyp, PtHyp -> JacHypPt
J ! [P, Q] : JacHyp, [PtHyp] -> JacHypPt
elt< J | P, Q > : JacHyp, PtHyp, PtHyp -> JacHypPt
For points P and Q on a hyperelliptic curve, this constructs the image of the divisor class [P - Q] as a point on its Jacobian J.
J ! [S, T] : JacHyp, [SeqEnum] -> JacHypPt
elt< J | S, T > : JacHyp, [PtHyp], [PtHyp] -> JacHypPt
Given two sequences S = [Pi] and T = [Qi] of points on the hyperelliptic curve with Jacobian J, each of length n, this returns the image of the divisor class ∑i [Pi] - ∑i [Qi] as a point on the Jacobian J.
JacobianPoint(J, D) : JacHyp, DivCrvElt -> JacHypPt
The point on the Jacobian J (of a hyperelliptic curve C) associated to the divisor D on C. If D does not have degree 0, then a suitable multiple of the divisor at infinity is subtracted. When the divisor at infinity on C has even degree, D is required to have even degree.

The function works for any divisor such that the corresponding point is definable in Magma. It is not implemented for characteristic 2.

J ! P : JacHyp, JacHypPt -> JacHypPt
Given a point P on a Jacobian J', construct the image of P on J, where J is a base extension of J'.
Points(J, a, d) : JacHyp, RngUPolElt, RngIntElt -> SetIndx
RationalPoints(J, a, d) : JacHyp, RngUPolElt, RngIntElt -> SetIndx
Find all points on the Jacobian J with first component a and degree d. Only implemented for genus 2 curves of the form y2 = f(x).

Example CrvHyp_point_creation_jacobian (H134E13)

Points on y2 = x6 - 3x - 1 and their images on the Jacobian.

> _<x> := PolynomialRing(Rationals());
> C := HyperellipticCurve(x^6-3*x-1);
> C;
Hyperelliptic Curve defined by y^2 = x^6 - 3*x - 1 over Rational Field
> J := Jacobian(C);
> J;
Jacobian of Hyperelliptic Curve defined by y^2 = x^6 - 3*x - 1 over Rational Field
Find some points on C and map them to J (using the first point to define the map C -> J):
> ptsC := Points(C : Bound := 100);
> ptsC;
{@ (1 : -1 : 0), (1 : 1 : 0), (-1 : -1 : 3), (-1 : 1 : 3) @}
> ptsJ := [ ptsC[i] - ptsC[1] : i in [2,3,4] ];
> ptsJ;
[ (1, x^3, 2), (x + 1/3, x^3, 2), (x + 1/3, x^3 + 2/27, 2) ]
We recreate the first of these, giving it in Magma's notation (which is read as "the divisor of degree 2 on C determined by z2 = y - x3 = 0").
> pt1 := elt< J | [1,x^3], 2 >;
> pt1 eq ptsJ[1];
true
The degree of the divisor must be specified here, otherwise it is assumed to be deg(a(x)) = 0:
> pt1 := J! [1,x^3];
> pt1; pt1 eq J!0;
(1, 0, 0)
true

Example CrvHyp_point_creation_jacobian2 (H134E14)

We define the nontrivial 2-torsion point on the Jacobian of C: y2 = (x2 + 1)(x6 + 7).

The divisor (y=0, x2 + 1=0) should define the only nontrivial 2-torsion point in J(Q).

> _<x> := PolynomialRing(Rationals());
> C := HyperellipticCurve( (x^2+1)*(x^6+7) );
> J := Jacobian(C);
Define the point corresponding to the divisor given by x2 + 1=0, y=0 on C. We can use the simpler syntax, not specifying that the degree of the divisor is 2, because this equals the degree of a(x)=x2 + 1.
> Ptors := J![x^2+1, 0];
> Ptors;
(x^2 + 1, 0, 2)
The alternative syntax would be as follows.
> Ptors1 := elt< J | [x^2+1, 0], 2 >;
> Ptors eq Ptors1;  // Are they the same?
true
Check that Ptors has order 2:
> 2*Ptors;
(1, 0, 0)
> $1 eq J!0;  // Is the previous result really the trivial point on J?
true
> Order(Ptors);  // Just to be absolutely sure ...
2
The other factor x6 + 7 will give the same point on J:
> Ptors2 := J![x^6+7,0];
(x^2 + 1, 0, 2)
Note that Magma returned the unique reduced triple representing this point, which means we can easily check whether or not it is the same point as Ptors.

For alternative ways to obtain 2-torsion points, see the section on torsion.

Example CrvHyp_point_creation_jacobian3 (H134E15)

A Jacobian with a point not coming from the curve.

Any curve containing the point (Sqrt(2), Sqrt(2)) has a Q-rational divisor D := (Sqrt(2), Sqrt(2)) + ( - Sqrt(2), - Sqrt(2)). This will give a nontrivial point in J(Q). For instance, take y2 = x6 - 6.

> _<x> := PolynomialRing(Rationals());
> C := HyperellipticCurve(x^6-6);
> J := Jacobian(C);
The divisor D has degree 2 and is given by x2 - 2 = y - x = 0.
> pt := J![x^2-2, x];
> pt;  // What is pt?
(x^2 - 2, x, 2)
> Parent(pt);  // Where does pt live?
Jacobian of Hyperelliptic Curve defined by y^2 = x^6 - 6 over Rational Field
> pt eq J!0;  // Is pt equal to 0 on J?
false
> Order(pt);
0
This means that P has infinite order in J(Q).

Alternatively we can construct the same point by first constructing the divisor, giving it as an ideal in the homogeneous coordinate ring in which C lives. (Note that Y has weight 3).

> P<X,Y,Z> := CoordinateRing(Ambient(C));
> D := Divisor(C, ideal<P | X^2-2*Z^2, Y-X*Z^2> );
> pt1 := J!D;
> pt eq pt1;
true

Random Points

Random(J) : JacHyp -> JacHypPt
A random point on a Jacobian J of a hyperelliptic curve defined over a finite field.

Booleans and Predicates for Points

For each of the following functions, P and Q are points on the Jacobian of a hyperelliptic curve.

P eq Q : JacHypPt, JacHypPt -> BoolElt
Returns true if and only if the points P and Q on the same Jacobian are equal.
P ne Q : JacHypPt, JacHypPt -> BoolElt
Returns false if and only if the two points P and Q on the same Jacobian are equal.
IsZero(P) : JacHypPt -> BoolElt
IsIdentity(P) : JacHypPt -> BoolElt
Returns true if and only if P is the zero element of the Jacobian.

Access Operations

P[i] : JacHypPt, RngIntElt -> RngElt
Given an integer 1 ≤i ≤2, this returns the i-th defining polynomial for the Jacobian point P. For i = 3, this returns the degree d of the associated reduced divisor.
Eltseq(P) : PtHyp -> SeqEnum, RngIntElt
ElementToSequence(P) : PtHyp -> SeqEnum, RngIntElt
Given a point P on the Jacobian J of a hyperelliptic curve, the function returns firstly, a sequence containing the two defining polynomials for the divisor associated to P and secondly, the degree of the divisor.

Arithmetic of Points

For each of the following functions, P and Q are points on the Jacobian of a hyperelliptic curve. The arithmetic operators on points can be called for curves over any base field, but the computations will be precise only for exact fields (i.e., not reals, complexes, p-adics or derivatives of these) so should really only be used for such fields.

- P : JacHypPt -> JacHypPt
The additive inverse of the point P on the Jacobian.
P + Q : JacHypPt, JacHypPt -> JacHypPt
The sum of the points P and Q on the same Jacobian.
P +:= Q : JacHypPt, JacHypPt ->
Given two points P and Q on the same Jacobian, set P equal to their sum.
P - Q : JacHypPt, JacHypPt -> JacHypPt
The difference of the points P and Q on the same Jacobian.
P -:= Q : JacHypPt, JacHypPt ->
Given two points P and Q on the same Jacobian, set P equal to the difference P - Q.
n * P : RngIntElt, JacHypPt -> JacHypPt
P * n : JacHypPt, RngIntElt -> JacHypPt
The n-th multiple of P in the group of rational points on the Jacobian.
P *:= n : JacHypPt, RngIntElt ->
Set P equal to the n-th multiple of itself.

Order of Points on the Jacobian

Order(P) : JacHypPt -> RngIntElt
Returns the order of the point P on the Jacobian J of a hyperelliptic curve defined over a finite field, the rationals or an algebraic number field, or 0 if P has infinite order. This first computes #J when J is defined over a finite field.
Order(P, l, u) : JacHypPt, RngIntElt, RngIntElt -> RngIntElt
    Alg: MonStg                         Default: "Shanks"
    UseInversion: BoolElt               Default: true
Returns the order of the point P where u and l are bounds for the order of P or for the order of J the parent of P. This does not compute #J.

If the parameter Alg is set to "Shanks" then the generic Shanks algorithm is used, otherwise, when Alg is "PollardRho", a Pollard-Rho variant is used (see [GH00]).

If UseInversion is true the search space is halved by using point negation.

Order(P, l, u, n, m) : JacHypPt, RngIntElt, RngIntElt ,RngIntElt, RngIntElt -> RngIntElt
    Alg: MonStg                         Default: "Shanks"
    UseInversion: BoolElt               Default: true
Returns the order of the point P where u and l are bounds for the order of P or for the order of J the parent of P and where n and m are such that the group order is n mod m. This does not compute #J.

The two parameters Alg and UseInversion have the same use as in the previous function.

HasOrder(P, n) : JacHypPt, RngIntElt -> BoolElt
Given a point P on the Jacobian J of a hyperelliptic curve and a positive integer n, this returns true if the order of the point is n.

Frobenius

Frobenius(P, k) : JacHypPt, FldFin -> JacHypPt
    Check: BoolElt                      Default: true
Given a point P that lies on the Jacobian J of a hyperelliptic curve that is defined over the finite field k=Fq, determine the image of P under the Frobenius map x -> xq. If Check is true, Magma verifies that the Jacobian of P is defined over k.

Weil Pairing

WeilPairing(P, Q, m) : JacHypPt, JacHypPt, RngIntElt -> RngElt
Computes the Weil pairing of P and Q, where P and Q are m-torsion points on the 2-dimensional Jacobian J defined over a finite field.

Example CrvHyp_Jac_WeilPairing (H134E16)

The following illustrates the use of the Weil Pairing in the MOV-reduction of the discrete logarithm problem on a Jacobian.

> PP<x>:=PolynomialRing(GF(2));
> h := PP!1;
> f := x^5 + x^4 + x^3 + 1;
> J := Jacobian(HyperellipticCurve(f,h));  // a supersingular curve
> Jext := BaseExtend(J, 41);
> Factorization(#Jext);
[ <7, 1>, <3887047, 1>, <177722253954175633, 1> ]
> m := 177722253954175633;                 // some big subgroup order
> cofact := 3887047*7;
> P := cofact*Random(Jext);
> Q := 876213876263897634*P;               // Q in <P>

Say we want to recompute the logarithm of Q in base P.

> Jext2 := BaseExtend(Jext, 6);            // go to an ext of deg 6
> NJ := #Jext2;
>
> R := Random(Jext2);
> R *:= NJ div m^Valuation(NJ, m);
>
> eP := WeilPairing(Jext2!P, R, m);
> eQ := WeilPairing(Jext2!Q, R, m);
> assert eP^876213876263897634 eq eQ;

So the discrete log problem on the Jacobian has been reduced to a discrete log problem in a finite field.

V2.28, 13 July 2023