Path: blob/master/src/sage/schemes/elliptic_curves/ell_number_field.py
8820 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 ell_torsion99from ell_generic import is_EllipticCurve100101from gp_simon import simon_two_descent102from constructor import EllipticCurve103from sage.rings.all import PolynomialRing, ZZ, RealField104import sage.misc.misc105from sage.misc.misc import verbose, forall106from sage.rings.integer import Integer107from sage.rings.arith import valuation108109import gal_reps_number_field110111class EllipticCurve_number_field(EllipticCurve_field):112r"""113Elliptic curve over a number field.114115EXAMPLES::116117sage: K.<i>=NumberField(x^2+1)118sage: EllipticCurve([i, i - 1, i + 1, 24*i + 15, 14*i + 35])119Elliptic 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 + 1120"""121def __init__(self, x, y=None):122r"""123Allow some ways to create an elliptic curve over a number124field in addition to the generic ones.125126INPUT:127128- ``x``, ``y`` -- see examples.129130EXAMPLES:131132A curve from the database of curves over `\QQ`, but over a larger field:133134sage: K.<i>=NumberField(x^2+1)135sage: EllipticCurve(K,'389a1')136Elliptic Curve defined by y^2 + y = x^3 + x^2 + (-2)*x over Number Field in i with defining polynomial x^2 + 1137138Making the field of definition explicitly larger::139140sage: EllipticCurve(K,[0,-1,1,0,0])141Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 over Number Field in i with defining polynomial x^2 + 1142143"""144if y is None:145if isinstance(x, list):146ainvs = x147field = ainvs[0].parent()148else:149if isinstance(y, str):150from sage.databases.cremona import CremonaDatabase151field = x152X = CremonaDatabase()[y]153ainvs = list(X.a_invariants())154else:155field = x156ainvs = y157if not (isinstance(field, Ring) and isinstance(ainvs,list)):158raise TypeError159160EllipticCurve_field.__init__(self, [field(x) for x in ainvs])161self._point = ell_point.EllipticCurvePoint_number_field162163def simon_two_descent(self, verbose=0, lim1=5, lim3=50, limtriv=10, maxprob=20, limbigprime=30):164r"""165Computes lower and upper bounds on the rank of the Mordell-Weil group,166and a list of independent points. Used internally by the :meth:`~rank`,167:meth:`~rank_bounds` and :meth:`~gens` methods.168169INPUT:170171- ``verbose`` -- 0, 1, 2, or 3 (default: 0), the verbosity level172173- ``lim1`` -- (default: 5) limit on trivial points on quartics174175- ``lim3`` -- (default: 50) limit on points on ELS quartics176177- ``limtriv`` -- (default: 10) limit on trivial points on elliptic curve178179- ``maxprob`` -- (default: 20)180181- ``limbigprime`` -- (default: 30) to distinguish between182small and large prime numbers. Use probabilistic tests for183large primes. If 0, don't use probabilistic tests.184185OUTPUT:186187``(lower, upper, list)`` where ``lower`` is a lower bound on188the rank, ``upper`` is an upper bound (the dimension of the1892-Selmer group) and190``list`` is a list of points of infinite order on the Weierstrass191model.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)])211sage: E.simon_two_descent(lim1=3, lim3=20, limtriv=5, maxprob=7, limbigprime=10)212(2, 2, [(-1 : 0 : 1), (-1/8*a + 5/8 : -3/16*a - 9/16 : 1)])213214::215216sage: K.<a> = NumberField(x^2 + 7, 'a')217sage: E = EllipticCurve(K, [0,0,0,1,a]); E218Elliptic Curve defined by y^2 = x^3 + x + a over Number Field in a with defining polynomial x^2 + 7219220sage: v = E.simon_two_descent(verbose=1); v221courbe elliptique : Y^2 = x^3 + Mod(1, y^2 + 7)*x + Mod(y, y^2 + 7)222points triviaux sur la courbe = [[1, 1, 0], [Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7), 1]]223#S(E/K)[2] = 2224#E(K)/2E(K) = 2225#III(E/K)[2] = 1226rang(E/K) = 1227listpointsmwr = [[Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7), 1]]228(1, 1, [(1/2*a + 3/2 : -a - 2 : 1)])229230sage: v = E.simon_two_descent(verbose=2) # random output231K = bnfinit(y^2 + 7);232a = Mod(y,K.pol);233bnfellrank(K, [0,0,0,1,a]);234courbe elliptique : Y^2 = x^3 + Mod(1, y^2 + 7)*x + Mod(y, y^2 + 7)235A = 0236B = Mod(1, y^2 + 7)237C = Mod(y, y^2 + 7)238LS2gen = [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))]239#LS2gen = 2240Recherche de points triviaux sur la courbe241points triviaux sur la courbe = [[1, 1, 0], [Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7), 1]]242zc = 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))243symbole de Hilbert (Mod(2, y^2 + 7),Mod(-5, y^2 + 7)) = -1244zc = 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))245symbole de Hilbert (Mod(-2*y + 2, y^2 + 7),Mod(1, y^2 + 7)) = 0246sol de Legendre = [1, 0, 1]~247zc*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))248quartique : (-1/2*y + 1/2)*Y^2 = x^4 + (-3*y - 15)*x^2 + (-8*y - 16)*x + (-11/2*y - 15/2)249reduite: 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)250non ELS en [2, [0, 1]~, 1, 1, [1, 1]~]251zc = 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))252vient du point trivial [Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7), 1]253m1 = 1254m2 = 1255#S(E/K)[2] = 2256#E(K)/2E(K) = 2257#III(E/K)[2] = 1258rang(E/K) = 1259listpointsmwr = [[Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7), 1]]260v = [1, 1, [[Mod(1/2*y + 3/2, y^2 + 7), Mod(-y - 2, y^2 + 7)]]]261sage: v262(1, 1, [(1/2*a + 3/2 : -a - 2 : 1)])263264265A curve with 2-torsion::266267sage: K.<a> = NumberField(x^2 + 7)268sage: E = EllipticCurve(K, '15a')269sage: E.simon_two_descent() # long time (3s on sage.math, 2013), points can vary270(1, 3, [...])271272A failure in the PARI/GP script ell.gp (VERSION 25/03/2009) is reported::273274sage: K = CyclotomicField(43).subfields(3)[0][0]275sage: E = EllipticCurve(K, '37')276sage: E.simon_two_descent() # long time (4s on sage.math, 2013)277Traceback (most recent call last):278...279RuntimeError:280*** at top-level: ans=bnfellrank(K,[0,0,1,281*** ^--------------------282*** in function bnfellrank: ...eqtheta,rnfeq,bbnf];rang=283*** bnfell2descent_gen(b284*** ^--------------------285*** in function bnfell2descent_gen: ...riv,r=nfsqrt(nf,norm(zc))286*** [1];if(DEBUGLEVEL_el287*** ^--------------------288*** array index (1) out of allowed range [none].289An error occurred while running Simon's 2-descent program290291"""292293try:294result = self._simon_two_descent_data[lim1,lim3,limtriv,maxprob,limbigprime]295if verbose == 0:296return result297except AttributeError:298self._simon_two_descent_data = {}299except KeyError:300pass301302t = simon_two_descent(self,303verbose=verbose, lim1=lim1, lim3=lim3, limtriv=limtriv,304maxprob=maxprob, limbigprime=limbigprime)305prob_rank = Integer(t[0])306two_selmer_rank = Integer(t[1])307prob_gens = [self(P) for P in t[2]]308self._simon_two_descent_data[lim1,lim3,limtriv,maxprob,limbigprime] = (prob_rank, two_selmer_rank, prob_gens)309return prob_rank, two_selmer_rank, prob_gens310311def division_field(self, p, names, map=False, **kwds):312"""313Given an elliptic curve over a number field `F` and a prime number `p`,314construct the field `F(E[p])`.315316INPUT:317318- ``p`` -- a prime number (an element of `\ZZ`)319320- ``names`` -- a variable name for the number field321322- ``map`` -- (default: ``False``) also return an embedding of323the :meth:`base_field` into the resulting field.324325- ``kwds`` -- additional keywords passed to326:func:`sage.rings.number_field.splitting_field.splitting_field`.327328OUTPUT:329330If ``map`` is ``False``, the division field as an absolute number331field. If ``map`` is ``True``, a tuple ``(K, phi)`` where ``phi``332is an embedding of the base field in the division field ``K``.333334.. WARNING::335336This takes a very long time when the degree of the division337field is large (e.g. when `p` is large or when the Galois338representation is surjective). The ``simplify`` flag also339has a big influence on the running time: sometimes340``simplify=False`` is faster, sometimes ``simplify=True``341(the default) is faster.342343EXAMPLES:344345The 2-division field is the same as the splitting field of346the 2-division polynomial (therefore, it has degree 1, 2, 3 or 6)::347348sage: E = EllipticCurve('15a1')349sage: K.<b> = E.division_field(2); K350Number Field in b with defining polynomial x351sage: E = EllipticCurve('14a1')352sage: K.<b> = E.division_field(2); K353Number Field in b with defining polynomial x^2 + 5*x + 92354sage: E = EllipticCurve('196b1')355sage: K.<b> = E.division_field(2); K356Number Field in b with defining polynomial x^3 + x^2 - 114*x - 127357sage: E = EllipticCurve('19a1')358sage: K.<b> = E.division_field(2); K359Number Field in b with defining polynomial x^6 + 10*x^5 + 24*x^4 - 212*x^3 + 1364*x^2 + 24072*x + 104292360361For odd primes `p`, the division field is either the splitting362field of the `p`-division polynomial, or a quadratic extension363of it. ::364365sage: E = EllipticCurve('50a1')366sage: F.<a> = E.division_polynomial(3).splitting_field(simplify_all=True); F367Number Field in a with defining polynomial x^6 - 3*x^5 + 4*x^4 - 3*x^3 - 2*x^2 + 3*x + 3368sage: K.<b> = E.division_field(3, simplify_all=True); K369Number Field in b with defining polynomial x^6 - 3*x^5 + 4*x^4 - 3*x^3 - 2*x^2 + 3*x + 3370371If we take any quadratic twist, the splitting field of the3723-division polynomial remains the same, but the 3-division field373becomes a quadratic extension::374375sage: E = E.quadratic_twist(5) # 50b3376sage: F.<a> = E.division_polynomial(3).splitting_field(simplify_all=True); F377Number Field in a with defining polynomial x^6 - 3*x^5 + 4*x^4 - 3*x^3 - 2*x^2 + 3*x + 3378sage: K.<b> = E.division_field(3, simplify_all=True); K379Number Field in b with defining polynomial x^12 - 3*x^11 + 8*x^10 - 15*x^9 + 30*x^8 - 63*x^7 + 109*x^6 - 144*x^5 + 150*x^4 - 120*x^3 + 68*x^2 - 24*x + 4380381Try another quadratic twist, this time over a subfield of `F`::382383sage: G.<c>,_,_ = F.subfields(3)[0]384sage: E = E.base_extend(G).quadratic_twist(c); E385Elliptic Curve defined by y^2 = x^3 + 5*a0*x^2 + (-200*a0^2)*x + (-42000*a0^2+42000*a0+126000) over Number Field in a0 with defining polynomial x^3 - 3*x^2 + 3*x + 9386sage: K.<b> = E.division_field(3, simplify_all=True); K387Number Field in b with defining polynomial x^12 - 10*x^10 + 55*x^8 - 60*x^6 + 75*x^4 + 1350*x^2 + 2025388389Some higher-degree examples::390391sage: E = EllipticCurve('11a1')392sage: K.<b> = E.division_field(2); K393Number Field in b with defining polynomial x^6 + 2*x^5 - 48*x^4 - 436*x^3 + 1668*x^2 + 28792*x + 73844394sage: K.<b> = E.division_field(3); K # long time (3s on sage.math, 2014)395Number Field in b with defining polynomial x^48 ...396sage: K.<b> = E.division_field(5); K397Number Field in b with defining polynomial x^4 - x^3 + x^2 - x + 1398sage: E.division_field(5, 'b', simplify=False)399Number Field in b with defining polynomial x^4 + x^3 + 11*x^2 + 41*x + 101400sage: E.base_extend(K).torsion_subgroup() # long time (2s on sage.math, 2014)401Torsion 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 b with defining polynomial x^4 - x^3 + x^2 - x + 1402403sage: E = EllipticCurve('27a1')404sage: K.<b> = E.division_field(3); K405Number Field in b with defining polynomial x^2 + 3*x + 9406sage: K.<b> = E.division_field(2); K407Number Field in b with defining polynomial x^6 + 6*x^5 + 24*x^4 - 52*x^3 - 228*x^2 + 744*x + 3844408sage: K.<b> = E.division_field(2, simplify_all=True); K409Number Field in b with defining polynomial x^6 - 3*x^5 + 5*x^3 - 3*x + 1410sage: K.<b> = E.division_field(5); K # long time (4s on sage.math, 2014)411Number Field in b with defining polynomial x^48 ...412sage: K.<b> = E.division_field(7); K # long time (8s on sage.math, 2014)413Number Field in b with defining polynomial x^72 ...414415Over a number field::416417sage: R.<x> = PolynomialRing(QQ)418sage: K.<i> = NumberField(x^2 + 1)419sage: E = EllipticCurve([0,0,0,0,i])420sage: L.<b> = E.division_field(2); L421Number Field in b with defining polynomial x^4 - x^2 + 1422sage: L.<b>, phi = E.division_field(2, map=True); phi423Ring morphism:424From: Number Field in i with defining polynomial x^2 + 1425To: Number Field in b with defining polynomial x^4 - x^2 + 1426Defn: i |--> -b^3427sage: L.<b>, phi = E.division_field(3, map=True)428sage: L429Number Field in b with defining polynomial x^24 - 6*x^22 - 12*x^21 - 21*x^20 + 216*x^19 + 48*x^18 + 804*x^17 + 1194*x^16 - 13488*x^15 + 21222*x^14 + 44196*x^13 - 47977*x^12 - 102888*x^11 + 173424*x^10 - 172308*x^9 + 302046*x^8 + 252864*x^7 - 931182*x^6 + 180300*x^5 + 879567*x^4 - 415896*x^3 + 1941012*x^2 + 650220*x + 443089430sage: phi431Ring morphism:432From: Number Field in i with defining polynomial x^2 + 1433To: Number Field in b with defining polynomial x^24 ...434Defn: i |--> -215621657062634529/183360797284413355040732*b^23 ...435436AUTHORS:437438- Jeroen Demeyer (2014-01-06): :trac:`11905`, use439``splitting_field`` method, moved from ``gal_reps.py``, make440it work over number fields.441"""442p = Integer(p)443if not p.is_prime():444raise ValueError("p must be a prime number")445446verbose("Adjoining X-coordinates of %s-torsion points"%p)447F = self.base_ring()448f = self.division_polynomial(p)449if p == 2:450# For p = 2, the division field is the splitting field of451# the division polynomial.452return f.splitting_field(names, map=map, **kwds)453454# Compute splitting field of X-coordinates.455# The Galois group of the division field is a subgroup of GL(2,p).456# The Galois group of the X-coordinates is a subgroup of GL(2,p)/{-1,+1}.457# We need the map to change the elliptic curve invariants to K.458deg_mult = F.degree()*p*(p+1)*(p-1)*(p-1)//2459K, F_to_K = f.splitting_field(names, degree_multiple=deg_mult, map=True, **kwds)460461verbose("Adjoining Y-coordinates of %s-torsion points"%p)462463# THEOREM (Cremona, http://trac.sagemath.org/ticket/11905#comment:21).464# Let K be a field, E an elliptic curve over K and p an odd465# prime number. Assume that K contains all roots of the466# p-division polynomial of E. Then either K contains all467# p-torsion points on E, or it doesn't contain any p-torsion468# point.469#470# PROOF. Let G be the absolute Galois group of K (every element471# in it fixes all elements of K). For any p-torsion point P472# over the algebraic closure and any sigma in G, we must have473# either sigma(P) = P or sigma(P) = -P (since K contains the474# X-coordinate of P). Now assume that K does not contain all475# p-torsion points. Then there exists a point P1 and a sigma in476# G such that sigma(P1) = -P1. Now take a different p-torsion477# point P2. Since sigma(P2) must be P2 or -P2 and478# sigma(P1+P2) = sigma(P1)+sigma(P2) = sigma(P1)-P2 must479# be P1+P2 or -(P1+P2), it follows that sigma(P2) = -sigma(P2).480# Therefore, K cannot contain any p-torsion point.481#482# This implies that it suffices to adjoin the Y-coordinate483# of just one point.484485# First factor f over F and then compute a root X of f over K.486g = prime_divisors(f)[0]487X = g.map_coefficients(F_to_K).roots(multiplicities=False)[0]488489# Polynomial defining the corresponding Y-coordinate490a1,a2,a3,a4,a6 = (F_to_K(ai) for ai in self.a_invariants())491rhs = X*(X*(X + a2) + a4) + a6492RK = PolynomialRing(K, 'x')493ypol = RK([-rhs, a1*X + a3, 1])494L = ypol.splitting_field(names, map=map, **kwds)495if map:496L, K_to_L = L497return L, F_to_K.post_compose(K_to_L)498else:499return L500501def height_pairing_matrix(self, points=None, precision=None):502r"""503Returns the height pairing matrix of the given points.504505INPUT:506507- points - either a list of points, which must be on this508curve, or (default) None, in which case self.gens() will be509used.510511- precision - number of bits of precision of result512(default: None, for default RealField precision)513514EXAMPLES::515516sage: E = EllipticCurve([0, 0, 1, -1, 0])517sage: E.height_pairing_matrix()518[0.0511114082399688]519520For rank 0 curves, the result is a valid 0x0 matrix::521522sage: EllipticCurve('11a').height_pairing_matrix()523[]524sage: E=EllipticCurve('5077a1')525sage: E.height_pairing_matrix([E.lift_x(x) for x in [-2,-7/4,1]], precision=100)526[ 1.3685725053539301120518194471 -1.3095767070865761992624519454 -0.63486715783715592064475542573]527[ -1.3095767070865761992624519454 2.7173593928122930896610589220 1.0998184305667292139777571432]528[-0.63486715783715592064475542573 1.0998184305667292139777571432 0.66820516565192793503314205089]529530sage: E = EllipticCurve('389a1')531sage: E = EllipticCurve('389a1')532sage: P,Q = E.point([-1,1,1]),E.point([0,-1,1])533sage: E.height_pairing_matrix([P,Q])534[0.686667083305587 0.268478098806726]535[0.268478098806726 0.327000773651605]536537Over a number field::538539sage: x = polygen(QQ)540sage: K.<t> = NumberField(x^2+47)541sage: EK = E.base_extend(K)542sage: EK.height_pairing_matrix([EK(P),EK(Q)])543[0.686667083305587 0.268478098806726]544[0.268478098806726 0.327000773651605]545546::547548sage: K.<i> = QuadraticField(-1)549sage: E = EllipticCurve([0,0,0,i,i])550sage: P = E(-9+4*i,-18-25*i)551sage: Q = E(i,-i)552sage: E.height_pairing_matrix([P,Q])553[ 2.16941934493768 -0.870059380421505]554[-0.870059380421505 0.424585837470709]555sage: E.regulator_of_points([P,Q])5560.164101403936070557"""558if points is None:559points = self.gens()560else:561for P in points:562assert P.curve() == self563564r = len(points)565if precision is None:566RR = RealField()567else:568RR = RealField(precision)569M = matrix.MatrixSpace(RR, r)570mat = M()571for j in range(r):572mat[j,j] = points[j].height(precision=precision)573for j in range(r):574for k in range(j+1,r):575mat[j,k]=((points[j]+points[k]).height(precision=precision) - mat[j,j] - mat[k,k])/2576mat[k,j]=mat[j,k]577return mat578579def regulator_of_points(self, points=[], precision=None):580"""581Returns the regulator of the given points on this curve.582583INPUT:584585- ``points`` -(default: empty list) a list of points on this curve586587- ``precision`` - int or None (default: None): the precision588in bits of the result (default real precision if None)589590EXAMPLES::591592sage: E = EllipticCurve('37a1')593sage: P = E(0,0)594sage: Q = E(1,0)595sage: E.regulator_of_points([P,Q])5960.000000000000000597sage: 2*P==Q598True599600::601602sage: E = EllipticCurve('5077a1')603sage: points = [E.lift_x(x) for x in [-2,-7/4,1]]604sage: E.regulator_of_points(points)6050.417143558758384606sage: E.regulator_of_points(points,precision=100)6070.41714355875838396981711954462608609::610611sage: E = EllipticCurve('389a')612sage: E.regulator_of_points()6131.00000000000000614sage: points = [P,Q] = [E(-1,1),E(0,-1)]615sage: E.regulator_of_points(points)6160.152460177943144617sage: E.regulator_of_points(points, precision=100)6180.15246017794314375162432475705619sage: E.regulator_of_points(points, precision=200)6200.15246017794314375162432475704945582324372707748663081784028621sage: E.regulator_of_points(points, precision=300)6220.152460177943143751624324757049455823243727077486630817840280980046053225683562463604114816623624Examples over number fields::625626sage: K.<a> = QuadraticField(97)627sage: E = EllipticCurve(K,[1,1])628sage: P = E(0,1)629sage: P.height()6300.476223106404866631sage: E.regulator_of_points([P])6320.476223106404866633634::635636sage: E = EllipticCurve('11a1')637sage: x = polygen(QQ)638sage: K.<t> = NumberField(x^2+47)639sage: EK = E.base_extend(K)640sage: T = EK(5,5)641sage: T.order()6425643sage: P = EK(-2, -1/2*t - 1/2)644sage: P.order()645+Infinity646sage: EK.regulator_of_points([P,T]) # random very small output647-1.23259516440783e-32648sage: EK.regulator_of_points([P,T]).abs() < 1e-30649True650651::652653sage: E = EllipticCurve('389a1')654sage: P,Q = E.gens()655sage: E.regulator_of_points([P,Q])6560.152460177943144657sage: K.<t> = NumberField(x^2+47)658sage: EK = E.base_extend(K)659sage: EK.regulator_of_points([EK(P),EK(Q)])6600.152460177943144661662::663664sage: K.<i> = QuadraticField(-1)665sage: E = EllipticCurve([0,0,0,i,i])666sage: P = E(-9+4*i,-18-25*i)667sage: Q = E(i,-i)668sage: E.height_pairing_matrix([P,Q])669[ 2.16941934493768 -0.870059380421505]670[-0.870059380421505 0.424585837470709]671sage: E.regulator_of_points([P,Q])6720.164101403936070673674"""675if points is None:676points = []677mat = self.height_pairing_matrix(points=points, precision=precision)678return mat.det(algorithm="hessenberg")679680681def is_local_integral_model(self,*P):682r"""683Tests if self is integral at the prime ideal `P`, or at all the684primes if `P` is a list or tuple.685686INPUT:687688- ``*P`` -- a prime ideal, or a list or tuple of primes.689690EXAMPLES::691692sage: K.<i> = NumberField(x^2+1)693sage: P1,P2 = K.primes_above(5)694sage: E = EllipticCurve([i/5,i/5,i/5,i/5,i/5])695sage: E.is_local_integral_model(P1,P2)696False697sage: Emin = E.local_integral_model(P1,P2)698sage: Emin.is_local_integral_model(P1,P2)699True700"""701if len(P)==1: P=P[0]702if isinstance(P,(tuple,list)):703return forall(P, lambda x : self.is_local_integral_model(x))[0]704return forall(self.ainvs(), lambda x : x.valuation(P) >= 0)[0]705706def local_integral_model(self,*P):707r"""708Return a model of self which is integral at the prime ideal709`P`.710711.. note::712713The integrality at other primes is not affected, even if714`P` is non-principal.715716INPUT:717718- ``*P`` -- a prime ideal, or a list or tuple of primes.719720EXAMPLES::721722sage: K.<i> = NumberField(x^2+1)723sage: P1,P2 = K.primes_above(5)724sage: E = EllipticCurve([i/5,i/5,i/5,i/5,i/5])725sage: E.local_integral_model((P1,P2))726Elliptic 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 + 1727"""728if len(P)==1: P=P[0]729if isinstance(P,(tuple,list)):730E=self731for Pi in P: E=E.local_integral_model(Pi)732return E733ai = self.a_invariants()734e = min([(ai[i].valuation(P)/[1,2,3,4,6][i]) for i in range(5)]).floor()735pi = self.base_field().uniformizer(P, 'negative')736return EllipticCurve([ai[i]/pi**(e*[1,2,3,4,6][i]) for i in range(5)])737738def is_global_integral_model(self):739r"""740Return true iff self is integral at all primes.741742EXAMPLES::743744sage: K.<i> = NumberField(x^2+1)745sage: E = EllipticCurve([i/5,i/5,i/5,i/5,i/5])746sage: P1,P2 = K.primes_above(5)747sage: Emin = E.global_integral_model()748sage: Emin.is_global_integral_model()749True750"""751return forall(self.a_invariants(), lambda x : x.is_integral())[0]752753def global_integral_model(self):754r"""755Return a model of self which is integral at all primes.756757EXAMPLES::758759sage: K.<i> = NumberField(x^2+1)760sage: E = EllipticCurve([i/5,i/5,i/5,i/5,i/5])761sage: P1,P2 = K.primes_above(5)762sage: E.global_integral_model()763Elliptic 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 + 1764765:trac:`7935`::766767sage: K.<a> = NumberField(x^2-38)768sage: E = EllipticCurve([a,1/2])769sage: E.global_integral_model()770Elliptic Curve defined by y^2 = x^3 + 1444*a*x + 27436 over Number Field in a with defining polynomial x^2 - 38771772:trac:`9266`::773774sage: K.<s> = NumberField(x^2-5)775sage: w = (1+s)/2776sage: E = EllipticCurve(K,[2,w])777sage: E.global_integral_model()778Elliptic Curve defined by y^2 = x^3 + 2*x + (1/2*s+1/2) over Number Field in s with defining polynomial x^2 - 5779780:trac:`12151`::781782sage: K.<v> = NumberField(x^2 + 161*x - 150)783sage: E = EllipticCurve([25105/216*v - 3839/36, 634768555/7776*v - 98002625/1296, 634768555/7776*v - 98002625/1296, 0, 0])784sage: E.global_integral_model()785Elliptic Curve defined by y^2 + (-502639783*v+465618899)*x*y + (-6603604211463489399460860*v+6117229527723443603191500)*y = x^3 + (1526887622075335620*v-1414427901517840500)*x^2 over Number Field in v with defining polynomial x^2 + 161*x - 150786787:trac:`14476`::788789sage: R.<t> = QQ[]790sage: K.<g> = NumberField(t^4 - t^3 - 3*t^2 - t + 1)791sage: E = EllipticCurve([ -43/625*g^3 + 14/625*g^2 - 4/625*g + 706/625, -4862/78125*g^3 - 4074/78125*g^2 - 711/78125*g + 10304/78125, -4862/78125*g^3 - 4074/78125*g^2 - 711/78125*g + 10304/78125, 0,0])792sage: E.global_integral_model()793Elliptic Curve defined by y^2 + (15*g^3-48*g-42)*x*y + (-111510*g^3-162162*g^2-44145*g+37638)*y = x^3 + (-954*g^3-1134*g^2+81*g+576)*x^2 over Number Field in g with defining polynomial t^4 - t^3 - 3*t^2 - t + 1794795"""796K = self.base_field()797ai = self.a_invariants()798Ps = [[ ff[0] for ff in a.denominator_ideal().factor() ] for a in ai if not a.is_integral() ]799Ps = sage.misc.misc.union(sage.misc.flatten.flatten(Ps))800for P in Ps:801pi = K.uniformizer(P,'positive')802e = min([(ai[i].valuation(P)/[1,2,3,4,6][i]) for i in range(5)]).floor()803if e < 0 :804ai = [ai[i]/pi**(e*[1,2,3,4,6][i]) for i in range(5)]805if all(a.is_integral() for a in ai):806break807for z in ai:808assert z.is_integral(), "bug in global_integral_model: %s" % list(ai)809return EllipticCurve(list(ai))810811integral_model = global_integral_model812813def _reduce_model(self):814r"""815816Transforms the elliptic curve to a model in which `a_1`,817`a_2`, `a_3` are reduced modulo 2, 3, 2 respectively.818819.. note::820821This only works on integral models, i.e. it requires that822`a_1`, `a_2` and `a_3` lie in the ring of integers of the base823field.824825EXAMPLES::826827sage: K.<a>=NumberField(x^2-38)828sage: E=EllipticCurve([a, -5*a + 19, -39*a + 237, 368258520200522046806318224*a - 2270097978636731786720858047, 8456608930180227786550494643437985949781*a - 52130038506835491453281450568107193773505])829sage: E.ainvs()830(a,831-5*a + 19,832-39*a + 237,833368258520200522046806318224*a - 2270097978636731786720858047,8348456608930180227786550494643437985949781*a - 52130038506835491453281450568107193773505)835sage: E._reduce_model().ainvs()836(a,837a + 1,838a + 1,839368258520200522046806318444*a - 2270097978636731786720859345,8408456608930173478039472018047583706316424*a - 52130038506793883217874390501829588391299)841sage: EllipticCurve([101,202,303,404,505])._reduce_model().ainvs()842(1, 1, 0, -2509254, 1528863051)843sage: EllipticCurve([-101,-202,-303,-404,-505])._reduce_model().ainvs()844(1, -1, 0, -1823195, 947995262)845846sage: E = EllipticCurve([a/4, 1])847sage: E._reduce_model()848Traceback (most recent call last):849...850TypeError: _reduce_model() requires an integral model.851"""852K = self.base_ring()853ZK = K.maximal_order()854try:855(a1, a2, a3, a4, a6) = [ZK(a) for a in self.a_invariants()]856except TypeError:857import sys858raise TypeError, "_reduce_model() requires an integral model.", sys.exc_info()[2]859860# N.B. Must define s, r, t in the right order.861if ZK.absolute_degree() == 1:862s = ((-a1)/2).round('up')863r = ((-a2 + s*a1 +s*s)/3).round()864t = ((-a3 - r*a1)/2).round('up')865else:866pariK = K._pari_()867s = K(pariK.nfeltdiveuc(-a1, 2))868r = K(pariK.nfeltdiveuc(-a2 + s*a1 + s*s, 3))869t = K(pariK.nfeltdiveuc(-a3 - r*a1, 2))870871return self.rst_transform(r, s, t)872873def local_data(self, P=None, proof=None, algorithm="pari", globally=False):874r"""875Local data for this elliptic curve at the prime `P`.876877INPUT:878879- ``P`` -- either None or a prime ideal of the base field of self.880881- ``proof`` -- whether to only use provably correct methods882(default controlled by global proof module). Note that the883proof module is number_field, not elliptic_curves, since the884functions that actually need the flag are in number fields.885886- ``algorithm`` (string, default: "pari") -- Ignored unless the887base field is `\QQ`. If "pari", use the PARI C-library888``ellglobalred`` implementation of Tate's algorithm over889`\QQ`. If "generic", use the general number field890implementation.891892- ``globally`` -- whether the local algorithm uses global generators893for the prime ideals. Default is False, which won't require any894information about the class group. If True, a generator for `P`895will be used if `P` is principal. Otherwise, or if ``globally``896is False, the minimal model returned will preserve integrality897at other primes, but not minimality.898899OUTPUT:900901If `P` is specified, returns the ``EllipticCurveLocalData``902object associated to the prime `P` for this curve. Otherwise,903returns a list of such objects, one for each prime `P` in the904support of the discriminant of this model.905906.. note::907908The model is not required to be integral on input.909910EXAMPLES::911912sage: K.<i> = NumberField(x^2+1)913sage: E = EllipticCurve([1 + i, 0, 1, 0, 0])914sage: E.local_data()915[Local data at Fractional ideal (i - 2):916Reduction type: bad non-split multiplicative917Local 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 + 1918Minimal discriminant valuation: 1919Conductor exponent: 1920Kodaira Symbol: I1921Tamagawa Number: 1,922Local data at Fractional ideal (-3*i - 2):923Reduction type: bad split multiplicative924Local 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 + 1925Minimal discriminant valuation: 2926Conductor exponent: 1927Kodaira Symbol: I2928Tamagawa Number: 2]929sage: E.local_data(K.ideal(3))930Local data at Fractional ideal (3):931Reduction type: good932Local 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 + 1933Minimal discriminant valuation: 0934Conductor exponent: 0935Kodaira Symbol: I0936Tamagawa Number: 1937938An example raised in :trac:`3897`::939940sage: E = EllipticCurve([1,1])941sage: E.local_data(3)942Local data at Principal ideal (3) of Integer Ring:943Reduction type: good944Local minimal model: Elliptic Curve defined by y^2 = x^3 + x + 1 over Rational Field945Minimal discriminant valuation: 0946Conductor exponent: 0947Kodaira Symbol: I0948Tamagawa Number: 1949"""950if proof is None:951import sage.structure.proof.proof952# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.953proof = sage.structure.proof.proof.get_flag(None, "number_field")954955if P is None:956primes = self.base_ring()(self.integral_model().discriminant()).support()957return [self._get_local_data(pr, proof) for pr in primes]958959from sage.schemes.elliptic_curves.ell_local_data import check_prime960P = check_prime(self.base_field(),P)961962return self._get_local_data(P,proof,algorithm,globally)963964def _get_local_data(self, P, proof, algorithm="pari", globally=False):965r"""966Internal function to create data for this elliptic curve at the prime `P`.967968This function handles the caching of local data. It is called969by local_data() which is the user interface and which parses970the input parameters `P` and proof.971972INPUT:973974- ``P`` -- either None or a prime ideal of the base field of self.975976- ``proof`` -- whether to only use provably correct methods977(default controlled by global proof module). Note that the978proof module is number_field, not elliptic_curves, since the979functions that actually need the flag are in number fields.980981- ``algorithm`` (string, default: "pari") -- Ignored unless the982base field is `\QQ`. If "pari", use the PARI C-library983``ellglobalred`` implementation of Tate's algorithm over984`\QQ`. If "generic", use the general number field985implementation.986987- ``globally`` -- whether the local algorithm uses global generators988for the prime ideals. Default is False, which won't require any989information about the class group. If True, a generator for `P`990will be used if `P` is principal. Otherwise, or if ``globally``991is False, the minimal model returned will preserve integrality992at other primes, but not minimality.993994EXAMPLES::995996sage: K.<i> = NumberField(x^2+1)997sage: E = EllipticCurve(K,[0,1,0,-160,308])998sage: p = K.ideal(i+1)999sage: E._get_local_data(p, False)1000Local data at Fractional ideal (i + 1):1001Reduction type: good1002Local 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 + 11003Minimal discriminant valuation: 01004Conductor exponent: 01005Kodaira Symbol: I01006Tamagawa Number: 110071008Verify that we cache based on the proof value and the algorithm choice::10091010sage: E._get_local_data(p, False) is E._get_local_data(p, True)1011False10121013sage: E._get_local_data(p, None, "pari") is E._get_local_data(p, None, "generic")1014False1015"""1016try:1017return self._local_data[P, proof, algorithm, globally]1018except AttributeError:1019self._local_data = {}1020except KeyError:1021pass1022from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData1023self._local_data[P, proof, algorithm, globally] = EllipticCurveLocalData(self, P, proof, algorithm, globally)1024return self._local_data[P, proof, algorithm, globally]10251026def local_minimal_model(self, P, proof = None, algorithm="pari"):1027r"""1028Returns a model which is integral at all primes and minimal at `P`.10291030INPUT:10311032- ``P`` -- either None or a prime ideal of the base field of self.10331034- ``proof`` -- whether to only use provably correct methods1035(default controlled by global proof module). Note that the1036proof module is number_field, not elliptic_curves, since the1037functions that actually need the flag are in number fields.10381039- ``algorithm`` (string, default: "pari") -- Ignored unless the1040base field is `\QQ`. If "pari", use the PARI C-library1041``ellglobalred`` implementation of Tate's algorithm over1042`\QQ`. If "generic", use the general number field1043implementation.10441045OUTPUT:10461047A model of the curve which is minimal (and integral) at `P`.10481049.. note::10501051The model is not required to be integral on input.10521053For principal `P`, a generator is used as a uniformizer,1054and integrality or minimality at other primes is not1055affected. For non-principal `P`, the minimal model1056returned will preserve integrality at other primes, but not1057minimality.10581059EXAMPLES::10601061sage: K.<a>=NumberField(x^2-5)1062sage: E=EllipticCurve([20, 225, 750, 1250*a + 6250, 62500*a + 15625])1063sage: P=K.ideal(a)1064sage: E.local_minimal_model(P).ainvs()1065(0, 1, 0, 2*a - 34, -4*a + 66)1066"""1067if proof is None:1068import sage.structure.proof.proof1069# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.1070proof = sage.structure.proof.proof.get_flag(None, "number_field")10711072return self.local_data(P, proof, algorithm).minimal_model()10731074def has_good_reduction(self, P):1075r"""1076Return True if this elliptic curve has good reduction at the prime `P`.10771078INPUT:10791080- ``P`` -- a prime ideal of the base field of self, or a field1081element generating such an ideal.10821083OUTPUT:10841085(bool) -- True if the curve has good reduction at `P`, else False.10861087.. note::10881089This requires determining a local integral minimal model;1090we do not just check that the discriminant of the current1091model has valuation zero.10921093EXAMPLES::10941095sage: E=EllipticCurve('14a1')1096sage: [(p,E.has_good_reduction(p)) for p in prime_range(15)]1097[(2, False), (3, True), (5, True), (7, False), (11, True), (13, True)]10981099sage: K.<a>=NumberField(x^3-2)1100sage: P17a, P17b = [P for P,e in K.factor(17)]1101sage: E = EllipticCurve([0,0,0,0,2*a+1])1102sage: [(p,E.has_good_reduction(p)) for p in [P17a,P17b]]1103[(Fractional ideal (4*a^2 - 2*a + 1), True),1104(Fractional ideal (2*a + 1), False)]1105"""1106return self.local_data(P).has_good_reduction()11071108def has_bad_reduction(self, P):1109r"""1110Return True if this elliptic curve has bad reduction at the prime `P`.11111112INPUT:11131114- ``P`` -- a prime ideal of the base field of self, or a field1115element generating such an ideal.11161117OUTPUT:11181119(bool) True if the curve has bad reduction at `P`, else False.11201121.. note::11221123This requires determining a local integral minimal model;1124we do not just check that the discriminant of the current1125model has valuation zero.11261127EXAMPLES::11281129sage: E=EllipticCurve('14a1')1130sage: [(p,E.has_bad_reduction(p)) for p in prime_range(15)]1131[(2, True), (3, False), (5, False), (7, True), (11, False), (13, False)]11321133sage: K.<a>=NumberField(x^3-2)1134sage: P17a, P17b = [P for P,e in K.factor(17)]1135sage: E = EllipticCurve([0,0,0,0,2*a+1])1136sage: [(p,E.has_bad_reduction(p)) for p in [P17a,P17b]]1137[(Fractional ideal (4*a^2 - 2*a + 1), False),1138(Fractional ideal (2*a + 1), True)]1139"""1140return self.local_data(P).has_bad_reduction()11411142def has_multiplicative_reduction(self, P):1143r"""1144Return True if this elliptic curve has (bad) multiplicative reduction at the prime `P`.11451146.. note::11471148See also ``has_split_multiplicative_reduction()`` and1149``has_nonsplit_multiplicative_reduction()``.11501151INPUT:11521153- ``P`` -- a prime ideal of the base field of self, or a field1154element generating such an ideal.11551156OUTPUT:11571158(bool) True if the curve has multiplicative reduction at `P`,1159else False.11601161EXAMPLES::11621163sage: E=EllipticCurve('14a1')1164sage: [(p,E.has_multiplicative_reduction(p)) for p in prime_range(15)]1165[(2, True), (3, False), (5, False), (7, True), (11, False), (13, False)]11661167sage: K.<a>=NumberField(x^3-2)1168sage: P17a, P17b = [P for P,e in K.factor(17)]1169sage: E = EllipticCurve([0,0,0,0,2*a+1])1170sage: [(p,E.has_multiplicative_reduction(p)) for p in [P17a,P17b]]1171[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]1172"""1173return self.local_data(P).has_multiplicative_reduction()11741175def has_split_multiplicative_reduction(self, P):1176r"""1177Return True if this elliptic curve has (bad) split multiplicative reduction at the prime `P`.11781179INPUT:11801181- ``P`` -- a prime ideal of the base field of self, or a field1182element generating such an ideal.11831184OUTPUT:11851186(bool) True if the curve has split multiplicative reduction at1187`P`, else False.11881189EXAMPLES::11901191sage: E=EllipticCurve('14a1')1192sage: [(p,E.has_split_multiplicative_reduction(p)) for p in prime_range(15)]1193[(2, False), (3, False), (5, False), (7, True), (11, False), (13, False)]11941195sage: K.<a>=NumberField(x^3-2)1196sage: P17a, P17b = [P for P,e in K.factor(17)]1197sage: E = EllipticCurve([0,0,0,0,2*a+1])1198sage: [(p,E.has_split_multiplicative_reduction(p)) for p in [P17a,P17b]]1199[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]1200"""1201return self.local_data(P).has_split_multiplicative_reduction()12021203def has_nonsplit_multiplicative_reduction(self, P):1204r"""1205Return True if this elliptic curve has (bad) non-split multiplicative reduction at the prime `P`.12061207INPUT:12081209- ``P`` -- a prime ideal of the base field of self, or a field1210element generating such an ideal.12111212OUTPUT:12131214(bool) True if the curve has non-split multiplicative1215reduction at `P`, else False.12161217EXAMPLES::12181219sage: E=EllipticCurve('14a1')1220sage: [(p,E.has_nonsplit_multiplicative_reduction(p)) for p in prime_range(15)]1221[(2, True), (3, False), (5, False), (7, False), (11, False), (13, False)]12221223sage: K.<a>=NumberField(x^3-2)1224sage: P17a, P17b = [P for P,e in K.factor(17)]1225sage: E = EllipticCurve([0,0,0,0,2*a+1])1226sage: [(p,E.has_nonsplit_multiplicative_reduction(p)) for p in [P17a,P17b]]1227[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]1228"""1229return self.local_data(P).has_nonsplit_multiplicative_reduction()12301231def has_additive_reduction(self, P):1232r"""1233Return True if this elliptic curve has (bad) additive reduction at the prime `P`.12341235INPUT:12361237- ``P`` -- a prime ideal of the base field of self, or a field1238element generating such an ideal.12391240OUTPUT:12411242(bool) True if the curve has additive reduction at `P`, else False.12431244EXAMPLES::12451246sage: E=EllipticCurve('27a1')1247sage: [(p,E.has_additive_reduction(p)) for p in prime_range(15)]1248[(2, False), (3, True), (5, False), (7, False), (11, False), (13, False)]12491250sage: K.<a>=NumberField(x^3-2)1251sage: P17a, P17b = [P for P,e in K.factor(17)]1252sage: E = EllipticCurve([0,0,0,0,2*a+1])1253sage: [(p,E.has_additive_reduction(p)) for p in [P17a,P17b]]1254[(Fractional ideal (4*a^2 - 2*a + 1), False),1255(Fractional ideal (2*a + 1), True)]1256"""1257return self.local_data(P).has_additive_reduction()12581259def tamagawa_number(self, P, proof = None):1260r"""1261Returns the Tamagawa number of this elliptic curve at the prime `P`.12621263INPUT:12641265- ``P`` -- either None or a prime ideal of the base field of self.12661267- ``proof`` -- whether to only use provably correct methods1268(default controlled by global proof module). Note that the1269proof module is number_field, not elliptic_curves, since the1270functions that actually need the flag are in number fields.12711272OUTPUT:12731274(positive integer) The Tamagawa number of the curve at `P`.12751276EXAMPLES::12771278sage: K.<a>=NumberField(x^2-5)1279sage: E=EllipticCurve([20, 225, 750, 625*a + 6875, 31250*a + 46875])1280sage: [E.tamagawa_number(P) for P in E.discriminant().support()]1281[1, 1, 1, 1]1282sage: K.<a> = QuadraticField(-11)1283sage: E = EllipticCurve('11a1').change_ring(K)1284sage: [E.tamagawa_number(P) for P in K(11).support()]1285[10]1286"""1287if proof is None:1288import sage.structure.proof.proof1289# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.1290proof = sage.structure.proof.proof.get_flag(None, "number_field")12911292return self.local_data(P, proof).tamagawa_number()12931294def tamagawa_numbers(self):1295"""1296Return a list of all Tamagawa numbers for all prime divisors of the1297conductor (in order).12981299EXAMPLES::13001301sage: e = EllipticCurve('30a1')1302sage: e.tamagawa_numbers()1303[2, 3, 1]1304sage: vector(e.tamagawa_numbers())1305(2, 3, 1)1306sage: K.<a>=NumberField(x^2+3)1307sage: eK = e.base_extend(K)1308sage: eK.tamagawa_numbers()1309[4, 6, 1]1310"""1311return [self.tamagawa_number(p) for p in prime_divisors(self.conductor())]13121313def tamagawa_exponent(self, P, proof = None):1314r"""1315Returns the Tamagawa index of this elliptic curve at the prime `P`.13161317INPUT:13181319- ``P`` -- either None or a prime ideal of the base field of self.13201321- ``proof`` -- whether to only use provably correct methods1322(default controlled by global proof module). Note that the1323proof module is number_field, not elliptic_curves, since the1324functions that actually need the flag are in number fields.13251326OUTPUT:13271328(positive integer) The Tamagawa index of the curve at P.13291330EXAMPLES::13311332sage: K.<a>=NumberField(x^2-5)1333sage: E=EllipticCurve([20, 225, 750, 625*a + 6875, 31250*a + 46875])1334sage: [E.tamagawa_exponent(P) for P in E.discriminant().support()]1335[1, 1, 1, 1]1336sage: K.<a> = QuadraticField(-11)1337sage: E = EllipticCurve('11a1').change_ring(K)1338sage: [E.tamagawa_exponent(P) for P in K(11).support()]1339[10]1340"""1341if proof is None:1342import sage.structure.proof.proof1343# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.1344proof = sage.structure.proof.proof.get_flag(None, "number_field")13451346return self.local_data(P, proof).tamagawa_exponent()13471348def tamagawa_product_bsd(self):1349r"""1350Given an elliptic curve `E` over a number field `K`, this function returns the1351integer `C(E/K)` that appears in the Birch and Swinnerton-Dyer conjecture accounting1352for the local information at finite places. If the model is a global minimal model then `C(E/K)` is1353simply 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.1354The definition is such that `C(E/K)` times the periods at the infinite places is invariant1355under change of the Weierstrass model. See [Ta2] and [Do] for details.13561357.. note::13581359This definition is slightly different from the definition of ``tamagawa_product``1360for curves defined over `\QQ`. Over the rational number it is always defined to be the product1361of the Tamagawa numbers, so the two definitions only agree when the model is global minimal.13621363OUTPUT:13641365A rational number13661367EXAMPLES::13681369sage: K.<i> = NumberField(x^2+1)1370sage: E = EllipticCurve([0,2+i])1371sage: E.tamagawa_product_bsd()1372113731374sage: E = EllipticCurve([(2*i+1)^2,i*(2*i+1)^7])1375sage: E.tamagawa_product_bsd()1376413771378An example where the Neron model changes over K::13791380sage: K.<t> = NumberField(x^5-10*x^3+5*x^2+10*x+1)1381sage: E = EllipticCurve(K,'75a1')1382sage: E.tamagawa_product_bsd()138351384sage: da = E.local_data()1385sage: [dav.tamagawa_number() for dav in da]1386[1, 1]13871388An example over `\mathbb{Q}` (:trac:`9413`)::13891390sage: E = EllipticCurve('30a')1391sage: E.tamagawa_product_bsd()1392613931394REFERENCES:13951396- [Ta2] Tate, John, On the conjectures of Birch and Swinnerton-Dyer and a geometric analog. Seminaire Bourbaki, Vol. 9, Exp. No. 306.13971398- [Do] Dokchitser, Tim and Vladimir, On the Birch-Swinnerton-Dyer quotients modulo squares, Annals of Math., 2010.13991400"""1401da = self.local_data()1402pr = 11403for dav in da:1404pp = dav.prime()1405cv = dav.tamagawa_number()1406# uu is the quotient of the Neron differential at pp divided by1407# the differential associated to this particular equation E1408uu = self.isomorphism_to(dav.minimal_model()).u1409if self.base_field().absolute_degree() == 1:1410p = pp.gens_reduced()[0]1411f = 11412v = valuation(ZZ(uu),p)1413else:1414p = pp.smallest_integer()1415f = pp.residue_class_degree()1416v = valuation(uu,pp)1417uu_abs_val = p**(f*v)1418pr *= cv * uu_abs_val1419return pr14201421def kodaira_symbol(self, P, proof = None):1422r"""1423Returns the Kodaira Symbol of this elliptic curve at the prime `P`.14241425INPUT:14261427- ``P`` -- either None or a prime ideal of the base field of self.14281429- ``proof`` -- whether to only use provably correct methods1430(default controlled by global proof module). Note that the1431proof module is number_field, not elliptic_curves, since the1432functions that actually need the flag are in number fields.14331434OUTPUT:14351436The Kodaira Symbol of the curve at P, represented as a string.14371438EXAMPLES::14391440sage: K.<a>=NumberField(x^2-5)1441sage: E=EllipticCurve([20, 225, 750, 625*a + 6875, 31250*a + 46875])1442sage: bad_primes = E.discriminant().support(); bad_primes1443[Fractional ideal (-a), Fractional ideal (7/2*a - 81/2), Fractional ideal (-a - 52), Fractional ideal (2)]1444sage: [E.kodaira_symbol(P) for P in bad_primes]1445[I0, I1, I1, II]1446sage: K.<a> = QuadraticField(-11)1447sage: E = EllipticCurve('11a1').change_ring(K)1448sage: [E.kodaira_symbol(P) for P in K(11).support()]1449[I10]1450"""1451if proof is None:1452import sage.structure.proof.proof1453# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.1454proof = sage.structure.proof.proof.get_flag(None, "number_field")14551456return self.local_data(P, proof).kodaira_symbol()145714581459def conductor(self):1460r"""1461Returns the conductor of this elliptic curve as a fractional1462ideal of the base field.14631464OUTPUT:14651466(fractional ideal) The conductor of the curve.14671468EXAMPLES::14691470sage: K.<i>=NumberField(x^2+1)1471sage: EllipticCurve([i, i - 1, i + 1, 24*i + 15, 14*i + 35]).conductor()1472Fractional ideal (21*i - 3)1473sage: K.<a>=NumberField(x^2-x+3)1474sage: EllipticCurve([1 + a , -1 + a , 1 + a , -11 + a , 5 -9*a ]).conductor()1475Fractional ideal (6*a)14761477A not so well known curve with everywhere good reduction::14781479sage: K.<a>=NumberField(x^2-38)1480sage: E=EllipticCurve([0,0,0, 21796814856932765568243810*a - 134364590724198567128296995, 121774567239345229314269094644186997594*a - 750668847495706904791115375024037711300])1481sage: E.conductor()1482Fractional ideal (1)14831484An example which used to fail (see :trac:`5307`)::14851486sage: K.<w>=NumberField(x^2+x+6)1487sage: E=EllipticCurve([w,-1,0,-w-6,0])1488sage: E.conductor()1489Fractional ideal (86304, w + 5898)14901491An example raised in :trac:`11346`::14921493sage: K.<g> = NumberField(x^2 - x - 1)1494sage: E1 = EllipticCurve(K,[0,0,0,-1/48,-161/864])1495sage: [(p.smallest_integer(),e) for p,e in E1.conductor().factor()]1496[(2, 4), (3, 1), (5, 1)]1497"""1498try:1499return self._conductor1500except AttributeError:1501pass15021503# Note: for number fields other than QQ we could initialize1504# N=K.ideal(1) or N=OK.ideal(1), which are the same, but for1505# K==QQ it has to be ZZ.ideal(1).1506OK = self.base_ring().ring_of_integers()1507self._conductor = prod([d.prime()**(d.conductor_valuation()) \1508for d in self.local_data()],\1509OK.ideal(1))1510return self._conductor15111512def global_minimal_model(self, proof = None):1513r"""1514Returns a model of self that is integral, minimal at all primes.15151516.. note::15171518This is only implemented for class number 1. In general,1519such a model may or may not exist.15201521INPUT:15221523- ``proof`` -- whether to only use provably correct methods1524(default controlled by global proof module). Note that the1525proof module is number_field, not elliptic_curves, since the1526functions that actually need the flag are in number fields.15271528OUTPUT:15291530A global integral and minimal model.15311532EXAMPLES::15331534sage: K.<a> = NumberField(x^2-38)1535sage: E = EllipticCurve([0,0,0, 21796814856932765568243810*a - 134364590724198567128296995, 121774567239345229314269094644186997594*a - 750668847495706904791115375024037711300])15361537sage: E2 = E.global_minimal_model()1538sage: E2 # random (the global minimal model is not unique)1539Elliptic 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 - 3815401541sage: E2.local_data()1542[]15431544See :trac:`11347`::15451546sage: K.<g> = NumberField(x^2 - x - 1)1547sage: E = EllipticCurve(K,[0,0,0,-1/48,161/864]).integral_model().global_minimal_model(); E1548Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 over Number Field in g with defining polynomial x^2 - x - 11549sage: [(p.norm(), e) for p, e in E.conductor().factor()]1550[(9, 1), (5, 1)]1551sage: [(p.norm(), e) for p, e in E.discriminant().factor()]1552[(-5, 2), (9, 1)]15531554See :trac:`14472`, this used not to work over a relative extension::15551556sage: K1.<w> = NumberField(x^2+x+1)1557sage: m = polygen(K1)1558sage: K2.<v> = K1.extension(m^2-w+1)1559sage: E = EllipticCurve([0*v,-432])1560sage: E.global_minimal_model()1561Elliptic Curve defined by y^2 + (v+w+1)*y = x^3 + ((6*w-10)*v+16*w+20) over Number Field in v with defining polynomial x^2 - w + 1 over its base field1562"""1563if proof is None:1564import sage.structure.proof.proof1565# We use the "number_field" flag because the actual proof dependence is in PARI's number field functions.1566proof = sage.structure.proof.proof.get_flag(None, "number_field")1567K = self.base_ring()1568if K.class_number() != 1:1569raise ValueError, "global minimal models only exist in general for class number 1"15701571E = self.global_integral_model()1572primes = E.base_ring()(E.discriminant()).support()1573for P in primes:1574E = E.local_data(P,proof, globally=True).minimal_model()1575return E._reduce_model()15761577def reduction(self,place):1578r"""1579Return the reduction of the elliptic curve at a place of good reduction.15801581INPUT:15821583- ``place`` -- a prime ideal in the base field of the curve15841585OUTPUT:15861587An elliptic curve over a finite field, the residue field of the place.15881589EXAMPLES::15901591sage: K.<i> = QuadraticField(-1)1592sage: EK = EllipticCurve([0,0,0,i,i+3])1593sage: v = K.fractional_ideal(2*i+3)1594sage: EK.reduction(v)1595Elliptic Curve defined by y^2 = x^3 + 5*x + 8 over Residue field of Fractional ideal (2*i + 3)1596sage: EK.reduction(K.ideal(1+i))1597Traceback (most recent call last):1598...1599ValueError: The curve must have good reduction at the place.1600sage: EK.reduction(K.ideal(2))1601Traceback (most recent call last):1602...1603ValueError: The ideal must be prime.1604sage: K=QQ.extension(x^2+x+1,"a")1605sage: E=EllipticCurve([1024*K.0,1024*K.0])1606sage: E.reduction(2*K)1607Elliptic Curve defined by y^2 + (abar+1)*y = x^3 over Residue field in abar of Fractional ideal (2)1608"""1609K = self.base_field()1610OK = K.ring_of_integers()1611try:1612place = K.ideal(place)1613except TypeError:1614raise TypeError, "The parameter must be an ideal of the base field of the elliptic curve"1615if not place.is_prime():1616raise ValueError, "The ideal must be prime."1617disc = self.discriminant()1618if not K.ideal(disc).valuation(place) == 0:1619local_data=self.local_data(place)1620if local_data.has_good_reduction():1621Fv = OK.residue_field(place)1622return local_data.minimal_model().change_ring(Fv)1623raise ValueError, "The curve must have good reduction at the place."1624Fv = OK.residue_field(place)1625return self.change_ring(Fv)16261627def _torsion_bound(self,number_of_places = 20):1628r"""1629An upper bound on the order of the torsion subgroup.16301631INPUT:16321633- ``number_of_places`` (positive integer, default = 20) -- the1634number of places that will be used to find the bound.16351636OUTPUT:16371638(integer) An upper bound on the torsion order.16391640ALGORITHM:16411642An upper bound on the order of the torsion.group of the1643elliptic curve is obtained by counting points modulo several1644primes of good reduction. Note that the upper bound returned1645by this function is a multiple of the order of the torsion1646group, and in general will be greater than the order.16471648EXAMPLES::16491650sage: CDB=CremonaDatabase()1651sage: [E._torsion_bound() for E in CDB.iter([14])]1652[6, 6, 6, 6, 6, 6]1653sage: [E.torsion_order() for E in CDB.iter([14])]1654[6, 6, 2, 6, 2, 6]1655"""1656E = self1657bound = 01658k = 01659K = E.base_field()1660OK = K.ring_of_integers()1661disc = E.discriminant()1662p = Integer(1)1663# runs through primes, decomposes them into prime ideals1664while k < number_of_places :1665p = p.next_prime()1666f = K.primes_above(p)1667# runs through prime ideals above p1668for qq in f:1669fqq = qq.residue_class_degree()1670charqq = qq.smallest_integer()1671# take only places with small residue field (so that the1672# number of points will be small)1673if fqq == 1 or charqq**fqq < 3*number_of_places:1674# check if the model is integral at the place1675if min([K.ideal(a).valuation(qq) for a in E.a_invariants()]) >= 0:1676eqq = qq.ramification_index()1677# check if the formal group at the place is torsion-free1678# if so the torsion injects into the reduction1679if eqq < charqq - 1 and disc.valuation(qq) == 0:1680Etilda = E.reduction(qq)1681Npp = Etilda.cardinality()1682bound = gcd(bound,Npp)1683if bound == 1:1684return bound1685k += 11686return bound16871688def torsion_subgroup(self):1689r"""1690Returns the torsion subgroup of this elliptic curve.16911692OUTPUT:16931694(``EllipticCurveTorsionSubgroup``) The1695``EllipticCurveTorsionSubgroup`` associated to this elliptic1696curve.16971698EXAMPLES::16991700sage: E = EllipticCurve('11a1')1701sage: K.<t>=NumberField(x^4 + x^3 + 11*x^2 + 41*x + 101)1702sage: EK = E.base_extend(K)1703sage: tor = EK.torsion_subgroup() # long time (2s on sage.math, 2014)1704sage: tor # long time1705Torsion 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 + 1011706sage: tor.gens() # long time1707((16 : 60 : 1), (t : 1/11*t^3 + 6/11*t^2 + 19/11*t + 48/11 : 1))17081709::17101711sage: E = EllipticCurve('15a1')1712sage: K.<t>=NumberField(x^2 + 2*x + 10)1713sage: EK=E.base_extend(K)1714sage: EK.torsion_subgroup()1715Torsion 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 + 1017161717::17181719sage: E = EllipticCurve('19a1')1720sage: 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)1721sage: EK=E.base_extend(K)1722sage: EK.torsion_subgroup()1723Torsion 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 + 117241725::17261727sage: K.<i> = QuadraticField(-1)1728sage: EK = EllipticCurve([0,0,0,i,i+3])1729sage: EK.torsion_subgroup ()1730Torsion 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 + 11731"""1732try:1733return self.__torsion_subgroup1734except AttributeError:1735self.__torsion_subgroup = ell_torsion.EllipticCurveTorsionSubgroup(self)1736return self.__torsion_subgroup17371738def torsion_order(self):1739r"""1740Returns the order of the torsion subgroup of this elliptic curve.17411742OUTPUT:17431744(integer) the order of the torsion subgroup of this elliptic curve.17451746EXAMPLES::17471748sage: E = EllipticCurve('11a1')1749sage: K.<t> = NumberField(x^4 + x^3 + 11*x^2 + 41*x + 101)1750sage: EK = E.base_extend(K)1751sage: EK.torsion_order() # long time (2s on sage.math, 2014)17522517531754::17551756sage: E = EllipticCurve('15a1')1757sage: K.<t> = NumberField(x^2 + 2*x + 10)1758sage: EK = E.base_extend(K)1759sage: EK.torsion_order()17601617611762::17631764sage: E = EllipticCurve('19a1')1765sage: 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)1766sage: EK = E.base_extend(K)1767sage: EK.torsion_order()1768917691770::17711772sage: K.<i> = QuadraticField(-1)1773sage: EK = EllipticCurve([0,0,0,i,i+3])1774sage: EK.torsion_order()177511776"""1777try:1778return self.__torsion_order1779except AttributeError:1780self.__torsion_order = self.torsion_subgroup().order()1781return self.__torsion_order17821783def torsion_points(self):1784r"""1785Returns a list of the torsion points of this elliptic curve.17861787OUTPUT:17881789(list) A sorted list of the torsion points.17901791EXAMPLES::17921793sage: E = EllipticCurve('11a1')1794sage: E.torsion_points()1795[(0 : 1 : 0), (5 : -6 : 1), (5 : 5 : 1), (16 : -61 : 1), (16 : 60 : 1)]1796sage: K.<t> = NumberField(x^4 + x^3 + 11*x^2 + 41*x + 101)1797sage: EK = E.base_extend(K)1798sage: EK.torsion_points() # long time (1s on sage.math, 2014)1799[(16 : 60 : 1),1800(5 : 5 : 1),1801(5 : -6 : 1),1802(16 : -61 : 1),1803(t : 1/11*t^3 + 6/11*t^2 + 19/11*t + 48/11 : 1),1804(-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),1805(-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),1806(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),1807(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),1808(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),1809(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),1810(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),1811(-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),1812(-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),1813(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),1814(-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),1815(-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),1816(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),1817(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),1818(t : -1/11*t^3 - 6/11*t^2 - 19/11*t - 59/11 : 1),1819(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),1820(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),1821(-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),1822(-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),1823(0 : 1 : 0)]18241825::18261827sage: E = EllipticCurve('15a1')1828sage: K.<t> = NumberField(x^2 + 2*x + 10)1829sage: EK = E.base_extend(K)1830sage: EK.torsion_points()1831[(-7 : -5*t - 2 : 1),1832(-7 : 5*t + 8 : 1),1833(-13/4 : 9/8 : 1),1834(-2 : -2 : 1),1835(-2 : 3 : 1),1836(-t - 2 : -t - 7 : 1),1837(-t - 2 : 2*t + 8 : 1),1838(-1 : 0 : 1),1839(t : t - 5 : 1),1840(t : -2*t + 4 : 1),1841(0 : 1 : 0),1842(1/2 : -5/4*t - 2 : 1),1843(1/2 : 5/4*t + 1/2 : 1),1844(3 : -2 : 1),1845(8 : -27 : 1),1846(8 : 18 : 1)]18471848::18491850sage: K.<i> = QuadraticField(-1)1851sage: EK = EllipticCurve(K,[0,0,0,0,-1])1852sage: EK.torsion_points ()1853[(-2 : -3*i : 1), (-2 : 3*i : 1), (0 : -i : 1), (0 : i : 1), (0 : 1 : 0), (1 : 0 : 1)]1854"""1855T = self.torsion_subgroup() # make sure it is cached1856return sorted(T.points()) # these are also cached in T18571858def rank_bounds(self,verbose=0, lim1=5, lim3=50, limtriv=10, maxprob=20, limbigprime=30):1859r"""1860Returns the lower and upper bounds using :meth:`~simon_two_descent`.1861The results of :meth:`~simon_two_descent` are cached.18621863.. NOTE::18641865The optional parameters control the Simon two descent algorithm;1866see the documentation of :meth:`~simon_two_descent` for more1867details.18681869INPUT:18701871- ``verbose`` -- 0, 1, 2, or 3 (default: 0), the verbosity level18721873- ``lim1`` -- (default: 5) limit on trivial points on quartics18741875- ``lim3`` -- (default: 50) limit on points on ELS quartics18761877- ``limtriv`` -- (default: 10) limit on trivial points on elliptic curve18781879- ``maxprob`` -- (default: 20)18801881- ``limbigprime`` -- (default: 30) to distinguish between1882small and large prime numbers. Use probabilistic tests for1883large primes. If 0, don't use probabilistic tests.18841885OUTPUT:18861887lower and upper bounds for the rank of the Mordell-Weil group188818891890.. NOTE::18911892For non-quadratic number fields, this code does1893return, but it takes a long time.18941895EXAMPLES::18961897sage: K.<a> = NumberField(x^2 + 23, 'a')1898sage: E = EllipticCurve(K, '37')1899sage: E == loads(dumps(E))1900True1901sage: E.rank_bounds()1902(2, 2)19031904Here is a curve with two-torsion, again the bounds coincide::19051906sage: Qrt5.<rt5>=NumberField(x^2-5)1907sage: E=EllipticCurve([0,5-rt5,0,rt5,0])1908sage: E.rank_bounds()1909(1, 1)19101911Finally an example with non-trivial 2-torsion in Sha. So the19122-descent will not be able to determine the rank, but can only1913give bounds::19141915sage: E = EllipticCurve("15a5")1916sage: K.<t> = NumberField(x^2-6)1917sage: EK = E.base_extend(K)1918sage: EK.rank_bounds(lim1=1,lim3=1,limtriv=1)1919(0, 2)19201921IMPLEMENTATION:19221923Uses Denis Simon's PARI/GP scripts from1924http://www.math.unicaen.fr/~simon/.19251926"""19271928lower, upper, gens = self.simon_two_descent(verbose=verbose,lim1=lim1,lim3=lim3,limtriv=limtriv,maxprob=maxprob,limbigprime=limbigprime)1929# this was corrected in trac 13593. upper is the dimension1930# of the 2-selmer group, so we can certainly remove the1931# 2-torsion of the Mordell-Weil group.1932upper -= self.two_torsion_rank()1933return lower, upper19341935def rank(self,verbose=0, lim1=5, lim3=50, limtriv=10, maxprob=20, limbigprime=30):1936r"""1937Return the rank of this elliptic curve, if it can be determined.19381939.. NOTE::19401941The optional parameters control the Simon two descent algorithm;1942see the documentation of :meth:`~simon_two_descent` for more1943details.19441945INPUT:19461947- ``verbose`` -- 0, 1, 2, or 3 (default: 0), the verbosity level19481949- ``lim1`` -- (default: 5) limit on trivial points on quartics19501951- ``lim3`` -- (default: 50) limit on points on ELS quartics19521953- ``limtriv`` -- (default: 10) limit on trivial points on elliptic curve19541955- ``maxprob`` -- (default: 20)19561957- ``limbigprime`` -- (default: 30) to distinguish between1958small and large prime numbers. Use probabilistic tests for1959large primes. If 0, don't use probabilistic tests.19601961OUTPUT:19621963If the upper and lower bounds given by Simon two-descent are1964the same, then the rank has been uniquely identified and we1965return this. Otherwise, we raise a ValueError with an error1966message specifying the upper and lower bounds.19671968.. NOTE::19691970For non-quadratic number fields, this code does return, but it takes1971a long time.19721973EXAMPLES::19741975sage: K.<a> = NumberField(x^2 + 23, 'a')1976sage: E = EllipticCurve(K, '37')1977sage: E == loads(dumps(E))1978True1979sage: E.rank()1980219811982Here is a curve with two-torsion in the Tate-Shafarevich group,1983so here the bounds given by the algorithm do not uniquely1984determine the rank::19851986sage: E = EllipticCurve("15a5")1987sage: K.<t> = NumberField(x^2-6)1988sage: EK = E.base_extend(K)1989sage: EK.rank(lim1=1, lim3=1, limtriv=1)1990Traceback (most recent call last):1991...1992ValueError: There is insufficient data to determine the rank -19932-descent gave lower bound 0 and upper bound 219941995IMPLEMENTATION:19961997Uses Denis Simon's PARI/GP scripts from1998http://www.math.unicaen.fr/~simon/.19992000"""20012002lower,upper = self.rank_bounds(verbose=verbose,lim1=lim1,lim3=lim3,limtriv=limtriv,maxprob=maxprob,limbigprime=limbigprime)2003if lower == upper:2004return lower2005else:2006raise ValueError, 'There is insufficient data to determine the rank - 2-descent gave lower bound %s and upper bound %s' % (lower, upper)20072008def gens(self,verbose=0, lim1=5, lim3=50, limtriv=10, maxprob=20, limbigprime=30):2009r"""2010Returns some points of infinite order on this elliptic curve.2011They are not necessarily linearly independent.20122013Check :meth:`~rank` or2014:meth:`~rank_bounds` to verify the number of generators.20152016.. NOTE::20172018The optional parameters control the Simon two descent algorithm;2019see the documentation of :meth:`~simon_two_descent` for more2020details.20212022INPUT:20232024- ``verbose`` -- 0, 1, 2, or 3 (default: 0), the verbosity level20252026- ``lim1`` -- (default: 5) limit on trivial points on quartics20272028- ``lim3`` -- (default: 50) limit on points on ELS quartics20292030- ``limtriv`` -- (default: 10) limit on trivial points on elliptic curve20312032- ``maxprob`` -- (default: 20)20332034- ``limbigprime`` -- (default: 30) to distinguish between2035small and large prime numbers. Use probabilistic tests for2036large primes. If 0, don't use probabilistic tests.20372038OUTPUT:20392040A set of points of infinite order given by the Simon two-descent.20412042.. NOTE::20432044For non-quadratic number fields, this code does return, but it takes2045a long time.20462047EXAMPLES::20482049sage: K.<a> = NumberField(x^2 + 23, 'a')2050sage: E = EllipticCurve(K, '37')2051sage: E == loads(dumps(E))2052True2053sage: E.gens()2054[(-1 : 0 : 1), (1/2*a - 5/2 : -1/2*a - 13/2 : 1)]205520562057Here is a curve of rank 2, yet the list contains many points::20582059sage: K.<t> = NumberField(x^2-17)2060sage: E = EllipticCurve(K,[-4,0])2061sage: E.gens()2062[(-1/2*t + 1/2 : -1/2*t + 1/2 : 1),2063(-2*t + 8 : -8*t + 32 : 1),2064(1/2*t + 3/2 : -1/2*t - 7/2 : 1),2065(-1/8*t - 7/8 : -1/16*t - 23/16 : 1),2066(1/8*t - 7/8 : -1/16*t + 23/16 : 1),2067(t + 3 : -2*t - 10 : 1),2068(2*t + 8 : -8*t - 32 : 1),2069(1/2*t + 1/2 : -1/2*t - 1/2 : 1),2070(-1/2*t + 3/2 : -1/2*t + 7/2 : 1),2071(t + 7 : -4*t - 20 : 1),2072(-t + 7 : -4*t + 20 : 1),2073(-t + 3 : -2*t + 10 : 1)]2074sage: E.rank()2075220762077Test that the points of finite order are not included :trac: `13593` ::20782079sage: E = EllipticCurve("17a3")2080sage: K.<t> = NumberField(x^2+3)2081sage: EK = E.base_extend(K)2082sage: EK.rank()208302084sage: EK.gens()2085[]208620872088IMPLEMENTATION:20892090Uses Denis Simon's PARI/GP scripts from2091http://www.math.unicaen.fr/~simon/.2092"""20932094lower,upper,gens = self.simon_two_descent(verbose=verbose,lim1=lim1,lim3=lim3,limtriv=limtriv,maxprob=maxprob,limbigprime=limbigprime)2095res = [P for P in gens if P.has_infinite_order()]2096return res209720982099def period_lattice(self, embedding):2100r"""2101Returns the period lattice of the elliptic curve for the given2102embedding of its base field with respect to the differential2103`dx/(2y + a_1x + a_3)`.21042105INPUT:21062107- ``embedding`` - an embedding of the base number field into `\RR` or `\CC`.21082109.. note::21102111The precision of the embedding is ignored: we only use the2112given embedding to determine which embedding into ``QQbar``2113to use. Once the lattice has been initialized, periods can2114be computed to arbitrary precision.211521162117EXAMPLES:21182119First define a field with two real embeddings::21202121sage: K.<a> = NumberField(x^3-2)2122sage: E=EllipticCurve([0,0,0,a,2])2123sage: embs=K.embeddings(CC); len(embs)2124321252126For each embedding we have a different period lattice::21272128sage: E.period_lattice(embs[0])2129Period 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:2130From: Number Field in a with defining polynomial x^3 - 22131To: Algebraic Field2132Defn: a |--> -0.6299605249474365? - 1.091123635971722?*I21332134sage: E.period_lattice(embs[1])2135Period 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:2136From: Number Field in a with defining polynomial x^3 - 22137To: Algebraic Field2138Defn: a |--> -0.6299605249474365? + 1.091123635971722?*I21392140sage: E.period_lattice(embs[2])2141Period 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:2142From: Number Field in a with defining polynomial x^3 - 22143To: Algebraic Field2144Defn: a |--> 1.259921049894873?21452146Although the original embeddings have only the default2147precision, we can obtain the basis with higher precision2148later::21492150sage: L=E.period_lattice(embs[0])2151sage: L.basis()2152(1.86405007647981 - 0.903761485143226*I, -0.149344633143919 - 2.06619546272945*I)21532154sage: L.basis(prec=100)2155(1.8640500764798108425920506200 - 0.90376148514322594749786960975*I, -0.14934463314391922099120107422 - 2.0661954627294548995621225062*I)2156"""2157from sage.schemes.elliptic_curves.period_lattice import PeriodLattice_ell2158return PeriodLattice_ell(self,embedding)21592160def is_isogenous(self, other, proof=True, maxnorm=100):2161"""2162Returns whether or not self is isogenous to other.21632164INPUT:21652166- ``other`` -- another elliptic curve.21672168- ``proof`` (default True) -- If ``False``, the function will2169return ``True`` whenever the two curves have the same2170conductor and are isogenous modulo `p` for all primes `p` of2171norm up to ``maxp``. If ``True``, the function returns2172False when the previous condition does not hold, and if it2173does hold we attempt to see if the curves are indeed2174isogenous. However, this has not been fully implemented2175(see examples below), so we may not be able to determine2176whether or not the curves are isogenous..21772178- ``maxnorm`` (integer, default 100) -- The maximum norm of2179primes `p` for which isogeny modulo `p` will be checked.21802181OUTPUT:21822183(bool) True if there is an isogeny from curve ``self`` to2184curve ``other``.21852186EXAMPLES::21872188sage: x = polygen(QQ, 'x')2189sage: F = NumberField(x^2 -2, 's'); F2190Number Field in s with defining polynomial x^2 - 22191sage: E1 = EllipticCurve(F, [7,8])2192sage: E2 = EllipticCurve(F, [0,5,0,1,0])2193sage: E3 = EllipticCurve(F, [0,-10,0,21,0])2194sage: E1.is_isogenous(E2)2195False2196sage: E1.is_isogenous(E1)2197True2198sage: E2.is_isogenous(E2)2199True2200sage: E2.is_isogenous(E1)2201False2202sage: E2.is_isogenous(E3)2203True22042205::22062207sage: x = polygen(QQ, 'x')2208sage: F = NumberField(x^2 -2, 's'); F2209Number Field in s with defining polynomial x^2 - 22210sage: E = EllipticCurve('14a1')2211sage: EE = EllipticCurve('14a2')2212sage: E1 = E.change_ring(F)2213sage: E2 = EE.change_ring(F)2214sage: E1.is_isogenous(E2)2215True22162217::22182219sage: x = polygen(QQ, 'x')2220sage: F = NumberField(x^2 -2, 's'); F2221Number Field in s with defining polynomial x^2 - 22222sage: k.<a> = NumberField(x^3+7)2223sage: E = EllipticCurve(F, [7,8])2224sage: EE = EllipticCurve(k, [2, 2])2225sage: E.is_isogenous(EE)2226Traceback (most recent call last):2227...2228ValueError: Second argument must be defined over the same number field.22292230Some examples from Cremona's 1981 tables::22312232sage: K.<i> = QuadraticField(-1)2233sage: E1 = EllipticCurve([i + 1, 0, 1, -240*i - 400, -2869*i - 2627])2234sage: E1.conductor()2235Fractional ideal (4*i + 7)2236sage: E2 = EllipticCurve([1+i,0,1,0,0])2237sage: E2.conductor()2238Fractional ideal (4*i + 7)2239sage: E1.is_isogenous(E2) # long time (2s on sage.math, 2014)2240Traceback (most recent call last):2241...2242NotImplementedError: Curves appear to be isogenous (same conductor, isogenous modulo all primes of norm up to 1000), but no isogeny has been constructed.2243sage: E1.is_isogenous(E2, proof=False)2244True22452246In this case E1 and E2 are in fact 9-isogenous, as may be2247deduced from the following::22482249sage: E3 = EllipticCurve([i + 1, 0, 1, -5*i - 5, -2*i - 5])2250sage: E3.is_isogenous(E1)2251True2252sage: E3.is_isogenous(E2)2253True2254sage: E1.isogeny_degree(E2)2255922562257"""2258if not is_EllipticCurve(other):2259raise ValueError, "Second argument is not an Elliptic Curve."2260if self.is_isomorphic(other):2261return True2262K = self.base_field()2263if K != other.base_field():2264raise ValueError, "Second argument must be defined over the same number field."22652266E1 = self.integral_model()2267E2 = other.integral_model()2268N = E1.conductor()2269if N != E2.conductor():2270return False22712272PI = K.primes_of_degree_one_iter()2273while True:2274P = PI.next()2275if P.norm() > maxnorm: break2276if not P.divides(N):2277OP = K.residue_field(P)2278if E1.change_ring(OP).cardinality() != E2.change_ring(OP).cardinality():2279return False22802281if not proof:2282return True22832284# We have not yet implemented isogenies of all possible2285# degrees, and do not know a bound on the possible degrees2286# over general number fields. But here we do at least try2287# some easy cases:22882289for l in [2,3,5,7,13]:2290if any([E2.is_isomorphic(f.codomain()) for f in E1.isogenies_prime_degree(l)]):2291return True22922293# Next we try looking modulo some more primes:22942295while True:2296if P.norm() > 10*maxnorm: break2297if not P.divides(N):2298OP = K.residue_field(P)2299if E1.change_ring(OP).cardinality() != E2.change_ring(OP).cardinality():2300return False2301P = PI.next()23022303# At this point is is highly likely that the curves are2304# isogenous, but we have not proved it.23052306raise 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)23072308def isogeny_degree(self, other):2309"""2310Returns the minimal degree of an isogeny between self and2311other, or 0 if no isogeny exists.23122313INPUT:23142315- ``other`` -- another elliptic curve.23162317OUTPUT:23182319(int) The degree of an isogeny from ``self`` to ``other``, or 0.23202321.. warning::23222323Not all isogenies over number fields are yet implemented.2324Currently the code only works if there is a chain of2325isogenies from ``self`` to ``other`` of degrees 2, 3, 5, 72326and 13.23272328EXAMPLES::23292330sage: x = QQ['x'].02331sage: F = NumberField(x^2 -2, 's'); F2332Number Field in s with defining polynomial x^2 - 22333sage: E = EllipticCurve('14a1')2334sage: EE = EllipticCurve('14a2')2335sage: E1 = E.change_ring(F)2336sage: E2 = EE.change_ring(F)2337sage: E1.isogeny_degree(E2)233822339sage: E2.isogeny_degree(E2)234012341sage: E5 = EllipticCurve('14a5').change_ring(F)2342sage: E1.isogeny_degree(E5)234362344"""2345if self.conductor() != other.conductor():2346return Integer(0)23472348if self.is_isomorphic(other):2349return Integer(1)23502351from sage.sets.set import Set23522353curves = [self]2354degrees = [Integer(1)]2355l_list = [l for l in Set([ZZ(f.degree()) for f in self.isogenies_prime_degree([2,3,5,7,13])])]23562357newcurves = []2358newdegs = []2359k = 02360while k<len(curves):2361newcurves.extend([f.codomain() for f in curves[k].isogenies_prime_degree(l_list)])2362newdegs.extend([degrees[k]*f.degree() for f in curves[k].isogenies_prime_degree(l_list)])2363newisogpairs = dict(zip(newcurves, newdegs))2364i = 02365while i<len(curves):2366j = 02367while j<len(newcurves):2368if curves[i].is_isomorphic(newcurves[j]):2369newdegs.remove(newisogpairs[newcurves[j]])2370newcurves.remove(newcurves[j])2371else:2372j = j+12373i = i+123742375m = 02376newisogpairs = dict(zip(newcurves, newdegs))2377while m<len(newcurves):2378if other.is_isomorphic(newcurves[m]):2379return newisogpairs[newcurves[m]]2380m = m+123812382curves.extend(newcurves)2383degrees.extend(newdegs)2384k = k+123852386raise NotImplementedError("Not all isogenies implemented over general number fields.")23872388def lll_reduce(self, points, height_matrix=None, precision=None):2389"""2390Returns an LLL-reduced basis from a given basis, with transform2391matrix.23922393INPUT:23942395- ``points`` - a list of points on this elliptic2396curve, which should be independent.23972398- ``height_matrix`` - the height-pairing matrix of2399the points, or ``None``. If ``None``, it will be computed.24002401- ``precision`` - number of bits of precision of intermediate2402computations (default: ``None``, for default RealField2403precision; ignored if ``height_matrix`` is supplied)24042405OUTPUT: A tuple (newpoints, U) where U is a unimodular integer2406matrix, new_points is the transform of points by U, such that2407new_points has LLL-reduced height pairing matrix24082409.. note::24102411If the input points are not independent, the output2412depends on the undocumented behaviour of PARI's2413``qflllgram()`` function when applied to a gram matrix which2414is not positive definite.24152416EXAMPLES:24172418Some examples over `\QQ`::24192420sage: E = EllipticCurve([0, 1, 1, -2, 42])2421sage: Pi = E.gens(); Pi2422[(-4 : 1 : 1), (-3 : 5 : 1), (-11/4 : 43/8 : 1), (-2 : 6 : 1)]2423sage: Qi, U = E.lll_reduce(Pi)2424sage: all(sum(U[i,j]*Pi[i] for i in range(4)) == Qi[j] for j in range(4))2425True2426sage: sorted(Qi)2427[(-4 : 1 : 1), (-3 : 5 : 1), (-2 : 6 : 1), (0 : 6 : 1)]2428sage: U.det()242912430sage: E.regulator_of_points(Pi)24314.590880369605732432sage: E.regulator_of_points(Qi)24334.5908803696057424342435::24362437sage: E = EllipticCurve([1,0,1,-120039822036992245303534619191166796374,504224992484910670010801799168082726759443756222911415116])2438sage: xi = [2005024558054813068,\2439-4690836759490453344,\24404700156326649806635,\24416785546256295273860,\24426823803569166584943,\24437788809602110240789,\244427385442304350994620556,\244554284682060285253719/4,\2446-94200235260395075139/25,\2447-3463661055331841724647/576,\2448-6684065934033506970637/676,\2449-956077386192640344198/2209,\2450-27067471797013364392578/2809,\2451-25538866857137199063309/3721,\2452-1026325011760259051894331/108241,\24539351361230729481250627334/1366561,\245410100878635879432897339615/1423249,\245511499655868211022625340735/17522596,\2456110352253665081002517811734/21353641,\2457414280096426033094143668538257/285204544,\245836101712290699828042930087436/4098432361,\245945442463408503524215460183165/5424617104,\2460983886013344700707678587482584/141566320009,\24611124614335716851053281176544216033/152487126016]2462sage: points = [E.lift_x(x) for x in xi]2463sage: newpoints, U = E.lll_reduce(points) # long time (35s on sage.math, 2011)2464sage: [P[0] for P in newpoints] # long time2465[6823803569166584943, 5949539878899294213, 2005024558054813068, 5864879778877955778, 23955263915878682727/4, 5922188321411938518, 5286988283823825378, 175620639884534615751/25, -11451575907286171572, 3502708072571012181, 1500143935183238709184/225, 27180522378120223419/4, -5811874164190604461581/625, 26807786527159569093, 7404442636649562303, 475656155255883588, 265757454726766017891/49, 7272142121019825303, 50628679173833693415/4, 6951643522366348968, 6842515151518070703, 111593750389650846885/16, 2607467890531740394315/9, -1829928525835506297]24662467An example to show the explicit use of the height pairing matrix::24682469sage: E = EllipticCurve([0, 1, 1, -2, 42])2470sage: Pi = E.gens()2471sage: H = E.height_pairing_matrix(Pi,3)2472sage: E.lll_reduce(Pi,height_matrix=H)2473(2474[1 0 0 1]2475[0 1 0 1]2476[0 0 0 1]2477[(-4 : 1 : 1), (-3 : 5 : 1), (-2 : 6 : 1), (1 : -7 : 1)], [0 0 1 1]2478)24792480Some examples over number fields (see :trac:`9411`)::24812482sage: K.<a> = QuadraticField(-23, 'a')2483sage: E = EllipticCurve(K, '37')2484sage: E.lll_reduce(E.gens())2485(2486[1 1]2487[(-1 : 0 : 1), (-2 : 1/2*a - 1/2 : 1)], [0 1]2488)24892490::24912492sage: K.<a> = QuadraticField(-5)2493sage: E = EllipticCurve(K,[0,a])2494sage: points = [E.point([-211/841*a - 6044/841,-209584/24389*a + 53634/24389]),E.point([-17/18*a - 1/9, -109/108*a - 277/108]) ]2495sage: E.lll_reduce(points)2496(2497[(-a + 4 : -3*a + 7 : 1), (-17/18*a - 1/9 : 109/108*a + 277/108 : 1)],2498[ 1 0]2499[ 1 -1]2500)2501"""2502r = len(points)2503if height_matrix is None:2504height_matrix = self.height_pairing_matrix(points, precision)2505U = height_matrix._pari_().lllgram().python()2506new_points = [sum([U[j, i]*points[j] for j in range(r)])2507for i in range(r)]2508return new_points, U25092510def galois_representation(self):2511r"""2512The compatible family of the Galois representation2513attached to this elliptic curve.25142515Given an elliptic curve `E` over a number field `K`2516and a rational prime number `p`, the `p^n`-torsion2517`E[p^n]` points of `E` is a representation of the2518absolute Galois group of `K`. As `n` varies2519we obtain the Tate module `T_p E` which is a2520a representation of `G_K` on a free `\ZZ_p`-module2521of rank `2`. As `p` varies the representations2522are compatible.25232524EXAMPLES::25252526sage: K = NumberField(x**2 + 1, 'a')2527sage: E = EllipticCurve('11a1').change_ring(K)2528sage: rho = E.galois_representation()2529sage: rho2530Compatible family of Galois representations associated to the Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 + (-10)*x + (-20) over Number Field in a with defining polynomial x^2 + 12531sage: rho.is_surjective(3)2532True2533sage: rho.is_surjective(5) # long time (4s on sage.math, 2014)2534False2535sage: rho.non_surjective()2536[5]2537"""2538return gal_reps_number_field.GaloisRepresentation(self)253925402541