Elliptic curves over finite fields#
AUTHORS:
William Stein (2005): Initial version
Robert Bradshaw et al….
John Cremona (2008-02): Point counting and group structure for non-prime fields, Frobenius endomorphism and order, elliptic logs
Mariah Lenox (2011-03): Added set_order method
- class sage.schemes.elliptic_curves.ell_finite_field.EllipticCurve_finite_field(K, ainvs)#
Bases:
EllipticCurve_field
,HyperellipticCurve_finite_field
Elliptic curve over a finite field.
EXAMPLES:
sage: EllipticCurve(GF(101),[2,3]) Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Finite Field of size 101 sage: F = GF(101^2, 'a') sage: EllipticCurve([F(2),F(3)]) Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Finite Field in a of size 101^2
Elliptic curves over \(\ZZ/N\ZZ\) with \(N\) prime are of type “elliptic curve over a finite field”:
sage: F = Zmod(101) sage: EllipticCurve(F, [2, 3]) Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Ring of integers modulo 101 sage: E = EllipticCurve([F(2), F(3)]) sage: type(E) <class 'sage.schemes.elliptic_curves.ell_finite_field.EllipticCurve_finite_field_with_category'> sage: E.category() Category of schemes over Ring of integers modulo 101
Elliptic curves over \(\ZZ/N\ZZ\) with \(N\) composite are of type “generic elliptic curve”:
sage: F = Zmod(95) sage: EllipticCurve(F, [2, 3]) Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Ring of integers modulo 95 sage: E = EllipticCurve([F(2), F(3)]) sage: type(E) <class 'sage.schemes.elliptic_curves.ell_generic.EllipticCurve_generic_with_category'> sage: E.category() Category of schemes over Ring of integers modulo 95 sage: TestSuite(E).run(skip=["_test_elements"])
- abelian_group()#
Return the abelian group structure of the group of points on this elliptic curve.
See also
If you do not need the complete abelian group structure but only generators of the group, use
gens()
which can be much faster in some cases.This method relies on
gens()
, which uses random points on the curve and hence the generators are likely to differ from one run to another. However, the group is cached, so the generators will not change in any one run of Sage.OUTPUT:
an
AdditiveAbelianGroupWrapper
object encapsulating the abelian group of rational points on this elliptic curve
ALGORITHM:
We first call
gens()
to obtain a generating set \((P,Q)\). Letting \(P\) denote the point of larger order \(n_1\), we extend \(P\) to a basis \((P,Q')\) by computing a scalar \(x\) such that \(Q'=Q-[x]P\) has order \(n_2=\#E/n_1\). Finding \(x\) involves a (typically easy) discrete-logarithm computation.The complexity of the algorithm is the cost of factoring the group order, plus \(\Theta(\sqrt{\ell})\) for each prime \(\ell\) such that the rational \(\ell^\infty\)-torsion of
self
is isomorphic to \(\ZZ/\ell^r\times\ZZ/\ell^s\) with \(r>s>0\), times a polynomial in the logarithm of the base-field size.AUTHORS:
John Cremona: original implementation
Lorenz Panny (2021): current implementation
EXAMPLES:
sage: E = EllipticCurve(GF(11),[2,5]) sage: E.abelian_group() Additive abelian group isomorphic to Z/10 embedded in Abelian group of points on Elliptic Curve defined by y^2 = x^3 + 2*x + 5 over Finite Field of size 11
sage: E = EllipticCurve(GF(41),[2,5]) sage: E.abelian_group() Additive abelian group isomorphic to Z/22 + Z/2 ...
sage: F.<a> = GF(3^6,'a') sage: E = EllipticCurve([a^4 + a^3 + 2*a^2 + 2*a, 2*a^5 + 2*a^3 + 2*a^2 + 1]) sage: E.abelian_group() Additive abelian group isomorphic to Z/26 + Z/26 ...
sage: F.<a> = GF(101^3,'a') sage: E = EllipticCurve([2*a^2 + 48*a + 27, 89*a^2 + 76*a + 24]) sage: E.abelian_group() Additive abelian group isomorphic to Z/1031352 ...
The group can be trivial:
sage: E = EllipticCurve(GF(2),[0,0,1,1,1]) sage: E.abelian_group() Trivial group embedded in Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 + x + 1 over Finite Field of size 2
Of course, there are plenty of points if we extend the field:
sage: E.cardinality(extension_degree=100) 1267650600228231653296516890625
This tests the patch for trac ticket #3111, using 10 primes randomly selected:
sage: E = EllipticCurve('389a') sage: for p in [5927, 2297, 1571, 1709, 3851, 127, 3253, 5783, 3499, 4817]: ....: G = E.change_ring(GF(p)).abelian_group() sage: for p in prime_range(10000): # long time (19s on sage.math, 2011) ....: if p != 389: ....: G = E.change_ring(GF(p)).abelian_group()
This tests that the bug reported in trac ticket #3926 has been fixed:
sage: K.<i> = QuadraticField(-1) sage: OK = K.ring_of_integers() sage: P=K.factor(10007)[0][0] sage: OKmodP = OK.residue_field(P) sage: E = EllipticCurve([0,0,0,i,i+3]) sage: Emod = E.change_ring(OKmodP); Emod Elliptic Curve defined by y^2 = x^3 + ibar*x + (ibar+3) over Residue field in ibar of Fractional ideal (10007) sage: Emod.abelian_group() #random generators (Multiplicative Abelian group isomorphic to C50067594 x C2, ((3152*ibar + 7679 : 7330*ibar + 7913 : 1), (8466*ibar + 1770 : 0 : 1)))
- cardinality(algorithm=None, extension_degree=1)#
Return the number of points on this elliptic curve.
INPUT:
algorithm
– (optional) string:'pari'
– use the PARI C-library functionellcard
.'bsgs'
– use the baby-step giant-step method asimplemented in Sage, with the Cremona-Sutherland version of Mestre’s trick.
'exhaustive'
– naive point counting.'subfield'
– reduce to a smaller field, provided that the j-invariant lies in a subfield.'all'
– compute cardinality with both'pari'
and'bsgs'
; return result if they agree or raise aAssertionError
if they do not
extension_degree
– an integer \(d\) (default: 1): if the base field is \(\GF{q}\), return the cardinality ofself
over the extension \(\GF{q^d}\) of degree \(d\).
OUTPUT:
The order of the group of rational points of
self
over its base field, or over an extension field of degree \(d\) as above. The result is cached.EXAMPLES:
sage: EllipticCurve(GF(4, 'a'), [1,2,3,4,5]).cardinality() 8 sage: k.<a> = GF(3^3) sage: l = [a^2 + 1, 2*a^2 + 2*a + 1, a^2 + a + 1, 2, 2*a] sage: EllipticCurve(k,l).cardinality() 29
sage: l = [1, 1, 0, 2, 0] sage: EllipticCurve(k, l).cardinality() 38
An even bigger extension (which we check against Magma):
sage: EllipticCurve(GF(3^100, 'a'), [1,2,3,4,5]).cardinality() 515377520732011331036459693969645888996929981504 sage: magma.eval("Order(EllipticCurve([GF(3^100)|1,2,3,4,5]))") # optional - magma '515377520732011331036459693969645888996929981504'
sage: EllipticCurve(GF(10007), [1,2,3,4,5]).cardinality() 10076 sage: EllipticCurve(GF(10007), [1,2,3,4,5]).cardinality(algorithm='pari') 10076 sage: EllipticCurve(GF(next_prime(10**20)), [1,2,3,4,5]).cardinality() 100000000011093199520
The cardinality is cached:
sage: E = EllipticCurve(GF(3^100, 'a'), [1,2,3,4,5]) sage: E.cardinality() is E.cardinality() True
The following is very fast since the curve is actually defined over the prime field:
sage: k.<a> = GF(11^100) sage: E1 = EllipticCurve(k, [3,3]) sage: N1 = E1.cardinality(algorithm="subfield"); N1 137806123398222701841183371720896367762643312000384671846835266941791510341065565176497846502742959856128 sage: E1.cardinality_pari() == N1 True sage: E2 = E1.quadratic_twist() sage: N2 = E2.cardinality(algorithm="subfield"); N2 137806123398222701841183371720896367762643312000384656816094284101308193849980588362304472492174093035876 sage: E2.cardinality_pari() == N2 True sage: N1 + N2 == 2*(k.cardinality() + 1) True
We can count points over curves defined as a reduction:
sage: x = polygen(QQ) sage: K.<w> = NumberField(x^2 + x + 1) sage: EK = EllipticCurve(K, [0, 0, w, 2, 1]) sage: E = EK.base_extend(K.residue_field(2)) sage: E Elliptic Curve defined by y^2 + wbar*y = x^3 + 1 over Residue field in wbar of Fractional ideal (2) sage: E.cardinality() 7 sage: E = EK.base_extend(K.residue_field(w - 1)) sage: E.abelian_group() Trivial group embedded in Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 + 2*x + 1 over Residue field of Fractional ideal (w - 1)
sage: R.<x> = GF(17)[] sage: pol = R.irreducible_element(5) sage: k.<a> = R.residue_field(pol) sage: E = EllipticCurve(R, [1, x]).base_extend(k) sage: E Elliptic Curve defined by y^2 = x^3 + x + a over Residue field in a of Principal ideal (x^5 + x + 14) of Univariate Polynomial Ring in x over Finite Field of size 17 sage: E.cardinality() 1421004
- cardinality_bsgs(verbose=False)#
Return the cardinality of
self
over the base field.ALGORITHM: A variant of “Mestre’s trick” extended to all finite fields by Cremona and Sutherland, 2008.
Note
The Mestre-Schoof-Cremona-Sutherland algorithm may fail for a small finite number of curves over \(F_q\) for \(q\) at most 49, so for \(q<50\) we use an exhaustive count.
Quadratic twists are not implemented in characteristic 2 when \(j=0 (=1728)\); but this case is treated separately.
EXAMPLES:
sage: p = next_prime(10^3) sage: E = EllipticCurve(GF(p),[3,4]) sage: E.cardinality_bsgs() 1020 sage: E = EllipticCurve(GF(3^4,'a'),[1,1]) sage: E.cardinality_bsgs() 64 sage: F.<a> = GF(101^3,'a') sage: E = EllipticCurve([2*a^2 + 48*a + 27, 89*a^2 + 76*a + 24]) sage: E.cardinality_bsgs() 1031352
- cardinality_exhaustive()#
Return the cardinality of
self
over the base field.This simply adds up the number of points with each x-coordinate: only used for small field sizes!
EXAMPLES:
sage: p = next_prime(10^3) sage: E = EllipticCurve(GF(p),[3,4]) sage: E.cardinality_exhaustive() 1020 sage: E = EllipticCurve(GF(3^4,'a'),[1,1]) sage: E.cardinality_exhaustive() 64
- cardinality_pari()#
Return the cardinality of
self
using PARI.This uses pari:ellcard.
EXAMPLES:
sage: p = next_prime(10^3) sage: E = EllipticCurve(GF(p),[3,4]) sage: E.cardinality_pari() 1020 sage: K = GF(next_prime(10^6)) sage: E = EllipticCurve(K,[1,0,0,1,1]) sage: E.cardinality_pari() 999945
Since trac ticket #16931, this now works over finite fields which are not prime fields:
sage: k.<a> = GF(7^3) sage: E = EllipticCurve_from_j(a) sage: E.cardinality_pari() 318 sage: K.<a> = GF(3^20) sage: E = EllipticCurve(K,[1,0,0,1,a]) sage: E.cardinality_pari() 3486794310
- count_points(n=1)#
Return the cardinality of this elliptic curve over the base field or extensions.
INPUT:
n
(int) – a positive integer
OUTPUT:
If \(n=1\), returns the cardinality of the curve over its base field.
If \(n>1\), returns a list \([c_1, c_2, ..., c_n]\) where \(c_d\) is the cardinality of the curve over the extension of degree \(d\) of its base field.
EXAMPLES:
sage: p = 101 sage: F = GF(p) sage: E = EllipticCurve(F, [2,3]) sage: E.count_points(1) 96 sage: E.count_points(5) [96, 10368, 1031904, 104053248, 10509895776]
sage: F.<a> = GF(p^2) sage: E = EllipticCurve(F, [a,a]) sage: E.cardinality() 10295 sage: E.count_points() 10295 sage: E.count_points(1) 10295 sage: E.count_points(5) [10295, 104072155, 1061518108880, 10828567126268595, 110462212555439192375]
- frobenius()#
Return the frobenius of
self
as an element of a quadratic order.Note
This computes the curve cardinality, which may be time-consuming.
Frobenius is only determined up to conjugacy.
EXAMPLES:
sage: E = EllipticCurve(GF(11),[3,3]) sage: E.frobenius() phi sage: E.frobenius().minpoly() x^2 - 4*x + 11
For some supersingular curves, Frobenius is in Z:
sage: E = EllipticCurve(GF(25,'a'),[0,0,0,0,1]) sage: E.frobenius() -5
- frobenius_endomorphism()#
Return the \(q\)-power Frobenius endomorphism of this elliptic curve, where \(q\) is the cardinality of the (finite) base field.
EXAMPLES:
sage: F.<t> = GF(11^4) sage: E = EllipticCurve([t,t]) sage: E.frobenius_endomorphism() Frobenius endomorphism of degree 14641 = 11^4: From: Elliptic Curve defined by y^2 = x^3 + t*x + t over Finite Field in t of size 11^4 To: Elliptic Curve defined by y^2 = x^3 + t*x + t over Finite Field in t of size 11^4 sage: E.frobenius_endomorphism() == E.frobenius_isogeny(4) True
See also
- frobenius_order()#
Return the quadratic order Z[phi] where phi is the Frobenius endomorphism of the elliptic curve.
Note
This computes the curve cardinality, which may be time-consuming.
EXAMPLES:
sage: E = EllipticCurve(GF(11),[3,3]) sage: E.frobenius_order() Order in Number Field in phi with defining polynomial x^2 - 4*x + 11
For some supersingular curves, Frobenius is in Z and the Frobenius order is Z:
sage: E = EllipticCurve(GF(25,'a'),[0,0,0,0,1]) sage: R = E.frobenius_order() sage: R Order in Number Field in phi with defining polynomial x + 5 sage: R.degree() 1
- frobenius_polynomial()#
Return the characteristic polynomial of Frobenius.
The Frobenius endomorphism of the elliptic curve has quadratic characteristic polynomial. In most cases this is irreducible and defines an imaginary quadratic order; for some supersingular curves, Frobenius is an integer a and the polynomial is \((x-a)^2\).
Note
This computes the curve cardinality, which may be time-consuming.
EXAMPLES:
sage: E = EllipticCurve(GF(11),[3,3]) sage: E.frobenius_polynomial() x^2 - 4*x + 11
For some supersingular curves, Frobenius is in Z and the polynomial is a square:
sage: E = EllipticCurve(GF(25,'a'),[0,0,0,0,1]) sage: E.frobenius_polynomial().factor() (x + 5)^2
- gens()#
Return points which generate the abelian group of points on this elliptic curve.
The algorithm involves factoring the group order of
self
, but is otherwise (randomized) polynomial-time.(The points returned by this function are not guaranteed to be the same each time, although they should remain fixed within a single run of Sage unless
abelian_group()
is called.)OUTPUT: a tuple of points on the curve.
if the group is trivial: an empty tuple.
if the group is cyclic: a tuple with 1 point, a generator.
if the group is not cyclic: a tuple with 2 points, where the order of the first point equals the exponent of the group.
Warning
In the case of 2 generators \(P\) and \(Q\), it is not guaranteed that the group is the cartesian product of the 2 cyclic groups \(\langle P \rangle\) and \(\langle Q \rangle\). In other words, the order of \(Q\) is not as small as possible. If you really need a basis (rather than just a generating set) of the group, use
abelian_group()
.EXAMPLES:
sage: E = EllipticCurve(GF(11),[2,5]) sage: P = E.gens()[0]; P # random (0 : 7 : 1) sage: E.cardinality(), P.order() (10, 10) sage: E = EllipticCurve(GF(41),[2,5]) sage: E.gens() # random ((20 : 38 : 1), (25 : 31 : 1)) sage: E.cardinality() 44
If the abelian group has been computed, return those generators instead:
sage: E.abelian_group() Additive abelian group isomorphic to Z/22 + Z/2 embedded in Abelian group of points on Elliptic Curve defined by y^2 = x^3 + 2*x + 5 over Finite Field of size 41 sage: ab_gens = E.abelian_group().gens() sage: ab_gens == E.gens() True sage: E.gens()[0].order() 22 sage: E.gens()[1].order() 2
Examples with 1 and 0 generators:
sage: F.<a> = GF(3^6) sage: E = EllipticCurve([a, a+1]) sage: pts = E.gens() sage: len(pts) 1 sage: pts[0].order() == E.cardinality() True sage: E = EllipticCurve(GF(2), [0,0,1,1,1]) sage: E.gens() ()
This works over larger finite fields where
abelian_group()
may be too expensive:sage: k.<a> = GF(5^60) sage: E = EllipticCurve([a, a]) sage: len(E.gens()) 2 sage: E.cardinality() 867361737988403547206134229616487867594472 sage: a = E.gens()[0].order(); a # random 433680868994201773603067114808243933797236 sage: b = E.gens()[1].order(); b # random 30977204928157269543076222486303138128374 sage: lcm(a,b) 433680868994201773603067114808243933797236
- is_isogenous(other, field=None, proof=True)#
Return whether or not self is isogenous to other.
INPUT:
other
– another elliptic curve.field
(default None) – a field containing the base fields of the two elliptic curves into which the two curves may be extended to test if they are isogenous over this field. By default is_isogenous will not try to find this field unless one of the curves can be extended into the base field of the other, in which case it will test over the larger base field.proof
(default True) – this parameter is here only to be consistent with versions for other types of elliptic curves.
OUTPUT:
(bool) True if there is an isogeny from curve
self
to curveother
defined overfield
.EXAMPLES:
sage: E1 = EllipticCurve(GF(11^2,'a'),[2,7]); E1 Elliptic Curve defined by y^2 = x^3 + 2*x + 7 over Finite Field in a of size 11^2 sage: E1.is_isogenous(5) Traceback (most recent call last): ... ValueError: Second argument is not an Elliptic Curve. sage: E1.is_isogenous(E1) True sage: E2 = EllipticCurve(GF(7^3,'b'),[3,1]); E2 Elliptic Curve defined by y^2 = x^3 + 3*x + 1 over Finite Field in b of size 7^3 sage: E1.is_isogenous(E2) Traceback (most recent call last): ... ValueError: The base fields must have the same characteristic. sage: E3 = EllipticCurve(GF(11^2,'c'),[4,3]); E3 Elliptic Curve defined by y^2 = x^3 + 4*x + 3 over Finite Field in c of size 11^2 sage: E1.is_isogenous(E3) False sage: E4 = EllipticCurve(GF(11^6,'d'),[6,5]); E4 Elliptic Curve defined by y^2 = x^3 + 6*x + 5 over Finite Field in d of size 11^6 sage: E1.is_isogenous(E4) True sage: E5 = EllipticCurve(GF(11^7,'e'),[4,2]); E5 Elliptic Curve defined by y^2 = x^3 + 4*x + 2 over Finite Field in e of size 11^7 sage: E1.is_isogenous(E5) Traceback (most recent call last): ... ValueError: Curves have different base fields: use the field parameter.
When the field is given:
sage: E1 = EllipticCurve(GF(13^2,’a’),[2,7]); E1 Elliptic Curve defined by y^2 = x^3 + 2*x + 7 over Finite Field in a of size 13^2 sage: E1.is_isogenous(5,GF(13^6,’f’)) Traceback (most recent call last): … ValueError: Second argument is not an Elliptic Curve. sage: E6 = EllipticCurve(GF(11^3,’g’),[9,3]); E6 Elliptic Curve defined by y^2 = x^3 + 9*x + 3 over Finite Field in g of size 11^3 sage: E1.is_isogenous(E6,QQ) Traceback (most recent call last): … ValueError: The base fields must have the same characteristic. sage: E7 = EllipticCurve(GF(13^5,’h’),[2,9]); E7 Elliptic Curve defined by y^2 = x^3 + 2*x + 9 over Finite Field in h of size 13^5 sage: E1.is_isogenous(E7,GF(13^4,’i’)) Traceback (most recent call last): … ValueError: Field must be an extension of the base fields of both curves sage: E1.is_isogenous(E7,GF(13^10,’j’)) False sage: E1.is_isogenous(E7,GF(13^30,’j’)) False
- is_ordinary(proof=True)#
Return True if this elliptic curve is ordinary, else False.
INPUT:
proof
(boolean, default True) – If True, returns a proved result. If False, then a return value of True is certain but a return value of False may be based on a probabilistic test. See the documentation of the functionis_j_supersingular()
for more details.
EXAMPLES:
sage: F = GF(101) sage: EllipticCurve(j=F(0)).is_ordinary() False sage: EllipticCurve(j=F(1728)).is_ordinary() True sage: EllipticCurve(j=F(66)).is_ordinary() False sage: EllipticCurve(j=F(99)).is_ordinary() True
- is_supersingular(proof=True)#
Return True if this elliptic curve is supersingular, else False.
INPUT:
proof
(boolean, default True) – If True, returns a proved result. If False, then a return value of False is certain but a return value of True may be based on a probabilistic test. See the documentation of the functionis_j_supersingular()
for more details.
EXAMPLES:
sage: F = GF(101) sage: EllipticCurve(j=F(0)).is_supersingular() True sage: EllipticCurve(j=F(1728)).is_supersingular() False sage: EllipticCurve(j=F(66)).is_supersingular() True sage: EllipticCurve(j=F(99)).is_supersingular() False
- order(algorithm=None, extension_degree=1)#
Return the number of points on this elliptic curve.
INPUT:
algorithm
– (optional) string:'pari'
– use the PARI C-library functionellcard
.'bsgs'
– use the baby-step giant-step method asimplemented in Sage, with the Cremona-Sutherland version of Mestre’s trick.
'exhaustive'
– naive point counting.'subfield'
– reduce to a smaller field, provided that the j-invariant lies in a subfield.'all'
– compute cardinality with both'pari'
and'bsgs'
; return result if they agree or raise aAssertionError
if they do not
extension_degree
– an integer \(d\) (default: 1): if the base field is \(\GF{q}\), return the cardinality ofself
over the extension \(\GF{q^d}\) of degree \(d\).
OUTPUT:
The order of the group of rational points of
self
over its base field, or over an extension field of degree \(d\) as above. The result is cached.EXAMPLES:
sage: EllipticCurve(GF(4, 'a'), [1,2,3,4,5]).cardinality() 8 sage: k.<a> = GF(3^3) sage: l = [a^2 + 1, 2*a^2 + 2*a + 1, a^2 + a + 1, 2, 2*a] sage: EllipticCurve(k,l).cardinality() 29
sage: l = [1, 1, 0, 2, 0] sage: EllipticCurve(k, l).cardinality() 38
An even bigger extension (which we check against Magma):
sage: EllipticCurve(GF(3^100, 'a'), [1,2,3,4,5]).cardinality() 515377520732011331036459693969645888996929981504 sage: magma.eval("Order(EllipticCurve([GF(3^100)|1,2,3,4,5]))") # optional - magma '515377520732011331036459693969645888996929981504'
sage: EllipticCurve(GF(10007), [1,2,3,4,5]).cardinality() 10076 sage: EllipticCurve(GF(10007), [1,2,3,4,5]).cardinality(algorithm='pari') 10076 sage: EllipticCurve(GF(next_prime(10**20)), [1,2,3,4,5]).cardinality() 100000000011093199520
The cardinality is cached:
sage: E = EllipticCurve(GF(3^100, 'a'), [1,2,3,4,5]) sage: E.cardinality() is E.cardinality() True
The following is very fast since the curve is actually defined over the prime field:
sage: k.<a> = GF(11^100) sage: E1 = EllipticCurve(k, [3,3]) sage: N1 = E1.cardinality(algorithm="subfield"); N1 137806123398222701841183371720896367762643312000384671846835266941791510341065565176497846502742959856128 sage: E1.cardinality_pari() == N1 True sage: E2 = E1.quadratic_twist() sage: N2 = E2.cardinality(algorithm="subfield"); N2 137806123398222701841183371720896367762643312000384656816094284101308193849980588362304472492174093035876 sage: E2.cardinality_pari() == N2 True sage: N1 + N2 == 2*(k.cardinality() + 1) True
We can count points over curves defined as a reduction:
sage: x = polygen(QQ) sage: K.<w> = NumberField(x^2 + x + 1) sage: EK = EllipticCurve(K, [0, 0, w, 2, 1]) sage: E = EK.base_extend(K.residue_field(2)) sage: E Elliptic Curve defined by y^2 + wbar*y = x^3 + 1 over Residue field in wbar of Fractional ideal (2) sage: E.cardinality() 7 sage: E = EK.base_extend(K.residue_field(w - 1)) sage: E.abelian_group() Trivial group embedded in Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 + 2*x + 1 over Residue field of Fractional ideal (w - 1)
sage: R.<x> = GF(17)[] sage: pol = R.irreducible_element(5) sage: k.<a> = R.residue_field(pol) sage: E = EllipticCurve(R, [1, x]).base_extend(k) sage: E Elliptic Curve defined by y^2 = x^3 + x + a over Residue field in a of Principal ideal (x^5 + x + 14) of Univariate Polynomial Ring in x over Finite Field of size 17 sage: E.cardinality() 1421004
- plot(*args, **kwds)#
Draw a graph of this elliptic curve over a prime finite field.
INPUT:
*args, **kwds
– all other options are passed to the circle graphing primitive.
EXAMPLES:
sage: E = EllipticCurve(FiniteField(17), [0,1]) sage: P = plot(E, rgbcolor=(0,0,1))
- points()#
All the points on this elliptic curve. The list of points is cached so subsequent calls are free.
EXAMPLES:
sage: p = 5 sage: F = GF(p) sage: E = EllipticCurve(F, [1, 3]) sage: a_sub_p = E.change_ring(QQ).ap(p); a_sub_p 2
sage: len(E.points()) 4 sage: p + 1 - a_sub_p 4 sage: E.points() [(0 : 1 : 0), (1 : 0 : 1), (4 : 1 : 1), (4 : 4 : 1)]
sage: K = GF((p, 2),'a') sage: E = E.change_ring(K) sage: len(E.points()) 32 sage: (p + 1)**2 - a_sub_p**2 32 sage: w = E.points(); w [(0 : 1 : 0), (0 : 2*a + 4 : 1), (0 : 3*a + 1 : 1), (1 : 0 : 1), (2 : 2*a + 4 : 1), (2 : 3*a + 1 : 1), (3 : 2*a + 4 : 1), (3 : 3*a + 1 : 1), (4 : 1 : 1), (4 : 4 : 1), (a : 1 : 1), (a : 4 : 1), (a + 2 : a + 1 : 1), (a + 2 : 4*a + 4 : 1), (a + 3 : a : 1), (a + 3 : 4*a : 1), (a + 4 : 0 : 1), (2*a : 2*a : 1), (2*a : 3*a : 1), (2*a + 4 : a + 1 : 1), (2*a + 4 : 4*a + 4 : 1), (3*a + 1 : a + 3 : 1), (3*a + 1 : 4*a + 2 : 1), (3*a + 2 : 2*a + 3 : 1), (3*a + 2 : 3*a + 2 : 1), (4*a : 0 : 1), (4*a + 1 : 1 : 1), (4*a + 1 : 4 : 1), (4*a + 3 : a + 3 : 1), (4*a + 3 : 4*a + 2 : 1), (4*a + 4 : a + 4 : 1), (4*a + 4 : 4*a + 1 : 1)]
Note that the returned list is an immutable sorted Sequence:
sage: w[0] = 9 Traceback (most recent call last): ... ValueError: object is immutable; please change a copy instead.
- random_element()#
Return a random point on this elliptic curve, uniformly chosen among all rational points.
ALGORITHM:
Choose the point at infinity with probability \(1/(2q + 1)\). Otherwise, take a random element from the field as x-coordinate and compute the possible y-coordinates. Return the i’th possible y-coordinate, where i is randomly chosen to be 0 or 1. If the i’th y-coordinate does not exist (either there is no point with the given x-coordinate or we hit a 2-torsion point with i == 1), try again.
This gives a uniform distribution because you can imagine \(2q + 1\) buckets, one for the point at infinity and 2 for each element of the field (representing the x-coordinates). This gives a 1-to-1 map of elliptic curve points into buckets. At every iteration, we simply choose a random bucket until we find a bucket containing a point.
AUTHORS:
Jeroen Demeyer (2014-09-09): choose points uniformly random, see trac ticket #16951.
EXAMPLES:
sage: k = GF(next_prime(7^5)) sage: E = EllipticCurve(k,[2,4]) sage: P = E.random_element(); P # random (16740 : 12486 : 1) sage: type(P) <class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'> sage: P in E True
sage: k.<a> = GF(7^5) sage: E = EllipticCurve(k,[2,4]) sage: P = E.random_element(); P # random (5*a^4 + 3*a^3 + 2*a^2 + a + 4 : 2*a^4 + 3*a^3 + 4*a^2 + a + 5 : 1) sage: type(P) <class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'> sage: P in E True
sage: k.<a> = GF(2^5) sage: E = EllipticCurve(k,[a^2,a,1,a+1,1]) sage: P = E.random_element();P # random (a^4 + a : a^4 + a^3 + a^2 : 1) sage: type(P) <class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'> sage: P in E True
Ensure that the entire point set is reachable:
sage: E = EllipticCurve(GF(11), [2,1]) sage: S = set() sage: while len(S) < E.cardinality(): ....: S.add(E.random_element())
- random_point()#
Return a random point on this elliptic curve, uniformly chosen among all rational points.
ALGORITHM:
Choose the point at infinity with probability \(1/(2q + 1)\). Otherwise, take a random element from the field as x-coordinate and compute the possible y-coordinates. Return the i’th possible y-coordinate, where i is randomly chosen to be 0 or 1. If the i’th y-coordinate does not exist (either there is no point with the given x-coordinate or we hit a 2-torsion point with i == 1), try again.
This gives a uniform distribution because you can imagine \(2q + 1\) buckets, one for the point at infinity and 2 for each element of the field (representing the x-coordinates). This gives a 1-to-1 map of elliptic curve points into buckets. At every iteration, we simply choose a random bucket until we find a bucket containing a point.
AUTHORS:
Jeroen Demeyer (2014-09-09): choose points uniformly random, see trac ticket #16951.
EXAMPLES:
sage: k = GF(next_prime(7^5)) sage: E = EllipticCurve(k,[2,4]) sage: P = E.random_element(); P # random (16740 : 12486 : 1) sage: type(P) <class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'> sage: P in E True
sage: k.<a> = GF(7^5) sage: E = EllipticCurve(k,[2,4]) sage: P = E.random_element(); P # random (5*a^4 + 3*a^3 + 2*a^2 + a + 4 : 2*a^4 + 3*a^3 + 4*a^2 + a + 5 : 1) sage: type(P) <class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'> sage: P in E True
sage: k.<a> = GF(2^5) sage: E = EllipticCurve(k,[a^2,a,1,a+1,1]) sage: P = E.random_element();P # random (a^4 + a : a^4 + a^3 + a^2 : 1) sage: type(P) <class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'> sage: P in E True
Ensure that the entire point set is reachable:
sage: E = EllipticCurve(GF(11), [2,1]) sage: S = set() sage: while len(S) < E.cardinality(): ....: S.add(E.random_element())
- rational_points()#
All the points on this elliptic curve. The list of points is cached so subsequent calls are free.
EXAMPLES:
sage: p = 5 sage: F = GF(p) sage: E = EllipticCurve(F, [1, 3]) sage: a_sub_p = E.change_ring(QQ).ap(p); a_sub_p 2
sage: len(E.points()) 4 sage: p + 1 - a_sub_p 4 sage: E.points() [(0 : 1 : 0), (1 : 0 : 1), (4 : 1 : 1), (4 : 4 : 1)]
sage: K = GF((p, 2),'a') sage: E = E.change_ring(K) sage: len(E.points()) 32 sage: (p + 1)**2 - a_sub_p**2 32 sage: w = E.points(); w [(0 : 1 : 0), (0 : 2*a + 4 : 1), (0 : 3*a + 1 : 1), (1 : 0 : 1), (2 : 2*a + 4 : 1), (2 : 3*a + 1 : 1), (3 : 2*a + 4 : 1), (3 : 3*a + 1 : 1), (4 : 1 : 1), (4 : 4 : 1), (a : 1 : 1), (a : 4 : 1), (a + 2 : a + 1 : 1), (a + 2 : 4*a + 4 : 1), (a + 3 : a : 1), (a + 3 : 4*a : 1), (a + 4 : 0 : 1), (2*a : 2*a : 1), (2*a : 3*a : 1), (2*a + 4 : a + 1 : 1), (2*a + 4 : 4*a + 4 : 1), (3*a + 1 : a + 3 : 1), (3*a + 1 : 4*a + 2 : 1), (3*a + 2 : 2*a + 3 : 1), (3*a + 2 : 3*a + 2 : 1), (4*a : 0 : 1), (4*a + 1 : 1 : 1), (4*a + 1 : 4 : 1), (4*a + 3 : a + 3 : 1), (4*a + 3 : 4*a + 2 : 1), (4*a + 4 : a + 4 : 1), (4*a + 4 : 4*a + 1 : 1)]
Note that the returned list is an immutable sorted Sequence:
sage: w[0] = 9 Traceback (most recent call last): ... ValueError: object is immutable; please change a copy instead.
- set_order(value, check, num_checks)#
Set the value of self._order to value.
Use this when you know a priori the order of the curve to avoid a potentially expensive order calculation.
INPUT:
value
– integer in the Hasse-Weil range for this curve.check
(boolean, default:True
) – whether or not to run sanity checks on the input.num_checks
(integer, default: 8) – ifcheck
isTrue
, the number of times to check whethervalue
times a random point on this curve equals the identity.
OUTPUT:
None
EXAMPLES:
This example illustrates basic usage.
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6 sage: E.set_order(6) sage: E.order() 6 sage: E.order() * E.random_point() (0 : 1 : 0)
We now give a more interesting case, the NIST-P521 curve. Its order is too big to calculate with Sage, and takes a long time using other packages, so it is very useful here.
sage: p = 2^521 - 1 sage: prev_proof_state = proof.arithmetic() sage: proof.arithmetic(False) # turn off primality checking sage: F = GF(p) sage: A = p - 3 sage: B = 1093849038073734274511112390766805569936207598951683748994586394495953116150735016013708737573759623248592132296706313309438452531591012912142327488478985984 sage: q = 6864797660130609714981900799081393217269435300143305409394463459185543183397655394245057746333217197532963996371363321113864768612440380340372808892707005449 sage: E = EllipticCurve([F(A), F(B)]) sage: E.set_order(q) sage: G = E.random_point() sage: G.order() * G # This takes practically no time. (0 : 1 : 0) sage: proof.arithmetic(prev_proof_state) # restore state
It is an error to pass a value which is not an integer in the Hasse-Weil range:
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6 sage: E.set_order("hi") Traceback (most recent call last): ... TypeError: unable to convert 'hi' to an integer sage: E.set_order(0) Traceback (most recent call last): ... ValueError: Value 0 illegal (not an integer in the Hasse range) sage: E.set_order(1000) Traceback (most recent call last): ... ValueError: Value 1000 illegal (not an integer in the Hasse range)
It is also very likely an error to pass a value which is not the actual order of this curve. How unlikely is determined by num_checks, the factorization of the actual order, and the actual group structure:
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6 sage: E.set_order(11) Traceback (most recent call last): ... ValueError: Value 11 illegal (multiple of random point not the identity)
However, set_order can be fooled, though it’s not likely in “real cases of interest”. For instance, the order can be set to a multiple of the actual order:
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6 sage: E.set_order(12) # 12 just fits in the Hasse range sage: E.order() 12
Or, the order can be set incorrectly along with num_checks set too small:
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6 sage: E.set_order(4, num_checks=0) sage: E.order() 4
The value of num_checks must be an integer. Negative values are interpreted as zero, which means don’t do any checking:
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6 sage: E.set_order(4, num_checks=-12) sage: E.order() 4
AUTHORS:
Mariah Lenox (2011-02-16)
- trace_of_frobenius()#
Return the trace of Frobenius acting on this elliptic curve.
Note
This computes the curve cardinality, which may be time-consuming.
EXAMPLES:
sage: E = EllipticCurve(GF(101),[2,3]) sage: E.trace_of_frobenius() 6 sage: E = EllipticCurve(GF(11^5,'a'),[2,5]) sage: E.trace_of_frobenius() 802
The following shows that the issue from trac ticket #2849 is fixed:
sage: E = EllipticCurve(GF(3^5,'a'),[-1,-1]) sage: E.trace_of_frobenius() -27
- sage.schemes.elliptic_curves.ell_finite_field.fill_ss_j_dict()#
Fill the global cache of supersingular j-_polynomials.
This function does nothing except the first time it is called, when it fills
supersingular_j_polynomials
with precomputed values for \(p<300\). Setting the values this way avoids start-up costs.
- sage.schemes.elliptic_curves.ell_finite_field.is_j_supersingular(j, proof=True)#
Return True if \(j\) is a supersingular \(j\)-invariant.
INPUT:
j
(finite field element) – an element of a finite fieldproof
(boolean, default True) – If True, returns a proved result. If False, then a return value of False is certain but a return value of True may be based on a probabilistic test. See the ALGORITHM section below for more details.
OUTPUT:
(boolean) True if \(j\) is supersingular, else False.
ALGORITHM:
For small characteristics \(p\) we check whether the \(j\)-invariant is in a precomputed list of supersingular values. Otherwise we next check the \(j\)-invariant. If \(j=0\), the curve is supersingular if and only if \(p=2\) or \(p\equiv3\pmod{4}\); if \(j=1728\), the curve is supersingular if and only if \(p=3\) or \(p\equiv2\pmod{3}\). Next, if the base field is the prime field \({\rm GF}(p)\), we check that \((p+1)P=0\) for several random points \(P\), returning False if any fail: supersingular curves over \({\rm GF}(p)\) have cardinality \(p+1\). If Proof is false we now return True. Otherwise we compute the cardinality and return True if and only if it is divisible by \(p\).
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_finite_field import is_j_supersingular, supersingular_j_polynomials sage: [(p,[j for j in GF(p) if is_j_supersingular(j)]) for p in prime_range(30)] [(2, [0]), (3, [0]), (5, [0]), (7, [6]), (11, [0, 1]), (13, [5]), (17, [0, 8]), (19, [7, 18]), (23, [0, 3, 19]), (29, [0, 2, 25])] sage: [j for j in GF(109) if is_j_supersingular(j)] [17, 41, 43] sage: PolynomialRing(GF(109),'j')(supersingular_j_polynomials[109]).roots() [(43, 1), (41, 1), (17, 1)] sage: [p for p in prime_range(100) if is_j_supersingular(GF(p)(0))] [2, 3, 5, 11, 17, 23, 29, 41, 47, 53, 59, 71, 83, 89] sage: [p for p in prime_range(100) if is_j_supersingular(GF(p)(1728))] [2, 3, 7, 11, 19, 23, 31, 43, 47, 59, 67, 71, 79, 83] sage: [p for p in prime_range(100) if is_j_supersingular(GF(p)(123456))] [2, 3, 59, 89]
- sage.schemes.elliptic_curves.ell_finite_field.supersingular_j_polynomial(p, use_cache=True)#
Return a polynomial whose roots are the supersingular \(j\)-invariants in characteristic \(p\), other than 0, 1728.
INPUT:
\(p\) (integer) – a prime number.
\(use_cache\) (boolean, default
True
) – use cached coefficients if they exist
ALGORITHM:
First compute H(X) whose roots are the Legendre \(\lambda\)-invariants of supersingular curves (Silverman V.4.1(b)) in characteristic \(p\). Then, using a resultant computation with the polynomial relating \(\lambda\) and \(j\) (Silverman III.1.7(b)), we recover the polynomial (in variable
j
) whose roots are the \(j\)-invariants. Factors of \(j\) and \(j-1728\) are removed if present.Note
The only point of the use_cache parameter is to allow checking the precomputed coefficients.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_finite_field import supersingular_j_polynomial sage: f = supersingular_j_polynomial(67); f j^5 + 53*j^4 + 4*j^3 + 47*j^2 + 36*j + 8 sage: f.factor() (j + 1) * (j^2 + 8*j + 45) * (j^2 + 44*j + 24)
sage: [supersingular_j_polynomial(p) for p in prime_range(30)] [1, 1, 1, 1, 1, j + 8, j + 9, j + 12, j + 4, j^2 + 2*j + 21]