Path: blob/master/sage/schemes/elliptic_curves/padic_lseries.py
4156 views
# -*- coding: utf-8 -*-1r"""2p-adic L-functions of elliptic curves34To an elliptic curve `E` over the rational numbers and a prime `p`, one5can associate a `p`-adic L-function; at least if `E` does not have additive6reduction at `p`. This function is defined by interpolation of L-values of `E`7at twists. Through the main conjecture of Iwasawa theory it should also be8equal to a characteristic series of a certain Selmer group.910If `E` is ordinary, then it is an element of the Iwasawa algebra11`\Lambda(\ZZ_p^\times) = \ZZ_p[\Delta][\![T]\!]`, where `\Delta` is the group12of `(p-1)`-st roots of unity in `\ZZ_p^\times`, and `T = [\gamma] - 1` where13`\gamma = 1 + p` is a generator of `1 + p\ZZ_p`. (There is a slightly different14description for `p = 2`.)1516One can decompose this algebra as the direct product of the subalgebras17corresponding to the characters of `\Delta`, which are simply the powers18`\tau^\eta` (`0 \le \eta \le p-2`) of the Teichmueller character `\tau: \Delta19\to \ZZ_p^\times`. Projecting the L-function into these components gives `p-1`20power series in `T`, each with coefficients in `\ZZ_p`.2122If `E` is supersingular, the series will have coefficients in a quadratic23extension of `\QQ_p`, and the coefficients will be unbounded. In this case we24have only implemented the series for `\eta = 0`. We have also implemented the25`p`-adic L-series as formulated by Perrin-Riou [BP], which has coefficients in26the Dieudonne module `D_pE = H^1_{dR}(E/\QQ_p)` of `E`. There is a different27description by Pollack [Po] which is not available here.2829According to the `p`-adic version of the Birch and Swinnerton-Dyer conjecture30[MTT], the order of vanishing of the `L`-function at the trivial character31(i.e. of the series for `\eta = 0` at `T = 0`) is just the rank of `E(\QQ)`, or32this rank plus one if the reduction at `p` is split multiplicative.3334See [SW] for more details.3536REFERENCES:3738- [MTT] B. Mazur, J. Tate, and J. Teitelbaum,39On `p`-adic analogues of the conjectures of Birch and40Swinnerton-Dyer, Inventiones mathematicae 84, (1986), 1-48.4142- [BP] Dominique Bernardi and Bernadette Perrin-Riou,43Variante `p`-adique de la conjecture de Birch et44Swinnerton-Dyer (le cas supersingulier), C. R. Acad. Sci. Paris,45Ser I. Math, 317 (1993), no 3, 227-232.4647- [Po] Robert Pollack, On the `p`-adic L-function of a modular form48at supersingular prime, Duke Math. J. 118 (2003), no 3, 523-558.4950- [SW] William Stein and Christian Wuthrich, Computations About Tate-Shafarevich Groups51using Iwasawa theory, preprint 2009.5253AUTHORS:5455- William Stein (2007-01-01): first version5657- Chris Wuthrich (22/05/2007): changed minor issues and added supersingular things5859- Chris Wuthrich (11/2008): added quadratic_twists6061- David Loeffler (01/2011): added nontrivial Teichmueller components6263"""6465######################################################################66# Copyright (C) 2007 William Stein <[email protected]>67#68# Distributed under the terms of the GNU General Public License (GPL)69#70# This code is distributed in the hope that it will be useful,71# but WITHOUT ANY WARRANTY; without even the implied warranty of72# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU73# General Public License for more details.74#75# The full text of the GPL is available at:76#77# http://www.gnu.org/licenses/78######################################################################798081from sage.rings.integer_ring import ZZ82from sage.rings.rational_field import QQ83from sage.rings.padics.factory import Qp84from sage.rings.infinity import infinity85from sage.rings.all import LaurentSeriesRing, PowerSeriesRing, PolynomialRing, Integers8687from sage.rings.integer import Integer88from sage.rings.arith import valuation, binomial, kronecker_symbol, gcd, prime_divisors, valuation8990from sage.structure.sage_object import SageObject9192from sage.misc.all import verbose, denominator, get_verbose93import sage.rings.arith as arith9495from sage.modules.free_module_element import vector96import sage.matrix.all as matrix97import monsky_washnitzer98# from sage.interfaces.all import gp99from sage.misc.functional import log100101from sage.misc.decorators import rename_keyword102103class pAdicLseries(SageObject):104r"""105The `p`-adic L-series of an elliptic curve.106107EXAMPLES:108An ordinary example::109110sage: e = EllipticCurve('389a')111sage: L = e.padic_lseries(5)112sage: L.series(0)113Traceback (most recent call last):114...115ValueError: n (=0) must be a positive integer116sage: L.series(1)117O(T^1)118sage: L.series(2)119O(5^4) + O(5)*T + (4 + O(5))*T^2 + (2 + O(5))*T^3 + (3 + O(5))*T^4 + O(T^5)120sage: L.series(3, prec=10)121O(5^5) + O(5^2)*T + (4 + 4*5 + O(5^2))*T^2 + (2 + 4*5 + O(5^2))*T^3 + (3 + O(5^2))*T^4 + (1 + O(5))*T^5 + O(5)*T^6 + (4 + O(5))*T^7 + (2 + O(5))*T^8 + O(5)*T^9 + O(T^10)122sage: L.series(2,quadratic_twist=-3)1232 + 4*5 + 4*5^2 + O(5^4) + O(5)*T + (1 + O(5))*T^2 + (4 + O(5))*T^3 + O(5)*T^4 + O(T^5)124125126A prime p such that E[p] is reducible::127128sage: L = EllipticCurve('11a').padic_lseries(5)129sage: L.series(1)1305 + O(5^2) + O(T)131sage: L.series(2)1325 + 4*5^2 + O(5^3) + O(5^0)*T + O(5^0)*T^2 + O(5^0)*T^3 + O(5^0)*T^4 + O(T^5)133sage: L.series(3)1345 + 4*5^2 + 4*5^3 + O(5^4) + O(5)*T + O(5)*T^2 + O(5)*T^3 + O(5)*T^4 + O(T^5)135136An example showing the calculation of nontrivial Teichmueller twists::137138sage: E=EllipticCurve('11a1')139sage: lp=E.padic_lseries(7)140sage: lp.series(4,eta=1)1416 + 2*7^3 + 5*7^4 + O(7^6) + (4*7 + 2*7^2 + O(7^3))*T + (2 + 3*7^2 + O(7^3))*T^2 + (1 + 2*7 + 2*7^2 + O(7^3))*T^3 + (1 + 3*7^2 + O(7^3))*T^4 + O(T^5)142sage: lp.series(4,eta=2)1435 + 6*7 + 4*7^2 + 2*7^3 + 3*7^4 + 2*7^5 + O(7^6) + (6 + 4*7 + 7^2 + O(7^3))*T + (3 + 2*7^2 + O(7^3))*T^2 + (1 + 4*7 + 7^2 + O(7^3))*T^3 + (6 + 6*7 + 6*7^2 + O(7^3))*T^4 + O(T^5)144sage: lp.series(4,eta=3)145O(7^6) + (3 + 2*7 + 5*7^2 + O(7^3))*T + (5 + 4*7 + 5*7^2 + O(7^3))*T^2 + (3*7 + 7^2 + O(7^3))*T^3 + (2*7 + 7^2 + O(7^3))*T^4 + O(T^5)146147(Note that the last series vanishes at `T = 0`, which is consistent with ::148149sage: E.quadratic_twist(-7).rank()1501151152This proves that `E` has rank 1 over `\QQ(\zeta_7)`.)153154the load-dumps test::155156sage: lp = EllipticCurve('11a').padic_lseries(5)157sage: lp == loads(dumps(lp))158True159"""160def __init__(self, E, p, use_eclib=False, normalize='L_ratio'):161r"""162INPUT:163164- ``E`` - an elliptic curve165- ``p`` - a prime of good reduction166- ``use_eclib`` - bool (default:True); whether or not to use167John Cremona's ``eclib`` for the computation of modular168symbols169- ``normalize`` - ``'L_ratio'`` (default), ``'period'`` or ``'none'``;170this is describes the way the modular symbols171are normalized. See ``modular_symbol`` of172an elliptic curve over Q for more details.173174EXAMPLES::175176sage: E = EllipticCurve('11a1')177sage: Lp = E.padic_lseries(3)178sage: Lp.series(2,prec=3)1792 + 3 + 3^2 + 2*3^3 + O(3^4) + (1 + O(3))*T + (1 + O(3))*T^2 + O(T^3)180"""181self._E = E182self._p = ZZ(p)183self._normalize = normalize184self._use_eclib = use_eclib185if not self._p.is_prime():186raise ValueError, "p (=%s) must be a prime"%p187if E.conductor() % (self._p)**2 == 0:188raise NotImplementedError, "p (=%s) must be a prime of semi-stable reduction"%p189190try :191crla = E.label()192except RuntimeError :193print "Warning : Curve outside Cremona's table. Computations of modular symbol space might take very long !"194195self._modular_symbol = E.modular_symbol(sign=+1, use_eclib = use_eclib, normalize=normalize)196197def __add_negative_space(self):198r"""199A helper function not designed for direct use.200201This function add the attribute ``_negative_modular_symbol`` to the class. This may take time202and will only be needed when twisting with negative fundamental discriminants.203204EXAMPLES::205206sage: E = EllipticCurve('11a1')207sage: lp = E.padic_lseries(5)208sage: lp.modular_symbol(1/7,sign=-1) #indirect doctest209-1210211"""212if self._use_eclib:213verbose('Currently there is no negative modular symbols in eclib, so we have to fall back on the implementation of modular symbols in sage')214# once there is a eclib implementation of -1, this should be changed.215self._negative_modular_symbol = self._E.modular_symbol(sign=-1, use_eclib = False, normalize=self._normalize)216else:217self._negative_modular_symbol = self._E.modular_symbol(sign=-1, use_eclib = False, normalize=self._normalize)218219def __cmp__(self,other):220r"""221Compare self and other.222223TESTS::224sage: lp1 = EllipticCurve('11a1').padic_lseries(5)225sage: lp2 = EllipticCurve('11a1').padic_lseries(7)226sage: lp3 = EllipticCurve('11a2').padic_lseries(5)227sage: lp1 == lp1228True229sage: lp1 == lp2230False231sage: lp1 == lp3232False233234"""235c = cmp(type(self), type(other))236if c:237return c238return cmp((self._E, self._p), (other._E, other._p))239240241def elliptic_curve(self):242r"""243Return the elliptic curve to which this `p`-adic L-series is associated.244245EXAMPLES::246247sage: L = EllipticCurve('11a').padic_lseries(5)248sage: L.elliptic_curve()249Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field250"""251return self._E252253def prime(self):254r"""255Returns the prime `p` as in 'p-adic L-function'.256257EXAMPLES::258259sage: L = EllipticCurve('11a').padic_lseries(5)260sage: L.prime()2615262"""263return self._p264265def _repr_(self):266r"""267Return print representation.268269EXAMPLES::270271sage: e = EllipticCurve('37a')272sage: e.padic_lseries(3)._repr_()273'3-adic L-series of Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field'274sage: e.padic_lseries(3,normalize='none')2753-adic L-series of Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field (not normalized)276sage: L = e.padic_lseries(3,normalize='none')277sage: L.rename('(factor)*L_3(T)')278sage: L279(factor)*L_3(T)280"""281s = "%s-adic L-series of %s"%(self._p, self._E)282if not self._normalize == 'L_ratio':283s += ' (not normalized)'284return s285286def modular_symbol(self, r, sign=+1, quadratic_twist= +1):287r"""288Return the modular symbol evaluated at `r`.289This is used to compute this `p`-adic290L-series.291292Note that the normalisation is not correct at this293stage: use ``_quotient_of periods_to_twist`` to correct.294295Note also that this function does not check if the condition296on the quadratic_twist=D is satisfied. So the result will only297be correct if for each prime `\ell` dividing `D`, we have298`ord_{\ell(N)}<= ord_{\ell}(D)`, where `N` is the conductor of the curve.299300INPUT:301302- ``r`` - a cusp given as either a rational number or oo303304- ``sign`` - +1 (default) or -1 (only implemented without twists)305306- ``quadratic_twist`` - a fundamental discriminant of a quadratic field or +1 (default)307308EXAMPLES::309310sage: E = EllipticCurve('11a1')311sage: lp = E.padic_lseries(5)312sage: [lp.modular_symbol(r) for r in [0,1/5,oo,1/11]]313[1/5, 6/5, 0, 0]314sage: [lp.modular_symbol(r,sign=-1) for r in [0,1/3,oo,1/7]]315[0, 1, 0, -1]316sage: [lp.modular_symbol(r,quadratic_twist=-20) for r in [0,1/5,oo,1/11]]317[2, 2, 0, 1]318319sage: lpt = E.quadratic_twist(-3).padic_lseries(5)320sage: et = E.padic_lseries(5)._quotient_of_periods_to_twist(-3)321sage: lpt.modular_symbol(0) == lp.modular_symbol(0,quadratic_twist=-3)/et322True323324"""325if quadratic_twist == +1 :326if sign == +1 :327return self._modular_symbol(r)328elif sign == -1:329try:330m = self._negative_modular_symbol331except (KeyError, AttributeError):332if not hasattr(self, '_modular_symbol_negative'):333self.__add_negative_space()334m = self._negative_modular_symbol335return m(r)336else :337D = quadratic_twist338if sign == -1:339raise NotImplementedError, "Quadratic twists for negative modular symbols are not yet implemented."340if D > 0:341m = self._modular_symbol342s = +1343else:344try:345m = self._negative_modular_symbol346except (KeyError, AttributeError):347if not hasattr(self, '_modular_symbol_negative'):348self.__add_negative_space()349m = self._negative_modular_symbol350s = -1351# without the ZZ here the u is treated as a 'int' and dividing by D gives 0352# this only happens when it is called from __init__ (?)353return s * sum([kronecker_symbol(D,u) * m(r+ZZ(u)/D) for u in range(1,abs(D))])354355356def measure(self, a, n, prec, quadratic_twist=+1, sign = +1):357r"""358Return the measure on `\ZZ_p^{\times}` defined by359360`\mu_{E,\alpha}^+ ( a + p^n \ZZ_p ) =361\frac{1}{\alpha^n} \left [\frac{a}{p^n}\right]^{+} -362\frac{1}{\alpha^{n+1}} \left[\frac{a}{p^{n-1}}\right]^{+}`363364where `[\cdot]^{+}` is the modular symbol. This is used to define365this `p`-adic L-function (at least when the reduction is good).366367The optional argument ``sign`` allows the minus symbol `[\cdot]^{-}` to368be substituted for the plus symbol.369370The optional argument ``quadratic_twist`` replaces `E` by the twist in371the above formula, but the twisted modular symbol is computed using a372sum over modular symbols of `E` rather then finding the modular symbols373for the twist. Quadratic twists are only implemented if the sign is374`+1`.375376Note that the normalisation is not correct at this377stage: use ``_quotient_of periods`` and ``_quotient_of periods_to_twist``378to correct.379380Note also that this function does not check if the condition381on the ``quadratic_twist=D`` is satisfied. So the result will only382be correct if for each prime `\ell` dividing `D`, we have383`ord_{\ell}(N)<= ord_{\ell}(D)`, where `N` is the conductor of the curve.384385INPUT:386387- ``a`` - an integer388389- ``n`` - a non-negative integer390391- ``prec`` - an integer392393- ``quadratic_twist`` (default = 1) - a fundamental discriminant of a quadratic field,394should be coprime to the conductor of `E`395396- ``sign`` (default = 1) - an integer, which should be `\pm 1`.397398EXAMPLES::399400sage: E = EllipticCurve('37a')401sage: L = E.padic_lseries(5)402sage: L.measure(1,2, prec=9)4032 + 3*5 + 4*5^3 + 2*5^4 + 3*5^5 + 3*5^6 + 4*5^7 + 4*5^8 + O(5^9)404sage: L.measure(1,2, quadratic_twist=8,prec=15)405O(5^15)406sage: L.measure(1,2, quadratic_twist=-4,prec=15)4074 + 4*5 + 4*5^2 + 3*5^3 + 2*5^4 + 5^5 + 3*5^6 + 5^8 + 2*5^9 + 3*5^12 + 2*5^13 + 4*5^14 + O(5^15)408409sage: E = EllipticCurve('11a1')410sage: a = E.quadratic_twist(-3).padic_lseries(5).measure(1,2,prec=15)411sage: b = E.padic_lseries(5).measure(1,2, quadratic_twist=-3,prec=15)412sage: a == b/E.padic_lseries(5)._quotient_of_periods_to_twist(-3)413True414415"""416s = ZZ(sign)417if s not in [1, -1]:418raise ValueError, "Sign must be +- 1"419if quadratic_twist != 1 and s != 1:420raise NotImplementedError, "Quadratic twists not implemented for sign -1"421422if quadratic_twist < 0:423s = -1424425try:426p, alpha, z, w, f = self.__measure_data[(n,prec,s)]427except (KeyError, AttributeError):428if not hasattr(self, '__measure_data'):429self.__measure_data = {}430p = self._p431alpha = self.alpha(prec=prec)432z = 1/(alpha**n)433w = p**(n-1)434if s == +1 :435f = self._modular_symbol436else :437try :438f = self._negative_modular_symbol439except (KeyError, AttributeError):440if not hasattr(self, '_modular_symbol_negative'):441self.__add_negative_space()442f = self._negative_modular_symbol443self.__measure_data[(n,prec,s)] = (p,alpha,z,w,f)444445if quadratic_twist == 1:446if self._E.conductor() % p == 0:447return z * f(a/(p*w))448return z * f(a/(p*w)) - (z/alpha) * f(a/w)449else:450D = quadratic_twist451chip = kronecker_symbol(D,p)452if self._E.conductor() % p == 0:453mu = chip**n * z * sum([kronecker_symbol(D,u) * f(a/(p*w)+ZZ(u)/D) for u in range(1,abs(D))])454else:455mu = chip**n * sum([kronecker_symbol(D,u) *(z * f(a/(p*w)+ZZ(u)/D) - chip *(z/alpha)* f(a/w+ZZ(u)/D)) for u in range(1,abs(D))])456return s*mu457458def alpha(self, prec=20):459r"""460Return a `p`-adic root `\alpha` of the polynomial `x^2 - a_p x461+ p` with `ord_p(\alpha) < 1`. In the ordinary case this is462just the unit root.463464INPUT:465- ``prec`` - positive integer, the `p`-adic precision of the root.466467EXAMPLES:468Consider the elliptic curve 37a::469470sage: E = EllipticCurve('37a')471472An ordinary prime::473474sage: L = E.padic_lseries(5)475sage: alpha = L.alpha(10); alpha4763 + 2*5 + 4*5^2 + 2*5^3 + 5^4 + 4*5^5 + 2*5^7 + 5^8 + 5^9 + O(5^10)477sage: alpha^2 - E.ap(5)*alpha + 5478O(5^10)479480A supersingular prime::481482sage: L = E.padic_lseries(3)483sage: alpha = L.alpha(10); alpha484(1 + O(3^10))*alpha485sage: alpha^2 - E.ap(3)*alpha + 3486(O(3^10))*alpha^2 + (O(3^11))*alpha + (O(3^11))487488A reducible prime::489490sage: L = EllipticCurve('11a').padic_lseries(5)491sage: L.alpha(5)4921 + 4*5 + 3*5^2 + 2*5^3 + 4*5^4 + O(5^5)493"""494try:495return self._alpha[prec]496except AttributeError:497self._alpha = {}498except KeyError:499pass500E = self._E501p = self._p502a_p = E.ap(p)503K = Qp(p, prec, print_mode='series')504505if E.conductor() % p == 0:506self._alpha[prec] = K(a_p)507return K(a_p)508509R = ZZ['x']510f = R([p, -a_p, 1])511if E.is_ordinary(p):512G = f.factor_padic(p, prec+5)513for pr, e in G:514a = -pr[0]515if a.valuation() < 1:516self._alpha[prec] = K(a)517return K(a)518raise ValueError, "bug in p-adic L-function alpha"519else: # supersingular case520f = f.change_ring(Qp(p, prec, print_mode='series'))521a = f.root_field('alpha', check_irreducible=False).gen()522self._alpha[prec] = a523return a524525def order_of_vanishing(self):526r"""527Return the order of vanishing of this `p`-adic L-series.528529The output of this function is provably correct, due to a530theorem of Kato [Ka].531532NOTE: currently `p` must be a prime of good ordinary reduction.533534REFERENCES:535536- [MTT] B. Mazur, J. Tate, and J. Teitelbaum,537On `p`-adic analogues of the conjectures of Birch and538Swinnerton-Dyer, Inventiones mathematicae 84, (1986), 1-48.539540- [Ka] Kayuza Kato, `p`-adic Hodge theory and values of zeta functions of modular541forms, Cohomologies `p`-adiques et applications arithmetiques III,542Asterisque vol 295, SMF, Paris, 2004.543544EXAMPLES::545546sage: L = EllipticCurve('11a').padic_lseries(3)547sage: L.order_of_vanishing()5480549sage: L = EllipticCurve('11a').padic_lseries(5)550sage: L.order_of_vanishing()5510552sage: L = EllipticCurve('37a').padic_lseries(5)553sage: L.order_of_vanishing()5541555sage: L = EllipticCurve('43a').padic_lseries(3)556sage: L.order_of_vanishing()5571558sage: L = EllipticCurve('37b').padic_lseries(3)559sage: L.order_of_vanishing()5600561sage: L = EllipticCurve('389a').padic_lseries(3)562sage: L.order_of_vanishing()5632564sage: L = EllipticCurve('389a').padic_lseries(5)565sage: L.order_of_vanishing()5662567sage: L = EllipticCurve('5077a').padic_lseries(5, use_eclib=True)568sage: L.order_of_vanishing()5693570"""571try:572return self.__ord573except AttributeError:574pass575576if not self.is_ordinary():577raise NotImplementedError578E = self.elliptic_curve()579if not E.is_good(self.prime()):580raise ValueError, "prime must be of good reduction"581r = E.rank()582n = 1583while True:584f = self.series(n)585v = f.valuation()586if v < r:587raise RuntimeError, "while computing p-adic order of vanishing, got a contradiction: the curve is %s, the curve has rank %s, but the p-adic L-series vanishes to order <= %s"%(E, r, v)588if v == r:589self.__ord = v590return v591n += 1592593594# def _c_bounds(self, n):595# raise NotImplementedError596597# def _prec_bounds(self, n,prec):598# raise NotImplementedError599600def teichmuller(self, prec):601r"""602Return Teichmuller lifts to the given precision.603604INPUT:605606- ``prec`` - a positive integer.607608OUTPUT:609610- a list of `p`-adic numbers, the cached Teichmuller lifts611612EXAMPLES::613614sage: L = EllipticCurve('11a').padic_lseries(7)615sage: L.teichmuller(1)616[0, 1, 2, 3, 4, 5, 6]617sage: L.teichmuller(2)618[0, 1, 30, 31, 18, 19, 48]619"""620p = self._p621K = Qp(p, prec, print_mode='series')622return [Integer(0)] + \623[a.residue(prec).lift() for a in K.teichmuller_system()]624625def _e_bounds(self, n, prec):626r"""627A helper function not designed for direct use.628629It computes the valuations of the coefficients of `\omega_n = (1+T)^{p^n}-1`.630631EXAMPLES::632633sage: E = EllipticCurve('11a1')634sage: Lp = E.padic_lseries(2)635sage: Lp._e_bounds(1,10)636[+Infinity, 1, 0, 0, 0, 0, 0, 0, 0, 0]637sage: Lp._e_bounds(2,10)638[+Infinity, 2, 1, 1, 0, 0, 0, 0, 0, 0]639sage: Lp._e_bounds(3,10)640[+Infinity, 3, 2, 2, 1, 1, 1, 1, 0, 0]641sage: Lp._e_bounds(4,10)642[+Infinity, 4, 3, 3, 2, 2, 2, 2, 1, 1]643"""644# trac 10280: replace with new corrected code, note that the sequence has to be decreasing.645pn = self._p**n646enj = infinity647res = [enj]648for j in range(1,prec):649bino = valuation(binomial(pn,j),self._p)650if bino < enj:651enj = bino652res.append(enj)653return res654655def _get_series_from_cache(self, n, prec, D, eta):656r"""657A helper function not designed for direct use.658659It picks up the series in the cache if it has been previously computed.660661EXAMPLES::662663sage: E = EllipticCurve('11a1')664sage: Lp = E.padic_lseries(5)665sage: Lp._get_series_from_cache(3,5,1,0)666sage: Lp.series(3,prec=5)6675 + 4*5^2 + 4*5^3 + O(5^4) + O(5)*T + O(5)*T^2 + O(5)*T^3 + O(5)*T^4 + O(T^5)668sage: Lp._get_series_from_cache(3,5,1,0)6695 + 4*5^2 + 4*5^3 + O(5^4) + O(5)*T + O(5)*T^2 + O(5)*T^3 + O(5)*T^4 + O(T^5)670"""671try:672return self.__series[(n,prec,D,eta)]673except AttributeError:674self.__series = {}675except KeyError:676for _n, _prec, _D, _eta in self.__series.keys():677if _n == n and _D == D and _eta == eta and _prec >= prec:678return self.__series[(_n,_prec,_D,_eta)].add_bigoh(prec)679return None680681def _set_series_in_cache(self, n, prec, D, eta, f):682r"""683A helper function not designed for direct use.684685It picks up the series in the cache if it has been previously computed.686687EXAMPLES::688689sage: E = EllipticCurve('11a1')690sage: Lp = E.padic_lseries(5)691sage: Lp.series(3,prec=5)6925 + 4*5^2 + 4*5^3 + O(5^4) + O(5)*T + O(5)*T^2 + O(5)*T^3 + O(5)*T^4 + O(T^5)693sage: Lp._set_series_in_cache(3,5,1,0,0)694sage: Lp.series(3,prec=5)6950696"""697self.__series[(n,prec,D, eta)] = f698699700def _quotient_of_periods_to_twist(self,D):701r"""702For a fundamental discriminant `D` of a quadratic number field this computes the constant `\eta` such that703`\sqrt{D}\cdot\Omega_{E_D}^{+} =\eta\cdot \Omega_E^{sign(D)}`. As in [MTT]_ page 40.704This is either 1 or 2 unless the condition on the twist is not satisfied, e.g. if we are 'twisting back'705to a semi-stable curve.706707REFERENCES:708709- [MTT] B. Mazur, J. Tate, and J. Teitelbaum,710On `p`-adic analogues of the conjectures of Birch and711Swinnerton-Dyer, Invertiones mathematicae 84, (1986), 1-48.712713.. note: No check on precision is made, so this may fail for huge `D`.714715EXAMPLES::716717sage: E = EllipticCurve('37b1')718sage: lp = E.padic_lseries(3)719sage: lp._quotient_of_periods_to_twist(-20)7201721sage: lp._quotient_of_periods_to_twist(-4)7221723sage: lp._quotient_of_periods_to_twist(-3)7241725sage: lp._quotient_of_periods_to_twist(-8)7262727sage: lp._quotient_of_periods_to_twist(8)7282729sage: lp._quotient_of_periods_to_twist(5)7301731sage: lp._quotient_of_periods_to_twist(12)7321733734sage: E = EllipticCurve('11a1')735sage: Et = E.quadratic_twist(-3)736sage: lpt = Et.padic_lseries(5)737sage: lpt._quotient_of_periods_to_twist(-3)7386739740"""741from sage.functions.all import sqrt742# This funciton does not depend on p and could be moved out of this file but it is needed only here743744# Note that the number of real components does not change by twisting.745if D == 1:746return 1747if D > 1:748Et = self._E.quadratic_twist(D)749qt = Et.period_lattice().basis()[0]/self._E.period_lattice().basis()[0]750qt *= sqrt(qt.parent()(D))751else:752Et = self._E.quadratic_twist(D)753qt = Et.period_lattice().basis()[0]/self._E.period_lattice().basis()[1].imag()754qt *= sqrt(qt.parent()(-D))755verbose('the real approximation is %s'%qt)756# we know from MTT that the result has a denominator 1757return QQ(int(round(8*qt)))/8758759760class pAdicLseriesOrdinary(pAdicLseries):761def series(self, n=2, quadratic_twist=+1, prec=5, eta=0):762r"""763Returns the `n`-th approximation to the `p`-adic L-series, in the764component corresponding to the `\eta`-th power of the Teichmueller765character, as a power series in `T` (corresponding to `\gamma-1` with766`\gamma=1+p` as a generator of `1+p\ZZ_p`). Each coefficient is a767`p`-adic number whose precision is provably correct.768769Here the normalization of the `p`-adic L-series is chosen770such that `L_p(E,1) = (1-1/\alpha)^2 L(E,1)/\Omega_E`771where `\alpha` is the unit root of the characteristic772polynomial of Frobenius on `T_pE` and `\Omega_E` is the773Neron period of `E`.774775INPUT:776777- ``n`` - (default: 2) a positive integer778- ``quadratic_twist`` - (default: +1) a fundamental discriminant of a779quadratic field, coprime to the conductor of the curve780- ``prec`` - (default: 5) maximal number of terms of the series to781compute; to compute as many as possible just give a very large782number for ``prec``; the result will still be correct.783- ``eta`` (default: 0) an integer (specifying the power of the784Teichmueller character on the group of roots of unity in785`\ZZ_p^\times`)786787ALIAS: power_series is identical to series.788789EXAMPLES:790We compute some `p`-adic L-functions associated to the elliptic791curve 11a::792793sage: E = EllipticCurve('11a')794sage: p = 3795sage: E.is_ordinary(p)796True797sage: L = E.padic_lseries(p)798sage: L.series(3)7992 + 3 + 3^2 + 2*3^3 + O(3^5) + (1 + 3 + O(3^2))*T + (1 + 2*3 + O(3^2))*T^2 + O(3)*T^3 + O(3)*T^4 + O(T^5)800801Another example at a prime of bad reduction, where the802`p`-adic L-function has an extra 0 (compared to the non803`p`-adic L-function)::804805sage: E = EllipticCurve('11a')806sage: p = 11807sage: E.is_ordinary(p)808True809sage: L = E.padic_lseries(p)810sage: L.series(2)811O(11^4) + (10 + O(11))*T + (6 + O(11))*T^2 + (2 + O(11))*T^3 + (5 + O(11))*T^4 + O(T^5)812813We compute a `p`-adic L-function that vanishes to order 2::814815sage: E = EllipticCurve('389a')816sage: p = 3817sage: E.is_ordinary(p)818True819sage: L = E.padic_lseries(p)820sage: L.series(1)821O(T^1)822sage: L.series(2)823O(3^4) + O(3)*T + (2 + O(3))*T^2 + O(T^3)824sage: L.series(3)825O(3^5) + O(3^2)*T + (2 + 2*3 + O(3^2))*T^2 + (2 + O(3))*T^3 + (1 + O(3))*T^4 + O(T^5)826827Checks if the precision can be changed (trac 5846)::828829sage: L.series(3,prec=4)830O(3^5) + O(3^2)*T + (2 + 2*3 + O(3^2))*T^2 + (2 + O(3))*T^3 + O(T^4)831sage: L.series(3,prec=6)832O(3^5) + O(3^2)*T + (2 + 2*3 + O(3^2))*T^2 + (2 + O(3))*T^3 + (1 + O(3))*T^4 + (1 + O(3))*T^5 + O(T^6)833834Rather than computing the `p`-adic L-function for the curve '15523a1', one can835compute it as a quadratic_twist::836837sage: E = EllipticCurve('43a1')838sage: lp = E.padic_lseries(3)839sage: lp.series(2,quadratic_twist=-19)8402 + 2*3 + 2*3^2 + O(3^4) + (1 + O(3))*T + (1 + O(3))*T^2 + O(T^3)841sage: E.quadratic_twist(-19).label() #optional --- conductor is greater than 10000842'15523a1'843844This proves that the rank of '15523a1' is zero, even if ``mwrank`` can not determine this.845846We calculate the `L`-series in the nontrivial Teichmueller components::847848sage: L = EllipticCurve('110a1').padic_lseries(5)849sage: for j in [0..3]: print L.series(4, eta=j)850O(5^6) + (2 + 2*5 + 2*5^2 + O(5^3))*T + (5 + 5^2 + O(5^3))*T^2 + (4 + 4*5 + 2*5^2 + O(5^3))*T^3 + (1 + 5 + 3*5^2 + O(5^3))*T^4 + O(T^5)8513 + 2*5 + 2*5^3 + 3*5^4 + O(5^6) + (2 + 5 + 4*5^2 + O(5^3))*T + (1 + 4*5 + 2*5^2 + O(5^3))*T^2 + (1 + 5 + 5^2 + O(5^3))*T^3 + (2 + 4*5 + 4*5^2 + O(5^3))*T^4 + O(T^5)8522 + O(5^6) + (1 + 5 + O(5^3))*T + (2 + 4*5 + 3*5^2 + O(5^3))*T^2 + (4 + 5 + 2*5^2 + O(5^3))*T^3 + (4 + O(5^3))*T^4 + O(T^5)8531 + 3*5 + 4*5^2 + 2*5^3 + 5^4 + 4*5^5 + O(5^6) + (2 + 4*5 + 3*5^2 + O(5^3))*T + (2 + 3*5 + 5^2 + O(5^3))*T^2 + (1 + O(5^3))*T^3 + (2*5 + 2*5^2 + O(5^3))*T^4 + O(T^5)854"""855n = ZZ(n)856if n < 1:857raise ValueError, "n (=%s) must be a positive integer"%n858eta = ZZ(eta) % (self._p - 1)859860# check if the conditions on quadratic_twist are satisfied861D = ZZ(quadratic_twist)862if D != 1:863if eta != 0: raise NotImplementedError, "quadratic twists only implemented for the 0th Teichmueller component"864if D % 4 == 0:865d = D//4866if not d.is_squarefree() or d % 4 == 1:867raise ValueError, "quadratic_twist (=%s) must be a fundamental discriminant of a quadratic field"%D868else:869if not D.is_squarefree() or D % 4 != 1:870raise ValueError, "quadratic_twist (=%s) must be a fundamental discriminant of a quadratic field"%D871if gcd(D,self._p) != 1:872raise ValueError, "quadratic twist (=%s) must be coprime to p (=%s) "%(D,self._p)873if gcd(D,self._E.conductor())!= 1:874for ell in prime_divisors(D):875if valuation(self._E.conductor(),ell) > valuation(D,ell) :876raise ValueError, "can not twist a curve of conductor (=%s) by the quadratic twist (=%s)."%(self._E.conductor(),D)877878879p = self._p880if p == 2 and self._normalize :881print 'Warning : For p=2 the normalization might not be correct !'882#verbose("computing L-series for p=%s, n=%s, and prec=%s"%(p,n,prec))883884bounds = self._prec_bounds(n,prec)885padic_prec = max(bounds[1:]) + 5886verbose("using p-adic precision of %s"%padic_prec)887888res_series_prec = min(p**(n-1), prec)889verbose("using series precision of %s"%res_series_prec)890891ans = self._get_series_from_cache(n, res_series_prec,D,eta)892if not ans is None:893verbose("found series in cache")894return ans895896K = QQ897gamma = K(1 + p)898R = PowerSeriesRing(K,'T',res_series_prec)899T = R(R.gen(),res_series_prec )900L = R(0)901one_plus_T_factor = R(1)902gamma_power = K(1)903teich = self.teichmuller(padic_prec)904p_power = p**(n-1)905906verbose("Now iterating over %s summands"%((p-1)*p_power))907verbose_level = get_verbose()908count_verb = 0909for j in range(p_power):910s = K(0)911if verbose_level >= 2 and j/p_power*100 > count_verb + 3:912verbose("%.2f percent done"%(float(j)/p_power*100))913count_verb += 3914for a in range(1,p):915b = teich[a] * gamma_power916s += teich[a]**eta * self.measure(b, n, padic_prec,quadratic_twist=D, sign = 1-2*(eta % 2)).lift()917L += s * one_plus_T_factor918one_plus_T_factor *= 1+T919gamma_power *= gamma920921verbose("the series before adjusting the precision is %s"%L)922# Now create series but with each coefficient truncated923# so it is proven correct:924K = Qp(p, padic_prec, print_mode='series')925R = PowerSeriesRing(K,'T',res_series_prec)926L = R(L,res_series_prec)927aj = L.list()928if len(aj) > 0:929aj = [aj[0].add_bigoh(padic_prec-2)] + [aj[j].add_bigoh(bounds[j]) for j in range(1,len(aj))]930L = R(aj,res_series_prec )931932L /= self._quotient_of_periods_to_twist(D)*self._E.real_components()933934self._set_series_in_cache(n, res_series_prec, D, eta, L)935936return L937938power_series = series939940941def is_ordinary(self):942r"""943Return True if the elliptic curve that this L-function is attached944to is ordinary.945946EXAMPLES::947948sage: L = EllipticCurve('11a').padic_lseries(5)949sage: L.is_ordinary()950True951"""952return True953954def is_supersingular(self):955r"""956Return True if the elliptic curve that this L function is attached957to is supersingular.958959EXAMPLES::960961sage: L = EllipticCurve('11a').padic_lseries(5)962sage: L.is_supersingular()963False964"""965return False966967def _c_bound(self):968r"""969A helper function not designed for direct use.970971It returns the maximal `p`-adic valuation of the possible denominators972of the modular symbols.973974EXAMPLES::975976sage: E = EllipticCurve('11a1')977sage: Lp = E.padic_lseries(5)978sage: Lp._c_bound()9791980sage: Lp = E.padic_lseries(17)981sage: Lp._c_bound()9820983984"""985try:986return self.__c_bound987except AttributeError:988pass989E = self._E990p = self._p991if E.galois_representation().is_irreducible(p):992ans = 0993else:994m = E.modular_symbol_space(sign=1)995b = m.boundary_map().codomain()996C = b._known_cusps() # all known, since computed the boundary map997ans = max([valuation(self.modular_symbol(a).denominator(), p) for a in C])998self.__c_bound = ans999return ans10001001def _prec_bounds(self, n, prec):1002r"""1003A helper function not designed for direct use.10041005It returns the `p`-adic precisions of the approximation1006to the `p`-adic L-function.10071008EXAMPLES::10091010sage: E = EllipticCurve('11a1')1011sage: Lp = E.padic_lseries(5)1012sage: Lp._prec_bounds(3,10)1013[+Infinity, 1, 1, 1, 1, 0, 0, 0, 0, 0]1014sage: Lp._prec_bounds(3,12)1015[+Infinity, 1, 1, 1, 1, 0, 0, 0, 0, 0, 0, 0]1016sage: Lp._prec_bounds(4,5)1017[+Infinity, 2, 2, 2, 2]1018sage: Lp._prec_bounds(15,10)1019[+Infinity, 13, 13, 13, 13, 12, 12, 12, 12, 12]10201021sage: Lp = E.padic_lseries(3)1022sage: Lp._prec_bounds(15,10)1023[+Infinity, 14, 14, 13, 13, 13, 13, 13, 13, 12]10241025"""1026p = self._p1027e = self._e_bounds(n-1, prec)1028c = self._c_bound()1029return [e[j] - c for j in range(len(e))]103010311032class pAdicLseriesSupersingular(pAdicLseries):1033def series(self, n=3, quadratic_twist = +1, prec=5, eta = 0):1034r"""1035Return the `n`-th approximation to the `p`-adic L-series as a1036power series in `T` (corresponding to `\gamma-1` with1037`\gamma=1+p` as a generator of `1+p\ZZ_p`). Each1038coefficient is an element of a quadratic extension of the `p`-adic1039number whose precision is probably correct.10401041Here the normalization of the `p`-adic L-series is chosen1042such that `L_p(E,1) = (1-1/\alpha)^2 L(E,1)/\Omega_E`1043where `\alpha` is the unit root of the characteristic1044polynomial of Frobenius on `T_pE` and `\Omega_E` is the1045Neron period of `E`.10461047INPUT:10481049- ``n`` - (default: 2) a positive integer1050- ``quadratic_twist`` - (default: +1) a fundamental discriminant of a1051quadratic field, coprime to the conductor of the curve1052- ``prec`` - (default: 5) maximal number of terms of the series to1053compute; to compute as many as possible just give a very large1054number for ``prec``; the result will still be correct.1055- ``eta`` (default: 0) an integer (specifying the power of the1056Teichmueller character on the group of roots of unity in1057`\ZZ_p^\times`)10581059ALIAS: power_series is identical to series.10601061EXAMPLES:1062A superingular example, where we must compute to higher precision to see anything::10631064sage: e = EllipticCurve('37a')1065sage: L = e.padic_lseries(3); L10663-adic L-series of Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field1067sage: L.series(2)1068O(T^3)1069sage: L.series(4) # takes a long time (several seconds)1070(O(3))*alpha + (O(3^2)) + ((O(3^-1))*alpha + (2*3^-1 + O(3^0)))*T + ((O(3^-1))*alpha + (2*3^-1 + O(3^0)))*T^2 + O(T^5)1071sage: L.alpha(2).parent()1072Univariate Quotient Polynomial Ring in alpha over 3-adic Field with capped1073relative precision 2 with modulus (1 + O(3^2))*x^2 + (3 + O(3^3))*x + (3 + O(3^3))1074"""1075n = ZZ(n)1076if n < 1:1077raise ValueError, "n (=%s) must be a positive integer"%n10781079# check if the conditions on quadratic_twist are satisfied1080D = ZZ(quadratic_twist)1081if D != 1:1082if eta != 0: raise NotImplementedError, "quadratic twists only implemented for the 0th Teichmueller component"1083if D % 4 == 0:1084d = D//41085if not d.is_squarefree() or d % 4 == 1:1086raise ValueError, "quadratic_twist (=%s) must be a fundamental discriminant of a quadratic field"%D1087else:1088if not D.is_squarefree() or D % 4 != 1:1089raise ValueError, "quadratic_twist (=%s) must be a fundamental discriminant of a quadratic field"%D1090if gcd(D,self._E.conductor())!= 1:1091for ell in prime_divisors(D):1092if valuation(self._E.conductor(),ell) > valuation(D,ell) :1093raise ValueError, "can not twist a curve of conductor (=%s) by the quadratic twist (=%s)."%(self._E.conductor(),D)10941095p = self._p1096if p == 2 and self._normalize :1097print 'Warning : for p == 2 the normalization might not be correct !'1098eta = ZZ(eta) % (p-1)10991100prec = min(p**(n-1), prec)1101bounds = self._prec_bounds(n,prec)1102padic_prec = max(sum(bounds[1:],[])) + 51103verbose("using p-adic precision of %s"%padic_prec)1104ans = self._get_series_from_cache(n, prec, quadratic_twist,eta)1105if not ans is None:1106verbose("found series in cache")1107return ans11081109alpha = self.alpha(prec=padic_prec)1110K = alpha.parent()1111gamma = 1 + p1112R = PowerSeriesRing(K,'T',prec)1113T = R(R.gen(), prec)1114L = R(0)1115one_plus_T_factor = R(1)1116gamma_power = 11117teich = self.teichmuller(padic_prec)11181119verbose("Now iterating over %s summands"%((p-1)*p**(n-1)))1120verbose_level = get_verbose()1121count_verb = 01122for j in range(p**(n-1)):1123s = K(0)1124if verbose_level >= 2 and j/p**(n-1)*100 > count_verb + 3:1125verbose("%.2f percent done"%(float(j)/p**(n-1)*100))1126count_verb += 31127for a in range(1,p):1128b = teich[a] * gamma_power1129s += teich[a]**eta * self.measure(b, n, padic_prec,quadratic_twist=D, sign=1-2*(eta % 2))1130L += s * one_plus_T_factor1131one_plus_T_factor *= 1+T1132gamma_power *= gamma11331134# Now create series but with each coefficient truncated1135# so it is proven correct:1136L = R(L,prec)1137aj = L.list()1138if len(aj) > 0:1139bj = [aj[0][0].add_bigoh(padic_prec-2) + alpha * aj[0][1].add_bigoh(padic_prec-2)]1140bj += [aj[j][0].add_bigoh(bounds[j][0]) + alpha * aj[j][1].add_bigoh(bounds[j][1]) for j in range(1,len(aj))]1141L = R(bj, prec)1142L /= self._quotient_of_periods_to_twist(D)*self._E.real_components()1143self._set_series_in_cache(n, prec, quadratic_twist, eta, L)1144return L11451146power_series = series11471148def is_ordinary(self):1149r"""1150Return True if the elliptic curve that this L-function is attached1151to is ordinary.11521153EXAMPLES::11541155sage: L = EllipticCurve('11a').padic_lseries(19)1156sage: L.is_ordinary()1157False1158"""1159return False11601161def is_supersingular(self):1162r"""1163Return True if the elliptic curve that this L function is attached1164to is supersingular.11651166EXAMPLES::11671168sage: L = EllipticCurve('11a').padic_lseries(19)1169sage: L.is_supersingular()1170True1171"""1172return True11731174def _prec_bounds(self, n,prec):1175r"""1176A helper function not designed for direct use.11771178It returns the `p`-adic precisions of the approximation1179to the `p`-adic L-function.11801181EXAMPLES::11821183sage: E = EllipticCurve('11a1')1184sage: Lp = E.padic_lseries(19)1185sage: Lp._prec_bounds(3,5)1186[[+Infinity, +Infinity], [-1, -1], [-1, -1], [-1, -1], [-1, -1]]1187sage: Lp._prec_bounds(2,5)1188[[+Infinity, +Infinity], [-1, -2], [-1, -2], [-1, -2], [-1, -2]]1189sage: Lp._prec_bounds(10,5)1190[[+Infinity, +Infinity], [3, 2], [3, 2], [3, 2], [3, 2]]1191"""1192p = self._p1193e = self._e_bounds(n-1,prec)1194c0 = ZZ(n+2)/21195c1 = ZZ(n+3)/21196return [[infinity,infinity]] + [[(e[j] - c0).floor(), (e[j] - c1).floor()] for j in range(1,len(e))]119711981199def Dp_valued_series(self, n=3, quadratic_twist = +1, prec=5):1200r"""1201Returns a vector of two components which are p-adic power series.1202The answer v is such that12031204`(1-\varphi)^{-2}\cdot L_p(E,T) =` ``v[1]`` `\cdot \omega +` ``v[2]`` `\cdot \varphi(\omega)`12051206as an element of the Dieudonne module `D_p(E) = H^1_{dR}(E/\QQ_p)` where1207`\omega` is the invariant differential and `\varphi` is the Frobenius on `D_p(E)`.1208According to the `p`-adic Birch and Swinnerton-Dyer1209conjecture [BP] this function has a zero of order1210rank of `E(\QQ)` and it's leading term is contains the order of1211the Tate-Shafarevich group, the Tamagawa numbers, the order of the1212torsion subgroup and the `D_p`-valued `p`-adic regulator.12131214INPUT:12151216- ``n`` - (default: 3) a positive integer1217- ``prec`` - (default: 5) a positive integer12181219REFERENCE:12201221- [BP] Dominique Bernardi and Bernadette Perrin-Riou,1222Variante `p`-adique de la conjecture de Birch et1223Swinnerton-Dyer (le cas supersingulier), C. R. Acad. Sci. Paris,1224Ser I. Math, 317 (1993), no 3, 227-232.12251226EXAMPLES::12271228sage: E = EllipticCurve('14a')1229sage: L = E.padic_lseries(5)1230sage: L.Dp_valued_series(4) # long time (9s on sage.math, 2011)1231(1 + 4*5 + 4*5^3 + O(5^4) + (4 + O(5))*T + (1 + O(5))*T^2 + (4 + O(5))*T^3 + (2 + O(5))*T^4 + O(T^5), O(5^4) + O(5)*T + O(5)*T^2 + O(5)*T^3 + (2 + O(5))*T^4 + O(T^5))1232"""1233E = self._E1234p = self._p1235lps = self.series(n, quadratic_twist=quadratic_twist, prec=prec)12361237# now split up the series in two lps = G + H * alpha1238R = lps.base_ring().base_ring() # Qp1239QpT , T = PowerSeriesRing(R,'T',prec).objgen()1240G = QpT([lps[n][0] for n in range(0,lps.prec())], prec)1241H = QpT([lps[n][1] for n in range(0,lps.prec())], prec)12421243# now compute phi1244phi = matrix.matrix([[0,-1/p],[1,E.ap(p)/p]])1245lpv = vector([G + (E.ap(p))*H , - R(p) * H ]) # this is L_p1246eps = (1-phi)**(-2)1247resu = lpv*eps.transpose()1248return resu124912501251@rename_keyword(deprecated='Sage version 4.6', method="algorithm")1252def frobenius(self, prec=20, algorithm = "mw"):1253r"""1254This returns a geometric Frobenius `\varphi` on the Diedonne module `D_p(E)`1255with respect to the basis `\omega`, the invariant differential, and `\eta=x\omega`.1256It satisfies `\varphi^2 - a_p/p\, \varphi + 1/p = 0`.12571258INPUT:12591260- ``prec`` - (default: 20) a positive integer12611262- ``algorithm`` - either 'mw' (default) for Monsky-Washintzer1263or 'approx' for the algorithm described by Bernardi and Perrin-Riou1264(much slower and not fully tested)126512661267EXAMPLES::12681269sage: E = EllipticCurve('14a')1270sage: L = E.padic_lseries(5)1271sage: phi = L.frobenius(5)1272sage: phi1273[ 2 + 5^2 + 5^4 + O(5^5) 3*5^-1 + 3 + 5 + 4*5^2 + 5^3 + O(5^4)]1274[ 3 + 3*5^2 + 4*5^3 + 3*5^4 + O(5^5) 3 + 4*5 + 3*5^2 + 4*5^3 + 3*5^4 + O(5^5)]1275sage: -phi^21276[5^-1 + O(5^4) O(5^4)]1277[ O(5^5) 5^-1 + O(5^4)]1278"""1279E = self._E1280p = self._p1281if algorithm != "mw" and algorithm !="approx":1282raise ValueError, "Unknown algorithm %s."%algorithm1283if algorithm == "approx":1284return self.__phi_bpr(prec=prec)1285if p < 4 and algorithm == "mw":1286print "Warning: If this fails try again using algorithm=\"approx\""1287Ew = E.integral_short_weierstrass_model()1288adjusted_prec = monsky_washnitzer.adjusted_prec(p, prec)1289modprecring = Integers(p**adjusted_prec)1290output_ring = Qp(p, prec)1291R, x = PolynomialRing(modprecring, 'x').objgen()1292Q = x**3 + modprecring(Ew.a4()) * x + modprecring(Ew.a6())1293trace = Ew.ap(p)1294fr = monsky_washnitzer.matrix_of_frobenius(Q, p, adjusted_prec, trace)1295fr = matrix.matrix(output_ring,2,2,fr)12961297# return a vector for PARI's ellchangecurve to pass from e1 to e21298def isom(e1,e2):1299if not e1.is_isomorphic(e2):1300raise ValueError, "Curves must be isomorphic."1301usq = (e1.discriminant()/e2.discriminant()).nth_root(6)1302u = usq.sqrt()1303s = (u * e2.a1() - e1.a1() )/ZZ(2)1304r = (usq * e2.a2() - e1.a2() + s**2 + e1.a1()*s)/ZZ(3)1305t = (u**3 * e2.a3() - e1.a3() - e1.a1()*r)/ZZ(2)1306return [u,r,s,t]13071308v = isom(E,Ew)1309u = v[0]1310r = v[1]13111312# change basis1313A = matrix.matrix([[u,-r/u],[0,1/u]])1314frn = A * fr * A**(-1)1315return 1/p*frn1316131713181319def __phi_bpr(self, prec=0):1320r"""1321This returns a geometric Frobenius `\varphi` on the Dieudonne module `D_p(E)`1322with respect to the basis `\omega`, the invariant differential, and `\eta=x\omega`.1323It satisfies `\varphi^2 - a_p/p\, \varphi + 1/p = 0`.13241325The algorithm used here is described in bernardi-perrin-riou on page 232.13261327.. note: Warning. This function has not been sufficiently tested. It is very slow.13281329EXAMPLES::13301331sage: E = EllipticCurve('11a1')1332sage: lp = E.padic_lseries(19)1333sage: lp.frobenius(prec=1,algorithm="approx") #indirect doctest1334[ O(19^0) 4*19^-1 + O(19^0)]1335[ 14 + O(19) O(19^0)]13361337sage: E = EllipticCurve('17a1')1338sage: lp = E.padic_lseries(3)1339sage: lp.frobenius(prec=3,algorithm="approx")1340[ O(3) 2*3^-1 + 2 + O(3)]1341[ 1 + O(3^2) O(3)]1342sage: lp.frobenius(prec=5,algorithm="approx")1343[ 3 + O(3^2) 2*3^-1 + 2 + 3 + O(3^2)]1344[ 1 + 2*3^2 + O(3^3) 2*3 + O(3^2)]134513461347"""1348E = self._E1349p = self._p1350if prec > 10:1351print "Warning: Very large value for the precision."1352if prec == 0:1353prec = floor((log(10000)/log(p)))1354verbose("prec set to %s"%prec)1355eh = E.formal()1356om = eh.differential(prec = p**prec+3)1357verbose("differential computed")1358xt = eh.x(prec=p**prec + 3)1359et = xt*om1360# c_(p^k) = cs[k] d...1361cs = [om[p**k-1] for k in range(0,prec+1)]1362ds = [et[p**k-1] for k in range(0,prec+1)]1363delta = 01364dpr = 01365gamma = 01366dga = 01367for k in range(1,prec+1):1368# this is the equation eq[0]*x+eq[1]*y+eq[2] == 01369# such that delta_ = delta + d^dpr*x ...1370eq = [(p**dpr*cs[k]) % p**k,(-p**dga*ds[k]) % p**k , (delta*cs[k]-gamma*ds[k]-cs[k-1]) % p**k ]1371verbose("valuations : %s"%([x.valuation(p) for x in eq]))1372v = min([x.valuation(p) for x in eq])1373if v == infinity:1374verbose("no new information at step k=%s"%k)1375else:1376eq = [ZZ(x/p**v) for x in eq]1377verbose("renormalised eq mod p^%s is now %s"%(k-v,eq))1378if eq[0].valuation(p) == 0:1379l = min(eq[1].valuation(p),k-v)1380if l == 0:1381verbose("not uniquely determined at step k=%s"%k)1382else:1383ainv = eq[0].inverse_mod(p**l)1384delta = delta - eq[2]*ainv*p**dpr1385dpr = dpr + l1386delta = delta % p**dpr1387verbose("delta_prec increased to %s\n delta is now %s"%(dpr,delta))1388elif eq[1].valuation(p) == 0:1389l = min(eq[0].valuation(p),k-v)1390ainv = eq[1].inverse_mod(p**l)1391gamma = gamma - eq[2]*ainv*p**dga1392dga = dga + l1393gamma = gamma % p**dga1394verbose("gamma_prec increased to %s\n gamma is now %s"%(dga,gamma))1395else:1396raise RuntimeError, "Bug: no delta or gamma can exist"13971398# end of approximation of delta and gamma1399R = Qp(p,max(dpr,dga)+1)1400delta = R(delta,absprec=dpr)1401gamma = R(gamma,absprec=dga)1402verbose("result delta = %s\n gamma = %s\n check : %s"%(delta,gamma, [Qp(p,k)(delta * cs[k] - gamma * ds[k] - cs[k-1]) for k in range(1,prec+1)] ))1403a = delta1404c = -gamma1405d = E.ap(p) - a1406b = (-1/p+a*d)/c1407phi = matrix.matrix([[a,b],[c,d]])1408return phi140914101411def bernardi_sigma_function(self, prec=20):1412r"""1413Return the `p`-adic sigma function of Bernardi in terms of `z = log(t)`.1414This is the same as ``padic_sigma`` with ``E2 = 0``.14151416EXAMPLES::14171418sage: E = EllipticCurve('14a')1419sage: L = E.padic_lseries(5)1420sage: L.bernardi_sigma_function(5) # Todo: some sort of consistency check!?1421z + 1/24*z^3 + 29/384*z^5 - 8399/322560*z^7 - 291743/92897280*z^9 - 4364831/5225472*z^10 + 2172371753/955514880*z^11 - 17875714529/6897623040*z^12 + 2839176621047/1605264998400*z^13 + 32012675789849/10042939146240*z^14 - 367444910151047/89894839910400*z^15 + 973773806885959/241030539509760*z^16 - 33997971208432501/17259809262796800*z^17 - 10331978660756704339/842918229599846400*z^18 + 18601407947897364480389/950670294194847744000*z^19 - 118837570440101901119321/8071784966648129126400*z^20 + O(z^21)1422"""1423E = self._E1424p = self._p14251426Eh = E.formal()1427lo = Eh.log(prec + 5)1428F = lo.reversion()14291430S = LaurentSeriesRing(QQ,'z')1431z = S.gen()1432F = F(z)1433xofF = Eh.x(prec + 2)(F)1434#r = ( E.a1()**2 + 4*E.a2() ) / ZZ(12)1435g = (1/z**2 - xofF ).power_series()1436h = g.integral().integral()1437sigma_of_z = z.power_series() * h.exp()14381439return sigma_of_z144014411442def Dp_valued_height(self,prec=20):1443r"""1444Returns the canonical `p`-adic height with values in the Dieudonne module `D_p(E)`.1445It is defined to be14461447`h_{\eta} \cdot \omega - h_{\omega} \cdot \eta`14481449where `h_{\eta}` is made out of the sigma function of Bernardi and1450`h_{\omega}` is `log_E^2`.1451The answer ``v`` is given as ``v[1]*omega + v[2]*eta``.1452The coordinates of ``v`` are dependent of the1453Weierstrass equation.14541455EXAMPLES::14561457sage: E = EllipticCurve('53a')1458sage: L = E.padic_lseries(5)1459sage: h = L.Dp_valued_height(7)1460sage: h(E.gens()[0])1461(3*5 + 5^2 + 2*5^3 + 3*5^4 + 4*5^5 + 5^6 + 5^7 + O(5^8), 5^2 + 4*5^4 + 2*5^7 + 3*5^8 + O(5^9))1462"""1463E = self._E1464p = self._p1465Ehat = E.formal()1466elog = Ehat.log(prec + Integer(3))14671468# we will have to do it properly with David Harvey's _multiply_point()1469n = arith.LCM(E.tamagawa_numbers())1470n = arith.LCM(n, E.Np(p)) # allowed here because E has good reduction at p14711472if p < 5:1473phi = self.frobenius(min(6,prec),algorithm="approx")1474else:1475phi = self.frobenius(prec+2,algorithm="mw")14761477def height(P,check=True):1478if P.is_finite_order():1479return Qp(p,prec)(0)1480if check:1481assert P.curve() == E, 'the point P must lie on the curve from which the height function was created'14821483Q = n * P1484tt = - Q[0]/Q[1]1485R = Qp(p,prec+5)1486tt = R(tt)1487zz = elog(tt)14881489homega = -zz**2/n**214901491eQ = denominator(Q[1])/denominator(Q[0])1492si = self.bernardi_sigma_function(prec=prec+4)1493heta = 2 * log(si(zz)/eQ) / n**214941495R = Qp(p,prec)14961497return vector([-R(heta),R(homega)])14981499return height1500150115021503def Dp_valued_regulator(self,prec=20,v1=0,v2=0):1504r"""1505Returns the canonical `p`-adic regulator with values in the Dieudonne module `D_p(E)`1506as defined by Perrin-Riou using the `p`-adic height with values in `D_p(E)`.1507The result is written in the basis `\omega`, `\varphi(\omega)`, and hence the1508coordinates of the result are independent of the chosen Weierstrass equation.15091510NOTE: The definition here is corrected with respect to Perrin-Riou's article [PR]. See1511[SW].151215131514REFERENCES:15151516- [PR] Perrin Riou, Arithmetique des courbes elliptiques a reduction supersinguliere en `p`,1517Experiment. Math. 12 (2003), no. 2, 155-186.15181519- [SW] William Stein and Christian Wuthrich, Computations About Tate-Shafarevich Groups1520using Iwasawa theory, preprint 2009.15211522EXAMPLES::15231524sage: E = EllipticCurve('43a')1525sage: L = E.padic_lseries(7)1526sage: L.Dp_valued_regulator(7)1527(5*7 + 6*7^2 + 4*7^3 + 4*7^4 + 7^5 + 4*7^7 + O(7^8), 4*7^2 + 2*7^3 + 3*7^4 + 7^5 + 6*7^6 + 4*7^7 + O(7^8))1528"""15291530p = self._p1531E = self._E15321533h = self.Dp_valued_height(prec=prec)15341535# this is the height_{v} (P) for a v in D_p1536def hv(vec,P):1537hP = h(P)1538return - vec[0]*hP[1] +vec[1]*hP[0]15391540# def hvpairing(vec,P,Q):1541# return (hv(vec, P+Q) - hv(vec,P)-hv(vec,Q))/21542K = Qp(p, prec)15431544if v1 ==0 and v2 ==0 :1545v1 = vector([K(0),K(1)]) # that is eta1546v2 = vector([K(-1),K(1)]) # and this is eta-omega.1547# the rest should not depend on this choice1548# as long as it is outside Q_p * omega15491550rk = E.rank()1551if rk == 0:1552return vector([K(1),K(0)])155315541555basis = E.gens()15561557def regv(vec):1558M = matrix.matrix(K,rk,rk,0)1559point_height = [hv(vec,P) for P in basis]1560for i in range(rk):1561for j in range(i+1, rk):1562M[i, j] = M[j, i] = (hv(vec,basis[i] + basis[j])- point_height[i] - point_height[j] )/21563for i in range(rk):1564M[i,i] = point_height[i]15651566return M.determinant()156715681569def Dp_pairing(vec1,vec2):1570return (vec1[0]*vec2[1]-vec1[1]*vec2[0])15711572omega_vec = vector([K(1),K(0)])15731574# note the correction here with respect to Perrin-Riou's definition.1575# only this way the result will be independent of the choice of v1 and v2.1576reg1 = regv(v1)/Dp_pairing(omega_vec,v1)**(rk-1)15771578reg2 = regv(v2)/Dp_pairing(omega_vec,v2)**(rk-1)157915801581# the regulator in the basis omega,eta1582reg_oe = (reg1 * v2 - reg2 * v1 ) / Dp_pairing(v2,v1)15831584if p < 5:1585phi = self.frobenius(min(6,prec),algorithm="approx")1586else:1587phi = self.frobenius(prec+2,algorithm="mw")15881589c = phi[1,0] # this is the 'period' [omega,phi(omega)]1590a = phi[0,0]15911592return vector([reg_oe[0] - a/c*reg_oe[1],reg_oe[1]/c])159315941595