Path: blob/master/sage/schemes/elliptic_curves/ell_number_field.py
4162 views
# -*- coding: utf-8 -*-1r"""2Elliptic curves over number fields34An elliptic curve `E` over a number field `K` can be given5by a Weierstrass equation whose coefficients lie in `K` or by6using ``base_extend`` on an elliptic curve defined over a subfield.78One major difference to elliptic curves over `\QQ` is that there9might not exist a global minimal equation over `K`, when `K` does10not have class number one.11Another difference is the lack of understanding of modularity for12general elliptic curves over general number fields.1314Currently Sage can obtain local information about `E/K_v` for finite places15`v`, it has an interface to Denis Simon's script for 2-descent, it can compute16the torsion subgroup of the Mordell-Weil group `E(K)`, and it can work with17isogenies defined over `K`.1819EXAMPLE::2021sage: K.<i> = NumberField(x^2+1)22sage: E = EllipticCurve([0,4+i])23sage: E.discriminant()24-3456*i - 648025sage: P= E([i,2])26sage: P+P27(-2*i + 9/16 : -9/4*i - 101/64 : 1)2829::3031sage: E.has_good_reduction(2+i)32True33sage: E.local_data(4+i)34Local data at Fractional ideal (i + 4):35Reduction type: bad additive36Local minimal model: Elliptic Curve defined by y^2 = x^3 + (i+4) over Number Field in i with defining polynomial x^2 + 137Minimal discriminant valuation: 238Conductor exponent: 239Kodaira Symbol: II40Tamagawa Number: 141sage: E.tamagawa_product_bsd()4214344::4546sage: E.simon_two_descent()47(1, 1, [(i : 2 : 1)])4849::5051sage: E.torsion_order()5215354::5556sage: E.isogenies_prime_degree(3)57[Isogeny of degree 3 from Elliptic Curve defined by y^2 = x^3 + (i+4) over Number Field in i with defining polynomial x^2 + 1 to Elliptic Curve defined by y^2 = x^3 + (-27*i-108) over Number Field in i with defining polynomial x^2 + 1]5859AUTHORS:6061- Robert Bradshaw 20076263- John Cremona6465- Chris Wuthrich6667REFERENCE:6869- [Sil] Silverman, Joseph H. The arithmetic of elliptic curves. Second edition. Graduate Texts in70Mathematics, 106. Springer, 2009.7172- [Sil2] Silverman, Joseph H. Advanced topics in the arithmetic of elliptic curves. Graduate Texts in73Mathematics, 151. Springer, 1994.74"""7576#*****************************************************************************77# Copyright (C) 2007 Robert Bradshaw <[email protected]>78# William Stein <[email protected]>79#80# Distributed under the terms of the GNU General Public License (GPL)81#82# This code is distributed in the hope that it will be useful,83# but WITHOUT ANY WARRANTY; without even the implied warranty of84# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU85# General Public License for more details.86#87# The full text of the GPL is available at:88#89# http://www.gnu.org/licenses/90#*****************************************************************************9192from ell_field import EllipticCurve_field93import ell_point94import sage.matrix.all as matrix95from sage.rings.ring import Ring96from sage.rings.arith import gcd, prime_divisors97from sage.misc.misc import prod98import sage.databases.cremona99import ell_torsion100from ell_generic import is_EllipticCurve101102from gp_simon import simon_two_descent103from constructor import EllipticCurve104from sage.rings.all import PolynomialRing, ZZ, RealField105from sage.misc.misc import verbose, forall106from sage.rings.integer import Integer107from sage.rings.infinity import Infinity # just for verbose output108from sage.rings.arith import valuation109110class EllipticCurve_number_field(EllipticCurve_field):111r"""112Elliptic curve over a number field.113114EXAMPLES::115116sage: K.<i>=NumberField(x^2+1)117sage: EllipticCurve([i, i - 1, i + 1, 24*i + 15, 14*i + 35])118Elliptic Curve defined by y^2 + i*x*y + (i+1)*y = x^3 + (i-1)*x^2 + (24*i+15)*x + (14*i+35) over Number Field in i with defining polynomial x^2 + 1119"""120def __init__(self, x, y=None):121r"""122Allow some ways to create an elliptic curve over a number123field in addition to the generic ones.124125INPUT:126127- ``x``, ``y`` -- see examples.128129EXAMPLES:130131A curve from the database of curves over `\QQ`, but over a larger field:132133sage: K.<i>=NumberField(x^2+1)134sage: EllipticCurve(K,'389a1')135Elliptic Curve defined by y^2 + y = x^3 + x^2 + (-2)*x over Number Field in i with defining polynomial x^2 + 1136137Making the field of definition explicitly larger::138139sage: EllipticCurve(K,[0,-1,1,0,0])140Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 over Number Field in i with defining polynomial x^2 + 1141142"""143if y is None:144if isinstance(x, list):145ainvs = x146field = ainvs[0].parent()147else:148if isinstance(y, str):149field = x150X = sage.databases.cremona.CremonaDatabase()[y]151ainvs = list(X.a_invariants())152else:153field = x154ainvs = y155if not (isinstance(field, Ring) and isinstance(ainvs,list)):156raise TypeError157158EllipticCurve_field.__init__(self, [field(x) for x in ainvs])159self._point = ell_point.EllipticCurvePoint_number_field160161def simon_two_descent(self, verbose=0, lim1=5, lim3=50, limtriv=10, maxprob=20, limbigprime=30):162r"""163Computes lower and upper bounds on the rank of the Mordell-Weil group,164and a list of independent points. Used internally by the :meth:`~rank`,165:meth:`~rank_bounds` and :meth:`~gens` methods.166167INPUT:168169- ``verbose`` -- 0, 1, 2, or 3 (default: 0), the verbosity level170171- ``lim1`` -- (default: 5) limit on trivial points on quartics172173- ``lim3`` -- (default: 50) limit on points on ELS quartics174175- ``limtriv`` -- (default: 10) limit on trivial points on elliptic curve176177- ``maxprob`` -- (default: 20)178179- ``limbigprime`` -- (default: 30) to distinguish between180small and large prime numbers. Use probabilistic tests for181large primes. If 0, don't use probabilistic tests.182183OUTPUT:184185``(lower, upper, list)`` where ``lower`` is a lower bound on186the rank, ``upper`` is an upper bound (the 2-Selmer rank) and187``list`` is a list of independent points on the Weierstrass188model. The length of ``list`` is equal to either ``lower``,189or ``lower-1``, since when ``lower`` is less than ``upper``190and of different parity, the value of ``lower`` is increased by1911.192193.. note::194195For non-quadratic number fields, this code does return, but196it takes a long time.197198ALGORITHM:199200Uses Denis Simon's PARI/GP scripts from201http://www.math.unicaen.fr/~simon/.202203EXAMPLES::204205sage: K.<a> = NumberField(x^2 + 23, 'a')206sage: E = EllipticCurve(K, '37')207sage: E == loads(dumps(E))208True209sage: E.simon_two_descent()210(2, 2, [(-1 : 0 : 1), (1/2*a - 5/2 : -1/2*a - 13/2 : 1)])211212::213214sage: K.<a> = NumberField(x^2 + 7, 'a')215sage: E = EllipticCurve(K, [0,0,0,1,a]); E216Elliptic Curve defined by y^2 = x^3 + x + a over Number Field in a with defining polynomial x^2 + 7217218sage: v = E.simon_two_descent(verbose=1); v219courbe elliptique : Y^2 = x^3 + Mod(1, y^2 + 7)*x + Mod(y, y^2 + 7)220points triviaux sur la courbe = [[1, 1, 0], [Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7), 1]]221#S(E/K)[2] = 2222#E(K)/2E(K) = 2223#III(E/K)[2] = 1224rang(E/K) = 1225listpointsmwr = [[Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7), 1]]226(1, 1, [(1/2*a + 3/2 : -a - 2 : 1)])227228sage: v = E.simon_two_descent(verbose=2) # random output229K = bnfinit(y^2 + 7);230a = Mod(y,K.pol);231bnfellrank(K, [0,0,0,1,a]);232courbe elliptique : Y^2 = x^3 + Mod(1, y^2 + 7)*x + Mod(y, y^2 + 7)233A = 0234B = Mod(1, y^2 + 7)235C = Mod(y, y^2 + 7)236LS2gen = [Mod(Mod(-5, y^2 + 7)*x^2 + Mod(-3*y, y^2 + 7)*x + Mod(8, y^2 + 7), x^3 + Mod(1, y^2 + 7)*x + Mod(y, y^2 + 7)), Mod(Mod(1, y^2 + 7)*x^2 + Mod(1/2*y - 1/2, y^2 + 7)*x - 1, x^3 + Mod(1, y^2 + 7)*x + Mod(y, y^2 + 7))]237#LS2gen = 2238Recherche de points triviaux sur la courbe239points triviaux sur la courbe = [[1, 1, 0], [Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7), 1]]240zc = Mod(Mod(-5, y^2 + 7)*x^2 + Mod(-3*y, y^2 + 7)*x + Mod(8, y^2 + 7), x^3 + Mod(1, y^2 + 7)*x + Mod(y, y^2 + 7))241symbole de Hilbert (Mod(2, y^2 + 7),Mod(-5, y^2 + 7)) = -1242zc = Mod(Mod(1, y^2 + 7)*x^2 + Mod(1/2*y - 1/2, y^2 + 7)*x + Mod(-1, y^2 + 7), x^3 + Mod(1, y^2 + 7)*x + Mod(y, y^2 + 7))243symbole de Hilbert (Mod(-2*y + 2, y^2 + 7),Mod(1, y^2 + 7)) = 0244sol de Legendre = [1, 0, 1]~245zc*z1^2 = Mod(Mod(2*y - 2, y^2 + 7)*x + Mod(2*y + 10, y^2 + 7), x^3 + Mod(1, y^2 + 7)*x + Mod(y, y^2 + 7))246quartique : (-1/2*y + 1/2)*Y^2 = x^4 + (-3*y - 15)*x^2 + (-8*y - 16)*x + (-11/2*y - 15/2)247reduite: Y^2 = (-1/2*y + 1/2)*x^4 - 4*x^3 + (-3*y + 3)*x^2 + (2*y - 2)*x + (1/2*y + 3/2)248non ELS en [2, [0, 1]~, 1, 1, [1, 1]~]249zc = Mod(Mod(1, y^2 + 7)*x^2 + Mod(1/2*y + 1/2, y^2 + 7)*x + Mod(-1, y^2 + 7), x^3 + Mod(1, y^2 + 7)*x + Mod(y, y^2 + 7))250vient du point trivial [Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7), 1]251m1 = 1252m2 = 1253#S(E/K)[2] = 2254#E(K)/2E(K) = 2255#III(E/K)[2] = 1256rang(E/K) = 1257listpointsmwr = [[Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7), 1]]258v = [1, 1, [[Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7)]]]259sage: v260(1, 1, [(1/2*a + 3/2 : -a - 2 : 1)])261262263A curve with 2-torsion::264265sage: K.<a> = NumberField(x^2 + 7, 'a')266sage: E = EllipticCurve(K, '15a')267sage: v = E.simon_two_descent(); v # long time (about 10 seconds), points can vary268(1, 3, [...])269"""270271try:272result = self._simon_two_descent_data[lim1,lim3,limtriv,maxprob,limbigprime]273if verbose == 0:274return result275except AttributeError:276self._simon_two_descent_data = {}277278x = PolynomialRing(self.base_ring(), 'x').gen(0)279t = simon_two_descent(self,280verbose=verbose, lim1=lim1, lim3=lim3, limtriv=limtriv,281maxprob=maxprob, limbigprime=limbigprime)282prob_rank = Integer(t[0])283two_selmer_rank = Integer(t[1])284prob_gens = [self(P) for P in t[2]]285self._simon_two_descent_data[lim1,lim3,limtriv,maxprob,limbigprime] = (prob_rank, two_selmer_rank, prob_gens)286return prob_rank, two_selmer_rank, prob_gens287288def height_pairing_matrix(self, points=None, precision=None):289r"""290Returns the height pairing matrix of the given points.291292INPUT:293294- points - either a list of points, which must be on this295curve, or (default) None, in which case self.gens() will be296used.297298- precision - number of bits of precision of result299(default: None, for default RealField precision)300301EXAMPLES::302303sage: E = EllipticCurve([0, 0, 1, -1, 0])304sage: E.height_pairing_matrix()305[0.0511114082399688]306307For rank 0 curves, the result is a valid 0x0 matrix::308309sage: EllipticCurve('11a').height_pairing_matrix()310[]311sage: E=EllipticCurve('5077a1')312sage: E.height_pairing_matrix([E.lift_x(x) for x in [-2,-7/4,1]], precision=100)313[ 1.3685725053539301120518194471 -1.3095767070865761992624519454 -0.63486715783715592064475542573]314[ -1.3095767070865761992624519454 2.7173593928122930896610589220 1.0998184305667292139777571432]315[-0.63486715783715592064475542573 1.0998184305667292139777571432 0.66820516565192793503314205089]316317sage: E = EllipticCurve('389a1')318sage: E = EllipticCurve('389a1')319sage: P,Q = E.point([-1,1,1]),E.point([0,-1,1])320sage: E.height_pairing_matrix([P,Q])321[0.686667083305587 0.268478098806726]322[0.268478098806726 0.327000773651605]323324Over a number field::325326sage: x = polygen(QQ)327sage: K.<t> = NumberField(x^2+47)328sage: EK = E.base_extend(K)329sage: EK.height_pairing_matrix([EK(P),EK(Q)])330[0.686667083305586 0.268478098806726]331[0.268478098806726 0.327000773651605]332333::334335sage: K.<i> = QuadraticField(-1)336sage: E = EllipticCurve([0,0,0,i,i])337sage: P = E(-9+4*i,-18-25*i)338sage: Q = E(i,-i)339sage: E.height_pairing_matrix([P,Q])340[ 2.16941934493768 -0.870059380421505]341[-0.870059380421505 0.424585837470709]342sage: E.regulator_of_points([P,Q])3430.164101403936070344"""345if points is None:346points = self.gens()347else:348for P in points:349assert P.curve() == self350351r = len(points)352if precision is None:353RR = RealField()354else:355RR = RealField(precision)356M = matrix.MatrixSpace(RR, r)357mat = M()358for j in range(r):359mat[j,j] = points[j].height(precision=precision)360for j in range(r):361for k in range(j+1,r):362mat[j,k]=((points[j]+points[k]).height(precision=precision) - mat[j,j] - mat[k,k])/2363mat[k,j]=mat[j,k]364return mat365366def regulator_of_points(self, points=[], precision=None):367"""368Returns the regulator of the given points on this curve.369370INPUT:371372- ``points`` -(default: empty list) a list of points on this curve373374- ``precision`` - int or None (default: None): the precision375in bits of the result (default real precision if None)376377EXAMPLES::378379sage: E = EllipticCurve('37a1')380sage: P = E(0,0)381sage: Q = E(1,0)382sage: E.regulator_of_points([P,Q])3830.000000000000000384sage: 2*P==Q385True386387::388389sage: E = EllipticCurve('5077a1')390sage: points = [E.lift_x(x) for x in [-2,-7/4,1]]391sage: E.regulator_of_points(points)3920.417143558758384393sage: E.regulator_of_points(points,precision=100)3940.41714355875838396981711954462395396::397398sage: E = EllipticCurve('389a')399sage: E.regulator_of_points()4001.00000000000000401sage: points = [P,Q] = [E(-1,1),E(0,-1)]402sage: E.regulator_of_points(points)4030.152460177943144404sage: E.regulator_of_points(points, precision=100)4050.15246017794314375162432475705406sage: E.regulator_of_points(points, precision=200)4070.15246017794314375162432475704945582324372707748663081784028408sage: E.regulator_of_points(points, precision=300)4090.152460177943143751624324757049455823243727077486630817840280980046053225683562463604114816410411Examples over number fields::412413sage: K.<a> = QuadraticField(97)414sage: E = EllipticCurve(K,[1,1])415sage: P = E(0,1)416sage: P.height()4170.476223106404866418sage: E.regulator_of_points([P])4190.476223106404866420421::422423sage: E = EllipticCurve('11a1')424sage: x = polygen(QQ)425sage: K.<t> = NumberField(x^2+47)426sage: EK = E.base_extend(K)427sage: T = EK(5,5)428sage: T.order()4295430sage: P = EK(-2, -1/2*t - 1/2)431sage: P.order()432+Infinity433sage: EK.regulator_of_points([P,T]) # random very small output434-1.23259516440783e-32435sage: EK.regulator_of_points([P,T]).abs() < 1e-30436True437438::439440sage: E = EllipticCurve('389a1')441sage: P,Q = E.gens()442sage: E.regulator_of_points([P,Q])4430.152460177943144444sage: K.<t> = NumberField(x^2+47)445sage: EK = E.base_extend(K)446sage: EK.regulator_of_points([EK(P),EK(Q)])4470.152460177943144448449::450451sage: K.<i> = QuadraticField(-1)452sage: E = EllipticCurve([0,0,0,i,i])453sage: P = E(-9+4*i,-18-25*i)454sage: Q = E(i,-i)455sage: E.height_pairing_matrix([P,Q])456[ 2.16941934493768 -0.870059380421505]457[-0.870059380421505 0.424585837470709]458sage: E.regulator_of_points([P,Q])4590.164101403936070460461"""462if points is None:463points = []464mat = self.height_pairing_matrix(points=points, precision=precision)465return mat.det(algorithm="hessenberg")466467468def is_local_integral_model(self,*P):469r"""470Tests if self is integral at the prime ideal `P`, or at all the471primes if `P` is a list or tuple.472473INPUT:474475- ``*P`` -- a prime ideal, or a list or tuple of primes.476477EXAMPLES::478479sage: K.<i> = NumberField(x^2+1)480sage: P1,P2 = K.primes_above(5)481sage: E = EllipticCurve([i/5,i/5,i/5,i/5,i/5])482sage: E.is_local_integral_model(P1,P2)483False484sage: Emin = E.local_integral_model(P1,P2)485sage: Emin.is_local_integral_model(P1,P2)486True487"""488if len(P)==1: P=P[0]489if isinstance(P,(tuple,list)):490return forall(P, lambda x : self.is_local_integral_model(x))[0]491return forall(self.ainvs(), lambda x : x.valuation(P) >= 0)[0]492493def local_integral_model(self,*P):494r"""495Return a model of self which is integral at the prime ideal496`P`.497498.. note::499500The integrality at other primes is not affected, even if501`P` is non-principal.502503INPUT:504505- ``*P`` -- a prime ideal, or a list or tuple of primes.506507EXAMPLES::508509sage: K.<i> = NumberField(x^2+1)510sage: P1,P2 = K.primes_above(5)511sage: E = EllipticCurve([i/5,i/5,i/5,i/5,i/5])512sage: E.local_integral_model((P1,P2))513Elliptic Curve defined by y^2 + (-i)*x*y + (-25*i)*y = x^3 + 5*i*x^2 + 125*i*x + 3125*i over Number Field in i with defining polynomial x^2 + 1514"""515if len(P)==1: P=P[0]516if isinstance(P,(tuple,list)):517E=self518for Pi in P: E=E.local_integral_model(Pi)519return E520ai = self.a_invariants()521e = min([(ai[i].valuation(P)/[1,2,3,4,6][i]) for i in range(5)]).floor()522pi = self.base_field().uniformizer(P, 'negative')523return EllipticCurve([ai[i]/pi**(e*[1,2,3,4,6][i]) for i in range(5)])524525def is_global_integral_model(self):526r"""527Return true iff self is integral at all primes.528529EXAMPLES::530531sage: K.<i> = NumberField(x^2+1)532sage: E = EllipticCurve([i/5,i/5,i/5,i/5,i/5])533sage: P1,P2 = K.primes_above(5)534sage: Emin = E.global_integral_model()535sage: Emin.is_global_integral_model()536True537"""538return forall(self.a_invariants(), lambda x : x.is_integral())[0]539540def global_integral_model(self):541r"""542Return a model of self which is integral at all primes.543544EXAMPLES::545546sage: K.<i> = NumberField(x^2+1)547sage: E = EllipticCurve([i/5,i/5,i/5,i/5,i/5])548sage: P1,P2 = K.primes_above(5)549sage: E.global_integral_model()550Elliptic Curve defined by y^2 + (-i)*x*y + (-25*i)*y = x^3 + 5*i*x^2 + 125*i*x + 3125*i over Number Field in i with defining polynomial x^2 + 1551552trac #7935::553554sage: K.<a> = NumberField(x^2-38)555sage: E = EllipticCurve([a,1/2])556sage: E.global_integral_model()557Elliptic Curve defined by y^2 = x^3 + 1444*a*x + 27436 over Number Field in a with defining polynomial x^2 - 38558559trac #9266::560561sage: K.<s> = NumberField(x^2-5)562sage: w = (1+s)/2563sage: E = EllipticCurve(K,[2,w])564sage: E.global_integral_model()565Elliptic Curve defined by y^2 = x^3 + 2*x + (1/2*s+1/2) over Number Field in s with defining polynomial x^2 - 5566567trac #12151::568569sage: K.<v> = NumberField(x^2 + 161*x - 150)570sage: E = EllipticCurve([25105/216*v - 3839/36, 634768555/7776*v - 98002625/1296, 634768555/7776*v - 98002625/1296, 0, 0])571sage: E.global_integral_model()572Elliptic Curve defined by y^2 + (33872485050625*v-31078224284250)*x*y + (2020602604156076340058146664245468750000*v-1871778534673615560803175189398437500000)*y = x^3 + (6933305282258321342920781250*v-6422644400723486559914062500)*x^2 over Number Field in v with defining polynomial x^2 + 161*x - 150573"""574K = self.base_field()575ai = self.a_invariants()576for a in ai:577if not a.is_integral():578for P, _ in a.denominator_ideal().factor():579pi = K.uniformizer(P,'positive')580e = min([(ai[i].valuation(P)/[1,2,3,4,6][i]) for i in range(5)]).floor()581ai = [ai[i]/pi**(e*[1,2,3,4,6][i]) for i in range(5)]582for z in ai:583assert z.is_integral(), "bug in global_integral_model: %s" % list(ai)584return EllipticCurve(list(ai))585586integral_model = global_integral_model587588def _reduce_model(self):589r"""590591Transforms the elliptic curve to a model in which `a_1`,592`a_2`, `a_3` are reduced modulo 2, 3, 2 respectively.593594.. note::595596This only works on integral models, i.e. it requires that597`a_1`, `a_2` and `a_3` lie in the ring of integers of the base598field.599600EXAMPLES::601602sage: K.<a>=NumberField(x^2-38)603sage: E=EllipticCurve([a, -5*a + 19, -39*a + 237, 368258520200522046806318224*a - 2270097978636731786720858047, 8456608930180227786550494643437985949781*a - 52130038506835491453281450568107193773505])604sage: E.ainvs()605(a,606-5*a + 19,607-39*a + 237,608368258520200522046806318224*a - 2270097978636731786720858047,6098456608930180227786550494643437985949781*a - 52130038506835491453281450568107193773505)610sage: E._reduce_model().ainvs()611(a,612a + 1,613a + 1,614368258520200522046806318444*a - 2270097978636731786720859345,6158456608930173478039472018047583706316424*a - 52130038506793883217874390501829588391299)616sage: EllipticCurve([101,202,303,404,505])._reduce_model().ainvs()617(1, 1, 0, -2509254, 1528863051)618sage: EllipticCurve([-101,-202,-303,-404,-505])._reduce_model().ainvs()619(1, -1, 0, -1823195, 947995262)620"""621ZK = self.base_ring().maximal_order()622try:623(a1, a2, a3, a4, a6) = [ZK(a) for a in self.a_invariants()]624except TypeError:625raise TypeError, "_reduce_model() requires an integral model."626# N.B. Must define s, r, t in the right order.627if ZK.degree() == 1:628s = ((-a1)/2).round('up')629r = ((-a2 + s*a1 +s*s)/3).round()630t = ((-a3 - r*a1)/2).round('up')631else:632s = ZK([(a/2).round('up') for a in (-a1).list()])633r = ZK([(a/3).round() for a in (-a2 + s*a1 +s*s).list()])634t = ZK([(a/2).round('up') for a in (-a3 - r*a1).list()])635636return self.rst_transform(r, s, t)637638def local_information(self, P=None, proof=None):639r"""640\code{local_information} has been renamed \code{local_data}641and is being deprecated.642"""643raise DeprecationWarning, "local_information is deprecated; use local_data instead"644return self.local_data(P,proof)645646def local_data(self, P=None, proof = None, algorithm="pari"):647r"""648Local data for this elliptic curve at the prime `P`.649650INPUT:651652- ``P`` -- either None or a prime ideal of the base field of self.653654- ``proof`` -- whether to only use provably correct methods655(default controlled by global proof module). Note that the656proof module is number_field, not elliptic_curves, since the657functions that actually need the flag are in number fields.658659- ``algorithm`` (string, default: "pari") -- Ignored unless the660base field is `\QQ`. If "pari", use the PARI C-library661``ellglobalred`` implementation of Tate's algorithm over662`\QQ`. If "generic", use the general number field663implementation.664665OUTPUT:666667If `P` is specified, returns the ``EllipticCurveLocalData``668object associated to the prime `P` for this curve. Otherwise,669returns a list of such objects, one for each prime `P` in the670support of the discriminant of this model.671672.. note::673674The model is not required to be integral on input.675676For principal `P`, a generator is used as a uniformizer,677and integrality or minimality at other primes is not678affected. For non-principal `P`, the minimal model679returned will preserve integrality at other primes, but not680minimality.681682EXAMPLES::683684sage: K.<i> = NumberField(x^2+1)685sage: E = EllipticCurve([1 + i, 0, 1, 0, 0])686sage: E.local_data()687[Local data at Fractional ideal (i - 2):688Reduction type: bad non-split multiplicative689Local minimal model: Elliptic Curve defined by y^2 + (i+1)*x*y + y = x^3 over Number Field in i with defining polynomial x^2 + 1690Minimal discriminant valuation: 1691Conductor exponent: 1692Kodaira Symbol: I1693Tamagawa Number: 1,694Local data at Fractional ideal (-3*i - 2):695Reduction type: bad split multiplicative696Local minimal model: Elliptic Curve defined by y^2 + (i+1)*x*y + y = x^3 over Number Field in i with defining polynomial x^2 + 1697Minimal discriminant valuation: 2698Conductor exponent: 1699Kodaira Symbol: I2700Tamagawa Number: 2]701sage: E.local_data(K.ideal(3))702Local data at Fractional ideal (3):703Reduction type: good704Local minimal model: Elliptic Curve defined by y^2 + (i+1)*x*y + y = x^3 over Number Field in i with defining polynomial x^2 + 1705Minimal discriminant valuation: 0706Conductor exponent: 0707Kodaira Symbol: I0708Tamagawa Number: 1709710An example raised in \#3897::711712sage: E = EllipticCurve([1,1])713sage: E.local_data(3)714Local data at Principal ideal (3) of Integer Ring:715Reduction type: good716Local minimal model: Elliptic Curve defined by y^2 = x^3 + x + 1 over Rational Field717Minimal discriminant valuation: 0718Conductor exponent: 0719Kodaira Symbol: I0720Tamagawa Number: 1721"""722if proof is None:723import sage.structure.proof.proof724# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.725proof = sage.structure.proof.proof.get_flag(None, "number_field")726727if P is None:728primes = self.base_ring()(self.integral_model().discriminant()).support()729return [self._get_local_data(pr, proof) for pr in primes]730731from sage.schemes.elliptic_curves.ell_local_data import check_prime732P = check_prime(self.base_field(),P)733734return self._get_local_data(P,proof,algorithm)735736def _get_local_data(self, P, proof, algorithm="pari"):737r"""738Internal function to create data for this elliptic curve at the prime `P`.739740This function handles the caching of local data. It is called741by local_data() which is the user interface and which parses742the input parameters `P` and proof.743744INPUT:745746- ``P`` -- either None or a prime ideal of the base field of self.747748- ``proof`` -- whether to only use provably correct methods749(default controlled by global proof module). Note that the750proof module is number_field, not elliptic_curves, since the751functions that actually need the flag are in number fields.752753- ``algorithm`` (string, default: "pari") -- Ignored unless the754base field is `\QQ`. If "pari", use the PARI C-library755``ellglobalred`` implementation of Tate's algorithm over756`\QQ`. If "generic", use the general number field757implementation.758759EXAMPLES::760761sage: K.<i> = NumberField(x^2+1)762sage: E = EllipticCurve(K,[0,1,0,-160,308])763sage: p = K.ideal(i+1)764sage: E._get_local_data(p, False)765Local data at Fractional ideal (i + 1):766Reduction type: good767Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + (-10)*x + (-10) over Number Field in i with defining polynomial x^2 + 1768Minimal discriminant valuation: 0769Conductor exponent: 0770Kodaira Symbol: I0771Tamagawa Number: 1772773Verify that we cache based on the proof value and the algorithm choice::774775sage: E._get_local_data(p, False) is E._get_local_data(p, True)776False777778sage: E._get_local_data(p, None, "pari") is E._get_local_data(p, None, "generic")779False780"""781try:782return self._local_data[P, proof, algorithm]783except AttributeError:784self._local_data = {}785except KeyError:786pass787from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData788self._local_data[P, proof, algorithm] = EllipticCurveLocalData(self, P, proof, algorithm)789return self._local_data[P, proof, algorithm]790791def local_minimal_model(self, P, proof = None, algorithm="pari"):792r"""793Returns a model which is integral at all primes and minimal at `P`.794795INPUT:796797- ``P`` -- either None or a prime ideal of the base field of self.798799- ``proof`` -- whether to only use provably correct methods800(default controlled by global proof module). Note that the801proof module is number_field, not elliptic_curves, since the802functions that actually need the flag are in number fields.803804- ``algorithm`` (string, default: "pari") -- Ignored unless the805base field is `\QQ`. If "pari", use the PARI C-library806``ellglobalred`` implementation of Tate's algorithm over807`\QQ`. If "generic", use the general number field808implementation.809810OUTPUT:811812A model of the curve which is minimal (and integral) at `P`.813814.. note::815816The model is not required to be integral on input.817818For principal `P`, a generator is used as a uniformizer,819and integrality or minimality at other primes is not820affected. For non-principal `P`, the minimal model821returned will preserve integrality at other primes, but not822minimality.823824EXAMPLES::825826sage: K.<a>=NumberField(x^2-5)827sage: E=EllipticCurve([20, 225, 750, 625*a + 6875, 31250*a + 46875])828sage: P=K.ideal(a)829sage: E.local_minimal_model(P).ainvs()830(0, 1, 0, a - 33, -2*a + 64)831"""832if proof is None:833import sage.structure.proof.proof834# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.835proof = sage.structure.proof.proof.get_flag(None, "number_field")836837return self.local_data(P, proof, algorithm).minimal_model()838839def has_good_reduction(self, P):840r"""841Return True if this elliptic curve has good reduction at the prime `P`.842843INPUT:844845- ``P`` -- a prime ideal of the base field of self, or a field846element generating such an ideal.847848OUTPUT:849850(bool) -- True if the curve has good reduction at `P`, else False.851852.. note::853854This requires determining a local integral minimal model;855we do not just check that the discriminant of the current856model has valuation zero.857858EXAMPLES::859860sage: E=EllipticCurve('14a1')861sage: [(p,E.has_good_reduction(p)) for p in prime_range(15)]862[(2, False), (3, True), (5, True), (7, False), (11, True), (13, True)]863864sage: K.<a>=NumberField(x^3-2)865sage: P17a, P17b = [P for P,e in K.factor(17)]866sage: E = EllipticCurve([0,0,0,0,2*a+1])867sage: [(p,E.has_good_reduction(p)) for p in [P17a,P17b]]868[(Fractional ideal (4*a^2 - 2*a + 1), True),869(Fractional ideal (2*a + 1), False)]870"""871return self.local_data(P).has_good_reduction()872873def has_bad_reduction(self, P):874r"""875Return True if this elliptic curve has bad reduction at the prime `P`.876877INPUT:878879- ``P`` -- a prime ideal of the base field of self, or a field880element generating such an ideal.881882OUTPUT:883884(bool) True if the curve has bad reduction at `P`, else False.885886.. note::887888This requires determining a local integral minimal model;889we do not just check that the discriminant of the current890model has valuation zero.891892EXAMPLES::893894sage: E=EllipticCurve('14a1')895sage: [(p,E.has_bad_reduction(p)) for p in prime_range(15)]896[(2, True), (3, False), (5, False), (7, True), (11, False), (13, False)]897898sage: K.<a>=NumberField(x^3-2)899sage: P17a, P17b = [P for P,e in K.factor(17)]900sage: E = EllipticCurve([0,0,0,0,2*a+1])901sage: [(p,E.has_bad_reduction(p)) for p in [P17a,P17b]]902[(Fractional ideal (4*a^2 - 2*a + 1), False),903(Fractional ideal (2*a + 1), True)]904"""905return self.local_data(P).has_bad_reduction()906907def has_multiplicative_reduction(self, P):908r"""909Return True if this elliptic curve has (bad) multiplicative reduction at the prime `P`.910911.. note::912913See also ``has_split_multiplicative_reduction()`` and914``has_nonsplit_multiplicative_reduction()``.915916INPUT:917918- ``P`` -- a prime ideal of the base field of self, or a field919element generating such an ideal.920921OUTPUT:922923(bool) True if the curve has multiplicative reduction at `P`,924else False.925926EXAMPLES::927928sage: E=EllipticCurve('14a1')929sage: [(p,E.has_multiplicative_reduction(p)) for p in prime_range(15)]930[(2, True), (3, False), (5, False), (7, True), (11, False), (13, False)]931932sage: K.<a>=NumberField(x^3-2)933sage: P17a, P17b = [P for P,e in K.factor(17)]934sage: E = EllipticCurve([0,0,0,0,2*a+1])935sage: [(p,E.has_multiplicative_reduction(p)) for p in [P17a,P17b]]936[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]937"""938return self.local_data(P).has_multiplicative_reduction()939940def has_split_multiplicative_reduction(self, P):941r"""942Return True if this elliptic curve has (bad) split multiplicative reduction at the prime `P`.943944INPUT:945946- ``P`` -- a prime ideal of the base field of self, or a field947element generating such an ideal.948949OUTPUT:950951(bool) True if the curve has split multiplicative reduction at952`P`, else False.953954EXAMPLES::955956sage: E=EllipticCurve('14a1')957sage: [(p,E.has_split_multiplicative_reduction(p)) for p in prime_range(15)]958[(2, False), (3, False), (5, False), (7, True), (11, False), (13, False)]959960sage: K.<a>=NumberField(x^3-2)961sage: P17a, P17b = [P for P,e in K.factor(17)]962sage: E = EllipticCurve([0,0,0,0,2*a+1])963sage: [(p,E.has_split_multiplicative_reduction(p)) for p in [P17a,P17b]]964[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]965"""966return self.local_data(P).has_split_multiplicative_reduction()967968def has_nonsplit_multiplicative_reduction(self, P):969r"""970Return True if this elliptic curve has (bad) non-split multiplicative reduction at the prime `P`.971972INPUT:973974- ``P`` -- a prime ideal of the base field of self, or a field975element generating such an ideal.976977OUTPUT:978979(bool) True if the curve has non-split multiplicative980reduction at `P`, else False.981982EXAMPLES::983984sage: E=EllipticCurve('14a1')985sage: [(p,E.has_nonsplit_multiplicative_reduction(p)) for p in prime_range(15)]986[(2, True), (3, False), (5, False), (7, False), (11, False), (13, False)]987988sage: K.<a>=NumberField(x^3-2)989sage: P17a, P17b = [P for P,e in K.factor(17)]990sage: E = EllipticCurve([0,0,0,0,2*a+1])991sage: [(p,E.has_nonsplit_multiplicative_reduction(p)) for p in [P17a,P17b]]992[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]993"""994return self.local_data(P).has_nonsplit_multiplicative_reduction()995996def has_additive_reduction(self, P):997r"""998Return True if this elliptic curve has (bad) additive reduction at the prime `P`.9991000INPUT:10011002- ``P`` -- a prime ideal of the base field of self, or a field1003element generating such an ideal.10041005OUTPUT:10061007(bool) True if the curve has additive reduction at `P`, else False.10081009EXAMPLES::10101011sage: E=EllipticCurve('27a1')1012sage: [(p,E.has_additive_reduction(p)) for p in prime_range(15)]1013[(2, False), (3, True), (5, False), (7, False), (11, False), (13, False)]10141015sage: K.<a>=NumberField(x^3-2)1016sage: P17a, P17b = [P for P,e in K.factor(17)]1017sage: E = EllipticCurve([0,0,0,0,2*a+1])1018sage: [(p,E.has_additive_reduction(p)) for p in [P17a,P17b]]1019[(Fractional ideal (4*a^2 - 2*a + 1), False),1020(Fractional ideal (2*a + 1), True)]1021"""1022return self.local_data(P).has_additive_reduction()10231024def tamagawa_number(self, P, proof = None):1025r"""1026Returns the Tamagawa number of this elliptic curve at the prime `P`.10271028INPUT:10291030- ``P`` -- either None or a prime ideal of the base field of self.10311032- ``proof`` -- whether to only use provably correct methods1033(default controlled by global proof module). Note that the1034proof module is number_field, not elliptic_curves, since the1035functions that actually need the flag are in number fields.10361037OUTPUT:10381039(positive integer) The Tamagawa number of the curve at `P`.10401041EXAMPLES::10421043sage: K.<a>=NumberField(x^2-5)1044sage: E=EllipticCurve([20, 225, 750, 625*a + 6875, 31250*a + 46875])1045sage: [E.tamagawa_number(P) for P in E.discriminant().support()]1046[1, 1, 1, 1]1047sage: K.<a> = QuadraticField(-11)1048sage: E = EllipticCurve('11a1').change_ring(K)1049sage: [E.tamagawa_number(P) for P in K(11).support()]1050[10]1051"""1052if proof is None:1053import sage.structure.proof.proof1054# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.1055proof = sage.structure.proof.proof.get_flag(None, "number_field")10561057return self.local_data(P, proof).tamagawa_number()10581059def tamagawa_numbers(self):1060"""1061Return a list of all Tamagawa numbers for all prime divisors of the1062conductor (in order).10631064EXAMPLES::10651066sage: e = EllipticCurve('30a1')1067sage: e.tamagawa_numbers()1068[2, 3, 1]1069sage: vector(e.tamagawa_numbers())1070(2, 3, 1)1071sage: K.<a>=NumberField(x^2+3)1072sage: eK = e.base_extend(K)1073sage: eK.tamagawa_numbers()1074[4, 6, 1]1075"""1076return [self.tamagawa_number(p) for p in prime_divisors(self.conductor())]10771078def tamagawa_exponent(self, P, proof = None):1079r"""1080Returns the Tamagawa index of this elliptic curve at the prime `P`.10811082INPUT:10831084- ``P`` -- either None or a prime ideal of the base field of self.10851086- ``proof`` -- whether to only use provably correct methods1087(default controlled by global proof module). Note that the1088proof module is number_field, not elliptic_curves, since the1089functions that actually need the flag are in number fields.10901091OUTPUT:10921093(positive integer) The Tamagawa index of the curve at P.10941095EXAMPLES::10961097sage: K.<a>=NumberField(x^2-5)1098sage: E=EllipticCurve([20, 225, 750, 625*a + 6875, 31250*a + 46875])1099sage: [E.tamagawa_exponent(P) for P in E.discriminant().support()]1100[1, 1, 1, 1]1101sage: K.<a> = QuadraticField(-11)1102sage: E = EllipticCurve('11a1').change_ring(K)1103sage: [E.tamagawa_exponent(P) for P in K(11).support()]1104[10]1105"""1106if proof is None:1107import sage.structure.proof.proof1108# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.1109proof = sage.structure.proof.proof.get_flag(None, "number_field")11101111return self.local_data(P, proof).tamagawa_exponent()11121113def tamagawa_product_bsd(self):1114r"""1115Given an elliptic curve `E` over a number field `K`, this function returns the1116integer `C(E/K)` that appears in the Birch and Swinnerton-Dyer conjecture accounting1117for the local information at finite places. If the model is a global minimal model then `C(E/K)` is1118simply the product of the Tamagawa numbers `c_v` where `v` runs over all prime ideals of `K`. Otherwise, if the model has to be changed at a place `v` a correction factor appears.1119The definition is such that `C(E/K)` times the periods at the infinite places is invariant1120under change of the Weierstrass model. See [Ta2] and [Do] for details.11211122.. note::11231124This definition is slightly different from the definition of ``tamagawa_product``1125for curves defined over `\QQ`. Over the rational number it is always defined to be the product1126of the Tamagawa numbers, so the two definitions only agree when the model is global minimal.11271128OUTPUT:11291130A rational number11311132EXAMPLES::11331134sage: K.<i> = NumberField(x^2+1)1135sage: E = EllipticCurve([0,2+i])1136sage: E.tamagawa_product_bsd()1137111381139sage: E = EllipticCurve([(2*i+1)^2,i*(2*i+1)^7])1140sage: E.tamagawa_product_bsd()1141411421143An example where the Neron model changes over K::11441145sage: K.<t> = NumberField(x^5-10*x^3+5*x^2+10*x+1)1146sage: E = EllipticCurve(K,'75a1')1147sage: E.tamagawa_product_bsd()114851149sage: da = E.local_data()1150sage: [dav.tamagawa_number() for dav in da]1151[1, 1]11521153An example over `\mathbb{Q}` (trac #9413)::11541155sage: E = EllipticCurve('30a')1156sage: E.tamagawa_product_bsd()1157611581159REFERENCES:11601161- [Ta2] Tate, John, On the conjectures of Birch and Swinnerton-Dyer and a geometric analog. Seminaire Bourbaki, Vol. 9, Exp. No. 306.11621163- [Do] Dokchitser, Tim and Vladimir, On the Birch-Swinnerton-Dyer quotients modulo squares, Annals of Math., 2010.11641165"""1166da = self.local_data()1167pr = 11168for dav in da:1169pp = dav.prime()1170cv = dav.tamagawa_number()1171# uu is the quotient of the Neron differential at pp divided by1172# the differential associated to this particular equation E1173uu = self.isomorphism_to(dav.minimal_model()).u1174if self.base_field().degree() == 1:1175p = pp.gens_reduced()[0]1176f = 11177v = valuation(ZZ(uu),p)1178else:1179p = pp.smallest_integer()1180f = pp.residue_class_degree()1181v = valuation(uu,pp)1182uu_abs_val = p**(f*v)1183pr *= cv * uu_abs_val1184return pr11851186def kodaira_symbol(self, P, proof = None):1187r"""1188Returns the Kodaira Symbol of this elliptic curve at the prime `P`.11891190INPUT:11911192- ``P`` -- either None or a prime ideal of the base field of self.11931194- ``proof`` -- whether to only use provably correct methods1195(default controlled by global proof module). Note that the1196proof module is number_field, not elliptic_curves, since the1197functions that actually need the flag are in number fields.11981199OUTPUT:12001201The Kodaira Symbol of the curve at P, represented as a string.12021203EXAMPLES::12041205sage: K.<a>=NumberField(x^2-5)1206sage: E=EllipticCurve([20, 225, 750, 625*a + 6875, 31250*a + 46875])1207sage: bad_primes = E.discriminant().support(); bad_primes1208[Fractional ideal (-a), Fractional ideal (7/2*a - 81/2), Fractional ideal (-a - 52), Fractional ideal (2)]1209sage: [E.kodaira_symbol(P) for P in bad_primes]1210[I0, I1, I1, II]1211sage: K.<a> = QuadraticField(-11)1212sage: E = EllipticCurve('11a1').change_ring(K)1213sage: [E.kodaira_symbol(P) for P in K(11).support()]1214[I10]1215"""1216if proof is None:1217import sage.structure.proof.proof1218# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.1219proof = sage.structure.proof.proof.get_flag(None, "number_field")12201221return self.local_data(P, proof).kodaira_symbol()122212231224def conductor(self):1225r"""1226Returns the conductor of this elliptic curve as a fractional1227ideal of the base field.12281229OUTPUT:12301231(fractional ideal) The conductor of the curve.12321233EXAMPLES::12341235sage: K.<i>=NumberField(x^2+1)1236sage: EllipticCurve([i, i - 1, i + 1, 24*i + 15, 14*i + 35]).conductor()1237Fractional ideal (21*i - 3)1238sage: K.<a>=NumberField(x^2-x+3)1239sage: EllipticCurve([1 + a , -1 + a , 1 + a , -11 + a , 5 -9*a ]).conductor()1240Fractional ideal (6*a)12411242A not so well known curve with everywhere good reduction::12431244sage: K.<a>=NumberField(x^2-38)1245sage: E=EllipticCurve([0,0,0, 21796814856932765568243810*a - 134364590724198567128296995, 121774567239345229314269094644186997594*a - 750668847495706904791115375024037711300])1246sage: E.conductor()1247Fractional ideal (1)12481249An example which used to fail (see trac #5307)::12501251sage: K.<w>=NumberField(x^2+x+6)1252sage: E=EllipticCurve([w,-1,0,-w-6,0])1253sage: E.conductor()1254Fractional ideal (86304, w + 5898)12551256An example raised in \#11346::12571258sage: K.<g> = NumberField(x^2 - x - 1)1259sage: E1 = EllipticCurve(K,[0,0,0,-1/48,-161/864])1260sage: [(p.smallest_integer(),e) for p,e in E1.conductor().factor()]1261[(2, 4), (3, 1), (5, 1)]1262"""1263try:1264return self._conductor1265except AttributeError:1266pass12671268# Note: for number fields other than QQ we could initialize1269# N=K.ideal(1) or N=OK.ideal(1), which are the same, but for1270# K==QQ it has to be ZZ.ideal(1).1271OK = self.base_ring().ring_of_integers()1272self._conductor = prod([d.prime()**(d.conductor_valuation()) \1273for d in self.local_data()],\1274OK.ideal(1))1275return self._conductor12761277def global_minimal_model(self, proof = None):1278r"""1279Returns a model of self that is integral, minimal at all primes.12801281.. note::12821283This is only implemented for class number 1. In general,1284such a model may or may not exist.12851286INPUT:12871288- ``proof`` -- whether to only use provably correct methods1289(default controlled by global proof module). Note that the1290proof module is number_field, not elliptic_curves, since the1291functions that actually need the flag are in number fields.12921293OUTPUT:12941295A global integral and minimal model.12961297EXAMPLES::12981299sage: K.<a> = NumberField(x^2-38)1300sage: E = EllipticCurve([0,0,0, 21796814856932765568243810*a - 134364590724198567128296995, 121774567239345229314269094644186997594*a - 750668847495706904791115375024037711300])13011302sage: E2 = E.global_minimal_model()1303sage: E2 # random (the global minimal model is not unique)1304Elliptic Curve defined by y^2 + a*x*y + (a+1)*y = x^3 + (a+1)*x^2 + (368258520200522046806318444*a-2270097978636731786720859345)*x + (8456608930173478039472018047583706316424*a-52130038506793883217874390501829588391299) over Number Field in a with defining polynomial x^2 - 3813051306sage: E2.local_data()1307[]13081309See trac \#11347::13101311sage: K.<g> = NumberField(x^2 - x - 1)1312sage: E = EllipticCurve(K,[0,0,0,-1/48,161/864]).integral_model().global_minimal_model(); E1313Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 over Number Field in g with defining polynomial x^2 - x - 11314sage: [(p.norm(), e) for p, e in E.conductor().factor()]1315[(9, 1), (5, 1)]1316sage: [(p.norm(), e) for p, e in E.discriminant().factor()]1317[(9, 1), (-5, 2)]1318"""1319if proof is None:1320import sage.structure.proof.proof1321# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.1322proof = sage.structure.proof.proof.get_flag(None, "number_field")1323K = self.base_ring()1324if K.class_number() != 1:1325raise ValueError, "global minimal models only exist in general for class number 1"13261327E = self.global_integral_model()1328primes = E.base_ring()(E.discriminant()).support()1329for P in primes:1330E = E.local_data(P,proof).minimal_model()1331return E._reduce_model()13321333def reduction(self,place):1334r"""1335Return the reduction of the elliptic curve at a place of good reduction.13361337INPUT:13381339- ``place`` -- a prime ideal in the base field of the curve13401341OUTPUT:13421343An elliptic curve over a finite field, the residue field of the place.13441345EXAMPLES::13461347sage: K.<i> = QuadraticField(-1)1348sage: EK = EllipticCurve([0,0,0,i,i+3])1349sage: v = K.fractional_ideal(2*i+3)1350sage: EK.reduction(v)1351Elliptic Curve defined by y^2 = x^3 + 5*x + 8 over Residue field of Fractional ideal (2*i + 3)1352sage: EK.reduction(K.ideal(1+i))1353Traceback (most recent call last):1354...1355ValueError: The curve must have good reduction at the place.1356sage: EK.reduction(K.ideal(2))1357Traceback (most recent call last):1358...1359ValueError: The ideal must be prime.1360sage: K=QQ.extension(x^2+x+1,"a")1361sage: E=EllipticCurve([1024*K.0,1024*K.0])1362sage: E.reduction(2*K)1363Elliptic Curve defined by y^2 + (abar+1)*y = x^3 over Residue field in abar of Fractional ideal (2)1364"""1365K = self.base_field()1366OK = K.ring_of_integers()1367try:1368place = K.ideal(place)1369except TypeError:1370raise TypeError, "The parameter must be an ideal of the base field of the elliptic curve"1371if not place.is_prime():1372raise ValueError, "The ideal must be prime."1373disc = self.discriminant()1374if not K.ideal(disc).valuation(place) == 0:1375local_data=self.local_data(place)1376if local_data.has_good_reduction():1377Fv = OK.residue_field(place)1378return local_data.minimal_model().change_ring(Fv)1379raise ValueError, "The curve must have good reduction at the place."1380Fv = OK.residue_field(place)1381return self.change_ring(Fv)13821383def _torsion_bound(self,number_of_places = 20):1384r"""1385An upper bound on the order of the torsion subgroup.13861387INPUT:13881389- ``number_of_places`` (positive integer, default = 20) -- the1390number of places that will be used to find the bound.13911392OUTPUT:13931394(integer) An upper bound on the torsion order.13951396ALGORITHM:13971398An upper bound on the order of the torsion.group of the1399elliptic curve is obtained by counting points modulo several1400primes of good reduction. Note that the upper bound returned1401by this function is a multiple of the order of the torsion1402group, and in general will be greater than the order.14031404EXAMPLES::14051406sage: CDB=CremonaDatabase()1407sage: [E._torsion_bound() for E in CDB.iter([14])]1408[6, 6, 6, 6, 6, 6]1409sage: [E.torsion_order() for E in CDB.iter([14])]1410[6, 6, 2, 6, 2, 6]1411"""1412E = self1413bound = 01414k = 01415K = E.base_field()1416OK = K.ring_of_integers()1417disc = E.discriminant()1418p = Integer(1)1419# runs through primes, decomposes them into prime ideals1420while k < number_of_places :1421p = p.next_prime()1422f = K.primes_above(p)1423# runs through prime ideals above p1424for qq in f:1425fqq = qq.residue_class_degree()1426charqq = qq.smallest_integer()1427# take only places with small residue field (so that the1428# number of points will be small)1429if fqq == 1 or charqq**fqq < 3*number_of_places:1430# check if the model is integral at the place1431if min([K.ideal(a).valuation(qq) for a in E.a_invariants()]) >= 0:1432eqq = qq.ramification_index()1433# check if the formal group at the place is torsion-free1434# if so the torsion injects into the reduction1435if eqq < charqq - 1 and disc.valuation(qq) == 0:1436Etilda = E.reduction(qq)1437Npp = Etilda.cardinality()1438bound = gcd(bound,Npp)1439if bound == 1:1440return bound1441k += 11442return bound14431444def torsion_subgroup(self):1445r"""1446Returns the torsion subgroup of this elliptic curve.14471448OUTPUT:14491450(``EllipticCurveTorsionSubgroup``) The1451``EllipticCurveTorsionSubgroup`` associated to this elliptic1452curve.14531454EXAMPLES::14551456sage: E = EllipticCurve('11a1')1457sage: K.<t>=NumberField(x^4 + x^3 + 11*x^2 + 41*x + 101)1458sage: EK=E.base_extend(K)1459sage: tor = EK.torsion_subgroup()1460sage: tor1461Torsion Subgroup isomorphic to Z/5 + Z/5 associated to the Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 + (-10)*x + (-20) over Number Field in t with defining polynomial x^4 + x^3 + 11*x^2 + 41*x + 1011462sage: tor.gens()1463((16 : 60 : 1), (t : 1/11*t^3 + 6/11*t^2 + 19/11*t + 48/11 : 1))14641465::14661467sage: E = EllipticCurve('15a1')1468sage: K.<t>=NumberField(x^2 + 2*x + 10)1469sage: EK=E.base_extend(K)1470sage: EK.torsion_subgroup()1471Torsion Subgroup isomorphic to Z/4 + Z/4 associated to the Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + (-10)*x + (-10) over Number Field in t with defining polynomial x^2 + 2*x + 1014721473::14741475sage: E = EllipticCurve('19a1')1476sage: K.<t>=NumberField(x^9-3*x^8-4*x^7+16*x^6-3*x^5-21*x^4+5*x^3+7*x^2-7*x+1)1477sage: EK=E.base_extend(K)1478sage: EK.torsion_subgroup()1479Torsion Subgroup isomorphic to Z/9 associated to the Elliptic Curve defined by y^2 + y = x^3 + x^2 + (-9)*x + (-15) over Number Field in t with defining polynomial x^9 - 3*x^8 - 4*x^7 + 16*x^6 - 3*x^5 - 21*x^4 + 5*x^3 + 7*x^2 - 7*x + 114801481::14821483sage: K.<i> = QuadraticField(-1)1484sage: EK = EllipticCurve([0,0,0,i,i+3])1485sage: EK.torsion_subgroup ()1486Torsion Subgroup isomorphic to Trivial group associated to the Elliptic Curve defined by y^2 = x^3 + i*x + (i+3) over Number Field in i with defining polynomial x^2 + 11487"""1488try:1489return self.__torsion_subgroup1490except AttributeError:1491self.__torsion_subgroup = ell_torsion.EllipticCurveTorsionSubgroup(self)1492return self.__torsion_subgroup14931494def torsion_order(self):1495r"""1496Returns the order of the torsion subgroup of this elliptic curve.14971498OUTPUT:14991500(integer) the order of the torsion subgroup of this elliptic curve.15011502EXAMPLES::15031504sage: E = EllipticCurve('11a1')1505sage: K.<t> = NumberField(x^4 + x^3 + 11*x^2 + 41*x + 101)1506sage: EK = E.base_extend(K)1507sage: EK.torsion_order()15082515091510::15111512sage: E = EllipticCurve('15a1')1513sage: K.<t> = NumberField(x^2 + 2*x + 10)1514sage: EK = E.base_extend(K)1515sage: EK.torsion_order()15161615171518::15191520sage: E = EllipticCurve('19a1')1521sage: K.<t> = NumberField(x^9-3*x^8-4*x^7+16*x^6-3*x^5-21*x^4+5*x^3+7*x^2-7*x+1)1522sage: EK = E.base_extend(K)1523sage: EK.torsion_order()1524915251526::15271528sage: K.<i> = QuadraticField(-1)1529sage: EK = EllipticCurve([0,0,0,i,i+3])1530sage: EK.torsion_order()153111532"""1533try:1534return self.__torsion_order1535except AttributeError:1536self.__torsion_order = self.torsion_subgroup().order()1537return self.__torsion_order15381539def torsion_points(self):1540r"""1541Returns a list of the torsion points of this elliptic curve.15421543OUTPUT:15441545(list) A sorted list of the torsion points.15461547EXAMPLES::15481549sage: E = EllipticCurve('11a1')1550sage: E.torsion_points()1551[(0 : 1 : 0), (5 : -6 : 1), (5 : 5 : 1), (16 : -61 : 1), (16 : 60 : 1)]1552sage: K.<t> = NumberField(x^4 + x^3 + 11*x^2 + 41*x + 101)1553sage: EK = E.base_extend(K)1554sage: EK.torsion_points()1555[(16 : 60 : 1),1556(5 : 5 : 1),1557(5 : -6 : 1),1558(16 : -61 : 1),1559(t : 1/11*t^3 + 6/11*t^2 + 19/11*t + 48/11 : 1),1560(-3/55*t^3 - 7/55*t^2 - 2/55*t - 133/55 : 6/55*t^3 + 3/55*t^2 + 25/11*t + 156/55 : 1),1561(-9/121*t^3 - 21/121*t^2 - 127/121*t - 377/121 : -7/121*t^3 + 24/121*t^2 + 197/121*t + 16/121 : 1),1562(5/121*t^3 - 14/121*t^2 - 158/121*t - 453/121 : -49/121*t^3 - 129/121*t^2 - 315/121*t - 207/121 : 1),1563(10/121*t^3 + 49/121*t^2 + 168/121*t + 73/121 : 32/121*t^3 + 60/121*t^2 - 261/121*t - 807/121 : 1),1564(1/11*t^3 - 5/11*t^2 + 19/11*t - 40/11 : -6/11*t^3 - 3/11*t^2 - 26/11*t - 321/11 : 1),1565(14/121*t^3 - 15/121*t^2 + 90/121*t + 232/121 : 16/121*t^3 - 69/121*t^2 + 293/121*t - 46/121 : 1),1566(3/55*t^3 + 7/55*t^2 + 2/55*t + 78/55 : 7/55*t^3 - 24/55*t^2 + 9/11*t + 17/55 : 1),1567(-5/121*t^3 + 36/121*t^2 - 84/121*t + 24/121 : 34/121*t^3 - 27/121*t^2 + 305/121*t + 708/121 : 1),1568(-26/121*t^3 + 20/121*t^2 - 219/121*t - 995/121 : 15/121*t^3 + 156/121*t^2 - 232/121*t + 2766/121 : 1),1569(1/11*t^3 - 5/11*t^2 + 19/11*t - 40/11 : 6/11*t^3 + 3/11*t^2 + 26/11*t + 310/11 : 1),1570(-26/121*t^3 + 20/121*t^2 - 219/121*t - 995/121 : -15/121*t^3 - 156/121*t^2 + 232/121*t - 2887/121 : 1),1571(-5/121*t^3 + 36/121*t^2 - 84/121*t + 24/121 : -34/121*t^3 + 27/121*t^2 - 305/121*t - 829/121 : 1),1572(3/55*t^3 + 7/55*t^2 + 2/55*t + 78/55 : -7/55*t^3 + 24/55*t^2 - 9/11*t - 72/55 : 1),1573(14/121*t^3 - 15/121*t^2 + 90/121*t + 232/121 : -16/121*t^3 + 69/121*t^2 - 293/121*t - 75/121 : 1),1574(t : -1/11*t^3 - 6/11*t^2 - 19/11*t - 59/11 : 1),1575(10/121*t^3 + 49/121*t^2 + 168/121*t + 73/121 : -32/121*t^3 - 60/121*t^2 + 261/121*t + 686/121 : 1),1576(5/121*t^3 - 14/121*t^2 - 158/121*t - 453/121 : 49/121*t^3 + 129/121*t^2 + 315/121*t + 86/121 : 1),1577(-9/121*t^3 - 21/121*t^2 - 127/121*t - 377/121 : 7/121*t^3 - 24/121*t^2 - 197/121*t - 137/121 : 1),1578(-3/55*t^3 - 7/55*t^2 - 2/55*t - 133/55 : -6/55*t^3 - 3/55*t^2 - 25/11*t - 211/55 : 1),1579(0 : 1 : 0)]15801581::15821583sage: E = EllipticCurve('15a1')1584sage: K.<t> = NumberField(x^2 + 2*x + 10)1585sage: EK = E.base_extend(K)1586sage: EK.torsion_points()1587[(8 : 18 : 1),1588(3 : -2 : 1),1589(8 : -27 : 1),1590(t : t - 5 : 1),1591(1/2 : 5/4*t + 1/2 : 1),1592(-t - 2 : 2*t + 8 : 1),1593(-7 : -5*t - 2 : 1),1594(-1 : 0 : 1),1595(-2 : 3 : 1),1596(-13/4 : 9/8 : 1),1597(-2 : -2 : 1),1598(t : -2*t + 4 : 1),1599(-7 : 5*t + 8 : 1),1600(-t - 2 : -t - 7 : 1),1601(1/2 : -5/4*t - 2 : 1),1602(0 : 1 : 0)]16031604::16051606sage: K.<i> = QuadraticField(-1)1607sage: EK = EllipticCurve(K,[0,0,0,0,-1])1608sage: EK.torsion_points ()1609[(-2 : -3*i : 1), (0 : -i : 1), (1 : 0 : 1), (0 : i : 1), (-2 : 3*i : 1), (0 : 1 : 0)]1610"""1611T = self.torsion_subgroup() # make sure it is cached1612return sorted(T.points()) # these are also cached in T16131614def rank_bounds(self,verbose=0, lim1=5, lim3=50, limtriv=10, maxprob=20, limbigprime=30):1615r"""1616Returns the lower and upper bounds using :meth:`~simon_two_descent`.1617The results of :meth:`~simon_two_descent` are cached.16181619.. NOTE::16201621The optional parameters control the Simon two descent algorithm;1622see the documentation of :meth:`~simon_two_descent` for more1623details.16241625INPUT:16261627- ``verbose`` -- 0, 1, 2, or 3 (default: 0), the verbosity level16281629- ``lim1`` -- (default: 5) limit on trivial points on quartics16301631- ``lim3`` -- (default: 50) limit on points on ELS quartics16321633- ``limtriv`` -- (default: 10) limit on trivial points on elliptic curve16341635- ``maxprob`` -- (default: 20)16361637- ``limbigprime`` -- (default: 30) to distinguish between1638small and large prime numbers. Use probabilistic tests for1639large primes. If 0, don't use probabilistic tests.16401641OUTPUT:16421643lower and upper bounds164416451646.. NOTE::16471648For non-quadratic number fields, this code does1649return, but it takes a long time.16501651EXAMPLES::16521653sage: K.<a> = NumberField(x^2 + 23, 'a')1654sage: E = EllipticCurve(K, '37')1655sage: E == loads(dumps(E))1656True1657sage: E.rank_bounds()1658(2, 2)16591660Here is a curve with two-torsion, so here the algorithm gives1661bounds on the rank::16621663sage: Qrt5.<rt5>=NumberField(x^2-5)1664sage: E=EllipticCurve([0,5-rt5,0,rt5,0])1665sage: E.rank_bounds()1666(1, 2)16671668IMPLEMENTATION:16691670Uses Denis Simon's PARI/GP scripts from1671http://www.math.unicaen.fr/~simon/.16721673"""16741675lower, upper, gens = self.simon_two_descent(verbose=verbose,lim1=lim1,lim3=lim3,limtriv=limtriv,maxprob=maxprob,limbigprime=limbigprime)1676return lower,upper16771678def rank(self,verbose=0, lim1=5, lim3=50, limtriv=10, maxprob=20, limbigprime=30):1679r"""1680Return the rank of this elliptic curve, if it can be determined.16811682.. NOTE::16831684The optional parameters control the Simon two descent algorithm;1685see the documentation of :meth:`~simon_two_descent` for more1686details.16871688INPUT:16891690- ``verbose`` -- 0, 1, 2, or 3 (default: 0), the verbosity level16911692- ``lim1`` -- (default: 5) limit on trivial points on quartics16931694- ``lim3`` -- (default: 50) limit on points on ELS quartics16951696- ``limtriv`` -- (default: 10) limit on trivial points on elliptic curve16971698- ``maxprob`` -- (default: 20)16991700- ``limbigprime`` -- (default: 30) to distinguish between1701small and large prime numbers. Use probabilistic tests for1702large primes. If 0, don't use probabilistic tests.17031704OUTPUT:17051706If the upper and lower bounds given by Simon two-descent are1707the same, then the rank has been uniquely identified and we1708return this. Otherwise, we raise a ValueError with an error1709message specifying the upper and lower bounds.17101711.. NOTE::17121713For non-quadratic number fields, this code does return, but it takes1714a long time.17151716EXAMPLES::17171718sage: K.<a> = NumberField(x^2 + 23, 'a')1719sage: E = EllipticCurve(K, '37')1720sage: E == loads(dumps(E))1721True1722sage: E.rank()1723217241725Here is a curve with two-torsion, so here the bounds given by the1726algorithm do not uniquely determine the rank::17271728sage: Qrt5.<rt5>=NumberField(x^2-5)1729sage: E=EllipticCurve([0,5-rt5,0,rt5,0])1730sage: E.rank()1731Traceback (most recent call last):1732...1733ValueError: There is insufficient data to determine the rank -17342-descent gave lower bound 1 and upper bound 217351736IMPLEMENTATION:17371738Uses Denis Simon's PARI/GP scripts from1739http://www.math.unicaen.fr/~simon/.17401741"""17421743lower,upper = self.rank_bounds(verbose=verbose,lim1=lim1,lim3=lim3,limtriv=limtriv,maxprob=maxprob,limbigprime=limbigprime)1744if lower == upper:1745return lower1746else:1747raise ValueError, 'There is insufficient data to determine the rank - 2-descent gave lower bound %s and upper bound %s' % (lower, upper)17481749def gens(self,verbose=0, lim1=5, lim3=50, limtriv=10, maxprob=20, limbigprime=30):1750r"""1751Returns some generators of this elliptic curve. Check :meth:`~rank` or1752:meth:`~rank_bounds` to verify the number of generators.17531754.. NOTE::17551756The optional parameters control the Simon two descent algorithm;1757see the documentation of :meth:`~simon_two_descent` for more1758details.17591760INPUT:17611762- ``verbose`` -- 0, 1, 2, or 3 (default: 0), the verbosity level17631764- ``lim1`` -- (default: 5) limit on trivial points on quartics17651766- ``lim3`` -- (default: 50) limit on points on ELS quartics17671768- ``limtriv`` -- (default: 10) limit on trivial points on elliptic curve17691770- ``maxprob`` -- (default: 20)17711772- ``limbigprime`` -- (default: 30) to distinguish between1773small and large prime numbers. Use probabilistic tests for1774large primes. If 0, don't use probabilistic tests.17751776OUTPUT:17771778The linearly independent elements given by the Simon two-descent.17791780.. NOTE::17811782For non-quadratic number fields, this code does return, but it takes1783a long time.17841785EXAMPLES::17861787sage: K.<a> = NumberField(x^2 + 23, 'a')1788sage: E = EllipticCurve(K, '37')1789sage: E == loads(dumps(E))1790True1791sage: E.gens()1792[(-1 : 0 : 1), (1/2*a - 5/2 : -1/2*a - 13/2 : 1)]179317941795Here is a curve with two-torsion, so here the algorithm does not1796uniquely determine the rank::17971798sage: Qrt5.<rt5>=NumberField(x^2-5)1799sage: E=EllipticCurve([0,5-rt5,0,rt5,0])1800sage: E.gens()1801[(3/2*rt5 + 5/2 : -9/2*rt5 - 15/2 : 1), (-1/2*rt5 + 3/2 : 3/2*rt5 - 9/2 : 1), (0 : 0 : 1)]18021803IMPLEMENTATION:18041805Uses Denis Simon's PARI/GP scripts from1806http://www.math.unicaen.fr/~simon/.1807"""18081809lower,upper,gens = self.simon_two_descent(verbose=verbose,lim1=lim1,lim3=lim3,limtriv=limtriv,maxprob=maxprob,limbigprime=limbigprime)1810return gens181118121813def period_lattice(self, embedding):1814r"""1815Returns the period lattice of the elliptic curve for the given1816embedding of its base field with respect to the differential1817`dx/(2y + a_1x + a_3)`.18181819INPUT:18201821- ``embedding`` - an embedding of the base number field into `\RR` or `\CC`.18221823.. note::18241825The precision of the embedding is ignored: we only use the1826given embedding to determine which embedding into ``QQbar``1827to use. Once the lattice has been initialized, periods can1828be computed to arbitrary precision.182918301831EXAMPLES:18321833First define a field with two real embeddings::18341835sage: K.<a> = NumberField(x^3-2)1836sage: E=EllipticCurve([0,0,0,a,2])1837sage: embs=K.embeddings(CC); len(embs)1838318391840For each embedding we have a different period lattice::18411842sage: E.period_lattice(embs[0])1843Period lattice associated to Elliptic Curve defined by y^2 = x^3 + a*x + 2 over Number Field in a with defining polynomial x^3 - 2 with respect to the embedding Ring morphism:1844From: Number Field in a with defining polynomial x^3 - 21845To: Algebraic Field1846Defn: a |--> -0.6299605249474365? - 1.091123635971722?*I18471848sage: E.period_lattice(embs[1])1849Period lattice associated to Elliptic Curve defined by y^2 = x^3 + a*x + 2 over Number Field in a with defining polynomial x^3 - 2 with respect to the embedding Ring morphism:1850From: Number Field in a with defining polynomial x^3 - 21851To: Algebraic Field1852Defn: a |--> -0.6299605249474365? + 1.091123635971722?*I18531854sage: E.period_lattice(embs[2])1855Period lattice associated to Elliptic Curve defined by y^2 = x^3 + a*x + 2 over Number Field in a with defining polynomial x^3 - 2 with respect to the embedding Ring morphism:1856From: Number Field in a with defining polynomial x^3 - 21857To: Algebraic Field1858Defn: a |--> 1.259921049894873?18591860Although the original embeddings have only the default1861precision, we can obtain the basis with higher precision1862later::18631864sage: L=E.period_lattice(embs[0])1865sage: L.basis()1866(1.86405007647981 - 0.903761485143226*I, -0.149344633143919 - 2.06619546272945*I)18671868sage: L.basis(prec=100)1869(1.8640500764798108425920506200 - 0.90376148514322594749786960975*I, -0.14934463314391922099120107422 - 2.0661954627294548995621225062*I)1870"""1871from sage.schemes.elliptic_curves.period_lattice import PeriodLattice_ell1872return PeriodLattice_ell(self,embedding)18731874def is_isogenous(self, other, proof=True, maxnorm=100):1875"""1876Returns whether or not self is isogenous to other.18771878INPUT:18791880- ``other`` -- another elliptic curve.18811882- ``proof`` (default True) -- If ``False``, the function will1883return ``True`` whenever the two curves have the same1884conductor and are isogenous modulo `p` for all primes `p` of1885norm up to ``maxp``. If ``True``, the function returns1886False when the previous condition does not hold, and if it1887does hold we attempt to see if the curves are indeed1888isogenous. However, this has not been fully implemented1889(see examples below), so we may not be able to determine1890whether or not the curves are isogenous..18911892- ``maxnorm`` (integer, default 100) -- The maximum norm of1893primes `p` for which isogeny modulo `p` will be checked.18941895OUTPUT:18961897(bool) True if there is an isogeny from curve ``self`` to1898curve ``other``.18991900EXAMPLES::19011902sage: x = polygen(QQ, 'x')1903sage: F = NumberField(x^2 -2, 's'); F1904Number Field in s with defining polynomial x^2 - 21905sage: E1 = EllipticCurve(F, [7,8])1906sage: E2 = EllipticCurve(F, [0,5,0,1,0])1907sage: E3 = EllipticCurve(F, [0,-10,0,21,0])1908sage: E1.is_isogenous(E2)1909False1910sage: E1.is_isogenous(E1)1911True1912sage: E2.is_isogenous(E2)1913True1914sage: E2.is_isogenous(E1)1915False1916sage: E2.is_isogenous(E3)1917True19181919::19201921sage: x = polygen(QQ, 'x')1922sage: F = NumberField(x^2 -2, 's'); F1923Number Field in s with defining polynomial x^2 - 21924sage: E = EllipticCurve('14a1')1925sage: EE = EllipticCurve('14a2')1926sage: E1 = E.change_ring(F)1927sage: E2 = EE.change_ring(F)1928sage: E1.is_isogenous(E2)1929True19301931::19321933sage: x = polygen(QQ, 'x')1934sage: F = NumberField(x^2 -2, 's'); F1935Number Field in s with defining polynomial x^2 - 21936sage: k.<a> = NumberField(x^3+7)1937sage: E = EllipticCurve(F, [7,8])1938sage: EE = EllipticCurve(k, [2, 2])1939sage: E.is_isogenous(EE)1940Traceback (most recent call last):1941...1942ValueError: Second argument must be defined over the same number field.19431944Some examples from Cremona's 1981 tables::19451946sage: K.<i> = QuadraticField(-1)1947sage: E1 = EllipticCurve([i + 1, 0, 1, -240*i - 400, -2869*i - 2627])1948sage: E1.conductor()1949Fractional ideal (4*i + 7)1950sage: E2 = EllipticCurve([1+i,0,1,0,0])1951sage: E2.conductor()1952Fractional ideal (4*i + 7)1953sage: E1.is_isogenous(E2)1954Traceback (most recent call last):1955...1956NotImplementedError: Curves appear to be isogenous (same conductor, isogenous modulo all primes of norm up to 1000), but no isogeny has been constructed.1957sage: E1.is_isogenous(E2, proof=False)1958True19591960In this case E1 and E2 are in fact 9-isogenous, as may be1961deduced from the following::19621963sage: E3 = EllipticCurve([i + 1, 0, 1, -5*i - 5, -2*i - 5])1964sage: E3.is_isogenous(E1)1965True1966sage: E3.is_isogenous(E2)1967True1968sage: E1.isogeny_degree(E2)1969919701971"""1972if not is_EllipticCurve(other):1973raise ValueError, "Second argument is not an Elliptic Curve."1974if self.is_isomorphic(other):1975return True1976K = self.base_field()1977if K != other.base_field():1978raise ValueError, "Second argument must be defined over the same number field."19791980E1 = self.integral_model()1981E2 = other.integral_model()1982N = E1.conductor()1983if N != E2.conductor():1984return False19851986PI = K.primes_of_degree_one_iter()1987while True:1988P = PI.next()1989if P.norm() > maxnorm: break1990if not P.divides(N):1991OP = K.residue_field(P)1992if E1.change_ring(OP).cardinality() != E2.change_ring(OP).cardinality():1993return False19941995if not proof:1996return True19971998# We have not yet implemented isogenies of all possible1999# degrees, and do not know a bound on the possible degrees2000# over general number fields. But here we do at least try2001# some easy cases:20022003for l in [2,3,5,7,13]:2004if any([E2.is_isomorphic(f.codomain()) for f in E1.isogenies_prime_degree(l)]):2005return True20062007# Next we try looking modulo some more primes:20082009while True:2010if P.norm() > 10*maxnorm: break2011if not P.divides(N):2012OP = K.residue_field(P)2013if E1.change_ring(OP).cardinality() != E2.change_ring(OP).cardinality():2014return False2015P = PI.next()20162017# At this point is is highly likely that the curves are2018# isogenous, but we have not proved it.20192020raise NotImplementedError, "Curves appear to be isogenous (same conductor, isogenous modulo all primes of norm up to %s), but no isogeny has been constructed." % (10*maxnorm)20212022def isogeny_degree(self, other):2023"""2024Returns the minimal degree of an isogeny between self and2025other, or 0 if no isogeny exists.20262027INPUT:20282029- ``other`` -- another elliptic curve.20302031OUTPUT:20322033(int) The degree of an isogeny from ``self`` to ``other``, or 0.20342035.. warning::20362037Not all isogenies over number fields are yet implemented.2038Currently the code only works if there is a chain of2039isogenies from ``self`` to ``other`` of degrees 2, 3, 5, 72040and 13.20412042EXAMPLES::20432044sage: x = QQ['x'].02045sage: F = NumberField(x^2 -2, 's'); F2046Number Field in s with defining polynomial x^2 - 22047sage: E = EllipticCurve('14a1')2048sage: EE = EllipticCurve('14a2')2049sage: E1 = E.change_ring(F)2050sage: E2 = EE.change_ring(F)2051sage: E1.isogeny_degree(E2)205222053sage: E2.isogeny_degree(E2)205412055sage: E5 = EllipticCurve('14a5').change_ring(F)2056sage: E1.isogeny_degree(E5)205762058"""2059if self.conductor() != other.conductor():2060return Integer(0)20612062if self.is_isomorphic(other):2063return Integer(1)20642065from sage.sets.set import Set20662067curves = [self]2068degrees = [Integer(1)]2069l_list = [l for l in Set([ZZ(f.degree()) for f in self.isogenies_prime_degree([2,3,5,7,13])])]20702071newcurves = []2072newdegs = []2073k = 02074while k<len(curves):2075newcurves.extend([f.codomain() for f in curves[k].isogenies_prime_degree(l_list)])2076newdegs.extend([degrees[k]*f.degree() for f in curves[k].isogenies_prime_degree(l_list)])2077newisogpairs = dict(zip(newcurves, newdegs))2078i = 02079while i<len(curves):2080j = 02081while j<len(newcurves):2082if curves[i].is_isomorphic(newcurves[j]):2083newdegs.remove(newisogpairs[newcurves[j]])2084newcurves.remove(newcurves[j])2085else:2086j = j+12087i = i+120882089m = 02090newisogpairs = dict(zip(newcurves, newdegs))2091while m<len(newcurves):2092if other.is_isomorphic(newcurves[m]):2093return newisogpairs[newcurves[m]]2094m = m+120952096curves.extend(newcurves)2097degrees.extend(newdegs)2098k = k+120992100raise NotImplementedError, "Not all isogenies implemented over general number fields."210121022103