Path: blob/master/src/sage/schemes/elliptic_curves/lseries_ell.py
8820 views
"""1Complex Elliptic Curve L-series23AUTHORS:45- Jeroen Demeyer (2013-10-17): compute L series with arbitrary precision6instead of floats.78- William Stein et al. (2005 and later)910"""11#*****************************************************************************12# Copyright (C) 2005 William Stein13# Copyright (C) 2013 Jeroen Demeyer14#15# Distributed under the terms of the GNU General Public License (GPL)16# as published by the Free Software Foundation; either version 2 of17# the License, or (at your option) any later version.18# http://www.gnu.org/licenses/19#*****************************************************************************2021from sage.structure.sage_object import SageObject22from sage.rings.all import RealField, RationalField23from math import sqrt, exp, log, ceil24import sage.functions.exp_integral as exp_integral25import sage.misc.all as misc2627class Lseries_ell(SageObject):28"""29An elliptic curve $L$-series.30"""31def __init__(self, E):32"""33Create an elliptic curve $L$-series.3435EXAMPLES::3637sage: EllipticCurve([1..5]).lseries()38Complex L-series of the Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field39"""40self.__E = E4142def elliptic_curve(self):43"""44Return the elliptic curve that this L-series is attached to.4546EXAMPLES::4748sage: E = EllipticCurve('389a')49sage: L = E.lseries()50sage: L.elliptic_curve ()51Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field52"""53return self.__E5455def taylor_series(self, a=1, prec=53, series_prec=6, var='z'):56"""57Return the Taylor series of this $L$-series about $a$ to58the given precision (in bits) and the number of terms.5960The output is a series in var, where you should view var as61equal to s-a. Thus this function returns the formal power62series whose coefficients are L^{(n)}(a)/n!.6364EXAMPLES::6566sage: E = EllipticCurve('389a')67sage: L = E.lseries()68sage: L.taylor_series(series_prec=3)69-1.28158145675273e-23 + (7.26268290541182e-24)*z + 0.759316500288427*z^2 + O(z^3) # 32-bit70-2.69129566562797e-23 + (1.52514901968783e-23)*z + 0.759316500288427*z^2 + O(z^3) # 64-bit71sage: L.taylor_series(series_prec=3)[2:]720.000000000000000 + 0.000000000000000*z + 0.759316500288427*z^2 + O(z^3)73"""74D = self.dokchitser(prec)75return D.taylor_series(a, series_prec, var)7677def _repr_(self):78"""79Return string representation of this L-series.8081EXAMPLES::8283sage: E = EllipticCurve('37a')84sage: L = E.lseries()85sage: L._repr_()86'Complex L-series of the Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field'87"""88return "Complex L-series of the %s"%self.__E8990def dokchitser(self, prec=53,91max_imaginary_part=0,92max_asymp_coeffs=40,93algorithm='gp'):94r"""95Return interface to Tim Dokchitser's program for computing96with the L-series of this elliptic curve; this provides a way97to compute Taylor expansions and higher derivatives of98$L$-series.99100INPUT:101prec -- integer (bits precision)102max_imaginary_part -- real number103max_asymp_coeffs -- integer104algorithm -- string: 'gp' or 'magma'105106\note{If algorithm='magma', then the precision is in digits rather107than bits and the object returned is a Magma L-series, which has108different functionality from the Sage L-series.}109110EXAMPLES::111112sage: E = EllipticCurve('37a')113sage: L = E.lseries().dokchitser()114sage: L(2)1150.381575408260711116sage: L = E.lseries().dokchitser(algorithm='magma') # optional - magma117sage: L.Evaluate(2) # optional - magma1180.38157540826071121129371040958008663667709753398892116119120If the curve has too large a conductor, it isn't possible to121compute with the L-series using this command. Instead a122RuntimeError is raised::123124sage: e = EllipticCurve([1,1,0,-63900,-1964465932632])125sage: L = e.lseries().dokchitser(15)126Traceback (most recent call last):127...128RuntimeError: Unable to create L-series, due to precision or other limits in PARI.129"""130if algorithm == 'magma':131from sage.interfaces.all import magma132return magma(self.__E).LSeries(Precision = prec)133134from sage.lfunctions.all import Dokchitser135key = (prec, max_imaginary_part, max_asymp_coeffs)136try:137return self.__dokchitser[key]138except KeyError:139pass140except AttributeError:141self.__dokchitser = {}142L = Dokchitser(conductor = self.__E.conductor(),143gammaV = [0,1],144weight = 2,145eps = self.__E.root_number(),146poles = [],147prec = prec)148gp = L.gp()149s = 'e = ellinit(%s);'%list(self.__E.minimal_model().a_invariants())150s += 'a(k) = ellak(e, k);'151L.init_coeffs('a(k)', 1, pari_precode = s,152max_imaginary_part=max_imaginary_part,153max_asymp_coeffs=max_asymp_coeffs)154L.rename('Dokchitser L-function associated to %s'%self.__E)155self.__dokchitser[key] = L156return L157158def sympow(self, n, prec):159r"""160Return $L(\Sym^{(n)}(E, \text{edge}))$ to prec digits161of precision.162163INPUT:164n -- integer165prec -- integer166167OUTPUT:168string -- real number to prec digits of precision as a string.169170\note{Before using this function for the first time for171a given $n$, you may have to type \code{sympow('-new_data <n>')},172where \code{<n>} is replaced by your value of $n$. This173command takes a long time to run.}174175EXAMPLES::176177sage: E = EllipticCurve('37a')178sage: a = E.lseries().sympow(2,16) # not tested - requires precomputing "sympow('-new_data 2')"179sage: a # not tested180'2.492262044273650E+00'181sage: RR(a) # not tested1822.49226204427365183"""184from sage.lfunctions.sympow import sympow185return sympow.L(self.__E, n, prec)186187def sympow_derivs(self, n, prec, d):188r"""189Return $0$th to $d$th derivatives of $L(\Sym^{(n)}(E,190\text{edge}))$ to prec digits of precision.191192INPUT:193n -- integer194prec -- integer195d -- integer196197OUTPUT:198a string, exactly as output by sympow199200\note{To use this function you may have to run a few commands201like \code{sympow('-new_data 1d2')}, each which takes a few202minutes. If this function fails it will indicate what203commands have to be run.}204205EXAMPLES::206207sage: E = EllipticCurve('37a')208sage: print E.lseries().sympow_derivs(1,16,2) # not tested -- requires precomputing "sympow('-new_data 2')"209sympow 1.018 RELEASE (c) Mark Watkins --- see README and COPYING for details210Minimal model of curve is [0,0,1,-1,0]211At 37: Inertia Group is C1 MULTIPLICATIVE REDUCTION212Conductor is 37213sp 1: Conductor at 37 is 1+0, root number is 1214sp 1: Euler factor at 37 is 1+1*x2151st sym power conductor is 37, global root number is -1216NT 1d0: 35217NT 1d1: 32218NT 1d2: 28219Maximal number of terms is 35220Done with small primes 1049221Computed: 1d0 1d1 1d2222Checked out: 1d12231n0: 3.837774351482055E-012241w0: 3.777214305638848E-012251n1: 3.059997738340522E-012261w1: 3.059997738340524E-012271n2: 1.519054910249753E-012281w2: 1.545605024269432E-01229"""230from sage.lfunctions.sympow import sympow231return sympow.Lderivs(self.__E, n, prec, d)232233def zeros(self, n):234"""235Return the imaginary parts of the first $n$ nontrivial zeros236on the critical line of the L-function in the upper half237plane, as 32-bit reals.238239EXAMPLES::240241sage: E = EllipticCurve('37a')242sage: E.lseries().zeros(2)243[0.000000000, 5.00317001]244245sage: a = E.lseries().zeros(20) # long time246sage: point([(1,x) for x in a]) # graph (long time)247248AUTHOR:249-- Uses Rubinstein's L-functions calculator.250"""251from sage.lfunctions.lcalc import lcalc252return lcalc.zeros(n, L=self.__E)253254def zeros_in_interval(self, x, y, stepsize):255r"""256Return the imaginary parts of (most of) the nontrivial zeros257on the critical line $\Re(s)=1$ with positive imaginary part258between $x$ and $y$, along with a technical quantity for each.259260INPUT:261x, y, stepsize -- positive floating point numbers262263OUTPUT:264list of pairs (zero, S(T)).265266Rubinstein writes: The first column outputs the imaginary part267of the zero, the second column a quantity related to S(T) (it268increases roughly by 2 whenever a sign change, i.e. pair of269zeros, is missed). Higher up the critical strip you should use270a smaller stepsize so as not to miss zeros.271272EXAMPLES::273274sage: E = EllipticCurve('37a')275sage: E.lseries().zeros_in_interval(6, 10, 0.1) # long time276[(6.87039122, 0.248922780), (8.01433081, -0.140168533), (9.93309835, -0.129943029)]277"""278from sage.lfunctions.lcalc import lcalc279return lcalc.zeros_in_interval(x, y, stepsize, L=self.__E)280281def values_along_line(self, s0, s1, number_samples):282"""283Return values of $L(E, s)$ at \code{number_samples}284equally-spaced sample points along the line from $s_0$ to285$s_1$ in the complex plane.286287\note{The L-series is normalized so that the center of the288critical strip is 1.}289290INPUT:291s0, s1 -- complex numbers292number_samples -- integer293294OUTPUT:295list -- list of pairs (s, zeta(s)), where the s are296equally spaced sampled points on the line from297s0 to s1.298299EXAMPLES::300301sage: E = EllipticCurve('37a')302sage: E.lseries().values_along_line(1, 0.5 + 20*I, 5)303[(0.500000000, ...),304(0.400000000 + 4.00000000*I, 3.31920245 - 2.60028054*I),305(0.300000000 + 8.00000000*I, -0.886341185 - 0.422640337*I),306(0.200000000 + 12.0000000*I, -3.50558936 - 0.108531690*I),307(0.100000000 + 16.0000000*I, -3.87043288 - 1.88049411*I)]308309"""310from sage.lfunctions.lcalc import lcalc311return lcalc.values_along_line(s0-RationalField()('1/2'),312s1-RationalField()('1/2'),313number_samples, L=self.__E)314315def twist_values(self, s, dmin, dmax):316r"""317Return values of $L(E, s, \chi_d)$ for each quadratic318character $\chi_d$ for $d_{\min} \leq d \leq d_{\max}$.319320\note{The L-series is normalized so that the center of the321critical strip is 1.}322323INPUT:324325- ``s`` -- complex numbers326- ``dmin`` -- integer327- ``dmax`` -- integer328329OUTPUT:330331list of pairs (d, L(E, s,chi_d))332333EXAMPLES::334335sage: E = EllipticCurve('37a')336sage: vals = E.lseries().twist_values(1, -12, -4)337sage: vals # abs tol 1e-17338[(-11, 1.47824342), (-8, 8.9590946e-18), (-7, 1.85307619), (-4, 2.45138938)]339sage: F = E.quadratic_twist(-8)340sage: F.rank()3411342sage: F = E.quadratic_twist(-7)343sage: F.rank()3440345"""346from sage.lfunctions.lcalc import lcalc347return lcalc.twist_values(s - RationalField()('1/2'), dmin, dmax, L=self.__E)348349def twist_zeros(self, n, dmin, dmax):350r"""351Return first $n$ real parts of nontrivial zeros of352$L(E,s,\chi_d)$ for each quadratic character $\chi_d$ with353$d_{\min} \leq d \leq d_{\max}$.354355\note{The L-series is normalized so that the center of the356critical strip is 1.}357358INPUT:359n -- integer360dmin -- integer361dmax -- integer362363OUTPUT:364dict -- keys are the discriminants $d$, and365values are list of corresponding zeros.366367EXAMPLES::368369sage: E = EllipticCurve('37a')370sage: E.lseries().twist_zeros(3, -4, -3) # long time371{-4: [1.60813783, 2.96144840, 3.89751747], -3: [2.06170900, 3.48216881, 4.45853219]}372"""373from sage.lfunctions.lcalc import lcalc374return lcalc.twist_zeros(n, dmin, dmax, L=self.__E)375376def at1(self, k=None, prec=None):377r"""378Compute `L(E,1)` using `k` terms of the series for `L(E,1)` as379explained in Section 7.5.3 of Henri Cohen's book "A Course in380Computational Algebraic Number Theory". If the argument `k`381is not specified, then it defaults to `\sqrt(N)`, where `N` is382the conductor.383384INPUT:385386- ``k`` -- number of terms of the series. If zero or ``None``,387use `k = \sqrt(N)`, where `N` is the conductor.388389- ``prec`` -- numerical precision in bits. If zero or ``None``,390use a reasonable automatic default.391392OUTPUT:393394A tuple of real numbers ``(L, err)`` where ``L`` is an395approximation for `L(E,1)` and ``err`` is a bound on the error396in the approximation.397398This function is disjoint from the PARI ``elllseries``399command, which is for a similar purpose. To use that command400(via the PARI C library), simply type401``E.pari_mincurve().elllseries(1)``.402403ALGORITHM:404405- Compute the root number eps. If it is -1, return 0.406407- Compute the Fourier coefficients `a_n`, for `n` up to and408including `k`.409410- Compute the sum411412.. MATH::4134142 * sum_{n=1}^{k} (a_n / n) * exp(-2*pi*n/Sqrt(N)),415416where `N` is the conductor of `E`.417418- Compute a bound on the tail end of the series, which is419420.. MATH::4214222 e^{-2 \pi (k+1) / \sqrt{N}} / (1 - e^{-2 \pi/\sqrt{N}}).423424For a proof see [Grigov-Jorza-Patrascu-Patrikis-Stein].425426EXAMPLES::427428sage: L, err = EllipticCurve('11a1').lseries().at1()429sage: L, err430(0.253804, 0.000181444)431sage: parent(L)432Real Field with 24 bits of precision433sage: E = EllipticCurve('37b')434sage: E.lseries().at1()435(0.7257177, 0.000800697)436sage: E.lseries().at1(100)437(0.7256810619361527823362055410263965487367603361763, 1.52469e-45)438sage: L,err = E.lseries().at1(100, prec=128)439sage: L4400.72568106193615278233620554102639654873441sage: parent(L)442Real Field with 128 bits of precision443sage: err4441.70693e-37445sage: parent(err)446Real Field with 24 bits of precision and rounding RNDU447448Rank 1 through 3 elliptic curves::449450sage: E = EllipticCurve('37a1')451sage: E.lseries().at1()452(0.0000000, 0.000000)453sage: E = EllipticCurve('389a1')454sage: E.lseries().at1()455(-0.001769566, 0.00911776)456sage: E = EllipticCurve('5077a1')457sage: E.lseries().at1()458(0.0000000, 0.000000)459"""460sqrtN = sqrt(self.__E.conductor())461if k:462k = int(k)463else:464k = int(ceil(sqrtN))465466if prec:467prec = int(prec)468else:469# Use the same precision as deriv_at1() below for470# consistency471prec = int(9.065*k/sqrtN + 1.443*log(k)) + 12472R = RealField(prec)473# Compute error term with bounded precision of 24 bits and474# round towards +infinity475Rerror = RealField(24, rnd='RNDU')476477if self.__E.root_number() == -1:478return (R.zero(), Rerror.zero())479480an = self.__E.anlist(k) # list of Sage Integers481pi = R.pi()482sqrtN = R(self.__E.conductor()).sqrt()483484z = (-2*pi/sqrtN).exp()485zpow = z486# Compute series sum and accumulate floating point errors487L = R.zero()488error = Rerror.zero()489490for n in xrange(1,k+1):491term = (zpow * an[n])/n492zpow *= z493L += term494# We express relative error in units of epsilon, where495# epsilon is a number divided by 2^precision.496# Instead of multiplying the error by 2 after the loop497# (to account for L *= 2), we already multiply it now.498#499# For multiplication and division, the relative error500# in epsilons is bounded by (1+e)^n - 1, where n is the501# number of operations (assuming exact inputs).502# exp(x) additionally multiplies this error by abs(x) and503# adds one epsilon. The inputs pi and sqrtN each contribute504# another epsilon.505# Assuming that 2*pi/sqrtN <= 2, the relative error for z is506# 7 epsilon. This implies a relative error of (8n-1) epsilon507# for zpow. We add 2 for the computation of term and 1/2 to508# compensate for the approximation (1+e)^n = 1+ne.509#510# The error of the addition is at most half an ulp of the511# result.512#513# Multiplying everything by two gives:514error += term.epsilon(Rerror)*(16*n + 3) + L.ulp(Rerror)515L *= 2516517# Add series error (we use (-2)/(z-1) instead of 2/(1-z)518# because this causes 1/(1-z) to be rounded up)519error += ((-2)*Rerror(zpow)) / Rerror(z - 1)520return (L, error)521522def deriv_at1(self, k=None, prec=None):523r"""524Compute `L'(E,1)` using `k` terms of the series for `L'(E,1)`,525under the assumption that `L(E,1) = 0`.526527The algorithm used is from Section 7.5.3 of Henri Cohen's book528``A Course in Computational Algebraic Number Theory.''529530INPUT:531532- ``k`` -- number of terms of the series. If zero or ``None``,533use `k = \sqrt(N)`, where `N` is the conductor.534535- ``prec`` -- numerical precision in bits. If zero or ``None``,536use a reasonable automatic default.537538OUTPUT:539540A tuple of real numbers ``(L1, err)`` where ``L1`` is an541approximation for `L'(E,1)` and ``err`` is a bound on the error542in the approximation.543544.. WARNING::545546This function only makes sense if `L(E)` has positive order547of vanishing at 1, or equivalently if `L(E,1) = 0`.548549ALGORITHM:550551- Compute the root number eps. If it is 1, return 0.552553- Compute the Fourier coefficients `a_n`, for `n` up to and554including `k`.555556- Compute the sum557558.. MATH::5595602 * \sum_{n=1}^{k} (a_n / n) * E_1(2 \pi n/\sqrt{N}),561562where `N` is the conductor of `E`, and `E_1` is the563exponential integral function.564565- Compute a bound on the tail end of the series, which is566567.. MATH::5685692 e^{-2 \pi (k+1) / \sqrt{N}} / (1 - e^{-2 \pi/\sqrt{N}}).570571For a proof see [Grigorov-Jorza-Patrascu-Patrikis-Stein]. This572is exactly the same as the bound for the approximation to573`L(E,1)` produced by :meth:`at1`.574575EXAMPLES::576577sage: E = EllipticCurve('37a')578sage: E.lseries().deriv_at1()579(0.3059866, 0.000801045)580sage: E.lseries().deriv_at1(100)581(0.3059997738340523018204836833216764744526377745903, 1.52493e-45)582sage: E.lseries().deriv_at1(1000)583(0.305999773834052301820483683321676474452637774590771998..., 2.75031e-449)584585With less numerical precision, the error is bounded by numerical accuracy::586587sage: L,err = E.lseries().deriv_at1(100, prec=64)588sage: L,err589(0.305999773834052302, 5.55318e-18)590sage: parent(L)591Real Field with 64 bits of precision592sage: parent(err)593Real Field with 24 bits of precision and rounding RNDU594595Rank 2 and rank 3 elliptic curves::596597sage: E = EllipticCurve('389a1')598sage: E.lseries().deriv_at1()599(0.0000000, 0.000000)600sage: E = EllipticCurve((1, 0, 1, -131, 558)) # curve 59450i1601sage: E.lseries().deriv_at1()602(-0.00010911444, 0.142428)603sage: E.lseries().deriv_at1(4000)604(6.9902290...e-50, 1.31318e-43)605"""606sqrtN = sqrt(self.__E.conductor())607if k:608k = int(k)609else:610k = int(ceil(sqrtN))611612if prec:613prec = int(prec)614else:615# Estimate number of bits for the computation, based on error616# estimate below (the denominator of that error is close enough617# to 1 that we can ignore it).618# 9.065 = 2*Pi/log(2)619# 1.443 = 1/log(2)620# 12 is an arbitrary extra number of bits (it is chosen621# such that the precision is 24 bits when the conductor622# equals 11 and k is the default value 4)623prec = int(9.065*k/sqrtN + 1.443*log(k)) + 12624R = RealField(prec)625# Compute error term with bounded precision of 24 bits and626# round towards +infinity627Rerror = RealField(24, rnd='RNDU')628629if self.__E.root_number() == 1:630# Order of vanishing at 1 of L(E) is even and assumed to be631# positive, so L'(E,1) = 0.632return (R.zero(), Rerror.zero())633634an = self.__E.anlist(k) # list of Sage Integers635pi = R.pi()636sqrtN = R(self.__E.conductor()).sqrt()637v = exp_integral.exponential_integral_1(2*pi/sqrtN, k)638639# Compute series sum and accumulate floating point errors640L = R.zero()641error = Rerror.zero()642# Sum of |an[n]|/n643sumann = Rerror.zero()644645for n in xrange(1,k+1):646term = (v[n-1] * an[n])/n647L += term648error += term.epsilon(Rerror)*5 + L.ulp(Rerror)649sumann += Rerror(an[n].abs())/n650L *= 2651652# Add error term for exponential_integral_1() errors.653# Absolute error for 2*v[i] is 4*max(1, v[0])*2^-prec654if v[0] > 1.0:655sumann *= Rerror(v[0])656error += (sumann >> (prec - 2))657658# Add series error (we use (-2)/(z-1) instead of 2/(1-z)659# because this causes 1/(1-z) to be rounded up)660z = (-2*pi/sqrtN).exp()661zpow = ((-2*(k+1))*pi/sqrtN).exp()662error += ((-2)*Rerror(zpow)) / Rerror(z - 1)663return (L, error)664665def __call__(self, s):666r"""667Returns the value of the L-series of the elliptic curve E at s, where s668must be a real number.669670.. NOTE::671672If the conductor of the curve is large, say `>10^{12}`,673then this function will take a very long time, since it674uses an `O(\sqrt{N})` algorithm.675676EXAMPLES::677678sage: E = EllipticCurve([1,2,3,4,5])679sage: L = E.lseries()680sage: L(1)6810.000000000000000682sage: L(1.1)6830.285491007678148684sage: L(1.1 + I)6850.174851377216615 + 0.816965038124457*I686"""687return self.dokchitser()(s)688689def L1_vanishes(self):690"""691Returns whether or not `L(E,1) = 0`. The result is provably692correct if the Manin constant of the associated optimal693quotient is <= 2. This hypothesis on the Manin constant694is true for all curves of conductor <= 40000 (by Cremona) and695all semistable curves (i.e., squarefree conductor).696697ALGORITHM: see :meth:`L_ratio`.698699EXAMPLES::700701sage: E = EllipticCurve([0, -1, 1, -10, -20]) # 11A = X_0(11)702sage: E.lseries().L1_vanishes()703False704sage: E = EllipticCurve([0, -1, 1, 0, 0]) # X_1(11)705sage: E.lseries().L1_vanishes()706False707sage: E = EllipticCurve([0, 0, 1, -1, 0]) # 37A (rank 1)708sage: E.lseries().L1_vanishes()709True710sage: E = EllipticCurve([0, 1, 1, -2, 0]) # 389A (rank 2)711sage: E.lseries().L1_vanishes()712True713sage: E = EllipticCurve([0, 0, 1, -38, 90]) # 361A (CM curve))714sage: E.lseries().L1_vanishes()715True716sage: E = EllipticCurve([0,-1,1,-2,-1]) # 141C (13-isogeny)717sage: E.lseries().L1_vanishes()718False719720AUTHOR: William Stein, 2005-04-20.721"""722return self.L_ratio() == 0723724def L_ratio(self):725r"""726Returns the ratio `L(E,1)/\Omega` as an exact rational727number. The result is \emph{provably} correct if the Manin728constant of the associated optimal quotient is `\leq 2`. This729hypothesis on the Manin constant is true for all semistable730curves (i.e., squarefree conductor), by a theorem of Mazur731from his \emph{Rational Isogenies of Prime Degree} paper.732733EXAMPLES::734735sage: E = EllipticCurve([0, -1, 1, -10, -20]) # 11A = X_0(11)736sage: E.lseries().L_ratio()7371/5738sage: E = EllipticCurve([0, -1, 1, 0, 0]) # X_1(11)739sage: E.lseries().L_ratio()7401/25741sage: E = EllipticCurve([0, 0, 1, -1, 0]) # 37A (rank 1)742sage: E.lseries().L_ratio()7430744sage: E = EllipticCurve([0, 1, 1, -2, 0]) # 389A (rank 2)745sage: E.lseries().L_ratio()7460747sage: E = EllipticCurve([0, 0, 1, -38, 90]) # 361A (CM curve))748sage: E.lseries().L_ratio()7490750sage: E = EllipticCurve([0,-1,1,-2,-1]) # 141C (13-isogeny)751sage: E.lseries().L_ratio()7521753sage: E = EllipticCurve(RationalField(), [1, 0, 0, 1/24624, 1/886464])754sage: E.lseries().L_ratio()7552756757See :trac:`3651` and :trac:`15299`::758759sage: EllipticCurve([0,0,0,-193^2,0]).sha().an()7604761sage: EllipticCurve([1, 0, 1, -131, 558]).sha().an() # long time7621.00000000000000763764ALGORITHM: Compute the root number. If it is -1 then L(E,s)765vanishes to odd order at 1, hence vanishes. If it is +1, use766a result about modular symbols and Mazur's "Rational Isogenies"767paper to determine a provably correct bound (assuming Manin768constant is <= 2) so that we can determine whether L(E,1) = 0.769770AUTHOR: William Stein, 2005-04-20.771"""772try:773return self.__lratio774except AttributeError:775pass776777if not self.__E.is_minimal():778self.__lratio = self.__E.minimal_model().lseries().L_ratio()779return self.__lratio780781QQ = RationalField()782if self.__E.root_number() == -1:783self.__lratio = QQ.zero()784return self.__lratio785786# Even root number. Decide if L(E,1) = 0. If E is a modular787# *OPTIMAL* quotient of J_0(N) elliptic curve, we know that T *788# L(E,1)/omega is an integer n, where T is the order of the789# image of the rational torsion point (0)-(oo) in E(Q), and790# omega is the least real Neron period. (This is proved in my791# Ph.D. thesis, but is probably well known.) We can easily792# compute omega to very high precision using AGM. So to prove793# that L(E,1) = 0 we compute T/omega * L(E,1) to sufficient794# precision to determine it as an integer. If eps is the795# error in computation of L(E,1), then the error in computing796# the product is (2T/Omega_E) * eps, and we need this to be797# less than 0.5, i.e.,798# (2T/Omega_E) * eps < 0.5,799# so800# eps < 0.5 * Omega_E / (2T) = Omega_E / (4*T).801#802# Since in general E need not be optimal, we have to choose803# eps = Omega_E/(8*t*B), where t is the exponent of E(Q)_tor,804# and is a multiple of the degree of an isogeny between E805# and the optimal curve.806#807# NOTES: We *do* have to worry about the Manin constant, since808# we are using the Neron model to compute omega, not the809# newform. My theorem replaces the omega above by omega/c,810# where c is the Manin constant, and the bound must be811# correspondingly smaller. If the level is square free, then812# the Manin constant is 1 or 2, so there's no problem (since813# we took 8 instead of 4 in the denominator). If the level814# is divisible by a square, then the Manin constant could815# be a divisible by an arbitrary power of that prime, except816# that Edixhoven claims the primes that appear are <= 7.817818t = self.__E.torsion_subgroup().order()819omega = self.__E.period_lattice().basis()[0]820d = self.__E._multiple_of_degree_of_isogeny_to_optimal_curve()821C = 8*d*t822eps = omega / C823824sqrtN = 2*int(sqrt(self.__E.conductor()))825k = sqrtN + 10826while True:827L1, error_bound = self.at1(k)828if error_bound < eps:829n = int(round(L1*C/omega))830quo = QQ((n,C))831self.__lratio = quo / self.__E.real_components()832return self.__lratio833k += sqrtN834misc.verbose("Increasing precision to %s terms."%k)835836837