NAME

libPARI - Technical Reference Guide for Elliptic curves and arithmetic geometry

This chapter is quite short, but is added as a placeholder, since we expect the library to expand in that direction.


Elliptic curves

Elliptic curves are represented in the Weierstrass model

   (E): y^2z + a_1xyz + a_3 yz = x^3 + a_2 x^2z + a_4 xz^2 + a_6z^3,

by the 5-tuple [a_1,a_2,a_3,a_4,a_6]. Points in the projective plane are represented as follows: the point at infinity (0:1:0) is coded as [0], a finite point (x:y:1) outside the projective line at infinity z = 0 is coded as [x,y]. Note that other points at infinity than (0:1:0) cannot be represented; this is harmless, since they do not belong to any of the elliptic curves E above.

Points on the curve are just projective points as described above, they are not tied to a curve in any way: the same point may be used in conjunction with different curves, provided it satisfies their equations (if it does not, the result is usually undefined). In particular, the point at infinity belongs to all elliptic curves.

As with factor for polynomial factorization, the 5-tuple [a_1,a_2,a_3,a_4,a_6] implicitly defines a base ring over which the curve is defined. Point coordinates must be operation-compatible with this base ring (gadd, gmul, gdiv involving them should not give errors).

Types of elliptic curves

We call a 5-tuble as above an ell5; most functions require an ell structure, as returned by ellinit, which contains additional data (usually dynamically computed as needed), depending on the base field.

GEN ellinit(GEN E, GEN D, long prec), returns an ell structure, associated to the elliptic curve E : either an ell5, a pair [a_4,a_6] or a t_STR in Cremona's notation, e.g. "11a1". The optional D (NULL to omit) describes the domain over which the curve is defined.

Type checking

void checkell(GEN e) raise an error unless e is a ell.

void checkell5(GEN e) raise an error unless e is an ell, a smallell or an ell5.

void checkellpt(GEN z) raise an error unless z is a point (either finite or at infinity).

long ell_get_type(GEN e) returns the domain type over which the curve is defined, one of

t_ELL_Q the field of rational numbers;

t_ELL_Qp the field of p-adic numbers, for some prime p;

t_ELL_Fp a prime finite field, base field elements are represented as Fp (t_INT reduced modulo p);

t_ELL_Fq a non-prime finite field (a prime finite field can also be represented by this subtype, but this is inefficient), base field elements are represented as t_FFELT;

t_ELL_Rg none of the above.

void checkell_Fq(GEN e) checks whether e is an ell, defined over a finite field (either prime or non-prime), raises pari_err_TYPE otherwise.

void checkell_Q(GEN e) checks whether e is an ell, defined over Q, raises pari_err_TYPE otherwise.

void checkell_Qp(GEN e) checks whether e is an ell, defined over some Q_p, raises pari_err_TYPE otherwise.

Extracting info from an ell structure

These functions expect an ell argument. If the required data is not part of the structure, it is computed then inserted, and the new value is returned.

All domains

GEN ell_get_a1(GEN e)

GEN ell_get_a2(GEN e)

GEN ell_get_a3(GEN e)

GEN ell_get_a4(GEN e)

GEN ell_get_a6(GEN e)

GEN ell_get_b2(GEN e)

GEN ell_get_b4(GEN e)

GEN ell_get_b6(GEN e)

GEN ell_get_b8(GEN e)

GEN ell_get_c4(GEN e)

GEN ell_get_c6(GEN e)

GEN ell_get_disc(GEN e)

GEN ell_get_j(GEN e)

Curves over Q

GEN ellQ_get_N(GEN e) returns the curve conductor

void ellQ_get_Nfa(GEN e, GEN *N, GEN *faN) sets N to the conductor and faN to its factorization

long ellrootno_global(GEN e) returns [c, [c_{p_1},...,c_{p_k}]], where the t_INT c belongs to {-1,1} is the global root number, and the c_{p_i} are the local root numbers at all prime divisors of the conductor, ordered as in faN above.

GEN elldatagenerators(GEN E) returns generators for E(Q) extracted from Cremona's table.

GEN ellanal_globalred(GEN e, GEN *v) takes an ell over Q and returns a global minimal model E (in ellinit form, over Q) for e suitable for analytic computations related to the curve L series: it contains ellglobalred data, as well as global and local root numbers. If v is not NULL, set *v to the needed change of variable: NULL if e was already the standard minimal model, such that E = ellchangecurve(e,v) otherwise. Compared to the direct use of ellchangecurve followed by ellrootno, this function avoids converting unneeded dynamic data and avoids potential memory leaks (the changed curve would have had to be deleted using obj_free). The original curve e is updated as well with the same information.

Curves over Q_p

GEN ellQp_get_p(GEN E) returns p

long ellQp_get_prec(GEN E) returns the default p-adic accuracy to which we must compute approximate results associated to E.

GEN ellQp_get_zero(GEN x) returns O(p^n), where n is the default p-adic accuracy as above.

The following functions are only defined when E has multiplicative reduction (Tate curves):

GEN ellQp_Tate_uniformization(GEN E, long prec) returns a t_VEC containing u^2, u, q, [a,b], at p-adic precision prec.

GEN ellQp_u(GEN E, long prec) returns u.

GEN ellQp_u2(GEN E, long prec) returns u^2.

GEN ellQp_q(GEN E, long prec) returns the Tate periode q.

GEN ellQp_ab(GEN E, long prec) returns [a,b].

GEN ellQp_root(GEN E, long prec) returns e_1.

Curves over a finite field F_q

GEN ellff_get_p(GEN E) returns the characteristic

GEN ellff_get_field(GEN E) returns p if F_q is a prime field, and a t_FFELT belonging to F_q otherwise.

GEN ellff_get_card(GEN E) returns #E(F_q)

GEN ellff_get_gens(GEN E) returns a minimal set of generators for E(F_q).

GEN ellff_get_group(GEN E) returns ellgroup(E).

GEN ellff_get_o(GEN E) returns [d, factor{d}], where d is the exponent of E(F_q).

GEN ellff_get_a4a6(GEN E) returns a canonical ``short model'' for E, and the corresponding change of variable [u,r,s,t]. For p != 2,3, this is [A_4,A_6,[u,r,s,t]], corresponding to y^2 = x^3 + A_4x + A_6, where A_4 = -27c_4, A_6 = -54c_6, [u,r,s,t] = [6, 3b_2,3a_1,108a_3].

@3* If p = 3 and the curve is ordinary (b_2 != 0), this is [[b_2], A_6, [1,v,-a_1,-a_3]], corresponding to

  y^2 = x^3 + b_2 x + A_6,

where v = b_4/b_2, A_6 = b_6 - v(b_4+v^2).

@3* If p = 3 and the curve is supersingular (b_2 = 0), this is [-b_4, b_6, [1,0,-a_1,-a_3]], corresponding to

  y^2 = x^3 + 2b_4 x + b_6.

@3* If p = 2 and the curve is ordinary (a_1 != 0), return [A_2,A_6,[a_1^{-1}, da_1^{-2}, 0, (a_4+d^2)a_1^{-1}]], corresponding to

   y^2 + xy = x^3 + A_2 x^2 + A_6,

where d = a_3/a_1, a_1^2 A_2 = (a_2 + d) and

   a_1^6 A_6 = d^3 + a_2 d^2 + a_4 d + a_6 + (a_4^2 + d^4)a_1^{-2}.

@3* If p = 2 and the curve is supersingular (a_1 = 0, a_3 != 0), return [[a_3, A_4, 1/a_3], A_6, [1,a_2,0,0]], corresponding to

   y^2 + a_3 y = x^3 + A_4 x + A_6,

where A_4 = a_2^2 + a_4, A_6 = a_2a_4 + a_6. The value 1/a_3 is included in the vector since it is frequently needed in computations.

Curves over C (This includes curves over Q!)

long ellR_get_prec(GEN E) returns the default accuracy to which we must compute approximate results associated to E.

GEN ellR_ab(GEN E, long prec) returns [a,b]

GEN ellR_omega(GEN x, long prec) returns periods [omega_1,omega_2].

GEN ellR_eta(GEN E, long prec) returns quasi-periods [eta_1,eta_2].

GEN ellR_roots(GEN E, long prec) returns [e_1,e_2,e_3]. If E is defined over R, then e_1 is real. If furthermore disc E > 0, then e_1 > e_2 > e_3.

long ellR_get_sign(GEN E) if E is defined over R returns the signe of its discriminant, otherwise return 0.

Points

int ell_is_inf(GEN z) tests whether the point z is the point at infinity.

GEN ellinf() returns the point at infinity [0].

Change of variables

GEN ellchangeinvert(GEN w) given a change of variables w = [u,r,s,t], returns the inverse change of variables w', such that if E' = ellchangecurve(E, w), then E = ellchangecurve(E, w').

Functions to handle elliptic curves over finite fields

Tolerant routines

GEN ellap(GEN E, GEN p) given a prime number p and an elliptic curve defined over Q or Q_p (assumed integral and minimal at p), computes the trace of Frobenius a_p = p+1 - #E(F_p). If E is defined over a non-prime finite field F_q, ignore p and return q+1 - #E(F_q). When p is implied (E defined over Q_p or a finite field), p can be omitted (set to NULL).

GEN ellsea(GEN E, GEN p, long s) available if the seadata package is installed. This function returns #E(F_p), using the Schoof-Elkies-Atkin algorithm; it is called by ellap: same conditions as above for E, except that t_ELL_Fq are not allowed. The extra flag s, if set to a non-zero value, causes the computation to return gen_0 (an impossible cardinality) if one of the small primes \ell > s divides the curve order. For cryptographic applications, where one is usually interested in curves of prime order, setting s = 1 efficiently weeds out most uninteresting curves; if curves of order a power of 2 times a prime are acceptable, set s = 2. There is no guarantee that the resulting cardinality is prime, only that it has no small prime divisor larger than s.

Curves defined a non-prime finite field

In this subsection, we assume that ell_get_type(E) is t_ELL_Fq. (As noted above, a curve defined over Z/pZ can be represented as a t_ELL_Fq.)

GEN FF_ellmul(GEN E, GEN P, GEN n) returns [n]P where n is an integer and P is a point on the curve E.

GEN FF_ellrandom(GEN E) returns a random point in E(F_q). This function never returns the point at infinity, unless this is the only point on the curve.

GEN FF_ellorder(GEN E, GEN P, GEN o) returns the order of the point P, where o is a multiple of the order of P, or its factorization.

GEN FF_ellcard(GEN E) returns #E(F_q).

GEN FF_ellgens(GEN E) returns the generators of the group E(F_q).

GEN FF_elllog(GEN E, GEN P, GEN G, GEN o) Let G be a point of order o, return e such that [e]P = G. If e does not exists, the result is undefined.

GEN FF_ellgroup(GEN E) returns the Abelian group E(F_q) in the form [h,cyc,gen].

GEN FF_ellweilpairing(GEN E, GEN P, GEN Q, GEN m) returns the Weil pairing of the points of m-torsion P and Q.

GEN FF_elltatepairing(GEN E, GEN P, GEN Q, GEN m) returns the Tate pairing of P and Q, where [m]P = 0.


Arithmetic on elliptic curve over a finite field in simple form

The functions in this section no longer operate on elliptic curve structures, as seen up to now. They are used to implement those higher-level functions without using cached information and thus require suitable explicitly enumerated data.

Helper functions

GEN elltrace_extension(GEN t, long n, GEN q) Let E some elliptic curve over F_q such that the trace of the Frobenius is t, returns the trace of the Frobenius over F_q^n.

Elliptic curves over F_p, p > 3

Let p a prime number and E the elliptic curve given by the equation E:y^2 = x^3+a_4 x+a_6, with a_4 and a_6 in F_p. A FpE is a point of E(F_p). Since an affine point and a_4 determine an unique a6, most functions do not take a_6 as an argument. A FpE is either the point at infinity (ellinf()) or a FpV whith two components. The parameters a_4 and a_6 are given as t_INTs when required.

GEN Fp_ellj(GEN a4, GEN a6, GEN p) returns the j-invariant of the curve E.

GEN Fp_ellcard(GEN a4, GEN a6, GEN p) returns the cardinal of the group E(F_p).

GEN Fp_ellcard_SEA(GEN a4, GEN a6, GEN p, long s) same as ellsea when only [a_4,a_6] are given.

GEN Fq_ellcard_SEA(GEN a4, GEN a6, GEN q, GEN T, GEN p, long s) same as ellsea when only [a_4,a_6] are given, over F_p[t]/(T). Assume p != 2,3.

GEN Fp_ffellcard(GEN a4, GEN a6, GEN q, long n, GEN p) returns the cardinal of the group E(F_q) where q = p^n.

GEN Fp_ellgroup(GEN a4, GEN a6, GEN N, GEN p, GEN *pt_m) returns the group structure D of the group E(F_p), which is assumed to be of order N and set *pt_m = m.

GEN Fp_ellgens(GEN a4, GEN a6, GEN ch, GEN D, GEN m, GEN p) returns generators of the group E(F_p) with the base change ch (see FpE_changepoint), where D and m are as returned by Fp_ellgroup.

GEN Fp_elldivpol(GEN a4, GEN a6, long n, GEN p) returns the n-division polynomial of the elliptic curve E.

FpE

GEN FpE_add(GEN P, GEN Q, GEN a4, GEN p) returns the sum P+Q in the group E(F_p), where E is defined by E:y^2 = x^3+a_4 x+a_6, for any value of a_6 compatible with the points given.

GEN FpE_sub(GEN P, GEN Q, GEN a4, GEN p) returns P-Q.

GEN FpE_dbl(GEN P, GEN a4, GEN p) returns 2.P.

GEN FpE_neg(GEN P, GEN p) returns -P.

GEN FpE_mul(GEN P, GEN n, GEN a4, GEN p) return n.P.

GEN FpE_changepoint(GEN P, GEN m, GEN a4, GEN p) returns the image Q of the point P on the curve E:y^2 = x^3+a_4 x+a_6 by the coordinate change m (which is a FpV).

GEN FpE_changepointinv(GEN P, GEN m, GEN a4, GEN p) returns the image Q on the curve E:y^2 = x^3+a_4 x+a_6 of the point P by the inverse of the coordinate change m (which is a FpV).

GEN random_FpE(GEN a4, GEN a6, GEN p) returns a random point on E(F_p), where E is defined by E:y^2 = x^3+a_4 x+a_6.

GEN FpE_order(GEN P, GEN o, GEN a4, GEN p) returns the order of P in the group E(F_p), where o is a multiple of the order of P, or its factorization.

GEN FpE_log(GEN P, GEN G, GEN o, GEN a4, GEN p) Let G be a point of order o, return e such that e.P = G. If e does not exists, the result is currently undefined.

GEN FpE_tatepairing(GEN P, GEN Q, GEN m, GEN a4, GEN p) returns the Tate pairing of the point of m-torsion P and the point Q.

GEN FpE_weilpairing(GEN P, GEN Q, GEN m, GEN a4, GEN p) returns the Weil pairing of the points of m-torsion P and Q.

GEN FpE_to_mod(GEN P, GEN p) returns P as a vector of t_INTMODs.

GEN RgE_to_FpE(GEN P, GEN p) returns the FpE obtained by applying Rg_to_Fp coefficientwise.

Fle

Let p be a prime ulong, and E the elliptic curve given by the equation E:y^2 = x^3+a_4 x+a_6, where a_4 and a_6 are ulong. A Fle is either the point at infinity (ellinf()), or a Flv with two components.

long Fl_elltrace(ulong a4, ulong a6, ulong p) returns the trace t of the Frobenius of E(F_p). The cardinal of E(F_p) is thus p+1-t, which might not fit in a ulong.

GEN Fle_add(GEN P, GEN Q, ulong a4, ulong p)

GEN Fle_dbl(GEN P, ulong a4, ulong p)

GEN Fle_sub(GEN P, GEN Q, ulong a4, ulong p)

GEN Fle_mul(GEN P, GEN n, ulong a4, ulong p)

GEN Fle_mulu(GEN P, ulong n, ulong a4, ulong p)

GEN Fle_order(GEN P, GEN o, ulong a4, ulong p)

GEN random_Fle(ulong a4, ulong a6, ulong p)

Elliptic curves over F_{2^n}

Let T be an irreducible F2x and E the elliptic curve given by either the equation E:y^2+x*y = x^3+a_2 x^2+a_6, where a_2, a_6 are F2x in F_2[X]/(T) (ordinary case) or E:y^2+a_3*y = x^3+a_4 x+a_6, where a_3, a_4, a_6 are F2x in F_2[X]/(T) (supersingular case).

A F2xqE is a point of E(F_2[X]/(T)). In the supersingular case, the parameter a2 is actually the t_VEC [a_3,a_4,a_3^{-1}].

GEN F2xq_ellcard(GEN a2, GEN a6, GEN T) Return the order of the group E(F_2[X]/(T)).

GEN F2xq_ellgroup(GEN a2, GEN a6, GEN N, GEN T, GEN *pt_m) Return the group structure D of the group E(F_2[X]/(T)), which is assumed to be of order N and set *pt_m = m.

GEN F2xq_ellgens(GEN a2, GEN a6, GEN ch, GEN D, GEN m, GEN T) Returns generators of the group E(F_2[X]/(T)) with the base change ch (see F2xqE_changepoint), where D and m are as returned by F2xq_ellgroup.

F2xqE

GEN F2xqE_changepoint(GEN P, GEN m, GEN a2, GEN T) returns the image Q of the point P on the curve E:y^2+x*y = x^3+a_2 x^2+a_6 by the coordinate change m (which is a F2xqV).

GEN F2xqE_changepointinv(GEN P, GEN m, GEN a2, GEN T) returns the image Q on the curve E:y^2 = x^3+a_4 x+a_6 of the point P by the inverse of the coordinate change m (which is a F2xqV).

GEN F2xqE_add(GEN P, GEN Q, GEN a2, GEN T)

GEN F2xqE_sub(GEN P, GEN Q, GEN a2, GEN T)

GEN F2xqE_dbl(GEN P, GEN a2, GEN T)

GEN F2xqE_neg(GEN P, GEN a2, GEN T)

GEN F2xqE_mul(GEN P, GEN n, GEN a2, GEN T)

GEN random_F2xqE(GEN a2, GEN a6, GEN T)

GEN F2xqE_order(GEN P, GEN o, GEN a2, GEN T) returns the order of P in the group E(F_2[X]/(T)), where o is a multiple of the order of P, or its factorization.

GEN F2xqE_log(GEN P, GEN G, GEN o, GEN a2, GEN T) Let G be a point of order o, return e such that e.P = G. If e does not exists, the result is currently undefined.

GEN F2xqE_tatepairing(GEN P, GEN Q, GEN m, GEN a2, GEN T) returns the Tate pairing of the point of m-torsion P and the point Q.

GEN F2xqE_weilpairing(GEN Q, GEN Q, GEN m, GEN a2, GEN T) returns the Weil pairing of the points of m-torsion P and Q.

GEN RgE_to_F2xqE(GEN P, GEN T) returns the F2xqE obtained by applying Rg_to_F2xq coefficientwise.

Elliptic curves over F_q, small characteristic p > 2

Let p be a prime ulong, T an irreducible Flx mod p, and E the elliptic curve given by the equation E:y^2 = x^3+a_4 x+a_6, where a_4 and a_6 are Flx in F_p[X]/(T). A FlxqE is a point of E(F_p[X]/(T)).

In the special case p = 3, ordinary elliptic curves (j(E) != 0) cannot be represented as above, but admit a model E:y^2 = x^3+a_2 x^2+a_6 with a_2 and a_6 being Flx in F_3[X]/(T). In that case, the parameter a2 is actually stored as a t_VEC, [a_2], to avoid ambiguities.

GEN Flxq_ellj(GEN a4, GEN a6, GEN T, ulong p) returns the j-invariant of the curve E.

GEN Flxq_ellcard(GEN a4, GEN a6, GEN T, ulong p) returns the order of E(F_p[X]/(T)).

GEN Flxq_ellgroup(GEN a4, GEN a6, GEN N, GEN T, ulong p, GEN *pt_m) returns the group structure D of the group E(F_p[X]/(T)), which is assumed to be of order N and set *pt_m = m.

GEN Flxq_ellgens(GEN a4, GEN a6, GEN ch, GEN D, GEN m, GEN T, ulong p) returns generators of the group E(F_p[X]/(T)) with the base change ch (see FlxqE_changepoint), where D and m are as returned by Flxq_ellgroup.

FlxqE

GEN FlxqE_changepoint(GEN P, GEN m, GEN a4, GEN T, ulong p) returns the image Q of the point P on the curve E:y^2 = x^3+a_4 x+a_6 by the coordinate change m (which is a FlxqV).

GEN FlxqE_changepointinv(GEN P, GEN m, GEN a4, GEN T, ulong p) returns the image Q on the curve E:y^2 = x^3+a_4 x+a_6 of the point P by the inverse of the coordinate change m (which is a FlxqV).

GEN FlxqE_add(GEN P, GEN Q, GEN a4, GEN T, ulong p)

GEN FlxqE_sub(GEN P, GEN Q, GEN a4, GEN T, ulong p)

GEN FlxqE_dbl(GEN P, GEN a4, GEN T, ulong p)

GEN FlxqE_neg(GEN P, GEN T, ulong p)

GEN FlxqE_mul(GEN P, GEN n, GEN a4, GEN T, ulong p)

GEN random_FlxqE(GEN a4, GEN a6, GEN T, ulong p)

GEN FlxqE_order(GEN P, GEN o, GEN a4, GEN T, ulong p) returns the order of P in the group E(F_p[X]/(T)), where o is a multiple of the order of P, or its factorization.

GEN FlxqE_log(GEN P, GEN G, GEN o, GEN a4, GEN T, ulong p) Let G be a point of order o, return e such that e.P = G. If e does not exists, the result is currently undefined.

GEN FlxqE_tatepairing(GEN P, GEN Q, GEN m, GEN a4, GEN T, ulong p) returns the Tate pairing of the point of m-torsion P and the point Q.

GEN FlxqE_weilpairing(GEN P, GEN Q, GEN m, GEN a4, GEN T, ulong p) returns the Weil pairing of the points of m-torsion P and Q.

GEN RgE_to_FlxqE(GEN P, GEN T, ulong p) returns the FlxqE obtained by applying Rg_to_Flxq coefficientwise.

Elliptic curves over F_q, large characteristic

Let p be a prime number, T an irreducible polynomial mod p, and E the elliptic curve given by the equation E:y^2 = x^3+a_4 x+a_6 with a_4 and a_6 in F_p[X]/(T). A FpXQE is a point of E(F_p[X]/(T)).

GEN FpXQ_ellj(GEN a4, GEN a6, GEN T, GEN p) returns the j-invariant of the curve E.

GEN FpXQ_ellcard(GEN a4, GEN a6, GEN T, GEN p) returns the order of E(F_p[X]/(T)).

GEN FpXQ_ellgroup(GEN a4, GEN a6, GEN N, GEN T, GEN p, GEN *pt_m) Return the group structure D of the group E(F_p[X]/(T)), which is assumed to be of order N and set *pt_m = m.

GEN FpXQ_ellgens(GEN a4, GEN a6, GEN ch, GEN D, GEN m, GEN T, GEN p) Returns generators of the group E(F_p[X]/(T)) with the base change ch (see FpXQE_changepoint), where D and m are as returned by FpXQ_ellgroup.

GEN FpXQ_elldivpol(GEN a4, GEN a6, long n, GEN T, GEN p) returns the n-division polynomial of the elliptic curve E.

GEN Fq_elldivpolmod(GEN a4,GEN a6, long n, GEN h, GEN T, GEN p) returns the n-division polynomial of the elliptic curve E modulo the polynomial h.

FpXQE

GEN FpXQE_changepoint(GEN P, GEN m, GEN a4, GEN T, GEN p) returns the image Q of the point P on the curve E:y^2 = x^3+a_4 x+a_6 by the coordinate change m (which is a FpXQV).

GEN FpXQE_changepointinv(GEN P, GEN m, GEN a4, GEN T, GEN p) returns the image Q on the curve E:y^2 = x^3+a_4 x+a_6 of the point P by the inverse of the coordinate change m (which is a FpXQV).

GEN FpXQE_add(GEN P, GEN Q, GEN a4, GEN T, GEN p)

GEN FpXQE_sub(GEN P, GEN Q, GEN a4, GEN T, GEN p)

GEN FpXQE_dbl(GEN P, GEN a4, GEN T, GEN p)

GEN FpXQE_neg(GEN P, GEN T, GEN p)

GEN FpXQE_mul(GEN P, GEN n, GEN a4, GEN T, GEN p)

GEN random_FpXQE(GEN a4, GEN a6, GEN T, GEN p)

GEN FpXQE_log(GEN P, GEN G, GEN o, GEN a4, GEN T, GEN p) Let G be a point of order o, return e such that e.P = G. If e does not exists, the result is currently undefined.

GEN FpXQE_order(GEN P, GEN o, GEN a4, GEN T, GEN p) returns the order of P in the group E(F_p[X]/(T)), where o is a multiple of the order of P, or its factorization.

GEN FpXQE_tatepairing(GEN P,GEN Q, GEN m, GEN a4, GEN T, GEN p) returns the Tate pairing of the point of m-torsion P and the point Q.

GEN FpXQE_weilpairing(GEN P,GEN Q, GEN m, GEN a4, GEN T, GEN p) returns the Weil pairing of the points of m-torsion P and Q.

GEN RgE_to_FpXQE(GEN P, GEN T, GEN p) returns the FpXQE obtained by applying Rg_to_FpXQ coefficientwise.


Other curves

The following functions deal with hyperelliptic curves in weighted projective space \P_{(1,d,1)}, with coordinates (x,y,z) and a model of the form y^2 = T(x,z), where T is homogeneous of degree 2d, and squarefree. Thus the curve is nonsingular of genus d-1.

long hyperell_locally_soluble(GEN T, GEN p) assumes that T belongs to Z[X] is integral. Returns 1 if the curve is locally soluble over Q_p, 0 otherwise.

long nf_hyperell_locally_soluble(GEN nf, GEN T, GEN pr) let K be a number field, associated to nf, pr a prid associated to some maximal ideal p; assumes that T belongs to Z_K[X] is integral. Returns 1 if the curve is locally soluble over K_{p}.

\newpage