Path: blob/master/sage/schemes/elliptic_curves/ell_finite_field.py
4156 views
"""1Elliptic curves over finite fields23AUTHORS:45- William Stein (2005): Initial version67- Robert Bradshaw et al....89- John Cremona (2008-02): Point counting and group structure for10non-prime fields, Frobenius endomorphism and order, elliptic logs1112- Mariah Lenox (2011-03): Added set_order method13"""1415#*****************************************************************************16# Copyright (C) 2005 William Stein <[email protected]>17#18# Distributed under the terms of the GNU General Public License (GPL)19#20# This code is distributed in the hope that it will be useful,21# but WITHOUT ANY WARRANTY; without even the implied warranty of22# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU23# General Public License for more details.24#25# The full text of the GPL is available at:26#27# http://www.gnu.org/licenses/28#*****************************************************************************2930from sage.misc.randstate import current_randstate3132from sage.schemes.plane_curves.projective_curve import Hasse_bounds33from ell_field import EllipticCurve_field34from constructor import EllipticCurve, EllipticCurve_from_j35from sage.schemes.hyperelliptic_curves.hyperelliptic_finite_field import HyperellipticCurve_finite_field36import sage.rings.ring as ring37from sage.rings.all import Integer, ZZ, PolynomialRing, GF, polygen38from sage.rings.finite_rings.all import is_FiniteFieldElement39import sage.groups.generic as generic40import ell_point41from sage.rings.arith import gcd, lcm42from sage.structure.sequence import Sequence4344import sage.plot.all as plot4546import sage.libs.pari47pari = sage.libs.pari.all.pari4849class EllipticCurve_finite_field(EllipticCurve_field, HyperellipticCurve_finite_field):50"""51Elliptic curve over a finite field.52"""53def __init__(self, x, y=None):54"""55Special constructor for elliptic curves over a finite field5657EXAMPLES::5859sage: EllipticCurve(GF(101),[2,3])60Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Finite Field of size 1016162::6364sage: F=GF(101^2, 'a')65sage: EllipticCurve([F(2),F(3)])66Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Finite Field in a of size 101^26768Elliptic curves over `\ZZ/N\ZZ` with `N` prime are of type69"elliptic curve over a finite field":7071sage: F = Zmod(101)72sage: EllipticCurve(F, [2, 3])73Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Ring of integers modulo 10174sage: E = EllipticCurve([F(2), F(3)])75sage: type(E)76<class 'sage.schemes.elliptic_curves.ell_finite_field.EllipticCurve_finite_field_with_category'>77sage: E.category()78Category of schemes over Ring of integers modulo 1017980Elliptic curves over `\ZZ/N\ZZ` with `N` composite are of type81"generic elliptic curve"::8283sage: F = Zmod(95)84sage: EllipticCurve(F, [2, 3])85Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Ring of integers modulo 9586sage: E = EllipticCurve([F(2), F(3)])87sage: type(E)88<class 'sage.schemes.elliptic_curves.ell_generic.EllipticCurve_generic_with_category'>89sage: E.category()90Category of schemes over Ring of integers modulo 9591sage: TestSuite(E).run(skip=["_test_elements"])92"""93if isinstance(x, list):94seq = Sequence(x)95else:96seq = Sequence(y, universe=x)97ainvs = list(seq)98field = seq.universe()99if not isinstance(field, ring.Ring):100raise TypeError101102EllipticCurve_field.__init__(self, ainvs)103104self._point = ell_point.EllipticCurvePoint_finite_field105106def plot(self, *args, **kwds):107"""108Draw a graph of this elliptic curve over a prime finite field.109110INPUT:111112113- ``*args, **kwds`` - all other options are passed114to the circle graphing primitive.115116117EXAMPLES::118119sage: E = EllipticCurve(FiniteField(17), [0,1])120sage: P = plot(E, rgbcolor=(0,0,1))121"""122R = self.base_ring()123if not R.is_prime_field():124raise NotImplementedError125126G = plot.Graphics()127G += plot.points([P[0:2] for P in self.points() if not P.is_zero()], *args, **kwds)128129return G130131def _points_via_group_structure(self):132"""133Return a list of all the points on the curve, for prime fields only134(see points() for the general case)135136EXAMPLES::137138sage: S=EllipticCurve(GF(97),[2,3])._points_via_group_structure()139sage: len(S)140100141142See trac #4687, where the following example did not work::143144sage: E=EllipticCurve(GF(2),[0, 0, 1, 1, 1])145sage: E.points()146[(0 : 1 : 0)]147148::149150sage: E=EllipticCurve(GF(2),[0, 0, 1, 0, 1])151sage: E.points()152[(0 : 1 : 0), (1 : 0 : 1), (1 : 1 : 1)]153154::155156sage: E=EllipticCurve(GF(4,'a'),[0, 0, 1, 0, 1])157sage: E.points()158[(0 : 1 : 0), (0 : a : 1), (0 : a + 1 : 1), (1 : 0 : 1), (1 : 1 : 1), (a : 0 : 1), (a : 1 : 1), (a + 1 : 0 : 1), (a + 1 : 1 : 1)]159"""160# TODO, eliminate when polynomial calling is fast161G = self.abelian_group()162pts = [x.element() for x in G.gens()]163164ni = G.generator_orders()165ngens = G.ngens()166167H0=[self(0)]168if ngens == 0: # trivial group169return H0170for m in range(1,ni[0]): H0.append(H0[-1]+pts[0])171if ngens == 1: # cyclic group172return H0173174# else noncyclic group175H1=[self(0)]176for m in range(1,ni[1]): H1.append(H1[-1]+pts[1])177return [P+Q for P in H0 for Q in H1]178179def points(self):180r"""181All the points on this elliptic curve. The list of points is cached182so subsequent calls are free.183184EXAMPLES::185186sage: p = 5187sage: F = GF(p)188sage: E = EllipticCurve(F, [1, 3])189sage: a_sub_p = E.change_ring(QQ).ap(p); a_sub_p1902191192::193194sage: len(E.points())1954196sage: p + 1 - a_sub_p1974198sage: E.points()199[(0 : 1 : 0), (1 : 0 : 1), (4 : 1 : 1), (4 : 4 : 1)]200201::202203sage: K = GF(p**2,'a')204sage: E = E.change_ring(K)205sage: len(E.points())20632207sage: (p + 1)**2 - a_sub_p**220832209sage: w = E.points(); w210[(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)]211212Note that the returned list is an immutable sorted Sequence::213214sage: w[0] = 9215Traceback (most recent call last):216...217ValueError: object is immutable; please change a copy instead.218"""219try:220return self.__points221except AttributeError: pass222223from sage.structure.sequence import Sequence224k = self.base_ring()225if k.is_prime_field() and k.order()>50:226v = self._points_via_group_structure()227else:228v =self._points_fast_sqrt()229v.sort()230self.__points = Sequence(v, immutable=True)231return self.__points232233rational_points = points234235def count_points(self, n=1):236"""237Returns the cardinality of this elliptic curve over the base field or extensions.238239INPUT:240241- ``n`` (int) -- a positive integer242243OUTPUT:244245If `n=1`, returns the cardinality of the curve over its base field.246247If `n>1`, returns a list `[c_1, c_2, ..., c_n]` where `c_d` is248the cardinality of the curve over the extension of degree `d`249of its base field.250251EXAMPLES::252253sage: p = 101254sage: F = GF(p)255sage: E = EllipticCurve(F, [2,3])256sage: E.count_points(1)25796258sage: E.count_points(5)259[96, 10368, 1031904, 104053248, 10509895776]260261::262263sage: F.<a> = GF(p^2)264sage: E = EllipticCurve(F, [a,a])265sage: E.cardinality()26610295267sage: E.count_points()26810295269sage: E.count_points(1)27010295271sage: E.count_points(5)272[10295, 104072155, 1061518108880, 10828567126268595, 110462212555439192375]273274"""275try:276n = Integer(n)277except TypeError:278raise TypeError, "n must be a positive integer"279280if n<1:281raise ValueError, "n must be a positive integer"282283if n==1:284return self.cardinality()285else:286return [ self.cardinality(extension_degree=i) for i in range(1,n+1)]287288def random_element(self):289"""290Returns a random point on this elliptic curve.291292If `q` is small, finds all points and returns one at random.293Otherwise, returns the point at infinity with probability294`1/(q+1)` where the base field has cardinality `q`, and then295picks random `x`-coordinates from the base field until one296gives a rational point.297298EXAMPLES::299300sage: k = GF(next_prime(7^5))301sage: E = EllipticCurve(k,[2,4])302sage: P = E.random_element(); P303(16740 : 12486 : 1)304sage: type(P)305<class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'>306sage: P in E307True308309::310311sage: k.<a> = GF(7^5)312sage: E = EllipticCurve(k,[2,4])313sage: P = E.random_element(); P314(2*a^4 + 3*a^3 + 5*a^2 + 6*a + 4 : 6*a^4 + 4*a^3 + a + 6 : 1)315sage: type(P)316<class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'>317sage: P in E318True319320::321322sage: k.<a> = GF(2^5)323sage: E = EllipticCurve(k,[a^2,a,1,a+1,1])324sage: P = E.random_element(); P325(a^4 + a^2 + 1 : a^3 + a : 1)326sage: type(P)327<class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'>328sage: P in E329True330331Ensure that the entire point set is reachable::332333sage: E = EllipticCurve(GF(11), [2,1])334sage: len(set(E.random_element() for _ in range(100)))33516336sage: E.cardinality()33716338339TESTS:340341See trac #8311::342343sage: E = EllipticCurve(GF(3), [0,0,0,2,2])344sage: E.random_element()345(0 : 1 : 0)346sage: E.cardinality()3471348349sage: E = EllipticCurve(GF(2), [0,0,1,1,1])350sage: E.random_point()351(0 : 1 : 0)352sage: E.cardinality()3531354355sage: F.<a> = GF(4)356sage: E = EllipticCurve(F, [0, 0, 1, 0, a])357sage: E.random_point()358(0 : 1 : 0)359sage: E.cardinality()3601361362"""363random = current_randstate().c_rand_double364k = self.base_field()365q = k.order()366367# For small fields we find all the rational points and pick368# one at random. Note that the group can be trivial for369# q=2,3,4 only (see #8311) so these cases need special370# treatment.371372if q < 5:373pts = self.points() # will be cached374return pts[ZZ.random_element(len(pts))]375376377# The following allows the origin self(0) to be picked378if random() <= 1/float(q+1):379return self(0)380381while True:382v = self.lift_x(k.random_element(), all=True)383if v:384return v[int(random() * len(v))]385386random_point = random_element387388389def trace_of_frobenius(self):390r"""391Return the trace of Frobenius acting on this elliptic curve.392393.. note::394395This computes the curve cardinality, which may be396time-consuming.397398EXAMPLES::399400sage: E=EllipticCurve(GF(101),[2,3])401sage: E.trace_of_frobenius()4026403sage: E=EllipticCurve(GF(11^5,'a'),[2,5])404sage: E.trace_of_frobenius()405802406407The following shows that the issue from trac #2849 is fixed::408409sage: E=EllipticCurve(GF(3^5,'a'),[-1,-1])410sage: E.trace_of_frobenius()411-27412"""413return 1 + self.base_field().order() - self.cardinality()414415def _cardinality_with_j_invariant_1728(self):416r"""417Special function to compute cardinality when j=1728.418419EXAMPLES: An example with q=p=1 (mod 4)420421::422423sage: F=GF(10009)424sage: [EllipticCurve(F,[0,0,0,11^i,0])._cardinality_with_j_invariant_1728() for i in range(4)]425[10016, 10210, 10004, 9810]426427An example with q=p=3 (mod 4)428429::430431sage: F=GF(10007)432sage: [EllipticCurve(F,[0,0,0,5^i,0])._cardinality_with_j_invariant_1728() for i in range(4)]433[10008, 10008, 10008, 10008]434435An example with `q=p^2`, p=3 (mod 4)436437::438439sage: F.<a>=GF(10007^2,'a')440sage: [EllipticCurve(F,[0,0,0,a^i,0])._cardinality_with_j_invariant_1728() for i in range(4)]441[100160064, 100140050, 100120036, 100140050]442443Examples with `q=2^d`, d odd (3 isomorphism classes)::444445sage: F.<a> = GF(2**15,'a')446sage: ais = [[0,0,1,0,0],[0,0,1,1,0],[0,0,1,1,1]]447sage: curves=[EllipticCurve(F,ai) for ai in ais]448sage: all([all([e1==e2 or not e1.is_isomorphic(e2) for e1 in curves]) for e2 in curves])449True450sage: [e._cardinality_with_j_invariant_1728() for e in curves]451[32769, 33025, 32513]452453Examples with `q=2^d`, d even (7 isomorphism classes)::454455sage: F.<a> = GF(2**16,'a')456sage: b = a^11 # trace 1457sage: ais = [[0,0,1,0,0],[0,0,1,0,b],[0,0,1,b,0],[0,0,a,0,0],[0,0,a,0,a^2*b],[0,0,a^2,0,0],[0,0,a^2,0,a^4*b]]458sage: curves=[EllipticCurve(F,ai) for ai in ais]459sage: all([all([e1==e2 or not e1.is_isomorphic(e2) for e1 in curves]) for e2 in curves])460True461sage: [e._cardinality_with_j_invariant_1728() for e in curves]462[65025, 66049, 65537, 65793, 65281, 65793, 65281]463464Examples with `q=3^d`, d odd (4 isomorphism classes)::465466sage: F.<a> = GF(3**15,'a')467sage: b=a^7 # has trace 1468sage: ais=[[0,0,0,1,0],[0,0,0,-1,0],[0,0,0,-1,b],[0,0,0,-1,-b]]469sage: curves=[EllipticCurve(F,ai) for ai in ais]470sage: all([all([e1==e2 or not e1.is_isomorphic(e2) for e1 in curves]) for e2 in curves])471True472sage: [e._cardinality_with_j_invariant_1728() for e in curves]473[14348908, 14348908, 14342347, 14355469]474475Examples with `q=3^d`, d even (6 isomorphism classes)::476477sage: F.<g>=GF(3^18,'g')478sage: i=F(-1).sqrt()479sage: a=g^8 # has trace 1480sage: ais= [[0,0,0,1,0],[0,0,0,1,i*a],[0,0,0,g,0],[0,0,0,g^3,0],[0,0,0,g^2,0], [0,0,0,g^2,i*a*g^3]]481sage: curves=[EllipticCurve(F,ai) for ai in ais]482sage: all([all([e1==e2 or not e1.is_isomorphic(e2) for e1 in curves]) for e2 in curves])483True484sage: [E._cardinality_with_j_invariant_1728() for E in curves]485[387459856, 387400807, 387420490, 387420490, 387381124, 387440173]486"""487try:488return self._order489except AttributeError:490pass491492k = self.base_ring()493assert self.j_invariant()==k(1728)494q = k.cardinality()495p = k.characteristic()496d = k.degree()497x=polygen(ZZ)498499# p=2, j=0=1728500#501# Number of isomorphism classes is 3 (in odd degree) or 7 (in even degree)502#503if p==2:504if d%2==1:505# The 3 classes are represented, independently of d,506# by [0,0,1,0,0], [0,0,1,1,0], [0,0,1,1,1]507E=EllipticCurve(k,[0,0,1,0,0])508if self.is_isomorphic(E):509N = q+1510else:511n = (d+1)//2512t = 2**n513n = n%4514if n==0 or n==1: t=-t515E=EllipticCurve(k,[0,0,1,1,1])516if self.is_isomorphic(E): t=-t517N = q+1-t518else:519# The 7 classes are represented by E1=[0,0,1,0,0],520# E2=[0,0,1,0,b], E3=[0,0,1,b,0], E4=[0,0,a,0,0],521# E4=[0,0,a,0,a^2*b], E6=[0,0,a^2,0,0],522# E7=[0,0,a^2,0,a^4*b], where a is a non-cube and b523# has trace 1. E1's Frobenius is pi=(-2)**(d//2); the524# Frobeniuses are then pi, -pi, 0; w*pi, -w*pi;525# w^2*pi, -w^2*pi where w is either cube root of526# unity, so the traces are 2*pi, -2*pi, 0, -pi, +pi;527# -pi, +pi.528delta = self.discriminant()529discube = (delta**((q-1)//3) == k(1))530pi = (-2)**(d//2)531if discube:532a = k.gen()533b = a534while b.trace()==0: b*=a535if self.is_isomorphic(EllipticCurve(k,[0,0,1,b,0])):536t = 0537else:538t = 2*pi539if not self.is_isomorphic(EllipticCurve(k,[0,0,1,0,0])):540t = -t541542else:543t = pi544if self.is_isomorphic(EllipticCurve(k,[0,0,delta,0,0])):545t = -t546N = q+1-t547548549# p=3, j=0=1728550#551# Number of isomorphism classes is 4 (odd degree) or 6 (even degree)552#553elif p==3:554if d%2==1:555# The 4 classes are represented by [0,0,0,1,0],556# [0,0,0,-1,0], [0,0,0,-1,a], [0,0,0,-1,-a] where a557# has trace 1558delta = self.discriminant()559if (-delta).is_square():560t = 0561else:562u = delta.sqrt()563if not u.is_square(): u=-u564tr = ((self.a3()**2+self.a6())/u).trace()565if tr==0:566t = 0567else:568d2 = (d+1)//2569t = 3**d2570if d2%2==1: t = -t571if tr==-1: t = -t572N = q+1-t573else:574# The 6 classes are represented by: [0,0,0,1,0],575# [0,0,0,1,i*a]; [0,0,0,g,0], [0,0,0,g^3,0];576# [0,0,0,g^2,0], [0,0,0,g^2,i*a*g^3]; where g577# generates the multiplicative group modulo 4th578# powers, and a has nonzero trace.579A4 = self.a4()-self.a1()*self.a3()580i = k(-1).sqrt()581t = 0582if A4.is_square():583u = A4.sqrt()584t = (-3)**(d//2)585if u.is_square():586A6 = (self.a3()**2+self.a6())/u587if (i*A6).trace()==0:588t = 2*t589else:590t = -t591else:592A6 = (self.a3()**2+self.a6())/(u*A4)593if (i*A6).trace()==0:594t = -2*t595N = q+1-t596597# p>3, j=1728598#599# Number of isomorphism classes is 4 if q=1 (mod 4), else 2600#601elif p%4==3:602if d%2==1:603t = 0604else:605t = (-p)**(d//2)606w = (self.c4()/k(48))**((q-1)//4)607if w==1: t = 2*t608elif w==-1: t = -2*t609else: t = 0610611N = q+1-t612613# p=1 (mod 4). First find Frobenius pi=a+b*i for [0,0,0,-1,0] over GF(p):614# N(pi)=p and N(pi-1)=0 (mod 8).615#616else: # p%4==1617R = ZZ.extension(x**2+1,'i')618i = R.gen(1)619pi = R.fraction_field().factor(p)[0][0].gens_reduced()[0]620a,b = pi.list()621if a%2==0:622a,b = -b,a623if (a+b+1)%4==0:624a,b = -a,-b625pi = a+b*i # Now pi=a+b*i with (a,b)=(1,0),(3,2) mod 4626627# Lift to Frobenius for [0,0,0,-1,0] over GF(p^d):628629if d>1:630pi = pi**d631a,b = pi.list()632633# Compute appropriate quartic twist:634635w = (self.c4()/k(48))**((q-1)//4)636if w==1:637t = 2*a638elif w==-1:639t = -2*a640elif k(b)==w*k(a):641t = 2*b642else:643t = -2*b644N = q+1-t645646self._order = Integer(N)647return self._order648649def _cardinality_with_j_invariant_0(self):650r"""651Special function to compute cardinality when j=0.652653EXAMPLES: An example with q=p=1 (mod 6)654655::656657sage: F=GF(1009)658sage: [EllipticCurve(F,[0,0,0,0,11^i])._cardinality_with_j_invariant_0() for i in range(6)]659[948, 967, 1029, 1072, 1053, 991]660661An example with q=p=5 (mod 6)662663::664665sage: F=GF(1013)666sage: [EllipticCurve(F,[0,0,0,0,3^i])._cardinality_with_j_invariant_0() for i in range(6)]667[1014, 1014, 1014, 1014, 1014, 1014]668669An example with `q=p^2`, p=5 (mod 6)670671::672673sage: F.<a>=GF(1013^2,'a')674sage: [EllipticCurve(F,[0,0,0,0,a^i])._cardinality_with_j_invariant_0() for i in range(6)]675[1028196, 1027183, 1025157, 1024144, 1025157, 1027183]676677For examples in characteristic 2 and 3, see the function678_cardinality_with_j_invariant_1728()679"""680681try:682return self._order683except AttributeError:684pass685686k = self.base_ring()687assert self.j_invariant()==k(0)688p = k.characteristic()689if p==2 or p==3: # then 0==1728690return self._cardinality_with_j_invariant_1728()691692q = k.cardinality()693d = k.degree()694x=polygen(ZZ)695696# p>3, j=0697#698# Number of isomorphism classes is 4 if q=1 (mod 4), else 2699#700if p%6==5:701if d%2==1:702t = 0703else:704t = (-p)**(d//2)705w = (self.c6()/k(-864))**((q-1)//6)706if w==1: t = 2*t707elif w==-1: t = -2*t708elif w**3==1: t = -t709710N = q+1-t711712# p=1 (mod 6). First find Frobenius pi=a+b*w for [0,0,0,0,1] over GF(p):713# N(pi)=p and N(pi-1)=0 (mod 12).714#715else: # p%6==1716R = ZZ.extension(x**2-x+1,'zeta6')717zeta6 = R.gen(1)718pi = R.fraction_field().factor(p)[0][0].gens_reduced()[0]719while (pi-1).norm()%12 !=0: pi*=zeta6720a,b = pi.list()721z = k(-b)/k(a) # a *specific* 6th root of unity in k722723# Now pi=a+b*zeta6 with N(pi-1)=0 (mod 12)724725# Lift to Frobenius for [0,0,0,0,1] over GF(p^d):726727if d>1:728pi = pi**d729a,b = pi.list()730731# Compute appropriate sextic twist:732733w = (self.c6()/k(-864))**((q-1)//6)734735if w==1: t = 2*a+b # = Trace(pi)736elif w==-1: t = -2*a-b # = Trace(-pi)737elif w==z: t = a-b # = Trace(pi*zeta6)738elif w==z**2: t = -a-2*b # = Trace(pi*zeta6**2)739elif w==z**4: t = b-a # = Trace(pi*zeta6**4)740elif w==z**5: t = a+2*b # = Trace(pi*zeta6**5)741742N = q+1-t743744self._order = Integer(N)745return self._order746747def cardinality(self, algorithm='heuristic', extension_degree=1):748r"""749Return the number of points on this elliptic curve over an750extension field (default: the base field).751752INPUT:753754755- ``algorithm`` - string (default: 'heuristic'), used756only for point counting over prime fields757758- ``'heuristic'`` - use a heuristic to choose between759``'pari'`` and ``'bsgs'``.760761- ``'pari'`` - use the baby step giant step or SEA methods762as implemented in PARI via the C-library function ellap.763764- ``'bsgs'`` - use the baby step giant step method as765implemented in Sage, with the Cremona-Sutherland version766of Mestre's trick.767768- ``'all'`` - (over prime fields only) compute cardinality769with all of PARI and bsgs; return result if they agree770or raise a RuntimeError if they do not.771772- ``extension_degree`` - int (default: 1); if the773base field is `k=GF(p^n)` and extension_degree=d, returns774the cardinality of `E(GF(p^{n d}))`.775776777OUTPUT: an integer778779The cardinality is cached.780781Over prime fields, one of the above algorithms is used. Over782non-prime fields, the serious point counting is done on a standard783curve with the same j-invariant over the field GF(p)(j), then784lifted to the base_field, and finally account is taken of twists.785786For j=0 and j=1728 special formulas are used instead.787788EXAMPLES::789790sage: EllipticCurve(GF(4,'a'),[1,2,3,4,5]).cardinality()7918792sage: k.<a> = GF(3^3)793sage: l = [a^2 + 1, 2*a^2 + 2*a + 1, a^2 + a + 1, 2, 2*a]794sage: EllipticCurve(k,l).cardinality()79529796797::798799sage: l = [1, 1, 0, 2, 0]800sage: EllipticCurve(k,l).cardinality()80138802803An even bigger extension (which we check against Magma)::804805sage: EllipticCurve(GF(3^100,'a'),[1,2,3,4,5]).cardinality()806515377520732011331036459693969645888996929981504807sage: magma.eval("Order(EllipticCurve([GF(3^100)|1,2,3,4,5]))") # optional - magma808'515377520732011331036459693969645888996929981504'809810::811812sage: EllipticCurve(GF(10007),[1,2,3,4,5]).cardinality()81310076814sage: EllipticCurve(GF(10007),[1,2,3,4,5]).cardinality(algorithm='pari')81510076816sage: EllipticCurve(GF(next_prime(10**20)),[1,2,3,4,5]).cardinality()817100000000011093199520818819The cardinality is cached::820821sage: E = EllipticCurve(GF(3^100,'a'),[1,2,3,4,5])822sage: E.cardinality() is E.cardinality()823True824sage: E=EllipticCurve(GF(11^2,'a'),[3,3])825sage: E.cardinality()826128827sage: EllipticCurve(GF(11^100,'a'),[3,3]).cardinality()828137806123398222701841183371720896367762643312000384671846835266941791510341065565176497846502742959856128829830TESTS::831832sage: EllipticCurve(GF(10007),[1,2,3,4,5]).cardinality(algorithm='foobar')833Traceback (most recent call last):834...835ValueError: Algorithm is not known836"""837if extension_degree>1:838# A recursive call to cardinality() with839# extension_degree=1, which will cache the cardinality, is840# made by the call to frobenius_order() here:841R=self.frobenius_order()842if R.degree()==1:843return (self.frobenius()**extension_degree-1)**2844else:845return (self.frobenius()**extension_degree-1).norm()846847# Now extension_degree==1848if algorithm != 'all':849try:850return self._order851except AttributeError:852pass853854k = self.base_ring()855q = k.cardinality()856857if q < 50:858return self.cardinality_exhaustive()859860# use special code for j=0, 1728 (for any field)861j = self.j_invariant()862if j==k(0):863return self._cardinality_with_j_invariant_0()864if j==k(1728):865return self._cardinality_with_j_invariant_1728()866867N = 0868p = k.characteristic()869d = k.degree()870871# Over prime fields, we have a variety of algorithms to choose from:872873if d == 1:874if algorithm == 'heuristic':875algorithm = 'pari'876if algorithm == 'pari':877N = self.cardinality_pari()878elif algorithm == 'sea':879N = self.cardinality_pari() # purely for backwards compatibility880elif algorithm == 'bsgs':881N = self.cardinality_bsgs()882elif algorithm == 'all':883N1 = self.cardinality_pari()884N2 = self.cardinality_bsgs()885if N1 == N2:886N = N1887else:888raise RuntimeError, "BUG! Cardinality with pari=%s but with bsgs=%s"%(N1, N2)889else:890raise ValueError, "Algorithm is not known"891self._order = Integer(N)892return self._order893894# now k is not a prime field and j is not 0, 1728895896# we count points on a standard curve with the same897# j-invariant, defined over the field it generates, then lift898# to the curve's own field, and finally allow for twists899900# Since j is not 0, 1728 the only twists are quadratic901902j_pol=j.minimal_polynomial()903j_deg=j_pol.degree()904905# if not possible to work over a smaller field:906if d==j_deg:907self._order = self.cardinality_bsgs()908return self._order909910kj=GF(p**j_deg,name='a',modulus=j_pol)911jkj=kj.gen() if j_deg>1 else j_pol.roots(multiplicities=False)[0]912913# recursive call which will do all the real work:914Ej = EllipticCurve_from_j(jkj)915N=Ej.cardinality(extension_degree=d//j_deg)916917# if curve ia a (quadratic) twist of the "standard" one:918if not self.is_isomorphic(EllipticCurve_from_j(j)): N=2*(q+1)-N919920self._order = N921return self._order922923order = cardinality # alias924925def frobenius_polynomial(self):926r"""927Return the characteristic polynomial of Frobenius.928929The Frobenius endomorphism of the elliptic curve has quadratic930characteristic polynomial. In most cases this is irreducible and931defines an imaginary quadratic order; for some supersingular932curves, Frobenius is an integer a and the polynomial is933`(x-a)^2`.934935.. note::936937This computes the curve cardinality, which may be938time-consuming.939940EXAMPLES::941942sage: E=EllipticCurve(GF(11),[3,3])943sage: E.frobenius_polynomial()944x^2 - 4*x + 11945946For some supersingular curves, Frobenius is in Z and the polynomial947is a square::948949sage: E=EllipticCurve(GF(25,'a'),[0,0,0,0,1])950sage: E.frobenius_polynomial().factor()951(x + 5)^2952"""953x=polygen(ZZ)954return x**2-self.trace_of_frobenius()*x+self.base_field().cardinality()955956def frobenius_order(self):957r"""958Return the quadratic order Z[phi] where phi is the Frobenius959endomorphism of the elliptic curve960961.. note::962963This computes the curve cardinality, which may be964time-consuming.965966EXAMPLES::967968sage: E=EllipticCurve(GF(11),[3,3])969sage: E.frobenius_order()970Order in Number Field in phi with defining polynomial x^2 - 4*x + 11971972For some supersingular curves, Frobenius is in Z and the Frobenius973order is Z::974975sage: E=EllipticCurve(GF(25,'a'),[0,0,0,0,1])976sage: R=E.frobenius_order()977sage: R978Order in Number Field in phi with defining polynomial x + 5979sage: R.degree()9801981"""982f = self.frobenius_polynomial().factor()[0][0]983return ZZ.extension(f,names='phi')984985def frobenius(self):986r"""987Return the frobenius of self as an element of a quadratic order988989.. note::990991This computes the curve cardinality, which may be992time-consuming.993994Frobenius is only determined up to conjugacy.995996EXAMPLES::997998sage: E=EllipticCurve(GF(11),[3,3])999sage: E.frobenius()1000phi1001sage: E.frobenius().minpoly()1002x^2 - 4*x + 1110031004For some supersingular curves, Frobenius is in Z::10051006sage: E=EllipticCurve(GF(25,'a'),[0,0,0,0,1])1007sage: E.frobenius()1008-51009"""1010R = self.frobenius_order()1011if R.degree()==1:1012return self.frobenius_polynomial().roots(multiplicities=False)[0]1013else:1014return R.gen(1)10151016def cardinality_exhaustive(self):1017r"""1018Return the cardinality of self over the base field. Simply adds up1019the number of points with each x-coordinate: only used for small1020field sizes!10211022EXAMPLES::10231024sage: p=next_prime(10^3)1025sage: E=EllipticCurve(GF(p),[3,4])1026sage: E.cardinality_exhaustive()102710201028sage: E=EllipticCurve(GF(3^4,'a'),[1,1])1029sage: E.cardinality_exhaustive()1030641031"""1032self._order = Integer(1+sum([len(self.lift_x(x,all=True)) for x in self.base_field()]))1033return self._order10341035def cardinality_pari(self):1036r"""1037Return the cardinality of self over the (prime) base field using PARI.10381039The result is not cached.10401041EXAMPLES::10421043sage: p=next_prime(10^3)1044sage: E=EllipticCurve(GF(p),[3,4])1045sage: E.cardinality_pari()104610201047sage: K=GF(next_prime(10^6))1048sage: E=EllipticCurve(K,[1,0,0,1,1])1049sage: E.cardinality_pari()105099994510511052TESTS::10531054sage: K.<a>=GF(3^20)1055sage: E=EllipticCurve(K,[1,0,0,1,a])1056sage: E.cardinality_pari()1057Traceback (most recent call last):1058...1059ValueError: cardinality_pari() only works over prime fields.1060sage: E.cardinality()1061348679431010621063"""1064k = self.base_ring()1065p = k.characteristic()1066if k.degree()==1:1067return ZZ(p + 1 - int(self._pari_().ellap(p)))1068else:1069raise ValueError, "cardinality_pari() only works over prime fields."10701071def cardinality_bsgs(self, verbose=False):1072r"""1073Return the cardinality of self over the base field. Will be called1074by user function cardinality only when necessary, i.e. when the1075j_invariant is not in the prime field.10761077ALGORITHM: A variant of "Mestre's trick" extended to all finite1078fields by Cremona and Sutherland, 2008.10791080.. note::108110821. The Mestre-Schoof-Cremona-Sutherland algorithm may fail for1083a small finite number of curves over `F_q` for `q` at most 49, so1084for `q<50` we use an exhaustive count.108510862. Quadratic twists are not implemented in characteristic 21087when `j=0 (=1728)`; but this case is treated separately.10881089EXAMPLES::10901091sage: p=next_prime(10^3)1092sage: E=EllipticCurve(GF(p),[3,4])1093sage: E.cardinality_bsgs()109410201095sage: E=EllipticCurve(GF(3^4,'a'),[1,1])1096sage: E.cardinality_bsgs()1097641098sage: F.<a>=GF(101^3,'a')1099sage: E=EllipticCurve([2*a^2 + 48*a + 27, 89*a^2 + 76*a + 24])1100sage: E.cardinality_bsgs()110110313521102"""1103E1 = self1104k = self.base_field()1105q = k.order()1106if q<50:1107if verbose:1108print "q=",q,"< 50 so using exhaustive count"1109return self.cardinality_exhaustive()11101111# Construct the quadratic twist:1112E2 = E1.quadratic_twist()1113if verbose:1114print "Quadratic twist is ",E2.ainvs()11151116bounds = Hasse_bounds(q)1117lower, upper = bounds1118B = upper-q-1 # = floor(2*sqrt(q))1119a = ZZ(0)1120N1 = N2 = M = ZZ(1)1121kmin = -B1122kmax = B1123q1 = q+11124# Throughout, we have #E=q+1-t where |t|<=B and t=a+k*M = a1125# (mod M) where kmin <= k <= kmax.11261127# M is the lcm of the orders of all the points found on E1 and1128# E2, which will eventually exceed 2*B, at which point1129# kmin=kmax.11301131if q > 2**10:1132N1 = ZZ(2)**sum([e for P,e in E1._p_primary_torsion_basis(2)])1133N2 = ZZ(2)**sum([e for P,e in E2._p_primary_torsion_basis(2)])1134if q > 2**20:1135N1 *= ZZ(3)**sum([e for P,e in E1._p_primary_torsion_basis(3)])1136N2 *= ZZ(3)**sum([e for P,e in E2._p_primary_torsion_basis(3)])1137if q > 2**40:1138N1 *= ZZ(5)**sum([e for P,e in E1._p_primary_torsion_basis(5)])1139N2 *= ZZ(5)**sum([e for P,e in E2._p_primary_torsion_basis(5)])1140# We now know that t=q+1 (mod N1) and t=-(q+1) (mod N2)1141a = q11142M = N11143g,u,v = M.xgcd(N2) # g==u*M+v*N21144if N2>g:1145a = (a*v*N2-q1*u*M)//g1146M *= (N2//g) # = lcm(M,N2)1147a = a%M1148if verbose:1149print "(a,M)=",(a,M)1150kmin = ((-B-a)/M).ceil()1151kmax = ((B-a)/M).floor()1152if kmin==kmax:1153self._order = q1-a-kmin*M1154if verbose: print "no random points were needed"1155return self._order1156if verbose: print "(2,3,5)-torsion subgroup gives M=",M11571158# N1, N2 are divisors of the orders of E1, E2 separately,1159# which are used to speed up the computation of the orders of1160# points on each curve. For large q it is worth initializing1161# these with the full order of the (2,3,5)-torsion which are1162# often non-trivial.11631164while kmax!=kmin:1165# Get a random point on E1 and find its order, using the1166# Hasse bounds and the fact that we know that the group1167# order is a multiple of N1:1168n = generic.order_from_bounds(E1.random_point(),bounds,N1,operation='+')1169if verbose: print "New point on E has order ",n1170# update N1 and M1171N1 = N1.lcm(n)1172g,u,v = M.xgcd(n) # g==u*M+v*n1173if n>g:1174# update congruence a (mod M) with q+1 (mod n)1175a = (a*v*n+q1*u*M)//g1176M *= (n//g) # = lcm(M,n)1177a = a%M1178if verbose: print "(a,M)=",(a,M)1179kmin = ((-B-a)/M).ceil()1180kmax = ((B-a)/M).floor()1181if kmin==kmax:1182self._order = q1-a-kmin*M1183return self._order1184if verbose: print "number of possibilities is now ",kmax-kmin+111851186# Get a random point on E2 and find its order, using the1187# Hasse bounds and the fact that we know that the group1188# order is a multiple of N2:1189n = generic.order_from_bounds(E2.random_point(),bounds,N2,operation='+')1190if verbose: print "New point on E' has order ",n1191# update N2 and M1192N2 = N2.lcm(n)1193g,u,v = M.xgcd(n) # g==u*M+v*n1194if n>g:1195# update congruence a (mod M) with -(q+1) (mod n)1196a = (a*v*n-q1*u*M)//g1197M *= (n//g) # = lcm(M,n)1198a = a%M1199if verbose: print "(a,M)=",(a,M)1200kmin = ((-B-a)/M).ceil()1201kmax = ((B-a)/M).floor()1202if kmin==kmax:1203self._order = q1-a-kmin*M1204return self._order1205if verbose: print "number of possibilities is now ",kmax-kmin+112061207def gens(self):1208"""1209Returns a tuple of length up to 2 of points which generate the1210abelian group of points on this elliptic curve. See1211abelian_group() for limitations.12121213The algorithm uses random points on the curve, and hence the1214generators are likely to differ from one run to another; but they1215are cached so will be consistent in any one run of Sage.12161217AUTHORS:12181219- John Cremona12201221EXAMPLES::12221223sage: E=EllipticCurve(GF(11),[2,5])1224sage: E.gens() # random output1225((0 : 7 : 1),)1226sage: EllipticCurve(GF(41),[2,5]).gens() # random output1227((21 : 1 : 1), (8 : 0 : 1))1228sage: F.<a>=GF(3^6,'a')1229sage: E=EllipticCurve([a,a+1])1230sage: pts=E.gens()1231sage: len(pts)123211233sage: pts[0].order()==E.cardinality()1234True1235"""1236try:1237G = self.abelian_group()1238return [x.element() for x in G.gens()]1239except AttributeError:1240pass12411242def __iter__(self):1243"""1244Return an iterator through the points of this elliptic curve.12451246EXAMPLES::12471248sage: E = EllipticCurve(GF(11), [1,2])1249sage: for P in E: print P, P.order()1250(0 : 1 : 0) 11251(1 : 2 : 1) 41252(1 : 9 : 1) 41253(2 : 1 : 1) 81254...1255(10 : 0 : 1) 21256"""1257for P in self.points():1258yield P12591260def __getitem__(self, n):1261"""1262Return the n'th point in self's __points list. This enables users1263to iterate over the curve's point set.12641265EXAMPLE::12661267sage: E=EllipticCurve(GF(97),[2,3])1268sage: S=E.points()1269sage: E[10]1270(10 : 76 : 1)1271sage: E[15]1272(17 : 10 : 1)1273sage: for P in E: print P.order()1274112755012765012775012785012795128051281501282...1283"""1284return self.points()[n]12851286def abelian_group(self, debug=False):1287r"""1288Returns the abelian group structure of the group of points on this1289elliptic curve.12901291.. warning::12921293The algorithm is definitely *not* intended for use with1294*large* finite fields! The factorization of the orders of1295elements must be feasible. Also, baby-step-giant-step1296methods are used which have space and time requirements1297which are `O(\sqrt{q})`.12981299Also, the algorithm uses random points on the curve and hence the1300generators are likely to differ from one run to another; but the1301group is cached so the generators will not change in any one run of1302Sage.13031304INPUT:130513061307- ``debug`` - (default: False): if True, print1308debugging messages130913101311OUTPUT:13121313- an abelian group13141315- tuple of images of each of the generators of the abelian1316group as points on this curve13171318AUTHORS:13191320- John Cremona13211322EXAMPLES::13231324sage: E=EllipticCurve(GF(11),[2,5])1325sage: E.abelian_group()1326Additive 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 1113271328::13291330sage: E=EllipticCurve(GF(41),[2,5])1331sage: E.abelian_group()1332Additive abelian group isomorphic to Z/2 + Z/22 ...13331334::13351336sage: F.<a>=GF(3^6,'a')1337sage: E=EllipticCurve([a^4 + a^3 + 2*a^2 + 2*a, 2*a^5 + 2*a^3 + 2*a^2 + 1])1338sage: E.abelian_group()1339Additive abelian group isomorphic to Z/26 + Z/26 ...13401341::13421343sage: F.<a>=GF(101^3,'a')1344sage: E=EllipticCurve([2*a^2 + 48*a + 27, 89*a^2 + 76*a + 24])1345sage: E.abelian_group()1346Additive abelian group isomorphic to Z/1031352 ...13471348The group can be trivial::13491350sage: E=EllipticCurve(GF(2),[0,0,1,1,1])1351sage: E.abelian_group()1352Trivial group embedded in Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 + x + 1 over Finite Field of size 213531354Of course, there are plenty of points if we extend the field::13551356sage: E.cardinality(extension_degree=100)1357126765060022823165329651689062513581359This tests the patch for trac #3111, using 10 primes randomly1360selected::13611362sage: E = EllipticCurve('389a')1363sage: for p in [5927, 2297, 1571, 1709, 3851, 127, 3253, 5783, 3499, 4817]:1364... G = E.change_ring(GF(p)).abelian_group()1365sage: for p in prime_range(10000): # long time (19s on sage.math, 2011)1366... if p != 389:1367... G = E.change_ring(GF(p)).abelian_group()13681369This tests that the bug reported in trac #3926 has been fixed::13701371sage: K.<i> = QuadraticField(-1)1372sage: OK = K.ring_of_integers()1373sage: P=K.factor(10007)[0][0]1374sage: OKmodP = OK.residue_field(P)1375sage: E = EllipticCurve([0,0,0,i,i+3])1376sage: Emod = E.change_ring(OKmodP); Emod1377Elliptic Curve defined by y^2 = x^3 + ibar*x + (ibar+3) over Residue field in ibar of Fractional ideal (10007)1378sage: Emod.abelian_group() #random generators1379(Multiplicative Abelian Group isomorphic to C50067594 x C2,1380((3152*ibar + 7679 : 7330*ibar + 7913 : 1), (8466*ibar + 1770 : 0 : 1)))1381"""1382if not debug:1383# if we're in debug mode, always recalculate1384try:1385return self.__abelian_group1386except AttributeError:1387pass13881389k = self.base_field()1390q = k.order()1391p = k.characteristic()1392d = k.degree()1393j = self.j_invariant()1394if d>1:1395d = j.minimal_polynomial().degree()139613971398# Before computing the group structure we compute the1399# cardinality. While this is not strictly necessary, it makes1400# the code simpler and also makes computation of orders of1401# points faster.14021403# j=0,172814041405if j==k(0):1406N = self._cardinality_with_j_invariant_0()1407if j==k(1728):1408N = self._cardinality_with_j_invariant_1728()14091410bounds = Hasse_bounds(q)1411lower, upper = bounds1412if debug:1413print "Lower and upper bounds on group order: [",lower,",",upper,"]"14141415try:1416N=self._order1417if debug:1418print "Group order already known to be ",N1419except:1420if (q<50):1421if debug:1422print "Computing group order naively"1423N=self.cardinality_exhaustive()1424elif d==1:1425if debug:1426print "Computing group order using PARI"1427N=self.cardinality()1428else:1429if debug:1430print "Computing group order using bsgs"1431N=self.cardinality_bsgs()1432if debug:1433print "... group order = ",N14341435self._order=N1436plist = N.prime_factors()1437P1=self(0)1438P2=self(0)1439n1= Integer(1)1440n2= Integer(1)1441P1._order=n11442P2._order=n21443npts = 014441445# At all stages the current subgroup is generated by P1, P2 with1446# orders n1,n2 which are disjoint. We stop when n1*n2 == N14471448if debug:1449"About to start generating random points"14501451while n1*n2 != N:1452if debug:1453"Getting a new random point"1454Q = self.random_point()1455while Q.is_zero(): Q = self.random_point()1456npts += 11457if debug:1458print "Q = ",Q,":",1459print " Order(Q) = ", Q.order()14601461Q1=n1*Q;14621463if Q1.is_zero() and npts>=10: # then P1,n1 will not change but we may increase n21464if debug: print "Case 2: n2 may increase"1465n1a = 1; n1b = n11466P1a = P11467n1a = n1.prime_to_m_part(N//n1)1468n1b = n1//n1a1469Q = n1a*Q # has order | n1b1470P1a = n1a*P1 # has order = n1b1471if debug: print "n1a=",n1a1472a = None1473for m in n1b.divisors():1474try:1475a = generic.bsgs(m*P1a,m*Q,(0,(n1b//m)-1),operation='+')1476break1477except ValueError:1478pass1479assert a != None1480a *= (m*n1a)1481if debug: print "linear relation gives m=",m,", a=",a1482if debug: assert m*Q==a*P11483if m>1: # else Q is in <P1>1484Q=Q-(a//m)*P1; # has order m and is disjoint from P11485if debug: assert Q.order()==m1486Q._order=m1487if n2==1: # this is our first nontrivial P21488P2=Q1489n2=m1490if debug:1491print "Adding second generator ",P2," of order ",n21492print "Subgroup order now ",n1*n2,"=",n1,"*",n21493else: # we must merge P2 and Q:1494oldn2=n2 # holds old value1495P2,n2=generic.merge_points((P2,n2),(Q,m),operation='+', check=debug)1496if debug: assert P2.order()==n21497P2._order=n21498if debug:1499if n2>oldn2:1500print "Replacing second generator by ",P2,1501print " of order ",n2, " gaining index ",n2//oldn21502print "Subgroup order now ",n1*n2,"=",n1,"*",n21503elif not Q1.is_zero(): # Q1 nonzero: n1 will increase1504if debug: print "Case 1: n1 may increase"1505oldn1=n11506if n2>1:1507P3=(n1//n2)*P1 # so P2,P3 are a basis for n2-torsion1508if debug: assert P3.order()==n21509P3._order=n21510if debug: print "storing generator ",P3," of ",n2,"-torsion"1511m = generic.order_from_multiple(Q,N,plist,operation='+', check=debug)1512P1,n1=generic.merge_points((P1,n1),(Q,m), check=debug)1513if debug: assert P1.order()==n11514P1._order=n11515if debug:1516print "Replacing first generator by ",P1," of order ",1517print n1,", gaining index ",n1//oldn11518print "Subgroup order now ",n1*n2,"=",n1,"*",n215191520# Now replace P2 by a point of order n2 s.t. it and1521# (n1//n2)*P1 are still a basis for n2-torsion:1522if n2>1:1523a,m = generic.linear_relation(P1,P3,operation='+')1524if debug: print "linear relation gives m=",m,", a=",a1525P3=P3-(a//m)*P11526if debug: assert P3.order()==m1527P3._order=m1528if debug: print "First P2 component =",P31529if m==n2:1530P2=P31531else:1532a,m = generic.linear_relation(P1,P2,operation='+')1533if debug: print "linear relation gives m=",m,", a=",a1534P2=P2-(a//m)*P1;1535if debug: assert P2.order()==m1536P2._order=m1537if debug: print "Second P2 component =",P21538P2,n2=generic.merge_points((P2,n2),(P3,m), check=debug)1539if debug: assert P2.order()==n21540P2._order=n21541if debug: print "Combined P2 component =",P215421543if debug:1544if P1.order()!=n1:1545print "Generator P1 = ",P1," has order ",P1.order(),1546print " and not ",n11547raise ValueError1548if P2.order()!=n2:1549print "Generator P2 = ",P2," has order ",P2.order()1550print " and not ",n21551raise ValueError1552if n2>1:1553if generic.linear_relation(P1,P2,operation='+')[1]!=n2:1554print "Generators not independent!"1555raise ValueError1556print "Generators: P1 = ",P1," of order ",n1,1557print ", P2 = ",P2," of order ",n21558print "Subgroup order is now ",n1*n2,"=",n1,"*",n215591560# Finished: record group order, structure and generators15611562from sage.groups.additive_abelian.additive_abelian_wrapper import AdditiveAbelianGroupWrapper1563self._order = n1*n21564if n1==1:1565self.__abelian_group = AdditiveAbelianGroupWrapper(self.point_homset(), [], [])1566else:1567if n2==1:1568self.__abelian_group = AdditiveAbelianGroupWrapper(self.point_homset(), [P1], [n1])1569else:1570self.__abelian_group = AdditiveAbelianGroupWrapper(self.point_homset(), [P1, P2], [n1, n2])1571return self.__abelian_group15721573def is_isogenous(self, other, field=None, proof=True):1574"""1575Returns whether or not self is isogenous to other15761577INPUT:15781579- ``other`` -- another elliptic curve.15801581- ``field`` (default None) -- a field containing the base1582fields of the two elliptic curves into which the two curves1583may be extended to test if they are isogenous over this1584field. By default is_isogenous will not try to find this1585field unless one of the curves can be extended into the base1586field of the other, in which case it will test over the1587larger base field.15881589- ``proof`` (default True) -- this parameter is here only to1590be consistent with versions for other types of elliptic1591curves.15921593OUTPUT:15941595(bool) True if there is an isogeny from curve ``self`` to1596curve ``other`` defined over ``field``.15971598EXAMPLES::15991600sage: E1 = EllipticCurve(GF(11^2,'a'),[2,7]); E11601Elliptic Curve defined by y^2 = x^3 + 2*x + 7 over Finite Field in a of size 11^21602sage: E1.is_isogenous(5)1603Traceback (most recent call last):1604...1605ValueError: Second argument is not an Elliptic Curve.1606sage: E1.is_isogenous(E1)1607True16081609sage: E2 = EllipticCurve(GF(7^3,'b'),[3,1]); E21610Elliptic Curve defined by y^2 = x^3 + 3*x + 1 over Finite Field in b of size 7^31611sage: E1.is_isogenous(E2)1612Traceback (most recent call last):1613...1614ValueError: The base fields must have the same characteristic.16151616sage: E3 = EllipticCurve(GF(11^2,'c'),[4,3]); E31617Elliptic Curve defined by y^2 = x^3 + 4*x + 3 over Finite Field in c of size 11^21618sage: E1.is_isogenous(E3)1619False16201621sage: E4 = EllipticCurve(GF(11^6,'d'),[6,5]); E41622Elliptic Curve defined by y^2 = x^3 + 6*x + 5 over Finite Field in d of size 11^61623sage: E1.is_isogenous(E4)1624True16251626sage: E5 = EllipticCurve(GF(11^7,'e'),[4,2]); E51627Elliptic Curve defined by y^2 = x^3 + 4*x + 2 over Finite Field in e of size 11^71628sage: E1.is_isogenous(E5)1629Traceback (most recent call last):1630...1631ValueError: Curves have different base fields: use the field parameter.16321633When the field is given:16341635sage: E1 = EllipticCurve(GF(13^2,'a'),[2,7]); E11636Elliptic Curve defined by y^2 = x^3 + 2*x + 7 over Finite Field in a of size 13^21637sage: E1.is_isogenous(5,GF(13^6,'f'))1638Traceback (most recent call last):1639...1640ValueError: Second argument is not an Elliptic Curve.1641sage: E6 = EllipticCurve(GF(11^3,'g'),[9,3]); E61642Elliptic Curve defined by y^2 = x^3 + 9*x + 3 over Finite Field in g of size 11^31643sage: E1.is_isogenous(E6,QQ)1644Traceback (most recent call last):1645...1646ValueError: The base fields must have the same characteristic.1647sage: E7 = EllipticCurve(GF(13^5,'h'),[2,9]); E71648Elliptic Curve defined by y^2 = x^3 + 2*x + 9 over Finite Field in h of size 13^51649sage: E1.is_isogenous(E7,GF(13^4,'i'))1650Traceback (most recent call last):1651...1652ValueError: Field must be an extension of the base fields of both curves1653sage: E1.is_isogenous(E7,GF(13^10,'j'))1654False1655sage: E1.is_isogenous(E7,GF(13^30,'j'))1656False1657"""1658from ell_generic import is_EllipticCurve1659if not is_EllipticCurve(other):1660raise ValueError, "Second argument is not an Elliptic Curve."1661if self.is_isomorphic(other):1662return True1663elif self.base_field().characteristic() != other.base_field().characteristic():1664raise ValueError, "The base fields must have the same characteristic."1665elif field==None:1666if self.base_field().degree() == other.base_field().degree():1667if self.cardinality() == other.cardinality():1668return True1669else:1670return False1671elif self.base_field().degree() == gcd(self.base_field().degree(),other.base_field().degree()):1672if self.cardinality(extension_degree=other.base_field().degree()//self.base_field().degree()) == other.cardinality():1673return True1674else:1675return False1676elif other.base_field().degree() == gcd(self.base_field().degree(),other.base_field().degree()):1677if other.cardinality(extension_degree=self.base_field().degree()//other.base_field().degree()) == self.cardinality():1678return True1679else:1680return False1681else:1682raise ValueError, "Curves have different base fields: use the field parameter."1683else:1684if not lcm(self.base_field().degree(), other.base_field().degree()).divides(field.degree()):1685raise ValueError, "Field must be an extension of the base fields of both curves"1686else:1687if \1688self.cardinality(extension_degree=field.degree()//self.base_field().degree())\1689== other.cardinality(extension_degree=field.degree()//other.base_field().degree()):1690return True1691else:1692return False16931694def is_supersingular(self, proof=True):1695r"""1696Return True if this elliptic curve is supersingular, else False.16971698INPUT:16991700- ``proof`` (boolean, default True) -- If True, returns a1701proved result. If False, then a return value of False is1702certain but a return value of True may be based on a1703probabilistic test. See the documentaion of the function1704:meth:`is_j_supersingular` for more details.17051706EXAMPLES::17071708sage: F = GF(101)1709sage: EllipticCurve(j=F(0)).is_supersingular()1710True1711sage: EllipticCurve(j=F(1728)).is_supersingular()1712False1713sage: EllipticCurve(j=F(66)).is_supersingular()1714True1715sage: EllipticCurve(j=F(99)).is_supersingular()1716False17171718TESTS::17191720sage: from sage.schemes.elliptic_curves.ell_finite_field import supersingular_j_polynomial, is_j_supersingular1721sage: F = GF(103)1722sage: ssjlist = [F(1728)] + supersingular_j_polynomial(103).roots(multiplicities=False)1723sage: Set([j for j in F if is_j_supersingular(j)]) == Set(ssjlist)1724True17251726"""1727return is_j_supersingular(self.j_invariant(), proof=proof)17281729def is_ordinary(self, proof=True):1730r"""1731Return True if this elliptic curve is ordinary, else False.17321733INPUT:17341735- ``proof`` (boolean, default True) -- If True, returns a1736proved result. If False, then a return value of True is1737certain but a return value of False may be based on a1738probabilistic test. See the documentaion of the function1739:meth:`is_j_supersingular` for more details.17401741EXAMPLES::17421743sage: F = GF(101)1744sage: EllipticCurve(j=F(0)).is_ordinary()1745False1746sage: EllipticCurve(j=F(1728)).is_ordinary()1747True1748sage: EllipticCurve(j=F(66)).is_ordinary()1749False1750sage: EllipticCurve(j=F(99)).is_ordinary()1751True17521753"""1754return not is_j_supersingular(self.j_invariant(), proof=proof)17551756def set_order(self, value, num_checks=8):1757r"""1758Set the value of self._order to value.17591760Use this when you know a priori the order of the curve to1761avoid a potentially expensive order calculation.17621763INPUT:17641765- ``value`` - Integer in the Hasse-Weil range for this1766curve.17671768- ``num_checks`` - Integer (default: 8) number of times to1769check whether value*(a random point on this curve) is1770equal to the identity.177117721773OUTPUT:17741775None17761777EXAMPLES:17781779This example illustrates basic usage.17801781::17821783sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61784sage: E.set_order(6)1785sage: E.order()178661787sage: E.order() * E.random_point()1788(0 : 1 : 0)17891790We now give a more interesting case, the NIST-P521 curve. Its1791order is too big to calculate with Sage, and takes a long time1792using other packages, so it is very useful here.17931794::17951796sage: p = 2^521 - 11797sage: prev_proof_state = proof.arithmetic()1798sage: proof.arithmetic(False) # turn off primality checking1799sage: F = GF(p)1800sage: A = p - 31801sage: B = 10938490380737342745111123907668055699362075989516837489945863944959531161507350160137087375737596232485921322967063133094384525315910129121423274884789859841802sage: q = 68647976601306097149819007990813932172694353001433054093944634591855431833976553942450577463332171975329639963713633211138647686124403803403728088927070054491803sage: E = EllipticCurve([F(A), F(B)])1804sage: E.set_order(q)1805sage: G = E.random_point()1806sage: E.order() * G # This takes practically no time.1807(0 : 1 : 0)1808sage: proof.arithmetic(prev_proof_state) # restore state18091810It is an error to pass a value which is not an integer in the1811Hasse-Weil range::18121813sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61814sage: E.set_order("hi")1815Traceback (most recent call last):1816...1817ValueError: Value hi illegal (not an integer in the Hasse range)1818sage: E.set_order(3.14159)1819Traceback (most recent call last):1820...1821ValueError: Value 3.14159000000000 illegal (not an integer in the Hasse range)1822sage: E.set_order(0)1823Traceback (most recent call last):1824...1825ValueError: Value 0 illegal (not an integer in the Hasse range)1826sage: E.set_order(1000)1827Traceback (most recent call last):1828...1829ValueError: Value 1000 illegal (not an integer in the Hasse range)18301831It is also very likely an error to pass a value which is not1832the actual order of this curve. How unlikely is determined by1833num_checks, the factorization of the actual order, and the1834actual group structure::18351836sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61837sage: E.set_order(11)1838Traceback (most recent call last):1839...1840ValueError: Value 11 illegal (multiple of random point not the identity)18411842However, set_order can be fooled, though it's not likely in1843"real cases of interest". For instance, the order can be set1844to a multiple of the actual order::18451846sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61847sage: E.set_order(12) # 12 just fits in the Hasse range1848sage: E.order()18491218501851Or, the order can be set incorrectly along with num_checks set1852too small::18531854sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61855sage: E.set_order(4, num_checks=0)1856WARNING: No checking done in set_order1857sage: E.order()1858418591860The value of num_checks must be an integer. Negative values1861are interpreted as zero, which means don't do any checking::18621863sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61864sage: E.set_order(4, num_checks=-12)1865WARNING: No checking done in set_order1866sage: E.order()1867418681869NOTES:18701871The implementation is based on the fact that orders of elliptic curves1872are cached in the (pseudo-private) _order slot.18731874AUTHORS:18751876- Mariah Lenox (2011-02-16)1877"""1878# Is value in the Hasse range?1879q = self.base_field().order()1880a,b = Hasse_bounds(q,1)1881#a = q + 1 - 2*q.isqrt()1882#b = q + 1 + 2*q.isqrt()1883if not value in ZZ:1884raise ValueError('Value %s illegal (not an integer in the Hasse range)'%value)1885if not a <= value <= b:1886raise ValueError('Value %s illegal (not an integer in the Hasse range)'%value)1887# Is value*random == identity?1888for i in range(num_checks):1889G = self.random_point()1890if value * G != self(0):1891raise ValueError('Value %s illegal (multiple of random point not the identity)'%value)1892if(num_checks <= 0):1893print 'WARNING: No checking done in set_order'1894self._order = value18951896def supersingular_j_polynomial(p):1897"""1898Return a polynomial whose roots are the supersingular `j`-invariants1899in characteristic `p`, other than 0, 1728.19001901INPUT:19021903- `p` (integer) -- a prime number.19041905ALGORITHM:19061907First compute H(X) whose roots are the Legendre1908`\lambda`-invariants of supersingular curves (Silverman V.4.1(b))1909in charactersitic `p`. Then, using a resultant computation with1910the polynomial relating `\lambda` and `j` (Silverman III.1.7(b)),1911we recover the polynomial (in variable ``j``) whose roots are the1912`j`-invariants. Factors of `j` and `j-1728` are removed if1913present.19141915EXAMPLES::19161917sage: from sage.schemes.elliptic_curves.ell_finite_field import supersingular_j_polynomial1918sage: f = supersingular_j_polynomial(67); f1919j^5 + 53*j^4 + 4*j^3 + 47*j^2 + 36*j + 81920sage: f.factor()1921(j + 1) * (j^2 + 8*j + 45) * (j^2 + 44*j + 24)19221923::19241925sage: [supersingular_j_polynomial(p) for p in prime_range(30)]1926[1, 1, 1, 1, 1, j + 8, j + 9, j + 12, j + 4, j^2 + 2*j + 21]19271928TESTS::19291930sage: supersingular_j_polynomial(6)1931Traceback (most recent call last):1932...1933ValueError: p (=6) should be a prime number19341935"""1936try:1937p = ZZ(p)1938except TypeError:1939raise ValueError, "p (=%s) should be a prime number"%p1940if not p.is_prime():1941raise ValueError, "p (=%s) should be a prime number"%p19421943J = polygen(GF(p),'j')1944if p<13:1945return J.parent().one()1946from sage.rings.all import binomial1947from sage.misc.all import prod1948m=(p-1)//21949X,T = PolynomialRing(GF(p),2,names=['X','T']).gens()1950H = sum([binomial(m,i)**2 * T**i for i in xrange(m+1)])1951F = T**2 * (T-1)**2 * X - 256*(T**2-T+1)**31952R = F.resultant(H,T)1953R = prod([fi for fi,e in R([J,0]).factor()])1954if R(0)==0:1955R = R//J1956if R(1728)==0:1957R = R//(J-1728)1958return R19591960# For p in [13..300] we have precomputed these polynomials and store1961# them (as lists of their coefficients in ZZ) in a dict:19621963supersingular_j_polynomials = dict()19641965supersingular_j_polynomials[13] = [8, 1]1966supersingular_j_polynomials[17] = [9, 1]1967supersingular_j_polynomials[19] = [12, 1]1968supersingular_j_polynomials[23] = [4, 1]1969supersingular_j_polynomials[29] = [21, 2, 1]1970supersingular_j_polynomials[31] = [8, 25, 1]1971supersingular_j_polynomials[37] = [11, 5, 23, 1]1972supersingular_j_polynomials[41] = [18, 10, 19, 1]1973supersingular_j_polynomials[43] = [32, 11, 21, 1]1974supersingular_j_polynomials[47] = [35, 33, 31, 1]1975supersingular_j_polynomials[53] = [24, 9, 30, 7, 1]1976supersingular_j_polynomials[59] = [39, 31, 35, 39, 1]1977supersingular_j_polynomials[61] = [60, 21, 27, 8, 60, 1]1978supersingular_j_polynomials[67] = [8, 36, 47, 4, 53, 1]1979supersingular_j_polynomials[71] = [18, 54, 28, 33, 1, 1]1980supersingular_j_polynomials[73] = [7, 39, 38, 9, 68, 60, 1]1981supersingular_j_polynomials[79] = [10, 25, 1, 63, 57, 55, 1]1982supersingular_j_polynomials[83] = [43, 72, 81, 81, 62, 11, 1]1983supersingular_j_polynomials[89] = [42, 79, 23, 22, 37, 86, 60, 1]1984supersingular_j_polynomials[97] = [19, 28, 3, 72, 2, 96, 10, 60, 1]1985supersingular_j_polynomials[101] = [9, 76, 45, 79, 1, 68, 87, 60, 1]1986supersingular_j_polynomials[103] = [64, 15, 24, 58, 70, 83, 84, 100, 1]1987supersingular_j_polynomials[107] = [6, 18, 72, 59, 43, 19, 17, 68, 1]1988supersingular_j_polynomials[109] = [107, 22, 39, 83, 30, 34, 108, 104, 60, 1]1989supersingular_j_polynomials[113] = [86, 71, 75, 6, 47, 97, 100, 4, 60, 1]1990supersingular_j_polynomials[127] = [32, 31, 5, 50, 115, 122, 114, 67, 38, 35, 1]1991supersingular_j_polynomials[131] = [65, 64, 10, 34, 129, 35, 94, 127, 7, 7, 1]1992supersingular_j_polynomials[137] = [104, 83, 3, 82, 112, 23, 77, 135, 18, 50, 60, 1]1993supersingular_j_polynomials[139] = [87, 79, 109, 21, 138, 9, 104, 130, 61, 118, 90, 1]1994supersingular_j_polynomials[149] = [135, 55, 80, 86, 87, 74, 32, 60, 130, 80, 146, 60, 1]1995supersingular_j_polynomials[151] = [94, 125, 8, 6, 93, 21, 114, 80, 107, 58, 42, 18, 1]1996supersingular_j_polynomials[157] = [14, 95, 22, 58, 110, 23, 71, 51, 47, 5, 147, 59, 60, 1]1997supersingular_j_polynomials[163] = [102, 26, 74, 95, 112, 151, 98, 107, 27, 37, 25, 111, 109, 1]1998supersingular_j_polynomials[167] = [14, 9, 27, 109, 97, 55, 51, 74, 145, 125, 36, 113, 89, 1]1999supersingular_j_polynomials[173] = [152, 73, 56, 12, 18, 96, 98, 49, 30, 43, 52, 79, 163, 60, 1]2000supersingular_j_polynomials[179] = [110, 51, 3, 94, 123, 90, 156, 90, 88, 119, 158, 27, 71, 29, 1]2001supersingular_j_polynomials[181] = [7, 65, 77, 29, 139, 34, 65, 84, 164, 73, 51, 136, 7, 141, 60, 1]2002supersingular_j_polynomials[191] = [173, 140, 144, 3, 135, 80, 182, 84, 93, 75, 83, 17, 22, 42, 160, 1]2003supersingular_j_polynomials[193] = [23, 48, 26, 15, 108, 141, 124, 44, 132, 49, 72, 173, 126, 101, 22, 60, 1]2004supersingular_j_polynomials[197] = [14, 111, 64, 170, 193, 32, 124, 91, 112, 163, 14, 112, 167, 191, 183, 60, 1]2005supersingular_j_polynomials[199] = [125, 72, 65, 30, 63, 45, 10, 177, 91, 102, 28, 27, 5, 150, 51, 128, 1]2006supersingular_j_polynomials[211] = [27, 137, 128, 90, 102, 141, 5, 77, 131, 144, 83, 108, 23, 105, 98, 13, 80, 1]2007supersingular_j_polynomials[223] = [56, 183, 46, 133, 191, 94, 20, 8, 92, 100, 57, 200, 166, 67, 59, 218, 28, 32, 1]2008supersingular_j_polynomials[227] = [79, 192, 142, 66, 11, 114, 100, 208, 57, 147, 32, 5, 144, 93, 185, 147, 92, 16, 1]2009supersingular_j_polynomials[229] = [22, 55, 182, 130, 228, 172, 63, 25, 108, 99, 100, 101, 220, 111, 205, 199, 91, 163, 60, 1]2010supersingular_j_polynomials[233] = [101, 148, 85, 113, 226, 68, 71, 103, 61, 44, 173, 175, 5, 225, 227, 99, 146, 170, 60, 1]2011supersingular_j_polynomials[239] = [225, 81, 47, 26, 133, 182, 238, 2, 144, 154, 234, 178, 165, 130, 35, 61, 144, 112, 207, 1]2012supersingular_j_polynomials[241] = [224, 51, 227, 139, 134, 186, 187, 152, 161, 175, 213, 59, 105, 88, 87, 124, 202, 40, 15, 60, 1]2013supersingular_j_polynomials[251] = [30, 183, 80, 127, 40, 56, 230, 168, 192, 48, 226, 61, 214, 54, 165, 147, 105, 88, 38, 171, 1]2014supersingular_j_polynomials[257] = [148, 201, 140, 146, 169, 147, 220, 4, 205, 224, 35, 42, 198, 97, 127, 7, 110, 229, 118, 202, 60, 1]2015supersingular_j_polynomials[263] = [245, 126, 72, 213, 14, 64, 152, 83, 169, 114, 9, 128, 138, 231, 103, 85, 114, 211, 173, 249, 135, 1]2016supersingular_j_polynomials[269] = [159, 32, 69, 95, 201, 266, 190, 176, 76, 151, 212, 21, 106, 49, 263, 105, 136, 194, 215, 181, 237, 60, 1]2017supersingular_j_polynomials[271] = [169, 87, 179, 109, 133, 101, 31, 167, 208, 99, 127, 120, 83, 62, 36, 23, 61, 50, 69, 263, 265, 111, 1]2018supersingular_j_polynomials[277] = [251, 254, 171, 72, 190, 237, 12, 231, 123, 217, 263, 151, 270, 183, 29, 228, 85, 4, 67, 101, 29, 169, 60, 1]2019supersingular_j_polynomials[281] = [230, 15, 146, 69, 41, 23, 142, 232, 18, 80, 58, 134, 270, 62, 272, 70, 247, 189, 118, 255, 274, 159, 60, 1]2020supersingular_j_polynomials[283] = [212, 4, 42, 155, 38, 1, 270, 175, 172, 256, 264, 232, 50, 82, 244, 127, 148, 46, 249, 72, 59, 124, 75, 1]2021supersingular_j_polynomials[293] = [264, 66, 165, 144, 243, 25, 163, 210, 18, 107, 160, 153, 70, 255, 91, 211, 22, 7, 256, 50, 150, 94, 225, 60, 1]202220232024def is_j_supersingular(j, proof=True):2025r"""2026Return True if `j` is a supersingular `j`-invariant.20272028INPUT:20292030- ``j`` (finite field element) -- an element of a finite field20312032- ``proof`` (boolean, default True) -- If True, returns a proved2033result. If False, then a return value of False is certain but a2034return value of True may be based on a probabilistic test. See2035the ALGORITHM section below for more details.20362037OUTPUT:20382039(boolean) True if `j` is supersingular, else False.20402041ALGORITHM:20422043For small characteristics `p` we check whether the `j`-invariant2044is in a precomputed list of supersingular values. Otherwise we2045next check the `j`-invariant. If `j=0`, the curve is2046supersingular if and only if `p=2` or `p\equiv3\pmod{4}`; if2047`j=1728`, the curve is supersingular if and only if `p=3` or2048`p\equiv2\pmod{3}`. Next, if the base field is the prime field2049`{\rm GF}(p)`, we check that `(p+1)P=0` for several random points2050`P`, returning False if any fail: supersingular curves over `{\rm2051GF}(p)` have cardinality `p+1`. If Proof is false we now return2052True. Otherwise we compute the cardinality and return True if and2053only if it is divisible by `p`.20542055EXAMPLES::20562057sage: from sage.schemes.elliptic_curves.ell_finite_field import is_j_supersingular, supersingular_j_polynomials2058sage: [(p,[j for j in GF(p) if is_j_supersingular(j)]) for p in prime_range(30)]2059[(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])]20602061sage: [j for j in GF(109) if is_j_supersingular(j)]2062[17, 41, 43]2063sage: PolynomialRing(GF(109),'j')(supersingular_j_polynomials[109]).roots()2064[(43, 1), (41, 1), (17, 1)]20652066sage: [p for p in prime_range(100) if is_j_supersingular(GF(p)(0))]2067[2, 3, 5, 11, 17, 23, 29, 41, 47, 53, 59, 71, 83, 89]2068sage: [p for p in prime_range(100) if is_j_supersingular(GF(p)(1728))]2069[2, 3, 7, 11, 19, 23, 31, 43, 47, 59, 67, 71, 79, 83]2070sage: [p for p in prime_range(100) if is_j_supersingular(GF(p)(123456))]2071[2, 3, 59, 89]20722073"""2074if not is_FiniteFieldElement(j):2075raise ValueError, "%s must be an element of a finite field"%j20762077F = j.parent()2078p = F.characteristic()2079d = F.degree()20802081if j.is_zero():2082return p==3 or p%3==220832084if (j-1728).is_zero():2085return p==2 or p%4==320862087# From now on we know that j != 0, 172820882089if p in (2,3,5,7,11):2090return False # since j=0, 1728 are the only s.s. invariants20912092# supersingular j-invariants have degree at most 2:20932094jpol = j.minimal_polynomial()2095degj = jpol.degree()2096if degj > 2:2097return False20982099# if p occurs in the precomputed list, use that:21002101try:2102coeffs = supersingular_j_polynomials[p]2103return PolynomialRing(F,'x')(coeffs)(j).is_zero()2104except KeyError:2105pass21062107# Over GF(p), supersingular elliptic curves have cardinality2108# exactly p+1, so we check some random points in order to detect2109# non-supersingularity. Over GF(p^2) (for p at least 5) the2110# cardinality is either (p-1)^2 or (p+1)^2, and the group has2111# exponent p+1 or p-1, so we can do a similar random check: unless2112# (p+1)*P=0 for all the random points, or (p-1)*P=0 for all of2113# them, we can certainly return False.21142115# First we replace j by an element of GF(p) or GF(p^2) (since F2116# might be a proper extension of these):21172118if degj==1:2119j = -jpol(0) # = j, but in GF(p)2120elif d>2:2121F = GF(p^2,'a')2122j = jpol.roots(F,multiplicities=False)[0] # j, but in GF(p^2)21232124E = EllipticCurve(j=j)2125if degj==1:2126for i in range(10):2127P = E.random_element()2128if not ((p+1)*P).is_zero():2129return False2130else:2131n = None # will hold either p+1 or p-1 later2132for i in range(10):2133P = E.random_element()2134# avoid 2-torsion; we know that a1=a3=0 and #E>4!2135while P[2].is_zero() or P[1].is_zero():2136P = E.random_element()21372138if n is None: # not yet decided between p+1 and p-12139pP = p*P2140if not pP[0]==P[0]: # i.e. pP is neither P nor -P2141return False2142if pP[1]==P[1]: # then p*P == P != -P2143n=p-12144else: # then p*P == -P != P2145n=p+12146else:2147if not (n*P).is_zero():2148return False214921502151# when proof is False we return True for any curve which passes2152# the probabilistic test:21532154if not proof:2155return True21562157# otherwise we check the trace of Frobenius (which could be2158# expensive since it involves counting the number of points on E):21592160return E.trace_of_frobenius() % p == 021612162216321642165