Path: blob/master/sage/schemes/elliptic_curves/lseries_ell.py
4156 views
"""1Complex Elliptic Curve L-series23"""45from sage.structure.sage_object import SageObject6from sage.rings.all import (7RealField,8RationalField,9ComplexField)10from math import sqrt, exp, ceil11import sage.functions.transcendental as transcendental12R = RealField()13Q = RationalField()14C = ComplexField()15import sage.misc.all as misc1617class Lseries_ell(SageObject):18"""19An elliptic curve $L$-series.2021EXAMPLES:2223"""24def __init__(self, E):25"""26Create an elliptic curve $L$-series.2728EXAMPLES:29sage: EllipticCurve([1..5]).lseries()30Complex L-series of the Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field31"""32self.__E = E3334def elliptic_curve(self):35"""36Return the elliptic curve that this L-series is attached to.3738EXAMPLES:39sage: E = EllipticCurve('389a')40sage: L = E.lseries()41sage: L.elliptic_curve ()42Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field43"""44return self.__E4546def taylor_series(self, a=1, prec=53, series_prec=6, var='z'):47"""48Return the Taylor series of this $L$-series about $a$ to49the given precision (in bits) and the number of terms.5051The output is a series in var, where you should view var as52equal to s-a. Thus this function returns the formal power53series whose coefficients are L^{(n)}(a)/n!.5455EXAMPLES:56sage: E = EllipticCurve('389a')57sage: L = E.lseries()58sage: L.taylor_series(series_prec=3) # random nearly 0 constant and linear terms59-2.69129566562797e-23 + (1.52514901968783e-23)*z + 0.759316500288427*z^2 + O(z^3)60sage: L.taylor_series(series_prec=3)[2:]610.000000000000000 + 0.000000000000000*z + 0.759316500288427*z^2 + O(z^3)62"""63D = self.dokchitser(prec)64return D.taylor_series(a, series_prec, var)6566def _repr_(self):67"""68Return string representation of this L-series.6970EXAMPLES:71sage: E = EllipticCurve('37a')72sage: L = E.lseries()73sage: L._repr_()74'Complex L-series of the Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field'75"""76return "Complex L-series of the %s"%self.__E7778def dokchitser(self, prec=53,79max_imaginary_part=0,80max_asymp_coeffs=40,81algorithm='gp'):82r"""83Return interface to Tim Dokchitser's program for computing84with the L-series of this elliptic curve; this provides a way85to compute Taylor expansions and higher derivatives of86$L$-series.8788INPUT:89prec -- integer (bits precision)90max_imaginary_part -- real number91max_asymp_coeffs -- integer92algorithm -- string: 'gp' or 'magma'9394\note{If algorithm='magma', then the precision is in digits rather95than bits and the object returned is a Magma L-series, which has96different functionality from the Sage L-series.}9798EXAMPLES:99sage: E = EllipticCurve('37a')100sage: L = E.lseries().dokchitser()101sage: L(2)1020.381575408260711103sage: L = E.lseries().dokchitser(algorithm='magma') # optional - magma104sage: L.Evaluate(2) # optional - magma1050.38157540826071121129371040958008663667709753398892116106107If the curve has too large a conductor, it isn't possible to108compute with the L-series using this command. Instead a109RuntimeError is raised:110sage: e = EllipticCurve([1,1,0,-63900,-1964465932632])111sage: L = e.lseries().dokchitser(15)112Traceback (most recent call last):113...114RuntimeError: Unable to create L-series, due to precision or other limits in PARI.115"""116if algorithm == 'magma':117from sage.interfaces.all import magma118return magma(self.__E).LSeries(Precision = prec)119120from sage.lfunctions.all import Dokchitser121key = (prec, max_imaginary_part, max_asymp_coeffs)122try:123return self.__dokchitser[key]124except KeyError:125pass126except AttributeError:127self.__dokchitser = {}128L = Dokchitser(conductor = self.__E.conductor(),129gammaV = [0,1],130weight = 2,131eps = self.__E.root_number(),132poles = [],133prec = prec)134gp = L.gp()135s = 'e = ellinit(%s);'%list(self.__E.minimal_model().a_invariants())136s += 'a(k) = ellak(e, k);'137L.init_coeffs('a(k)', 1, pari_precode = s,138max_imaginary_part=max_imaginary_part,139max_asymp_coeffs=max_asymp_coeffs)140L.rename('Dokchitser L-function associated to %s'%self.__E)141self.__dokchitser[key] = L142return L143144def sympow(self, n, prec):145r"""146Return $L(\Sym^{(n)}(E, \text{edge}))$ to prec digits147of precision.148149INPUT:150n -- integer151prec -- integer152153OUTPUT:154string -- real number to prec digits of precision as a string.155156\note{Before using this function for the first time for157a given $n$, you may have to type \code{sympow('-new_data <n>')},158where \code{<n>} is replaced by your value of $n$. This159command takes a long time to run.}160161EXAMPLES:162sage: E = EllipticCurve('37a')163sage: a = E.lseries().sympow(2,16) # optional - requires precomputing "sympow('-new_data 2')"164sage: a # optional165'2.492262044273650E+00'166sage: RR(a) # optional1672.49226204427365168"""169from sage.lfunctions.sympow import sympow170return sympow.L(self.__E, n, prec)171172def sympow_derivs(self, n, prec, d):173r"""174Return $0$th to $d$th derivatives of $L(\Sym^{(n)}(E,175\text{edge}))$ to prec digits of precision.176177INPUT:178n -- integer179prec -- integer180d -- integer181182OUTPUT:183a string, exactly as output by sympow184185\note{To use this function you may have to run a few commands186like \code{sympow('-new_data 1d2')}, each which takes a few187minutes. If this function fails it will indicate what188commands have to be run.}189190EXAMPLES:191sage: E = EllipticCurve('37a')192sage: print E.lseries().sympow_derivs(1,16,2) # optional -- requires precomputing "sympow('-new_data 2')"193sympow 1.018 RELEASE (c) Mark Watkins --- see README and COPYING for details194Minimal model of curve is [0,0,1,-1,0]195At 37: Inertia Group is C1 MULTIPLICATIVE REDUCTION196Conductor is 37197sp 1: Conductor at 37 is 1+0, root number is 1198sp 1: Euler factor at 37 is 1+1*x1991st sym power conductor is 37, global root number is -1200NT 1d0: 35201NT 1d1: 32202NT 1d2: 28203Maximal number of terms is 35204Done with small primes 1049205Computed: 1d0 1d1 1d2206Checked out: 1d12071n0: 3.837774351482055E-012081w0: 3.777214305638848E-012091n1: 3.059997738340522E-012101w1: 3.059997738340524E-012111n2: 1.519054910249753E-012121w2: 1.545605024269432E-01213"""214from sage.lfunctions.sympow import sympow215return sympow.Lderivs(self.__E, n, prec, d)216217def zeros(self, n):218"""219Return the imaginary parts of the first $n$ nontrivial zeros220on the critical line of the L-function in the upper half221plane, as 32-bit reals.222223EXAMPLES:224sage: E = EllipticCurve('37a')225sage: E.lseries().zeros(2)226[0.000000000, 5.00317001]227228sage: a = E.lseries().zeros(20) # long time229sage: point([(1,x) for x in a]) # graph (long time)230231AUTHOR:232-- Uses Rubinstein's L-functions calculator.233"""234from sage.lfunctions.lcalc import lcalc235return lcalc.zeros(n, L=self.__E)236237def zeros_in_interval(self, x, y, stepsize):238r"""239Return the imaginary parts of (most of) the nontrivial zeros240on the critical line $\Re(s)=1$ with positive imaginary part241between $x$ and $y$, along with a technical quantity for each.242243INPUT:244x, y, stepsize -- positive floating point numbers245246OUTPUT:247list of pairs (zero, S(T)).248249Rubinstein writes: The first column outputs the imaginary part250of the zero, the second column a quantity related to S(T) (it251increases roughly by 2 whenever a sign change, i.e. pair of252zeros, is missed). Higher up the critical strip you should use253a smaller stepsize so as not to miss zeros.254255EXAMPLES:256sage: E = EllipticCurve('37a')257sage: E.lseries().zeros_in_interval(6, 10, 0.1) # long time258[(6.87039122, 0.248922780), (8.01433081, -0.140168533), (9.93309835, -0.129943029)]259"""260from sage.lfunctions.lcalc import lcalc261return lcalc.zeros_in_interval(x, y, stepsize, L=self.__E)262263def values_along_line(self, s0, s1, number_samples):264"""265Return values of $L(E, s)$ at \code{number_samples}266equally-spaced sample points along the line from $s_0$ to267$s_1$ in the complex plane.268269\note{The L-series is normalized so that the center of the270critical strip is 1.}271272INPUT:273s0, s1 -- complex numbers274number_samples -- integer275276OUTPUT:277list -- list of pairs (s, zeta(s)), where the s are278equally spaced sampled points on the line from279s0 to s1.280281EXAMPLES:282sage: I = CC.0283sage: E = EllipticCurve('37a')284sage: E.lseries().values_along_line(1, 0.5+20*I, 5) # long time and slightly random output285[(0.500000000, 0), (0.400000000 + 4.00000000*I, 3.31920245 - 2.60028054*I), (0.300000000 + 8.00000000*I, -0.886341185 - 0.422640337*I), (0.200000000 + 12.0000000*I, -3.50558936 - 0.108531690*I), (0.100000000 + 16.0000000*I, -3.87043288 - 1.88049411*I)]286"""287from sage.lfunctions.lcalc import lcalc288return lcalc.values_along_line(s0-RationalField()('1/2'),289s1-RationalField()('1/2'),290number_samples, L=self.__E)291292def twist_values(self, s, dmin, dmax):293r"""294Return values of $L(E, s, \chi_d)$ for each quadratic295character $\chi_d$ for $d_{\min} \leq d \leq d_{\max}$.296297\note{The L-series is normalized so that the center of the298critical strip is 1.}299300INPUT:301s -- complex numbers302dmin -- integer303dmax -- integer304305OUTPUT:306list -- list of pairs (d, L(E, s,chi_d))307308EXAMPLES:309sage: E = EllipticCurve('37a')310sage: E.lseries().twist_values(1, -12, -4) # slightly random output depending on architecture311[(-11, 1.4782434171), (-8, 0), (-7, 1.8530761916), (-4, 2.4513893817)]312sage: F = E.quadratic_twist(-8)313sage: F.rank()3141315sage: F = E.quadratic_twist(-7)316sage: F.rank()3170318"""319from sage.lfunctions.lcalc import lcalc320return lcalc.twist_values(s - RationalField()('1/2'), dmin, dmax, L=self.__E)321322def twist_zeros(self, n, dmin, dmax):323r"""324Return first $n$ real parts of nontrivial zeros of325$L(E,s,\chi_d)$ for each quadratic character $\chi_d$ with326$d_{\min} \leq d \leq d_{\max}$.327328\note{The L-series is normalized so that the center of the329critical strip is 1.}330331INPUT:332n -- integer333dmin -- integer334dmax -- integer335336OUTPUT:337dict -- keys are the discriminants $d$, and338values are list of corresponding zeros.339340EXAMPLES:341sage: E = EllipticCurve('37a')342sage: E.lseries().twist_zeros(3, -4, -3) # long time343{-4: [1.60813783, 2.96144840, 3.89751747], -3: [2.06170900, 3.48216881, 4.45853219]}344"""345from sage.lfunctions.lcalc import lcalc346return lcalc.twist_zeros(n, dmin, dmax, L=self.__E)347348def at1(self, k=0):349r"""350Compute $L(E,1)$ using $k$ terms of the series for $L(E,1)$ as351explained on page 406 of Henri Cohen's book"A Course in Computational352Algebraic Number Theory". If the argument $k$ is not specified,353then it defaults to $\sqrt(N)$, where $N$ is the conductor.354355The real precision used in each step of the computation is the356precision of machine floats.357358INPUT:359k -- (optional) an integer, defaults to sqrt(N).360361OUTPUT:362float -- L(E,1)363float -- a bound on the error in the approximation; this364is a provably correct upper bound on the sum365of the tail end of the series used to compute L(E,1).366367This function is disjoint from the PARI \code{elllseries}368command, which is for a similar purpose. To use that command369(via the PARI C library), simply type370\code{E.pari_mincurve().elllseries(1)}371372ALGORITHM:373\begin{enumerate}374\item Compute the root number eps. If it is -1, return 0.375376\item Compute the Fourier coefficients a_n, for n up to and377including k.378379\item Compute the sum380$$3812 * sum_{n=1}^{k} (a_n / n) * exp(-2*pi*n/Sqrt(N)),382$$383where N is the conductor of E.384385\item Compute a bound on the tail end of the series, which is386$$3872 * e^(-2 * pi * (k+1) / sqrt(N)) / (1 - e^(-2*pi/sqrt(N))).388$$389For a proof see [Grigov-Jorza-Patrascu-Patrikis-Stein].390\end{enumerate}391392EXAMPLES:393sage: E = EllipticCurve('37b')394sage: E.lseries().at1(100)395(0.725681061936153, 1.52437502288743e-45)396"""397if self.__E.root_number() == -1:398return 0399sqrtN = float(self.__E.conductor().sqrt())400k = int(k)401if k == 0: k = int(ceil(sqrtN))402an = self.__E.anlist(k) # list of Sage ints403# Compute z = e^(-2pi/sqrt(N))404pi = 3.14159265358979323846405z = exp(-2*pi/sqrtN)406zpow = z407s = 0.0408for n in xrange(1,k+1):409s += (zpow * float(an[n]))/n410zpow *= z411412error = 2*zpow / (1 - z)413414return R(2*s), R(error)415416def deriv_at1(self, k=0):417r"""418Compute $L'(E,1)$ using$ k$ terms of the series for $L'(E,1)$.419420The algorithm used is from page 406 of Henri Cohen's book ``A421Course in Computational Algebraic Number Theory.''422423The real precision of the computation is the precision of424Python floats.425426INPUT:427k -- int; number of terms of the series428429OUTPUT:430real number -- an approximation for L'(E,1)431real number -- a bound on the error in the approximation432433ALGORITHM:434\begin{enumerate}435\item Compute the root number eps. If it is 1, return 0.436437\item Compute the Fourier coefficients $a_n$, for $n$ up to and438including $k$.439440\item Compute the sum441$$4422 * \sum_{n=1}^{k} (a_n / n) * E_1(2 \pi n/\sqrt{N}),443$$444where $N$ is the conductor of $E$, and $E_1$ is the445exponential integral function.446447\item Compute a bound on the tail end of the series, which is448$$4492 * e^{-2 \pi (k+1) / \sqrt{N}} / (1 - e^{-2 \ pi/\sqrt{N}}).450$$451For a proof see [Grigorov-Jorza-Patrascu-Patrikis-Stein]. This452is exactly the same as the bound for the approximation to453$L(E,1)$ produced by \code{E.lseries().at1}.454\end{enumerate}455456EXAMPLES::457458sage: E = EllipticCurve('37a')459sage: E.lseries().deriv_at1()460(0.305986660898516, 0.000800351433106958)461sage: E.lseries().deriv_at1(100)462(0.305999773834052, 1.52437502288740e-45)463sage: E.lseries().deriv_at1(1000)464(0.305999773834052, 0.000000000000000)465"""466if self.__E.root_number() == 1: return 0467k = int(k)468sqrtN = float(self.__E.conductor().sqrt())469if k == 0: k = int(ceil(sqrtN))470an = self.__E.anlist(k) # list of Sage Integers471# Compute z = e^(-2pi/sqrt(N))472pi = 3.14159265358979323846473v = transcendental.exponential_integral_1(2*pi/sqrtN, k)474L = 2*float(sum([ (v[n-1] * an[n])/n for n in xrange(1,k+1)]))475error = 2*exp(-2*pi*(k+1)/sqrtN)/(1-exp(-2*pi/sqrtN))476return R(L), R(error)477478def __call__(self, s):479r"""480Returns the value of the L-series of the elliptic curve E at s, where s481must be a real number.482483Use self.extended for s complex.484485\note{If the conductor of the curve is large, say $>10^{12}$,486then this function will take a very long time, since it uses487an $O(\sqrt{N})$ algorithm.}488489EXAMPLES:490sage: E = EllipticCurve([1,2,3,4,5])491sage: L = E.lseries()492sage: L(1)4930.000000000000000494sage: L(1.1)4950.285491007678148496sage: L(1.1 + I)4970.174851377216615 + 0.816965038124457*I498"""499return self.dokchitser()(s)500501#def extended(self, s, prec):502# r"""503# Returns the value of the L-series of the elliptic curve E at s504# can be any complex number using prec terms of the power series505# expansion.506#507#508# WARNING: This may be slow. Consider using \code{dokchitser()}509# instead.510#511# INPUT:512# s -- complex number513# prec -- integer514#515# EXAMPLES:516# sage: E = EllipticCurve('389a')517# sage: E.lseries().extended(1 + I, 50)518# -0.638409959099589 + 0.715495262192901*I519# sage: E.lseries().extended(1 + 0.1*I, 50)520# -0.00761216538818315 + 0.000434885704670107*I521#522# NOTE: You might also want to use Tim Dokchitser's523# L-function calculator, which is available by typing524# L = E.lseries().dokchitser(), then evaluating L. It525# gives the same information but is sometimes much faster.526#527# """528# try:529# s = C(s)530# except TypeError:531# raise TypeError, "Input argument %s must be coercible to a complex number"%s532# prec = int(prec)533# if abs(s.imag()) < R(0.0000000000001):534# return self(s.real())535# N = self.__E.conductor()536# from sage.symbolic.constants import pi537# pi = R(pi)538# Gamma = transcendental.gamma539# Gamma_inc = transcendental.gamma_inc540# a = self.__E.anlist(prec)541# eps = self.__E.root_number()542# sqrtN = float(N.sqrt())543# def F(n, t):544# return Gamma_inc(t+1, 2*pi*n/sqrtN) * C(sqrtN/(2*pi*n))**(t+1)545# return C(N)**(-s/2) * C(2*pi)**s * Gamma(s)**(-1)\546# * sum([a[n]*(F(n,s-1) + eps*F(n,1-s)) for n in xrange(1,prec+1)])547548def L1_vanishes(self):549"""550Returns whether or not L(E,1) = 0. The result is provably551correct if the Manin constant of the associated optimal552quotient is <= 2. This hypothesis on the Manin constant553is true for all curves of conductor <= 40000 (by Cremona) and554all semistable curves (i.e., squarefree conductor).555556EXAMPLES:557sage: E = EllipticCurve([0, -1, 1, -10, -20]) # 11A = X_0(11)558sage: E.lseries().L1_vanishes()559False560sage: E = EllipticCurve([0, -1, 1, 0, 0]) # X_1(11)561sage: E.lseries().L1_vanishes()562False563sage: E = EllipticCurve([0, 0, 1, -1, 0]) # 37A (rank 1)564sage: E.lseries().L1_vanishes()565True566sage: E = EllipticCurve([0, 1, 1, -2, 0]) # 389A (rank 2)567sage: E.lseries().L1_vanishes()568True569sage: E = EllipticCurve([0, 0, 1, -38, 90]) # 361A (CM curve))570sage: E.lseries().L1_vanishes()571True572sage: E = EllipticCurve([0,-1,1,-2,-1]) # 141C (13-isogeny)573sage: E.lseries().L1_vanishes()574False575576WARNING: It's conceivable that machine floats are not large577enough precision for the computation; if this could be the578case a RuntimeError is raised. The curve's real period would579have to be very small for this to occur.580581ALGORITHM: Compute the root number. If it is -1 then L(E,s)582vanishes to odd order at 1, hence vanishes. If it is +1, use583a result about modular symbols and Mazur's "Rational Isogenies"584paper to determine a provably correct bound (assuming Manin585constant is <= 2) so that we can determine whether L(E,1) = 0.586587AUTHOR: William Stein, 2005-04-20.588"""589return self.L_ratio() == 0590591def L_ratio(self):592r"""593Returns the ratio $L(E,1)/\Omega$ as an exact rational594number. The result is \emph{provably} correct if the Manin595constant of the associated optimal quotient is $\leq 2$. This596hypothesis on the Manin constant is true for all semistable597curves (i.e., squarefree conductor), by a theorem of Mazur598from his \emph{Rational Isogenies of Prime Degree} paper.599600EXAMPLES:601sage: E = EllipticCurve([0, -1, 1, -10, -20]) # 11A = X_0(11)602sage: E.lseries().L_ratio()6031/5604sage: E = EllipticCurve([0, -1, 1, 0, 0]) # X_1(11)605sage: E.lseries().L_ratio()6061/25607sage: E = EllipticCurve([0, 0, 1, -1, 0]) # 37A (rank 1)608sage: E.lseries().L_ratio()6090610sage: E = EllipticCurve([0, 1, 1, -2, 0]) # 389A (rank 2)611sage: E.lseries().L_ratio()6120613sage: E = EllipticCurve([0, 0, 1, -38, 90]) # 361A (CM curve))614sage: E.lseries().L_ratio()6150616sage: E = EllipticCurve([0,-1,1,-2,-1]) # 141C (13-isogeny)617sage: E.lseries().L_ratio()6181619sage: E = EllipticCurve(RationalField(), [1, 0, 0, 1/24624, 1/886464])620sage: E.lseries().L_ratio()6212622623# See trac #3651:624sage: EllipticCurve([0,0,0,-193^2,0]).sha().an()6254626627WARNING: It's conceivable that machine floats are not large628enough precision for the computation; if this could be the629case a RuntimeError is raised. The curve's real period would630have to be very small for this to occur.631632ALGORITHM: Compute the root number. If it is -1 then L(E,s)633vanishes to odd order at 1, hence vanishes. If it is +1, use634a result about modular symbols and Mazur's "Rational Isogenies"635paper to determine a provably correct bound (assuming Manin636constant is <= 2) so that we can determine whether L(E,1) = 0.637638AUTHOR: William Stein, 2005-04-20.639"""640try:641return self.__lratio642except AttributeError:643pass644645if not self.__E.is_minimal():646self.__lratio = self.__E.minimal_model().lseries().L_ratio()647return self.__lratio648649if self.__E.root_number() == -1:650self.__lratio = Q(0)651return self.__lratio652653# Even root number. Decide if L(E,1) = 0. If E is a modular654# *OPTIMAL* quotient of J_0(N) elliptic curve, we know that T *655# L(E,1)/omega is an integer n, where T is the order of the656# image of the rational torsion point (0)-(oo) in E(Q), and657# omega is the least real Neron period. (This is proved in my658# Ph.D. thesis, but is probably well known.) We can easily659# compute omega to very high precision using AGM. So to prove660# that L(E,1) = 0 we compute T/omega * L(E,1) to sufficient661# precision to determine it as an integer. If eps is the662# error in computation of L(E,1), then the error in computing663# the product is (2T/Omega_E) * eps, and we need this to be664# less than 0.5, i.e.,665# (2T/Omega_E) * eps < 0.5,666# so667# eps < 0.5 * Omega_E / (2T) = Omega_E / (4*T).668#669# Since in general E need not be optimal, we have to choose670# eps = Omega_E/(8*t*B), where t is the exponent of E(Q)_tor,671# and is a multiple of the degree of an isogeny between E672# and the optimal curve.673#674# NOTES: We *do* have to worry about the Manin constant, since675# we are using the Neron model to compute omega, not the676# newform. My theorem replaces the omega above by omega/c,677# where c is the Manin constant, and the bound must be678# correspondingly smaller. If the level is square free, then679# the Manin constant is 1 or 2, so there's no problem (since680# we took 8 instead of 4 in the denominator). If the level681# is divisible by a square, then the Manin constant could682# be a divisible by an arbitrary power of that prime, except683# that Edixhoven claims the primes that appear are <= 7.684685t = self.__E.torsion_subgroup().order()686omega = self.__E.period_lattice().basis()[0]687d = self.__E._multiple_of_degree_of_isogeny_to_optimal_curve()688C = 8*d*t689eps = omega / C690# coercion of 10**(-15) to our real field is needed to691# make unambiguous comparison692if eps < R(10**(-15)): # liberal bound on precision of float693raise RuntimeError, "Insufficient machine precision (=%s) for computation."%eps694sqrtN = 2*int(sqrt(self.__E.conductor()))695k = sqrtN + 10696while True:697L1, error_bound = self.at1(k)698if error_bound < eps:699n = int(round(L1*C/omega))700quo = Q(n) / Q(C)701self.__lratio = quo / self.__E.real_components()702return self.__lratio703k += sqrtN704misc.verbose("Increasing precision to %s terms."%k)705706707