def sigma(self, prec=10): """ EXAMPLES:: sage: E = EllipticCurve('14a') sage: F = E.formal_group() sage: F.sigma(5) t + 1/2*t^2 + (1/2*c + 1/3)*t^3 + (3/4*c + 3/4)*t^4 + O(t^5) """ a1,a2,a3,a4,a6 = self.curve().ainvs() k = self.curve().base_ring() fl = self.log(prec) R = rings.PolynomialRing(k,'c'); c = R.gen() F = fl.reverse() S = rings.LaurentSeriesRing(R,'z') c = S(c) z = S.gen() F = F(z + O(z**prec)) wp = self.x()(F) e2 = 12*c - a1**2 - 4*a2 g = (1/z**2 - wp + e2/12).power_series() h = g.integral().integral() sigma_of_z = z.power_series() * h.exp() T = rings.PowerSeriesRing(R,'t') fl = fl(T.gen()+O(T.gen()**prec)) sigma_of_t = sigma_of_z(fl) return sigma_of_t
def local_coordinates_at_infinity(self, prec=20, name='t'): """ For the genus `g` hyperelliptic curve `y^2 = f(x)`, return `(x(t), y(t))` such that `(y(t))^2 = f(x(t))`, where `t = x^g/y` is the local parameter at infinity INPUT: - ``prec`` -- desired precision of the local coordinates - ``name`` -- generator of the power series ring (default: ``t``) OUTPUT: `(x(t),y(t))` such that `y(t)^2 = f(x(t))` and `t = x^g/y` is the local parameter at infinity EXAMPLES:: sage: R.<x> = QQ['x'] sage: H = HyperellipticCurve(x^5-5*x^2+1) sage: x,y = H.local_coordinates_at_infinity(10) sage: x t^-2 + 5*t^4 - t^8 - 50*t^10 + O(t^12) sage: y t^-5 + 10*t - 2*t^5 - 75*t^7 + 50*t^11 + O(t^12) :: sage: R.<x> = QQ['x'] sage: H = HyperellipticCurve(x^3-x+1) sage: x,y = H.local_coordinates_at_infinity(10) sage: x t^-2 + t^2 - t^4 - t^6 + 3*t^8 + O(t^12) sage: y t^-3 + t - t^3 - t^5 + 3*t^7 - 10*t^11 + O(t^12) Note: if even degree model, just returns local coordinate above one point AUTHOR: - Jennifer Balakrishnan (2007-12) """ g = self.genus() pol = self.hyperelliptic_polynomials()[0] K = LaurentSeriesRing(self.base_ring(), name) t = K.gen() K.set_default_prec(prec + 2) L = PolynomialRing(K, 'x') x = L.gen() i = 0 w = (x**g / t)**2 - pol wprime = w.derivative(x) if pol.degree() == 2 * g + 1: x = t**-2 else: x = t**-1 for i in range((RR(log(prec + 2) / log(2))).ceil()): x = x - w(x) / wprime(x) y = x**g / t return x + O(t**(prec + 2)), y + O(t**(prec + 2))
def local_coordinates_at_nonweierstrass(self, P, prec=20, name='t'): """ For a non-Weierstrass point `P = (a,b)` on the hyperelliptic curve `y^2 = f(x)`, return `(x(t), y(t))` such that `(y(t))^2 = f(x(t))`, where `t = x - a` is the local parameter. INPUT: - ``P = (a, b)`` -- a non-Weierstrass point on self - ``prec`` -- desired precision of the local coordinates - ``name`` -- gen of the power series ring (default: ``t``) OUTPUT: `(x(t),y(t))` such that `y(t)^2 = f(x(t))` and `t = x - a` is the local parameter at `P` EXAMPLES:: sage: R.<x> = QQ['x'] sage: H = HyperellipticCurve(x^5-23*x^3+18*x^2+40*x) sage: P = H(1,6) sage: x,y = H.local_coordinates_at_nonweierstrass(P,prec=5) sage: x 1 + t + O(t^5) sage: y 6 + t - 7/2*t^2 - 1/2*t^3 - 25/48*t^4 + O(t^5) sage: Q = H(-2,12) sage: x,y = H.local_coordinates_at_nonweierstrass(Q,prec=5) sage: x -2 + t + O(t^5) sage: y 12 - 19/2*t - 19/32*t^2 + 61/256*t^3 - 5965/24576*t^4 + O(t^5) AUTHOR: - Jennifer Balakrishnan (2007-12) """ d = P[1] if d == 0: raise TypeError( "P = %s is a Weierstrass point. Use local_coordinates_at_weierstrass instead!" % P) pol = self.hyperelliptic_polynomials()[0] L = PowerSeriesRing(self.base_ring(), name) t = L.gen() L.set_default_prec(prec) K = PowerSeriesRing(L, 'x') pol = K(pol) x = K.gen() b = P[0] f = pol(t + b) for i in range((RR(log(prec) / log(2))).ceil()): d = (d + f / d) / 2 return t + b + O(t**(prec)), d + O(t**(prec))
def local_coordinates_at_weierstrass(self, P, prec=20, name='t'): """ For a finite Weierstrass point on the hyperelliptic curve y^2 = f(x), returns (x(t), y(t)) such that (y(t))^2 = f(x(t)), where t = y is the local parameter. INPUT: - P a finite Weierstrass point on self - prec: desired precision of the local coordinates - name: gen of the power series ring (default: 't') OUTPUT: (x(t),y(t)) such that y(t)^2 = f(x(t)) and t = y is the local parameter at P EXAMPLES: sage: R.<x> = QQ['x'] sage: H = HyperellipticCurve(x^5-23*x^3+18*x^2+40*x) sage: A = H(4,0) sage: x,y = H.local_coordinates_at_weierstrass(A,prec =5) sage: x 4 + 1/360*t^2 - 191/23328000*t^4 + 7579/188956800000*t^6 + O(t^7) sage: y t + O(t^7) sage: B = H(-5,0) sage: x,y = H.local_coordinates_at_weierstrass(B,prec = 5) sage: x -5 + 1/1260*t^2 + 887/2000376000*t^4 + 643759/1587898468800000*t^6 + O(t^7) sage: y t + O(t^7) AUTHOR: - Jennifer Balakrishnan (2007-12) """ if P[1] != 0: raise TypeError, "P = %s is not a finite Weierstrass point. Use local_coordinates_at_nonweierstrass instead!" % P pol = self.hyperelliptic_polynomials()[0] L = PowerSeriesRing(self.base_ring(), name) t = L.gen() L.set_default_prec(prec + 2) K = PowerSeriesRing(L, 'x') pol = K(pol) x = K.gen() b = P[0] g = pol / (x - b) c = b + 1 / g(b) * t**2 f = pol - t**2 fprime = f.derivative() for i in range((RR(log(prec + 2) / log(2))).ceil()): c = c - f(c) / fprime(c) return c + O(t**(prec + 2)), t + O(t**(prec + 2))
def sigma(self, prec=10): r""" Return the Weierstrass sigma function as a formal power series solution to the differential equation .. MATH:: \frac{d^2 \log \sigma}{dz^2} = - \wp(z) with initial conditions `\sigma(O)=0` and `\sigma'(O)=1`, expressed in the variable `t=\log_E(z)` of the formal group. INPUT: - ``prec`` - integer (default 10) OUTPUT: a power series with given precision Other solutions can be obtained by multiplication with a function of the form `\exp(c z^2)`. If the curve has good ordinary reduction at a prime `p` then there is a canonical choice of `c` that produces the canonical `p`-adic sigma function. To obtain that, please use ``E.padic_sigma(p)`` instead. See :meth:`~sage.schemes.elliptic_curves.ell_rational_field.EllipticCurve_rational_field.padic_sigma` EXAMPLES:: sage: E = EllipticCurve('14a') sage: F = E.formal_group() sage: F.sigma(5) t + 1/2*t^2 + 1/3*t^3 + 3/4*t^4 + O(t^5) """ a1, a2, a3, a4, a6 = self.curve().ainvs() k = self.curve().base_ring() fl = self.log(prec) F = fl.reverse() S = rings.LaurentSeriesRing(k, 'z') z = S.gen() F = F(z + O(z**prec)) wp = self.x()(F) + (a1**2 + 4 * a2) / 12 g = (1 / z**2 - wp).power_series() h = g.integral().integral() sigma_of_z = z.power_series() * h.exp() T = rings.PowerSeriesRing(k, 't') fl = fl(T.gen() + O(T.gen()**prec)) sigma_of_t = sigma_of_z(fl) return sigma_of_t
def local_coordinates_at_infinity(self, prec=20, name='t'): """ For the genus g hyperelliptic curve y^2 = f(x), returns (x(t), y(t)) such that (y(t))^2 = f(x(t)), where t = x^g/y is the local parameter at infinity INPUT: - prec: desired precision of the local coordinates - name: gen of the power series ring (default: 't') OUTPUT: (x(t),y(t)) such that y(t)^2 = f(x(t)) and t = x^g/y is the local parameter at infinity EXAMPLES: sage: R.<x> = QQ['x'] sage: H = HyperellipticCurve(x^5-5*x^2+1) sage: x,y = H.local_coordinates_at_infinity(10) sage: x t^-2 + 5*t^4 - t^8 - 50*t^10 + O(t^12) sage: y t^-5 + 10*t - 2*t^5 - 75*t^7 + 50*t^11 + O(t^12) sage: R.<x> = QQ['x'] sage: H = HyperellipticCurve(x^3-x+1) sage: x,y = H.local_coordinates_at_infinity(10) sage: x t^-2 + t^2 - t^4 - t^6 + 3*t^8 + O(t^12) sage: y t^-3 + t - t^3 - t^5 + 3*t^7 - 10*t^11 + O(t^12) AUTHOR: - Jennifer Balakrishnan (2007-12) """ g = self.genus() pol = self.hyperelliptic_polynomials()[0] K = LaurentSeriesRing(self.base_ring(), name) t = K.gen() K.set_default_prec(prec + 2) L = PolynomialRing(self.base_ring(), 'x') x = L.gen() i = 0 w = (x**g / t)**2 - pol wprime = w.derivative(x) x = t**-2 for i in range((RR(log(prec + 2) / log(2))).ceil()): x = x - w(x) / wprime(x) y = x**g / t return x + O(t**(prec + 2)), y + O(t**(prec + 2))
def y(self, prec=20): r""" Return the formal series `y(t) = -1/w(t)` in terms of the local parameter `t = -x/y` at infinity. INPUT: - ``prec`` - integer (default 20) OUTPUT: a Laurent series with given precision DETAILS: Return the formal series .. MATH:: y(t) = - t^{-3} + a_1 t^{-2} + a_2 t + a_3 + \cdots to precision `O(t^{prec})` of page 113 of [Silverman AEC1]. The result is cached, and a cached version is returned if possible. .. warning:: The resulting series will have precision prec, but its parent PowerSeriesRing will have default precision 20 (or whatever the default default is). EXAMPLES:: sage: EllipticCurve([0, 0, 1, -1, 0]).formal_group().y(10) -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) """ prec = max(prec,0) try: pr, y = self.__y except AttributeError: pr = -1 if prec <= pr: t = y.parent().gen() return y + O(t**prec) w = self.w(prec+6) # XXX why 6? t = w.parent().gen() y = -(w**(-1)) + O(t**prec) self.__y = (prec, y) return self.__y[1]
def inverse(self, prec=20): r""" The formal group inverse law i(t), which satisfies F(t, i(t)) = 0. INPUT: - ``prec`` - integer (default 20) OUTPUT: a power series with given precision DETAILS: Return the formal power series .. MATH:: i(t) = - t + a_1 t^2 + \cdots to precision `O(t^{prec})` of page 114 of [Silverman AEC1]. The result is cached, and a cached version is returned if possible. .. warning:: The resulting power series will have precision prec, but its parent PowerSeriesRing will have default precision 20 (or whatever the default default is). EXAMPLES:: sage: P.<a1, a2, a3, a4, a6> = ZZ[] sage: E = EllipticCurve(list(P.gens())) sage: i = E.formal_group().inverse(6); i -t - a1*t^2 - a1^2*t^3 + (-a1^3 - a3)*t^4 + (-a1^4 - 3*a1*a3)*t^5 + O(t^6) sage: F = E.formal_group().group_law(6) sage: F(i.parent().gen(), i) O(t^6) """ prec = max(prec,0) try: pr, inv = self.__inverse except AttributeError: pr = -1 if prec <= pr: t = inv.parent().gen() return inv + O(t**prec) x = self.x(prec) y = self.y(prec) a1, _, a3, _, _ = self.curve().ainvs() inv = x / ( y + a1*x + a3) # page 114 of Silverman, AEC I inv = inv.power_series().add_bigoh(prec) self.__inverse = (prec, inv) return inv
def newton_sqrt(self,f,x0, prec): r""" Takes the square root of the power series `f` by Newton's method NOTE: this function should eventually be moved to `p`-adic power series ring INPUT: - f power series wtih coefficients in `\QQ_p` or an extension - x0 seeds the Newton iteration - prec precision OUTPUT: the square root of `f` EXAMPLES:: sage: R.<x> = QQ['x'] sage: H = HyperellipticCurve(x^5-23*x^3+18*x^2+40*x) sage: Q = H(0,0) sage: u,v = H.local_coord(Q,prec=100) sage: K = Qp(11,5) sage: HK = H.change_ring(K) sage: L.<a> = K.extension(x^20-11) sage: HL = H.change_ring(L) sage: S = HL(u(a),v(a)) sage: f = H.hyperelliptic_polynomials()[0] sage: y = HK.newton_sqrt( f(u(a)^11), a^11,5) sage: y^2 - f(u(a)^11) O(a^122) AUTHOR: - Jennifer Balakrishnan """ z = x0 try: x = f.parent().variable_name() if x!='a' : #this is to distinguish between extensions of Qp that are finite vs. not S = f.base_ring()[[x]] x = S.gen() except ValueError: pass z = x0 loop_prec = (log(RR(prec))/log(RR(2))).ceil() for i in range(loop_prec): z = (z+f/z)/2 try: return z + O(x**prec) except (NameError,ArithmeticError,TypeError): return z
def inverse(self, prec=20): r""" The formal group inverse law i(t), which satisfies F(t, i(t)) = 0. INPUT: - ``prec`` - integer (default 20) OUTPUT: a power series with given precision DETAILS: Return the formal power series .. math:: i(t) = - t + a_1 t^2 + \cdots to precision `O(t^{prec})` of page 114 of [Silverman AEC1]. The result is cached, and a cached version is returned if possible. .. warning:: The resulting power series will have precision prec, but its parent PowerSeriesRing will have default precision 20 (or whatever the default default is). EXAMPLES:: sage: e = EllipticCurve([1, 2]) sage: F = e.formal_group().group_law(5) sage: i = e.formal_group().inverse(5) sage: Fx = F.base_extend(F.base_ring().base_extend(i.parent())) sage: Fx (i) (i.parent().gen()) O(t^5) """ prec = max(prec, 0) try: pr, inv = self.__inverse except AttributeError: pr = -1 if prec <= pr: t = inv.parent().gen() return inv + O(t**prec) x = self.x(prec) y = self.y(prec) a1, _, a3, _, _ = self.curve().ainvs() inv = x / (y + a1 * x + a3) # page 114 of Silverman, AEC I inv = inv.power_series().add_bigoh(prec) self.__inverse = (prec, inv) return inv
def x(self, prec=20): r""" Return the formal series `x(t) = t/w(t)` in terms of the local parameter `t = -x/y` at infinity. INPUT: - ``prec`` - integer (default 20) OUTPUT: a Laurent series with given precision DETAILS: Return the formal series .. MATH:: x(t) = t^{-2} - a_1 t^{-1} - a_2 - a_3 t - \cdots to precision `O(t^{prec})` of page 113 of [Silverman AEC1]. .. warning:: The resulting series will have precision prec, but its parent PowerSeriesRing will have default precision 20 (or whatever the default default is). EXAMPLES:: sage: EllipticCurve([0, 0, 1, -1, 0]).formal_group().x(10) t^-2 - t + t^2 - t^4 + 2*t^5 - t^6 - 2*t^7 + 6*t^8 - 6*t^9 + O(t^10) """ prec = max(prec, 0) y = self.y(prec) t = y.parent().gen() return -t*y + O(t**prec)
def mult_by_n(self, n, prec=10): """ The formal 'multiplication by n' endomorphism `[n]`. INPUT: - ``prec`` - integer (default 10) OUTPUT: a power series with given precision DETAILS: Return the formal power series .. math:: [n](t) = n t + \cdots to precision `O(t^{prec})` of Proposition 2.3 of [Silverman AEC1]. .. warning:: The resulting power series will have precision prec, but its parent PowerSeriesRing will have default precision 20 (or whatever the default default is). AUTHORS: - Nick Alexander: minor fixes, docstring - David Harvey (2007-03): faster algorithm for char 0 field case - Hamish Ivey-Law (2009-06): double-and-add algorithm for non char 0 field case. - Tom Boothby (2009-06): slight improvement to double-and-add EXAMPLES:: sage: e = EllipticCurve([1, 2, 3, 4, 6]) sage: e.formal_group().mult_by_n(0, 5) O(t^5) sage: e.formal_group().mult_by_n(1, 5) t + O(t^5) We verify an identity of low degree:: sage: none = e.formal_group().mult_by_n(-1, 5) sage: two = e.formal_group().mult_by_n(2, 5) sage: ntwo = e.formal_group().mult_by_n(-2, 5) sage: ntwo - none(two) O(t^5) sage: ntwo - two(none) O(t^5) It's quite fast:: sage: E = EllipticCurve("37a"); F = E.formal_group() sage: F.mult_by_n(100, 20) 100*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) TESTS:: sage: F = EllipticCurve(GF(17), [1, 1]).formal_group() sage: F.mult_by_n(10, 50) # long time (13s on sage.math, 2011) 10*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) sage: F = EllipticCurve(GF(101), [1, 1]).formal_group() sage: F.mult_by_n(100, 20) 100*t + O(t^20) sage: P.<a1, a2, a3, a4, a6> = PolynomialRing(ZZ, 5) sage: E = EllipticCurve(list(P.gens())) sage: E.formal().mult_by_n(2,prec=5) 2*t - a1*t^2 - 2*a2*t^3 + (a1*a2 - 7*a3)*t^4 + O(t^5) sage: E = EllipticCurve(QQ, [1,2,3,4,6]) sage: E.formal().mult_by_n(2,prec=5) 2*t - t^2 - 4*t^3 - 19*t^4 + O(t^5) """ if self.curve().base_ring().is_field() and self.curve().base_ring( ).characteristic() == 0 and n != 0: # The following algorithm only works over a field of # characteristic zero. I don't know whether something similar # can be done over a general ring. It would be nice if it did, # since it's much faster than using the formal group law. # -- dmharvey # Create a "formal point" on the original curve E. # Our answer only needs prec-1 coefficients (since lowest term # is t^1), and x(t) = t^(-2) + ... and y(t) = t^(-3) + ..., # so we only need x(t) mod t^(prec-3) and y(t) mod t^(prec-4) x = self.x(prec - 3) y = self.y(prec - 4) R = x.parent() # the Laurent series ring over the base ring X = self.curve().change_ring(R) P = X(x, y) # and multiply it by n, using the group law on E Q = n * P # express it in terms of the formal parameter return -Q[0] / Q[1] # Now the general case, not necessarily over a field. # For unknown reasons, this seems to lose one place of precision: # the coefficient of t**(prec-1) seems off. So we apply the easy fix. orig_prec = max(prec, 0) prec = orig_prec + 1 R = rings.PowerSeriesRing(self.curve().base_ring(), "t") t = R.gen() if n == 1: return t + O(t**orig_prec) if n == 0: return R(0) + O(t**orig_prec) if n == -1: return R(self.inverse(orig_prec)) if n < 0: F = self.inverse(prec)(self.mult_by_n(-n, orig_prec)) return R(F.add_bigoh(orig_prec)) F = self.group_law(prec) g = F.parent().base_ring().gen() # Double and add is faster than the naive method when n >= 4. if n < 4: for m in range(2, n + 1): g = F(g) return R(g.add_bigoh(orig_prec)) if n & 1: result = g else: result = F.parent().base_ring()(0) n = n >> 1 while n > 0: g = F(g, g) if n & 1: result = F(result, g) n = n >> 1 return R(result.add_bigoh(orig_prec))
def group_law(self, prec=10): r""" The formal group law. INPUT: - ``prec`` - integer (default 10) OUTPUT: a power series with given precision in ZZ[[ ZZ[['t1']],'t2']] DETAILS: Return the formal power series .. math:: F(t_1, t_2) = t_1 + t_2 - a_1 t_1 t_2 - \cdots to precision `O(t^{prec})` of page 115 of [Silverman AEC1]. The result is cached, and a cached version is returned if possible. .. warning:: The resulting power series will have precision prec, but its parent PowerSeriesRing will have default precision 20 (or whatever the default default is). AUTHORS: - Nick Alexander: minor fixes, docstring EXAMPLES:: sage: e = EllipticCurve([1, 2]) sage: e.formal_group().group_law(5) t1 + O(t1^5) + (1 - 2*t1^4 + O(t1^5))*t2 + (-4*t1^3 + O(t1^5))*t2^2 + (-4*t1^2 - 30*t1^4 + O(t1^5))*t2^3 + (-2*t1 - 30*t1^3 + O(t1^5))*t2^4 + O(t2^5) sage: e = EllipticCurve('14a1') sage: ehat = e.formal() sage: ehat.group_law(3) t1 + O(t1^3) + (1 - t1 - 4*t1^2 + O(t1^3))*t2 + (-4*t1 + O(t1^3))*t2^2 + O(t2^3) sage: ehat.group_law(5) t1 + O(t1^5) + (1 - t1 - 2*t1^3 - 32*t1^4 + O(t1^5))*t2 + (-3*t1^2 - 59*t1^3 - 120*t1^4 + O(t1^5))*t2^2 + (-2*t1 - 59*t1^2 - 141*t1^3 - 347*t1^4 + O(t1^5))*t2^3 + (-32*t1 - 120*t1^2 - 347*t1^3 - 951*t1^4 + O(t1^5))*t2^4 + O(t2^5) sage: e = EllipticCurve(GF(7),[3,4]) sage: ehat = e.formal() sage: ehat.group_law(3) t1 + O(t1^3) + (1 + O(t1^3))*t2 + O(t2^3) sage: F = ehat.group_law(7); F t1 + O(t1^7) + (1 + t1^4 + 2*t1^6 + O(t1^7))*t2 + (2*t1^3 + 6*t1^5 + O(t1^7))*t2^2 + (2*t1^2 + 3*t1^4 + 2*t1^6 + O(t1^7))*t2^3 + (t1 + 3*t1^3 + 4*t1^5 + O(t1^7))*t2^4 + (6*t1^2 + 4*t1^4 + O(t1^7))*t2^5 + (2*t1 + 2*t1^3 + O(t1^7))*t2^6 + O(t2^7) TESTS:: sage: i = e.formal_group().inverse(7) sage: Fx = F.base_extend(F.base_ring().base_extend(i.parent())) sage: Fx (i.parent().gen()) (i) O(t^7) Let's ensure caching with changed precision is working:: sage: e.formal_group().group_law(4) t1 + O(t1^4) + (1 + O(t1^4))*t2 + (2*t1^3 + O(t1^4))*t2^2 + (2*t1^2 + O(t1^4))*t2^3 + O(t2^4) Test for trac ticket 9646:: sage: P.<a1, a2, a3, a4, a6> = PolynomialRing(ZZ, 5) sage: E = EllipticCurve(list(P.gens())) sage: F = E.formal().group_law(prec = 4) sage: t2 = F.parent().gen() sage: t1 = F.parent().base_ring().gen() sage: F(t1, 0) t1 sage: F(0, t2) t2 sage: F[2][1] -a2 """ prec = max(prec, 0) if prec <= 0: raise ValueError, "The precision must be positive." R1 = rings.PowerSeriesRing(self.curve().base_ring(), "t1") R2 = rings.PowerSeriesRing(R1, "t2") t1 = R1.gen().add_bigoh(prec) t2 = R2.gen().add_bigoh(prec) if prec == 1: return R2(O(t2)) elif prec == 2: return R2(t1 + t2 - self.curve().a1() * t1 * t2) fix_prec = lambda F, final_prec: R2([c + O(t1**final_prec) for c in F]) + O(t2**final_prec) try: pr, F = self.__group_law except AttributeError: pr = -1 if prec <= pr: # we have to 'fix up' coefficient precisions return fix_prec(F, prec) w = self.w(prec + 1) tsum = lambda n: sum([t2**m * t1**(n - m - 1) for m in range(n)]) lam = sum([tsum(n) * w[n] for n in range(3, prec + 1)]) w1 = R1(w, prec + 1) nu = w1 - lam * t1 + O(t1**prec) a1, a2, a3, a4, a6 = self.curve().ainvs() lam2 = lam * lam lam3 = lam2 * lam # note that the following formula differs from the one in Silverman page 119. See trac ticket 9646 for the explanation and justification. t3 = -t1 - t2 - \ (a1*lam + a3*lam2 + a2*nu + 2*a4*lam*nu + 3*a6*lam2*nu)/ \ (1 + a2*lam + a4*lam2 + a6*lam3) inv = self.inverse(prec) # we have to 'fix up' precision issues F = fix_prec(inv(t3), prec) self.__group_law = (prec, F) return F