Path: blob/master/src/sage/schemes/elliptic_curves/formal_group.py
8821 views
# -*- coding: utf-8 -*-1r"""2Formal groups of elliptic curves34AUTHORS:56- William Stein: original implementations78- David Harvey: improved asymptotics of some methods910- Nick Alexander: separation from ell_generic.py, bugfixes and11docstrings12"""1314from sage.structure.sage_object import SageObject1516import sage.misc.misc as misc17import sage.rings.all as rings18from sage.rings.all import O192021class EllipticCurveFormalGroup(SageObject):22r"""23The formal group associated to an elliptic curve.24"""25def __init__(self, E):26"""27EXAMPLE::2829sage: E = EllipticCurve('11a')30sage: F = E.formal_group(); F31Formal Group associated to the Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field32sage: F == loads(dumps(F))33True34"""35self.__E = E3637def __cmp__(self, other):38"""39Compare self and other.4041TESTS::4243sage: E = EllipticCurve('35a')44sage: F1 = E.formal_group()45sage: F2 = E.formal_group()46sage: F1 == F247True48"""49c = cmp(type(self), type(other))50if c: return c51return cmp(self.__E, other.__E)5253def _repr_(self):54"""55EXAMPLE::5657sage: E = EllipticCurve('43a')58sage: F = E.formal_group()59sage: F._repr_()60'Formal Group associated to the Elliptic Curve defined by y^2 + y = x^3 + x^2 over Rational Field'61"""62return "Formal Group associated to the %s"%self.__E6364def curve(self):65r"""66The elliptic curve this formal group is associated to.6768EXAMPLES::6970sage: E = EllipticCurve("37a")71sage: F = E.formal_group()72sage: F.curve()73Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field74"""75return self.__E7677def w(self, prec=20):78r"""79The formal group power series w.8081INPUT:828384- ``prec`` - integer (default 20)858687OUTPUT: a power series with given precision8889DETAILS: Return the formal power series9091.. math::9293w(t) = t^3 + a_1 t^4 + (a_2 + a_1^2) t^5 + \cdots949596to precision `O(t^{prec})` of Proposition IV.1.1 of97[Silverman AEC1]. This is the formal expansion of98`w = -1/y` about the formal parameter `t = -x/y` at99`\\infty`.100101The result is cached, and a cached version is returned if102possible.103104.. warning::105106The resulting power series will have precision prec, but107its parent PowerSeriesRing will have default precision 20108(or whatever the default default is).109110ALGORITHM: Uses Newton's method to solve the elliptic curve111equation at the origin. Complexity is roughly `O(M(n))`112where `n` is the precision and `M(n)` is the time113required to multiply polynomials of length `n` over the114coefficient ring of `E`.115116AUTHOR:117118- David Harvey (2006-09-09): modified to use Newton's119method instead of a recurrence formula.120121EXAMPLES::122123sage: e = EllipticCurve([0, 0, 1, -1, 0])124sage: e.formal_group().w(10)125t^3 + t^6 - t^7 + 2*t^9 + O(t^10)126127Check that caching works::128129sage: e = EllipticCurve([3, 2, -4, -2, 5])130sage: e.formal_group().w(20)131t^3 + 3*t^4 + 11*t^5 + 35*t^6 + 101*t^7 + 237*t^8 + 312*t^9 - 949*t^10 - 10389*t^11 - 57087*t^12 - 244092*t^13 - 865333*t^14 - 2455206*t^15 - 4366196*t^16 + 6136610*t^17 + 109938783*t^18 + 688672497*t^19 + O(t^20)132sage: e.formal_group().w(7)133t^3 + 3*t^4 + 11*t^5 + 35*t^6 + O(t^7)134sage: e.formal_group().w(35)135t^3 + 3*t^4 + 11*t^5 + 35*t^6 + 101*t^7 + 237*t^8 + 312*t^9 - 949*t^10 - 10389*t^11 - 57087*t^12 - 244092*t^13 - 865333*t^14 - 2455206*t^15 - 4366196*t^16 + 6136610*t^17 + 109938783*t^18 + 688672497*t^19 + 3219525807*t^20 + 12337076504*t^21 + 38106669615*t^22 + 79452618700*t^23 - 33430470002*t^24 - 1522228110356*t^25 - 10561222329021*t^26 - 52449326572178*t^27 - 211701726058446*t^28 - 693522772940043*t^29 - 1613471639599050*t^30 - 421817906421378*t^31 + 23651687753515182*t^32 + 181817896829144595*t^33 + 950887648021211163*t^34 + O(t^35)136"""137prec = max(prec, 0)138k = self.curve().base_ring()139140try:141# Try cached version142w = self.__w143cached_prec = w.prec()144R = w.parent()145except AttributeError:146# No cached version available147R = rings.PowerSeriesRing(k, "t")148w = R([k(0), k(0), k(0), k(1)], 4)149cached_prec = 4150self.__w = w151152if prec < cached_prec:153return R(w, prec)154155# We use the following iteration, which doubles the precision156# at each step:157#158# z^3 - a_3 w^2 - a_4 z w^2 - 2 a_6 w^3159# w' = -----------------------------------------------------160# 1 - a_1 z - a_2 z^2 - 2 a_3 w - 2 a_4 z w - 3 a_6 w^2161162a1, a2, a3, a4, a6 = self.curve().ainvs()163current_prec = cached_prec164w = w.truncate() # work with polynomials instead of power series165166numerator_const = w.parent()([0, 0, 0, 1]) # z^3167denominator_const = w.parent()([1, -a1, -a2]) # 1 - a_1 z - a_2 z^2168169last_prec = 0170for next_prec in misc.newton_method_sizes(prec):171if next_prec > current_prec:172if w.degree() - 1 > last_prec:173# Here it's best to throw away some precision to get us174# in sync with the sizes recommended by175# newton_method_sizes(). This is especially counter-176# intuitive when we throw away almost half of our177# cached data!178179# todo: this might not actually be true, depending on180# the overhead of truncate(), which is currently very181# high e.g. for NTL based polynomials (but this is182# slated to be fixed...)183184w = w.truncate(last_prec)185186w_squared = w.square()187w_cubed = (w_squared * w).truncate(next_prec)188189numerator = numerator_const \190- a3 * w_squared \191- a4 * w_squared.shift(1) \192- (2*a6) * w_cubed193194denominator = denominator_const \195- (2*a3) * w \196- (2*a4) * w.shift(1) \197- (3*a6) * w_squared198199# todo: this is quite inefficient, because it gets200# converted to a power series, then the power series201# inversion works with polynomials again, and then202# it gets converted *back* to a power series, and203# then we convert it to a polynomial again! That's four204# useless conversions!!!205206inverse = ~R(denominator, prec=next_prec)207inverse = inverse.truncate(next_prec)208209w = (numerator * inverse).truncate(next_prec)210211last_prec = next_prec212213# convert back to power series214w = R(w, prec)215self.__w = (prec, w)216return self.__w[1]217218def x(self, prec=20):219r"""220Return the formal series `x(t) = t/w(t)` in terms of the221local parameter `t = -x/y` at infinity.222223INPUT:224225226- ``prec`` - integer (default 20)227228229OUTPUT: a Laurent series with given precision230231DETAILS: Return the formal series232233.. math::234235x(t) = t^{-2} - a_1 t^{-1} - a_2 - a_3 t - \cdots236237238to precision `O(t^{prec})` of page 113 of [Silverman AEC1].239240.. warning::241242The resulting series will have precision prec, but its243parent PowerSeriesRing will have default precision 20 (or244whatever the default default is).245246EXAMPLES::247248sage: EllipticCurve([0, 0, 1, -1, 0]).formal_group().x(10)249t^-2 - t + t^2 - t^4 + 2*t^5 - t^6 - 2*t^7 + 6*t^8 - 6*t^9 + O(t^10)250"""251prec = max(prec, 0)252y = self.y(prec)253t = y.parent().gen()254return -t*y + O(t**prec)255256def y(self, prec=20):257r"""258Return the formal series `y(t) = -1/w(t)` in terms of the259local parameter `t = -x/y` at infinity.260261INPUT:262263264- ``prec`` - integer (default 20)265266267OUTPUT: a Laurent series with given precision268269DETAILS: Return the formal series270271.. math::272273y(t) = - t^{-3} + a_1 t^{-2} + a_2 t + a_3 + \cdots274275276to precision `O(t^{prec})` of page 113 of [Silverman AEC1].277278The result is cached, and a cached version is returned if279possible.280281.. warning::282283The resulting series will have precision prec, but its284parent PowerSeriesRing will have default precision 20 (or285whatever the default default is).286287EXAMPLES::288289sage: EllipticCurve([0, 0, 1, -1, 0]).formal_group().y(10)290-t^-3 + 1 - t + t^3 - 2*t^4 + t^5 + 2*t^6 - 6*t^7 + 6*t^8 + 3*t^9 + O(t^10)291"""292prec = max(prec,0)293try:294pr, y = self.__y295except AttributeError:296pr = -1297if prec <= pr:298t = y.parent().gen()299return y + O(t**prec)300w = self.w(prec+6) # XXX why 6?301t = w.parent().gen()302y = -(w**(-1)) + O(t**prec)303self.__y = (prec, y)304return self.__y[1]305306def differential(self, prec=20):307r"""308Returns the power series `f(t) = 1 + \cdots` such that309`f(t) dt` is the usual invariant differential310`dx/(2y + a_1 x + a_3)`.311312INPUT:313314315- ``prec`` - nonnegative integer (default 20), answer316will be returned `O(t^{\mathrm{prec}})`317318319OUTPUT: a power series with given precision320321DETAILS: Return the formal series322323.. math::324325f(t) = 1 + a_1 t + ({a_1}^2 + a_2) t^2 + \cdots326327328to precision `O(t^{prec})` of page 113 of [Silverman AEC1].329330The result is cached, and a cached version is returned if331possible.332333.. warning::334335The resulting series will have precision prec, but its336parent PowerSeriesRing will have default precision 20 (or337whatever the default default is).338339EXAMPLES::340341sage: EllipticCurve([-1, 1/4]).formal_group().differential(15)3421 - 2*t^4 + 3/4*t^6 + 6*t^8 - 5*t^10 - 305/16*t^12 + 105/4*t^14 + O(t^15)343sage: EllipticCurve(Integers(53), [-1, 1/4]).formal_group().differential(15)3441 + 51*t^4 + 14*t^6 + 6*t^8 + 48*t^10 + 24*t^12 + 13*t^14 + O(t^15)345346AUTHOR:347348- David Harvey (2006-09-10): factored out of log349"""350prec = max(prec,0)351try:352cached_prec, omega = self.__omega353except AttributeError:354cached_prec = -1355if prec <= cached_prec:356return omega.add_bigoh(prec)357358a = self.curve().ainvs()359x = self.x(prec+1)360y = self.y(prec+1)361xprime = x.derivative()362g = xprime / (2*y + a[0]*x + a[2])363self.__omega = (prec, g.power_series().add_bigoh(prec))364return self.__omega[1]365366def log(self, prec=20):367r"""368Returns the power series `f(t) = t + \cdots` which is an369isomorphism to the additive formal group.370371Generally this only makes sense in characteristic zero, although372the terms before `t^p` may work in characteristic373`p`.374375INPUT:376377378- ``prec`` - nonnegative integer (default 20)379380381OUTPUT: a power series with given precision382383EXAMPLES::384385sage: EllipticCurve([-1, 1/4]).formal_group().log(15)386t - 2/5*t^5 + 3/28*t^7 + 2/3*t^9 - 5/11*t^11 - 305/208*t^13 + O(t^15)387388AUTHORS:389390- David Harvey (2006-09-10): rewrote to use differential391"""392return self.differential(prec-1).integral().add_bigoh(prec)393394def inverse(self, prec=20):395r"""396The formal group inverse law i(t), which satisfies F(t, i(t)) = 0.397398INPUT:399400401- ``prec`` - integer (default 20)402403404OUTPUT: a power series with given precision405406DETAILS: Return the formal power series407408.. math::409410i(t) = - t + a_1 t^2 + \cdots411412413to precision `O(t^{prec})` of page 114 of [Silverman AEC1].414415The result is cached, and a cached version is returned if416possible.417418.. warning::419420The resulting power series will have precision prec, but421its parent PowerSeriesRing will have default precision 20422(or whatever the default default is).423424EXAMPLES::425426sage: P.<a1, a2, a3, a4, a6> = ZZ[]427sage: E = EllipticCurve(list(P.gens()))428sage: i = E.formal_group().inverse(6); i429-t - a1*t^2 - a1^2*t^3 + (-a1^3 - a3)*t^4 + (-a1^4 - 3*a1*a3)*t^5 + O(t^6)430sage: F = E.formal_group().group_law(6)431sage: F(i.parent().gen(), i)432O(t^6)433"""434prec = max(prec,0)435try:436pr, inv = self.__inverse437except AttributeError:438pr = -1439if prec <= pr:440t = inv.parent().gen()441return inv + O(t**prec)442x = self.x(prec)443y = self.y(prec)444a1, _, a3, _, _ = self.curve().ainvs()445inv = x / ( y + a1*x + a3) # page 114 of Silverman, AEC I446inv = inv.power_series().add_bigoh(prec)447self.__inverse = (prec, inv)448return inv449450def group_law(self, prec=10):451r"""452The formal group law.453454INPUT:455456457- ``prec`` - integer (default 10)458459460OUTPUT: a power series with given precision in R[['t1','t2']], where461the curve is defined over R.462463DETAILS: Return the formal power series464465.. math::466467F(t_1, t_2) = t_1 + t_2 - a_1 t_1 t_2 - \cdots468469470to precision `O(t1,t2)^{prec}` of page 115 of [Silverman AEC1].471472The result is cached, and a cached version is returned if possible.473474475AUTHORS:476477- Nick Alexander: minor fixes, docstring478479- Francis Clarke (2012-08): modified to use two-variable power series ring480481EXAMPLES::482483sage: e = EllipticCurve([1, 2])484sage: e.formal_group().group_law(6)485t1 + t2 - 2*t1^4*t2 - 4*t1^3*t2^2 - 4*t1^2*t2^3 - 2*t1*t2^4 + O(t1, t2)^6486487sage: e = EllipticCurve('14a1')488sage: ehat = e.formal()489sage: ehat.group_law(3)490t1 + t2 - t1*t2 + O(t1, t2)^3491sage: ehat.group_law(5)492t1 + t2 - t1*t2 - 2*t1^3*t2 - 3*t1^2*t2^2 - 2*t1*t2^3 + O(t1, t2)^5493494sage: e = EllipticCurve(GF(7), [3, 4])495sage: ehat = e.formal()496sage: ehat.group_law(3)497t1 + t2 + O(t1, t2)^3498sage: F = ehat.group_law(7); F499t1 + t2 + t1^4*t2 + 2*t1^3*t2^2 + 2*t1^2*t2^3 + t1*t2^4 + O(t1, t2)^7500501TESTS::502503sage: R.<x,y,z> = GF(7)[[]]504sage: F(x, ehat.inverse()(x))5050 + O(x, y, z)^7506sage: F(x, y) == F(y, x)507True508sage: F(x, F(y, z)) == F(F(x, y), z)509True510511Let's ensure caching with changed precision is working::512513sage: e.formal_group().group_law(4)514t1 + t2 + O(t1, t2)^4515516Test for :trac:`9646`::517518sage: P.<a1, a2, a3, a4, a6> = PolynomialRing(ZZ, 5)519sage: E = EllipticCurve(list(P.gens()))520sage: F = E.formal().group_law(prec=5)521sage: t1, t2 = F.parent().gens()522sage: F(t1, 0)523t1 + O(t1, t2)^5524sage: F(0, t2)525t2 + O(t1, t2)^5526sage: F.coefficients()[t1*t2^2]527-a2528529"""530prec = max(prec,0)531if prec <= 0:532raise ValueError, "The precision must be positive."533534R = rings.PowerSeriesRing(self.curve().base_ring(), 2, 't1,t2')535t1, t2 = R.gens()536537if prec == 1:538return R(0)539elif prec == 2:540return t1 + t2 - self.curve().a1()*t1*t2541542try:543pr, F = self.__group_law544if prec <= pr:545return F.add_bigoh(prec)546except AttributeError:547pass548549w = self.w(prec+1)550lam = sum([w[n]*sum(t2**m * t1**(n-m-1) for m in range(n)) for n in range(3, prec+1)])551lam = lam.add_bigoh(prec)552nu = w(t1) - lam*t1553a1, a2, a3, a4, a6 = self.curve().ainvs()554lam2 = lam*lam555lam3 = lam2*lam556# note that the following formula differs from the one in Silverman page 119.557# See trac ticket 9646 for the explanation and justification.558t3 = -t1 - t2 - \559(a1*lam + a3*lam2 + a2*nu + 2*a4*lam*nu + 3*a6*lam2*nu)/ \560(1 + a2*lam + a4*lam2 + a6*lam3)561inv = self.inverse(prec)562563F = inv(t3).add_bigoh(prec)564self.__group_law = (prec, F)565return F566567def mult_by_n(self, n, prec=10):568"""569The formal 'multiplication by n' endomorphism `[n]`.570571INPUT:572573574- ``prec`` - integer (default 10)575576577OUTPUT: a power series with given precision578579DETAILS: Return the formal power series580581.. math::582583[n](t) = n t + \cdots584585586to precision `O(t^{prec})` of Proposition 2.3 of [Silverman587AEC1].588589.. warning::590591The resulting power series will have precision prec, but592its parent PowerSeriesRing will have default precision 20593(or whatever the default default is).594595AUTHORS:596597- Nick Alexander: minor fixes, docstring598599- David Harvey (2007-03): faster algorithm for char 0 field600case601602- Hamish Ivey-Law (2009-06): double-and-add algorithm for603non char 0 field case.604605- Tom Boothby (2009-06): slight improvement to double-and-add606607- Francis Clarke (2012-08): adjustments and simplifications using group_law608code as modified to yield a two-variable power series.609610EXAMPLES::611612sage: e = EllipticCurve([1, 2, 3, 4, 6])613sage: e.formal_group().mult_by_n(0, 5)614O(t^5)615sage: e.formal_group().mult_by_n(1, 5)616t + O(t^5)617618We verify an identity of low degree::619620sage: none = e.formal_group().mult_by_n(-1, 5)621sage: two = e.formal_group().mult_by_n(2, 5)622sage: ntwo = e.formal_group().mult_by_n(-2, 5)623sage: ntwo - none(two)624O(t^5)625sage: ntwo - two(none)626O(t^5)627628It's quite fast::629630sage: E = EllipticCurve("37a"); F = E.formal_group()631sage: F.mult_by_n(100, 20)632100*t - 49999950*t^4 + 3999999960*t^5 + 14285614285800*t^7 - 2999989920000150*t^8 + 133333325333333400*t^9 - 3571378571674999800*t^10 + 1402585362624965454000*t^11 - 146666057066712847999500*t^12 + 5336978000014213190385000*t^13 - 519472790950932256570002000*t^14 + 93851927683683567270392002800*t^15 - 6673787211563812368630730325175*t^16 + 320129060335050875009191524993000*t^17 - 45670288869783478472872833214986000*t^18 + 5302464956134111125466184947310391600*t^19 + O(t^20)633634TESTS::635636sage: F = EllipticCurve(GF(17), [1, 1]).formal_group()637sage: F.mult_by_n(10, 50) # long time (13s on sage.math, 2011)63810*t + 5*t^5 + 7*t^7 + 13*t^9 + t^11 + 16*t^13 + 13*t^15 + 9*t^17 + 16*t^19 + 15*t^23 + 15*t^25 + 2*t^27 + 10*t^29 + 8*t^31 + 15*t^33 + 6*t^35 + 7*t^37 + 9*t^39 + 10*t^41 + 5*t^43 + 4*t^45 + 6*t^47 + 13*t^49 + O(t^50)639640sage: F = EllipticCurve(GF(101), [1, 1]).formal_group()641sage: F.mult_by_n(100, 20)642100*t + O(t^20)643644sage: P.<a1, a2, a3, a4, a6> = PolynomialRing(ZZ, 5)645sage: E = EllipticCurve(list(P.gens()))646sage: E.formal().mult_by_n(2,prec=5)6472*t - a1*t^2 - 2*a2*t^3 + (a1*a2 - 7*a3)*t^4 + O(t^5)648649sage: E = EllipticCurve(QQ, [1,2,3,4,6])650sage: E.formal().mult_by_n(2,prec=5)6512*t - t^2 - 4*t^3 - 19*t^4 + O(t^5)652653654"""655if self.curve().base_ring().is_field() and self.curve().base_ring().characteristic() == 0 and n != 0:656# The following algorithm only works over a field of657# characteristic zero. I don't know whether something similar658# can be done over a general ring. It would be nice if it did,659# since it's much faster than using the formal group law.660# -- dmharvey661662# Create a "formal point" on the original curve E.663# Our answer only needs prec-1 coefficients (since lowest term664# is t^1), and x(t) = t^(-2) + ... and y(t) = t^(-3) + ...,665# so we only need x(t) mod t^(prec-3) and y(t) mod t^(prec-4)666x = self.x(prec-3)667y = self.y(prec-4)668R = x.parent() # the Laurent series ring over the base ring669X = self.curve().change_ring(R)670P = X(x, y)671672# and multiply it by n, using the group law on E673Q = n*P674675# express it in terms of the formal parameter676return -Q[0] / Q[1]677678679# Now the general case, not necessarily over a field.680681R = rings.PowerSeriesRing(self.curve().base_ring(), "t")682t = R.gen()683684if n == 1:685return t.add_bigoh(prec)686687if n == 0:688return R(0).add_bigoh(prec)689690if n == -1:691return R(self.inverse(prec))692693if n < 0:694return self.inverse(prec)(self.mult_by_n(-n, prec))695696F = self.group_law(prec)697698result = t699if n < 4:700for m in range(n - 1):701result = F(result, t)702return result703704# Double and add is faster than the naive method when n >= 4.705g = t706if n & 1:707result = g708else:709result = 0710n = n >> 1711712while n > 0:713g = F(g, g)714if n & 1:715result = F(result, g)716n = n >> 1717718return result719720def sigma(self, prec=10):721"""722EXAMPLE::723724sage: E = EllipticCurve('14a')725sage: F = E.formal_group()726sage: F.sigma(5)727t + 1/2*t^2 + (1/2*c + 1/3)*t^3 + (3/4*c + 3/4)*t^4 + O(t^5)728"""729a1,a2,a3,a4,a6 = self.curve().ainvs()730731k = self.curve().base_ring()732fl = self.log(prec)733R = rings.PolynomialRing(k,'c'); c = R.gen()734F = fl.reversion()735736S = rings.LaurentSeriesRing(R,'z')737c = S(c)738z = S.gen()739F = F(z + O(z**prec))740wp = self.x()(F)741e2 = 12*c - a1**2 - 4*a2742g = (1/z**2 - wp + e2/12).power_series()743h = g.integral().integral()744sigma_of_z = z.power_series() * h.exp()745746T = rings.PowerSeriesRing(R,'t')747fl = fl(T.gen()+O(T.gen()**prec))748sigma_of_t = sigma_of_z(fl)749return sigma_of_t750751752