Path: blob/master/src/sage/schemes/elliptic_curves/monsky_washnitzer.py
8820 views
r"""1Computation of Frobenius matrix on Monsky-Washnitzer cohomology23The most interesting functions to be exported here are4:func:`matrix_of_frobenius` and :func:`adjusted_prec`.56Currently this code is limited to the case `p \geq 5` (no7`GF(p^n)` for `n > 1`), and only handles the8elliptic curve case (not more general hyperelliptic curves).910REFERENCES:1112.. [Ked2001] Kedlaya, K., "Counting points on hyperelliptic curves using13Monsky-Washnitzer cohomology", J. Ramanujan Math. Soc. 16 (2001) no144, 323-3381516.. [Edix] Edixhoven, B., "Point counting after Kedlaya", EIDMA-Stieltjes17graduate course, Lieden (lecture notes?).1819AUTHORS:2021- David Harvey and Robert Bradshaw: initial code developed at the 200622MSRI graduate workshop, working with Jennifer Balakrishnan and Liang23Xiao2425- David Harvey (2006-08): cleaned up, rewrote some chunks, lots more26documentation, added Newton iteration method, added more complete27'trace trick', integrated better into Sage.2829- David Harvey (2007-02): added algorithm with sqrt(p) complexity30(removed in May 2007 due to better C++ implementation)3132- Robert Bradshaw (2007-03): keep track of exact form in reduction33algorithms3435- Robert Bradshaw (2007-04): generalization to hyperelliptic curves36"""3738#*****************************************************************************39# Copyright (C) 2006 William Stein <[email protected]>40# 2006 Robert Bradshaw <[email protected]>41# 2006 David Harvey <[email protected]>42#43# Distributed under the terms of the GNU General Public License (GPL)44# http://www.gnu.org/licenses/45#*****************************************************************************464748from sage.rings.all import Integers, Integer, PolynomialRing, PowerSeriesRing, Rationals, Rational, LaurentSeriesRing4950from sage.rings.polynomial.polynomial_element import is_Polynomial5152from sage.modules.module import Module53from sage.structure.element import ModuleElement54from sage.matrix.all import matrix55from sage.modules.all import vector56from sage.rings.ring import CommutativeAlgebra57from sage.structure.element import CommutativeAlgebraElement58from sage.rings.infinity import Infinity5960from sage.rings.arith import binomial, integer_ceil as ceil61from sage.misc.functional import log62from sage.misc.misc import newton_method_sizes6364from ell_generic import is_EllipticCurve65from constructor import EllipticCurve666768class SpecialCubicQuotientRing(CommutativeAlgebra):69r"""70Specialised class for representing the quotient ring71`R[x,T]/(T - x^3 - ax - b)`, where `R` is an72arbitrary commutative base ring (in which 2 and 3 are invertible),73`a` and `b` are elements of that ring.7475Polynomials are represented internally in the form76`p_0 + p_1 x + p_2 x^2` where the `p_i` are77polynomials in `T`. Multiplication of polynomials always78reduces high powers of `x` (i.e. beyond `x^2`) to79powers of `T`.8081Hopefully this ring is faster than a general quotient ring because82it uses the special structure of this ring to speed multiplication83(which is the dominant operation in the frobenius matrix84calculation). I haven't actually tested this theory though...8586.. TODO::8788Eventually we will want to run this in characteristic 3, so we89need to: (a) Allow `Q(x)` to contain an `x^2` term, and (b) Remove90the requirement that 3 be invertible. Currently this is used in91the Toom-Cook algorithm to speed multiplication.9293EXAMPLES::9495sage: B.<t> = PolynomialRing(Integers(125))96sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))97sage: R98SpecialCubicQuotientRing over Ring of integers modulo 125 with polynomial T = x^3 + 124*x + 9499100Get generators::101102sage: x, T = R.gens()103sage: x104(0) + (1)*x + (0)*x^2105sage: T106(T) + (0)*x + (0)*x^2107108Coercions::109110sage: R(7)111(7) + (0)*x + (0)*x^2112113Create elements directly from polynomials::114115sage: A, z = R.poly_ring().objgen()116sage: A117Univariate Polynomial Ring in T over Ring of integers modulo 125118sage: R.create_element(z^2, z+1, 3)119(T^2) + (T + 1)*x + (3)*x^2120121Some arithmetic::122123sage: x^3124(T + 31) + (1)*x + (0)*x^2125sage: 3 * x**15 * T**2 + x - T126(3*T^7 + 90*T^6 + 110*T^5 + 20*T^4 + 58*T^3 + 26*T^2 + 124*T) + (15*T^6 + 110*T^5 + 35*T^4 + 63*T^2 + 1)*x + (30*T^5 + 40*T^4 + 8*T^3 + 38*T^2)*x^2127128Retrieve coefficients (output is zero-padded)::129130sage: x^10131(3*T^2 + 61*T + 8) + (T^3 + 93*T^2 + 12*T + 40)*x + (3*T^2 + 61*T + 9)*x^2132sage: (x^10).coeffs()133[[8, 61, 3, 0], [40, 12, 93, 1], [9, 61, 3, 0]]134135.. TODO::136137write an example checking multiplication of these polynomials138against Sage's ordinary quotient ring arithmetic. I can't seem139to get the quotient ring stuff happening right now...140"""141def __init__(self, Q, laurent_series=False):142"""143Constructor.144145INPUT:146147- ``Q`` -- a polynomial of the form148`Q(x) = x^3 + ax + b`, where `a`, `b` belong to a ring in which1492, 3 are invertible.150151- ``laurent_series`` -- whether or not to allow152negative powers of `T` (default=False)153154EXAMPLES::155156sage: B.<t> = PolynomialRing(Integers(125))157sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))158sage: R159SpecialCubicQuotientRing over Ring of integers modulo 125 with polynomial T = x^3 + 124*x + 94160161::162163sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 + 2*t^2 - t + B(1/4))164Traceback (most recent call last):165...166ValueError: Q (=t^3 + 2*t^2 + 124*t + 94) must be of the form x^3 + ax + b167168::169170sage: B.<t> = PolynomialRing(Integers(10))171sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + 1)172Traceback (most recent call last):173...174ArithmeticError: 2 and 3 must be invertible in the coefficient ring (=Ring of integers modulo 10) of Q175"""176if not is_Polynomial(Q):177raise TypeError("Q (=%s) must be a polynomial" % Q)178179if Q.degree() != 3 or not Q[2].is_zero():180raise ValueError("Q (=%s) must be of the form x^3 + ax + b" % Q)181182base_ring = Q.parent().base_ring()183184if not base_ring(6).is_unit():185raise ArithmeticError("2 and 3 must be invertible in the "186"coefficient ring (=%s) of Q" % base_ring)187188# CommutativeAlgebra.__init__ tries to establish a coercion189# from the base ring, by trac ticket #9138. The corresponding190# hom set is cached. In order to use self as cache key, its191# string representation is used. In otder to get the string192# representation, we need to know the attributes _a and193# _b. Hence, in #9138, we have to move CommutativeAlgebra.__init__194# further down:195self._a = Q[1]196self._b = Q[0]197if laurent_series:198self._poly_ring = LaurentSeriesRing(base_ring, 'T') # R[T]199else:200self._poly_ring = PolynomialRing(base_ring, 'T') # R[T]201self._poly_generator = self._poly_ring.gen(0) # the generator T202CommutativeAlgebra.__init__(self, base_ring)203204# Precompute a matrix that is used in the Toom-Cook multiplication.205# This is where we need 2 and 3 invertible.206207# (a good description of Toom-Cook is online at:208# http://www.gnu.org/software/gmp/manual/html_node/Toom-Cook-3-Way-Multiplication.html)209210self._speedup_matrix = (matrix(Integers(), 3, 3, [2, 4, 8,2111, 1, 1,2128, 4, 2])**(-1)).change_ring(base_ring).list()213214def __repr__(self):215"""216String representation217218EXAMPLES::219220sage: B.<t> = PolynomialRing(Integers(125))221sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))222sage: print R223SpecialCubicQuotientRing over Ring of integers modulo 125 with polynomial T = x^3 + 124*x + 94224"""225return "SpecialCubicQuotientRing over %s with polynomial T = %s" % \226(self.base_ring(), PolynomialRing(self.base_ring(), 'x')(227[self._b, self._a, 0, 1]))228229def poly_ring(self):230"""231Return the underlying polynomial ring in `T`.232233EXAMPLES::234235sage: B.<t> = PolynomialRing(Integers(125))236sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))237sage: R.poly_ring()238Univariate Polynomial Ring in T over Ring of integers modulo 125239"""240return self._poly_ring241242def gens(self):243"""244Return a list [x, T] where x and T are the generators of the ring245(as element *of this ring*).246247.. note::248249I have no idea if this is compatible with the usual Sage250'gens' interface.251252EXAMPLES::253254sage: B.<t> = PolynomialRing(Integers(125))255sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))256sage: x, T = R.gens()257sage: x258(0) + (1)*x + (0)*x^2259sage: T260(T) + (0)*x + (0)*x^2261"""262return [SpecialCubicQuotientRingElement(self, self._poly_ring(0),263self._poly_ring(1),264self._poly_ring(0),265check=False),266SpecialCubicQuotientRingElement(self, self._poly_generator,267self._poly_ring(0),268self._poly_ring(0),269check=False)]270271def create_element(self, p0, p1, p2, check=True):272"""273Creates the element `p_0 + p_1*x + p_2*x^2`, where the `p_i`274are polynomials in `T`.275276INPUT:277278- ``p0, p1, p2`` -- coefficients; must be coercible279into poly_ring()280281- ``check`` -- bool (default True): whether to carry282out coercion283284EXAMPLES::285286sage: B.<t> = PolynomialRing(Integers(125))287sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))288sage: A, z = R.poly_ring().objgen()289sage: R.create_element(z^2, z+1, 3)290(T^2) + (T + 1)*x + (3)*x^2291"""292return SpecialCubicQuotientRingElement(self, p0, p1, p2, check)293294def __call__(self, value):295"""296EXAMPLES::297298sage: B.<t> = PolynomialRing(Integers(125))299sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))300sage: R(3)301(3) + (0)*x + (0)*x^2302"""303return self._coerce_(value)304305def _coerce_impl(self, value):306"""307EXAMPLES::308309sage: B.<t> = PolynomialRing(Integers(125))310sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))311sage: R._coerce_impl(3)312(3) + (0)*x + (0)*x^2313"""314# coerce to underlying polynomial ring (possibly via base ring):315value = self._poly_ring._coerce_(value)316317return SpecialCubicQuotientRingElement(self, value, self._poly_ring(0),318self._poly_ring(0), check=False)319320321class SpecialCubicQuotientRingElement(CommutativeAlgebraElement):322"""323An element of a SpecialCubicQuotientRing.324"""325def __init__(self, parent, p0, p1, p2, check=True):326"""327Constructs the element `p_0 + p_1*x + p_2*x^2`, where328the `p_i` are polynomials in `T`.329330INPUT:331332- ``parent`` -- a SpecialCubicQuotientRing333334- ``p0, p1, p2`` -- coefficients; must be coercible335into parent.poly_ring()336337- ``check`` -- bool (default True): whether to carry338out coercion339340EXAMPLES::341342sage: B.<t> = PolynomialRing(Integers(125))343sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))344sage: from sage.schemes.elliptic_curves.monsky_washnitzer import SpecialCubicQuotientRingElement345sage: SpecialCubicQuotientRingElement(R, 2, 3, 4)346(2) + (3)*x + (4)*x^2347"""348if not isinstance(parent, SpecialCubicQuotientRing):349raise TypeError("parent (=%s) must be a SpecialCubicQuotientRing" % parent)350351CommutativeAlgebraElement.__init__(self, parent)352353if check:354poly_ring = parent.poly_ring()355p0 = poly_ring(p0)356p1 = poly_ring(p1)357p2 = poly_ring(p2)358359self._triple = (p0, p1, p2)360361def coeffs(self):362"""363Returns list of three lists of coefficients, corresponding to the364`x^0`, `x^1`, `x^2` coefficients. The lists365are zero padded to the same length. The list entries belong to the366base ring.367368EXAMPLES::369370sage: B.<t> = PolynomialRing(Integers(125))371sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))372sage: p = R.create_element(t, t^2 - 2, 3)373sage: p.coeffs()374[[0, 1, 0], [123, 0, 1], [3, 0, 0]]375"""376coeffs = [column.coeffs() for column in self._triple]377degree = max([len(x) for x in coeffs])378base_ring = self.parent().base_ring()379for column in coeffs:380column.extend([base_ring(0)] * (degree - len(column)))381return coeffs382383def __nonzero__(self):384"""385EXAMPLES::386387sage: B.<t> = PolynomialRing(Integers(125))388sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))389sage: x, T = R.gens()390sage: not x391False392sage: not T393False394sage: not R.create_element(0, 0, 0)395True396"""397return not not self._triple[0] or not not self._triple[1] or not not self._triple[2]398399def __cmp__(self, other):400"""401EXAMPLES::402403sage: B.<t> = PolynomialRing(Integers(125))404sage: x, t = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4)).gens()405sage: x == t406False407sage: x == x408True409sage: x == x + x - x410True411"""412return cmp(self._triple, other._triple)413414def _repr_(self):415"""416EXAMPLES::417418sage: B.<t> = PolynomialRing(Integers(125))419sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))420sage: x, T = R.gens()421sage: x + T*x - 2*T^2422(123*T^2) + (T + 1)*x + (0)*x^2423"""424return "(%s) + (%s)*x + (%s)*x^2" % self._triple425426def _latex_(self):427"""428EXAMPLES::429430sage: B.<t> = PolynomialRing(Integers(125))431sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))432sage: x, T = R.gens()433sage: f = x + T*x - 2*T^2434sage: latex(f)435(123 T^{2}) + (T + 1)x + (0)x^2436"""437return "(%s) + (%s)x + (%s)x^2" % \438tuple([column._latex_() for column in self._triple])439440def _add_(self, other):441"""442EXAMPLES::443444sage: B.<t> = PolynomialRing(Integers(125))445sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))446sage: f = R.create_element(2, t, t^2 - 3)447sage: g = R.create_element(3 + t, -t, t)448sage: f + g449(T + 5) + (0)*x + (T^2 + T + 122)*x^2450"""451return SpecialCubicQuotientRingElement(self.parent(),452self._triple[0] + other._triple[0],453self._triple[1] + other._triple[1],454self._triple[2] + other._triple[2],455check=False)456457def _sub_(self, other):458"""459EXAMPLES::460461sage: B.<t> = PolynomialRing(Integers(125))462sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))463sage: f = R.create_element(2, t, t^2 - 3)464sage: g = R.create_element(3 + t, -t, t)465sage: f - g466(124*T + 124) + (2*T)*x + (T^2 + 124*T + 122)*x^2467"""468return SpecialCubicQuotientRingElement(self.parent(),469self._triple[0] - other._triple[0],470self._triple[1] - other._triple[1],471self._triple[2] - other._triple[2],472check=False)473474def shift(self, n):475"""476Returns this element multiplied by `T^n`.477478EXAMPLES::479480sage: B.<t> = PolynomialRing(Integers(125))481sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))482sage: f = R.create_element(2, t, t^2 - 3)483sage: f484(2) + (T)*x + (T^2 + 122)*x^2485sage: f.shift(1)486(2*T) + (T^2)*x + (T^3 + 122*T)*x^2487sage: f.shift(2)488(2*T^2) + (T^3)*x + (T^4 + 122*T^2)*x^2489"""490return SpecialCubicQuotientRingElement(self.parent(),491self._triple[0].shift(n),492self._triple[1].shift(n),493self._triple[2].shift(n),494check=False)495496def scalar_multiply(self, scalar):497"""498Multiplies this element by a scalar, i.e. just multiply each499coefficient of `x^j` by the scalar.500501INPUT:502503- ``scalar`` -- either an element of base_ring, or an504element of poly_ring.505506EXAMPLES::507508sage: B.<t> = PolynomialRing(Integers(125))509sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))510sage: x, T = R.gens()511sage: f = R.create_element(2, t, t^2 - 3)512sage: f513(2) + (T)*x + (T^2 + 122)*x^2514sage: f.scalar_multiply(2)515(4) + (2*T)*x + (2*T^2 + 119)*x^2516sage: f.scalar_multiply(t)517(2*T) + (T^2)*x + (T^3 + 122*T)*x^2518"""519scalar = self.parent()._poly_ring(scalar)520return SpecialCubicQuotientRingElement(self.parent(),521scalar * self._triple[0],522scalar * self._triple[1],523scalar * self._triple[2],524check=False)525526def square(self):527"""528EXAMPLES::529530sage: B.<t> = PolynomialRing(Integers(125))531sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))532sage: x, T = R.gens()533534::535536sage: f = R.create_element(1 + 2*t + 3*t^2, 4 + 7*t + 9*t^2, 3 + 5*t + 11*t^2)537sage: f.square()538(73*T^5 + 16*T^4 + 38*T^3 + 39*T^2 + 70*T + 120) + (121*T^5 + 113*T^4 + 73*T^3 + 8*T^2 + 51*T + 61)*x + (18*T^4 + 60*T^3 + 22*T^2 + 108*T + 31)*x^2539"""540return self * self541542def _mul_(self, other):543"""544EXAMPLES::545546sage: B.<t> = PolynomialRing(Integers(125))547sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))548sage: x, T = R.gens()549550::551552sage: f = R.create_element(1 + 2*t + 3*t^2, 4 + 7*t + 9*t^2, 3 + 5*t + 11*t^2)553sage: g = R.create_element(4 + 3*t + 7*t^2, 2 + 3*t + t^2, 8 + 4*t + 6*t^2)554sage: f * g555(65*T^5 + 27*T^4 + 33*T^3 + 75*T^2 + 120*T + 57) + (66*T^5 + T^4 + 123*T^3 + 95*T^2 + 24*T + 50)*x + (45*T^4 + 75*T^3 + 37*T^2 + 2*T + 52)*x^2556"""557if not isinstance(other, SpecialCubicQuotientRingElement):558return self.scalar_multiply(other)559560# Here we do Toom-Cook three-way multiplication, which reduces the561# naive 9 polynomial multiplications to only 5 polynomial multiplications.562563a0, a1, a2 = self._triple564b0, b1, b2 = other._triple565M = self.parent()._speedup_matrix566567if self is other:568# faster method if we're squaring569p0 = a0 * a0570temp = a0 + 2*a1 + 4*a2571p1 = temp * temp572temp = a0 + a1 + a2573p2 = temp * temp574temp = 4*a0 + 2*a1 + a2575p3 = temp * temp576p4 = a2 * a2577578else:579p0 = a0 * b0580p1 = (a0 + 2*a1 + 4*a2) * (b0 + 2*b1 + 4*b2)581p2 = (a0 + a1 + a2) * (b0 + b1 + b2)582p3 = (4*a0 + 2*a1 + a2) * (4*b0 + 2*b1 + b2)583p4 = a2 * b2584585q1 = p1 - p0 - 16*p4586q2 = p2 - p0 - p4587q3 = p3 - 16*p0 - p4588589c0 = p0590c1 = M[0]*q1 + M[1]*q2 + M[2]*q3591c2 = M[3]*q1 + M[4]*q2 + M[5]*q3592c3 = M[6]*q1 + M[7]*q2 + M[8]*q3593c4 = p4594595# Now the product is c0 + c1 x + c2 x^2 + c3 x^3 + c4 x^4.596# We need to reduce mod y = x^3 + ax + b and return result.597598parent = self.parent()599T = parent._poly_generator600b = parent._b601a = parent._a602603# todo: These lines are necessary to get binop stuff working604# for certain base rings, e.g. when we compute b*c3 in the605# final line. They shouldn't be necessary. Need to fix this606# somewhere else in Sage.607a = parent._poly_ring(a)608b = parent._poly_ring(b)609610return SpecialCubicQuotientRingElement(parent,611-b*c3 + c0 + c3*T,612-b*c4 - a*c3 + c1 + c4*T,613-a*c4 + c2,614check=False)615616617def transpose_list(input):618"""619INPUT:620621- ``input`` -- a list of lists, each list of the same622length623624OUTPUT:625626- ``output`` -- a list of lists such that output[i][j]627= input[j][i]628629EXAMPLES::630631sage: from sage.schemes.elliptic_curves.monsky_washnitzer import transpose_list632sage: L = [[1, 2], [3, 4], [5, 6]]633sage: transpose_list(L)634[[1, 3, 5], [2, 4, 6]]635"""636h = len(input)637w = len(input[0])638639output = []640for i in range(w):641row = []642for j in range(h):643row.append(input[j][i])644output.append(row)645return output646647648def helper_matrix(Q):649"""650Computes the (constant) matrix used to calculate the linear651combinations of the `d(x^i y^j)` needed to eliminate the652negative powers of `y` in the cohomology (i.e. in653reduce_negative()).654655INPUT:656657- ``Q`` -- cubic polynomial658659EXAMPLES::660661sage: t = polygen(QQ,'t')662sage: from sage.schemes.elliptic_curves.monsky_washnitzer import helper_matrix663sage: helper_matrix(t**3-4*t-691)664[ 64/12891731 -16584/12891731 4297329/12891731]665[ 6219/12891731 -32/12891731 8292/12891731]666[ -24/12891731 6219/12891731 -32/12891731]667"""668a = Q[1]669b = Q[0]670671# Discriminant (should be invertible for a curve of good reduction)672D = 4*a**3 + 27*b**2673Dinv = D**(-1) # NB do not use 1/D674675# This is the inverse of the matrix676# [ a, -3b, 0 ]677# [ 0, -2a, -3b ]678# [ 3, 0, -2a ]679680return Dinv * matrix([[4*a**2, -6*b*a, 9*b**2],681[-9*b, -2*a**2, 3*b*a],682[6*a, -9*b, -2*a**2]])683684685def lift(x):686r"""687Tries to call x.lift(), presumably from the `p`-adics to ZZ.688689If this fails, it assumes the input is a power series, and tries to690lift it to a power series over QQ.691692This function is just a very kludgy solution to the problem of693trying to make the reduction code (below) work over both Zp and694Zp[[t]].695696EXAMPLES::697698sage: from sage.schemes.elliptic_curves.monsky_washnitzer import lift699sage: l = lift(Qp(13)(131)); l700131701sage: l.parent()702Integer Ring703704sage: x=PowerSeriesRing(Qp(17),'x').gen()705sage: l = lift(4+5*x+17*x**6); l7064 + 5*t + 17*t^6707sage: l.parent()708Power Series Ring in t over Rational Field709"""710try:711return x.lift()712except AttributeError:713return PowerSeriesRing(Rationals(), "t")(x.list(), x.prec())714715716def reduce_negative(Q, p, coeffs, offset, exact_form=None):717"""718Applies cohomology relations to incorporate negative powers of719`y` into the `y^0` term.720721INPUT:722723- ``p`` -- prime724725- ``Q`` -- cubic polynomial726727- ``coeffs`` -- list of length 3 lists. The728`i^{th}` list [a, b, c] represents729`y^{2(i - offset)} (a + bx + cx^2) dx/y`.730731- ``offset`` -- nonnegative integer732733OUTPUT: The reduction is performed in-place. The output is placed734in coeffs[offset]. Note that coeffs[i] will be meaningless for i735offset after this function is finished.736737EXAMPLE::738739sage: R.<x> = Integers(5^3)['x']740sage: Q = x^3 - x + R(1/4)741sage: coeffs = [[10, 15, 20], [1, 2, 3], [4, 5, 6], [7, 8, 9]]742sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]743sage: monsky_washnitzer.reduce_negative(Q, 5, coeffs, 3)744sage: coeffs[3]745[28, 52, 9]746747::748749sage: R.<x> = Integers(7^3)['x']750sage: Q = x^3 - x + R(1/4)751sage: coeffs = [[7, 14, 21], [1, 2, 3], [4, 5, 6], [7, 8, 9]]752sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]753sage: monsky_washnitzer.reduce_negative(Q, 7, coeffs, 3)754sage: coeffs[3]755[245, 332, 9]756"""757758m = helper_matrix(Q).list()759base_ring = Q.base_ring()760next_a = coeffs[0]761762if exact_form is not None:763x = exact_form.parent().gen(0)764y = exact_form.parent()(exact_form.parent().base_ring().gen(0))765766try:767three_j_plus_5 = 5 - base_ring(6*offset)768three_j_plus_7 = 7 - base_ring(6*offset)769six = base_ring(6)770771for i in range(0, offset):772773j = 2*(i-offset)774a = next_a775next_a = coeffs[i+1]776777# todo: the following divisions will sometimes involve778# a division by (a power of) p. In all cases, we know (from779# Kedlaya's estimates) that the answer should be p-integral.780# However, since we're working over $Z/p^k Z$, we're not allowed781# to "divide by p". So currently we lift to Q, divide, and coerce782# back. Eventually, when pAdicInteger is implemented, and plays783# nicely with pAdicField, we should reimplement this stuff784# using pAdicInteger.785786if (p.divides(j+1)):787# need to lift here to perform the division788a[0] = base_ring(lift(a[0]) / (j+1))789a[1] = base_ring(lift(a[1]) / (j+1))790a[2] = base_ring(lift(a[2]) / (j+1))791else:792j_plus_1_inv = ~base_ring(j+1)793a[0] = a[0] * j_plus_1_inv794a[1] = a[1] * j_plus_1_inv795a[2] = a[2] * j_plus_1_inv796797c1 = m[3]*a[0] + m[4]*a[1] + m[5]*a[2]798c2 = m[6]*a[0] + m[7]*a[1] + m[8]*a[2]799next_a[0] = next_a[0] - three_j_plus_5 * c1800next_a[1] = next_a[1] - three_j_plus_7 * c2801802three_j_plus_7 = three_j_plus_7 + six803three_j_plus_5 = three_j_plus_5 + six804805if exact_form is not None:806c0 = m[0]*a[0] + m[1]*a[1] + m[2]*a[2]807exact_form += (c0 + c1*x + c2 * x**2) * y**(j+1)808809except NotImplementedError:810raise NotImplementedError("It looks like you've found a "811"non-integral matrix of Frobenius! "812"(Q=%s, p=%s)\nTime to write a paper." % (Q, p))813814coeffs[int(offset)] = next_a815816return exact_form817818819def reduce_positive(Q, p, coeffs, offset, exact_form=None):820"""821Applies cohomology relations to incorporate positive powers of822`y` into the `y^0` term.823824INPUT:825826- ``Q`` -- cubic polynomial827828- ``coeffs`` -- list of length 3 lists. The829`i^{th}` list [a, b, c] represents830`y^{2(i - offset)} (a + bx + cx^2) dx/y`.831832- ``offset`` -- nonnegative integer833834OUTPUT: The reduction is performed in-place. The output is placed835in coeffs[offset]. Note that coeffs[i] will be meaningless for i836offset after this function is finished.837838EXAMPLE::839840sage: R.<x> = Integers(5^3)['x']841sage: Q = x^3 - x + R(1/4)842843::844845sage: coeffs = [[1, 2, 3], [10, 15, 20]]846sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]847sage: monsky_washnitzer.reduce_positive(Q, 5, coeffs, 0)848sage: coeffs[0]849[16, 102, 88]850851::852853sage: coeffs = [[9, 8, 7], [10, 15, 20]]854sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]855sage: monsky_washnitzer.reduce_positive(Q, 5, coeffs, 0)856sage: coeffs[0]857[24, 108, 92]858"""859860base_ring = Q.base_ring()861next_a = coeffs[len(coeffs) - 1]862863Qa = Q[1]864Qb = Q[0]865866A = 2*Qa867B = 3*Qb868869offset = Integer(offset)870871if exact_form is not None:872x = exact_form.parent().gen(0)873y = exact_form.parent().base_ring().gen(0)874# y = exact_form.parent()(exact_form.parent().base_ring().gen(0))875876for i in range(len(coeffs)-1, offset, -1):877j = 2*(i-offset) - 2878a = next_a879next_a = coeffs[i-1]880881a[0] = a[0] - Qa*a[2]/3 # subtract d(y^j + 3)882if exact_form is not None:883exact_form += Q.base_ring()(a[2].lift() / (3*j+9)) * y**(j+3)884885# todo: see comments about pAdicInteger in reduceNegative()886887# subtract off c1 of d(x y^j + 1), and888if p.divides(3*j + 5):889c1 = base_ring(lift(a[0]) / (3*j + 5))890else:891c1 = a[0] / (3*j + 5)892893# subtract off c2 of d(x^2 y^j + 1)894if p.divides(3*j + 7):895c2 = base_ring(lift(a[1]) / (3*j + 7))896else:897c2 = a[1] / (3*j + 7)898899next_a[0] = next_a[0] + B*c1*(j+1)900next_a[1] = next_a[1] + A*c1*(j+1) + B*c2*(j+1)901next_a[2] = next_a[2] + A*c2*(j+1)902903if exact_form is not None:904exact_form += (c1*x + c2 * x**2) * y**(j+1)905906coeffs[int(offset)] = next_a907908return exact_form909910911def reduce_zero(Q, coeffs, offset, exact_form=None):912"""913Applies cohomology relation to incorporate `x^2 y^0` term914into `x^0 y^0` and `x^1 y^0` terms.915916INPUT:917918- ``Q`` -- cubic polynomial919920- ``coeffs`` -- list of length 3 lists. The921`i^{th}` list [a, b, c] represents922`y^{2(i - offset)} (a + bx + cx^2) dx/y`.923924- ``offset`` -- nonnegative integer925926OUTPUT: The reduction is performed in-place. The output is placed927in coeffs[offset]. This method completely ignores coeffs[i] for i928!= offset.929930EXAMPLE::931932sage: R.<x> = Integers(5^3)['x']933sage: Q = x^3 - x + R(1/4)934sage: coeffs = [[1, 2, 3], [4, 5, 6], [7, 8, 9]]935sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]936sage: monsky_washnitzer.reduce_zero(Q, coeffs, 1)937sage: coeffs[1]938[6, 5, 0]939"""940941a = coeffs[int(offset)]942if a[2] == 0:943return exact_form944945Qa = Q[1]946947a[0] = a[0] - a[2]*Qa/3 # $3x^2 dx/y = -a dx/y$948949coeffs[int(offset)] = a950951if exact_form is not None:952y = exact_form.parent()(exact_form.parent().base_ring().gen(0))953exact_form += Q.base_ring()(a[2] / 3) * y954955a[2] = 0956957coeffs[int(offset)] = a958return exact_form959960961def reduce_all(Q, p, coeffs, offset, compute_exact_form=False):962"""963Applies cohomology relations to reduce all terms to a linear964combination of `dx/y` and `x dx/y`.965966INPUT:967968- ``Q`` -- cubic polynomial969970- ``coeffs`` -- list of length 3 lists. The971`i^{th}` list [a, b, c] represents972`y^{2(i - offset)} (a + bx + cx^2) dx/y`.973974- ``offset`` -- nonnegative integer975976977OUTPUT:978979980- ``A, B`` - pair such that the input differential is981cohomologous to (A + Bx) dx/y.982983984.. note::985986The algorithm operates in-place, so the data in coeffs is987destroyed.988989EXAMPLE::990991sage: R.<x> = Integers(5^3)['x']992sage: Q = x^3 - x + R(1/4)993sage: coeffs = [[1, 2, 3], [4, 5, 6], [7, 8, 9]]994sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]995sage: monsky_washnitzer.reduce_all(Q, 5, coeffs, 1)996(21, 106)997"""998999R = Q.base_ring()10001001if compute_exact_form:1002# exact_form = SpecialCubicQuotientRing(Q, laurent_series=True)(0)1003exact_form = PolynomialRing(LaurentSeriesRing(Q.base_ring(), 'y'), 'x')(0)1004# t = (Q.base_ring().order().factor())[0]1005# from sage.rings.padics.qp import pAdicField1006# exact_form = PolynomialRing(LaurentSeriesRing(pAdicField(p, t[1]), 'y'), 'x')(0)1007else:1008exact_form = None10091010while len(coeffs) <= offset:1011coeffs.append([R(0), R(0), R(0)])10121013exact_form = reduce_negative(Q, p, coeffs, offset, exact_form)1014exact_form = reduce_positive(Q, p, coeffs, offset, exact_form)1015exact_form = reduce_zero(Q, coeffs, offset, exact_form)10161017if exact_form is None:1018return coeffs[int(offset)][0], coeffs[int(offset)][1]1019else:1020return (coeffs[int(offset)][0], coeffs[int(offset)][1]), exact_form102110221023def frobenius_expansion_by_newton(Q, p, M):1024r"""1025Computes the action of Frobenius on `dx/y` and on1026`x dx/y`, using Newton's method (as suggested in Kedlaya's1027paper [Ked2001]_).10281029(This function does *not* yet use the cohomology relations - that1030happens afterwards in the "reduction" step.)10311032More specifically, it finds `F_0` and `F_1` in1033the quotient ring `R[x, T]/(T - Q(x))`, such that10341035.. math::10361037F( dx/y) = T^{-r} F0 dx/y, \text{\ and\ } F(x dx/y) = T^{-r} F1 dx/y10381039where10401041.. math::10421043r = ( (2M-3)p - 1 )/2.104410451046(Here `T` is `y^2 = z^{-2}`, and `R` is the1047coefficient ring of `Q`.)10481049`F_0` and `F_1` are computed in the1050SpecialCubicQuotientRing associated to `Q`, so all powers1051of `x^j` for `j \geq 3` are reduced to powers of1052`T`.10531054INPUT:10551056- ``Q`` -- cubic polynomial of the form1057`Q(x) = x^3 + ax + b`, whose coefficient ring is a1058`Z/(p^M)Z`-algebra10591060- ``p`` -- residue characteristic of the p-adic field10611062- ``M`` -- p-adic precision of the coefficient ring1063(this will be used to determine the number of Newton iterations)10641065OUTPUT:10661067- ``F0, F1`` - elements of1068SpecialCubicQuotientRing(Q), as described above10691070- ``r`` - non-negative integer, as described above10711072EXAMPLES::10731074sage: from sage.schemes.elliptic_curves.monsky_washnitzer import frobenius_expansion_by_newton1075sage: R.<x> = Integers(5^3)['x']1076sage: Q = x^3 - x + R(1/4)1077sage: frobenius_expansion_by_newton(Q,5,3)1078((25*T^5 + 75*T^3 + 100*T^2 + 100*T + 100) + (5*T^6 + 80*T^5 + 100*T^31079+ 25*T + 50)*x + (55*T^5 + 50*T^4 + 75*T^3 + 25*T^2 + 25*T + 25)*x^2,1080(5*T^8 + 15*T^7 + 95*T^6 + 10*T^5 + 25*T^4 + 25*T^3 + 100*T^2 + 50)1081+ (65*T^7 + 55*T^6 + 70*T^5 + 100*T^4 + 25*T^2 + 100*T)*x1082+ (15*T^6 + 115*T^5 + 75*T^4 + 100*T^3 + 50*T^2 + 75*T + 75)*x^2, 7)1083"""10841085S = SpecialCubicQuotientRing(Q)1086x, _ = S.gens() # T = y^21087base_ring = S.base_ring()10881089# When we compute Frob(1/y) we actually only need precision M-1, since1090# we're going to multiply by p at the end anyway.1091M = float(M - 1)10921093# Kedlaya sets s = Q(x^p)/T^p = 1 + p T^{-p} E, where1094# E = (Q(x^p) - Q(x)^p) / p (has integral coefficients).1095# Then he computes s^{-1/2} in S, using Newton's method to find1096# successive approximations. We follow this plan, but we normalise our1097# approximations so that we only ever need positive powers of T.10981099# Start by setting r = Q(x^p)/2 = 1/2 T^p s.1100# (The 1/2 is for convenience later on.)1101x_to_p_less_one = x**(p-1)1102x_to_p = x_to_p_less_one * x1103x_to_p_cubed = x_to_p.square() * x_to_p1104r = (base_ring(1) / base_ring(2)) * (x_to_p_cubed + Q[1]*x_to_p + S(Q[0]))11051106# todo: this next loop would be clearer if it used the newton_method_sizes()1107# function11081109# We will start with a hard-coded initial approximation, which we provide1110# up to precision 3. First work out what precision is best to start with.1111if M <= 3:1112initial_precision = M1113elif ceil(log(M/2, 2)) == ceil(log(M/3, 2)):1114# In this case there is no advantage to starting with precision three,1115# because we'll overshoot at the end. E.g. suppose the final precision1116# is 8. If we start with precision 2, we need two iterations to get us1117# to 8. If we start at precision 3, we will still need two iterations,1118# but we do more work along the way. So may as well start with only 2.1119initial_precision = 21120else:1121initial_precision = 311221123# Now compute the first approximation. In the main loop below, X is the1124# normalised approximation, and k is the precision. More specifically,1125# X = T^{p(k-1)} x_i, where x_i is an approximation to s^{-1/2}, and the1126# approximation is correct mod p^k.1127if initial_precision == 1:1128k = 11129X = S(1)1130elif initial_precision == 2:1131# approximation is 3/2 - 1/2 s1132k = 21133X = S(base_ring(3) / base_ring(2)).shift(p) - r1134elif initial_precision == 3:1135# approximation is (15 - 10 s + 3 s^2) / 81136k = 31137X = (base_ring(1) / base_ring(8)) * (S(15).shift(2*p)1138- (base_ring(20) * r).shift(p) +1139(base_ring(12) * r.square()))1140# The key to the following calculation is that the T^{-m} coefficient1141# of every x_i is divisible by p^(ceil(m/p)) (for m >= 0). Therefore if1142# we are only expecting an answer correct mod p^k, we can truncate1143# beyond the T^{-(k-1)p} term without any problems.11441145# todo: what would be really nice is to be able to work in a lower1146# precision *coefficient ring* when we start the iteration, and move up to1147# higher precision rings as the iteration proceeds. This would be feasible1148# over Integers(p**n), but quite complicated (maybe impossible) over a more1149# general base ring. This might give a decent constant factor speedup;1150# or it might not, depending on how much the last iteration dominates the1151# whole runtime. My guess is that it isn't worth the effort.11521153three_halves = base_ring(3) / base_ring(2)11541155# Newton iteration loop1156while k < M:1157# target_k = k' = precision we want our answer to be after this iteration1158target_k = 2*k11591160# This prevents us overshooting. For example if the current precision1161# is 3 and we want to get to 10, we're better off going up to 51162# instead of 6, because it is less work to get from 5 to 10 than it1163# is to get from 6 to 10.1164if ceil(log(M/target_k, 2)) == ceil(log(M/(target_k-1), 2)):1165target_k -= 111661167# temp = T^{p(3k-2)} 1/2 s x_i^31168temp = X.square() * (X * r)11691170# We know that the final result is only going to be correct mod1171# p^(target_k), so we might as well truncate the extraneous terms now.1172# temp = T^{p(k'-1)} 1/2 s x_i^31173temp = temp.shift(-p*(3*k - target_k - 1))11741175# X = T^{p(k'-1)} (3/2 x_i - 1/2 s x_i^3)1176# = T^{p(k'-1)} x_{i+1}1177X = (three_halves * X).shift(p*(target_k - k)) - temp11781179k = target_k11801181# Now k should equal M, since we're up to the correct precision1182assert k == M, "Oops, something went wrong in the iteration"11831184# We should have s^{-1/2} correct to precision M.1185# The following line can be uncommented to verify this.1186# (It is a slow verification though, can double the whole computation time.)11871188#assert (p * X.square() * r * base_ring(2)).coeffs() == \1189# R(p).shift(p*(2*M - 1)).coeffs()11901191# Finally incorporate frobenius of dx and x dx, and choose offset that1192# compensates for our normalisations by powers of T.1193F0 = base_ring(p) * x_to_p_less_one * X1194F1 = F0 * x_to_p1195offset = ((2*k-1)*p - 1)/211961197return F0, F1, offset119811991200def frobenius_expansion_by_series(Q, p, M):1201r"""1202Computes the action of Frobenius on `dx/y` and on `x dx/y`, using a1203series expansion.12041205(This function computes the same thing as1206frobenius_expansion_by_newton(), using a different method.1207Theoretically the Newton method should be asymptotically faster,1208when the precision gets large. However, in practice, this functions1209seems to be marginally faster for moderate precision, so I'm1210keeping it here until I figure out exactly why it is faster.)12111212(This function does *not* yet use the cohomology relations - that1213happens afterwards in the "reduction" step.)12141215More specifically, it finds F0 and F1 in the quotient ring1216`R[x, T]/(T - Q(x))`, such that1217`F( dx/y) = T^{-r} F0 dx/y`, and1218`F(x dx/y) = T^{-r} F1 dx/y` where1219`r = ( (2M-3)p - 1 )/2`. (Here `T` is `y^2 = z^{-2}`,1220and `R` is the coefficient ring of `Q`.)12211222`F_0` and `F_1` are computed in the1223SpecialCubicQuotientRing associated to `Q`, so all powers1224of `x^j` for `j \geq 3` are reduced to powers of1225`T`.12261227It uses the sum12281229.. math::12301231F0 = \sum_{k=0}^{M-2} {-1/2 \choose k} p x^{p-1} E^k T^{(M-2-k)p}12321233and12341235.. math::12361237F1 = x^p F0,12381239where `E = Q(x^p) - Q(x)^p`.12401241INPUT:12421243- ``Q`` -- cubic polynomial of the form1244`Q(x) = x^3 + ax + b`, whose coefficient ring is a1245`\ZZ/(p^M)\ZZ` -algebra12461247- ``p`` -- residue characteristic of the `p`-adic field12481249- ``M`` -- `p`-adic precision of the coefficient ring1250(this will be used to determine the number of terms in the1251series)12521253OUTPUT:12541255- ``F0, F1`` - elements of1256SpecialCubicQuotientRing(Q), as described above12571258- ``r`` - non-negative integer, as described above12591260EXAMPLES::12611262sage: from sage.schemes.elliptic_curves.monsky_washnitzer import frobenius_expansion_by_series1263sage: R.<x> = Integers(5^3)['x']1264sage: Q = x^3 - x + R(1/4)1265sage: frobenius_expansion_by_series(Q,5,3)1266((25*T^5 + 75*T^3 + 100*T^2 + 100*T + 100) + (5*T^6 + 80*T^5 + 100*T^31267+ 25*T + 50)*x + (55*T^5 + 50*T^4 + 75*T^3 + 25*T^2 + 25*T + 25)*x^2,1268(5*T^8 + 15*T^7 + 95*T^6 + 10*T^5 + 25*T^4 + 25*T^3 + 100*T^2 + 50)1269+ (65*T^7 + 55*T^6 + 70*T^5 + 100*T^4 + 25*T^2 + 100*T)*x1270+ (15*T^6 + 115*T^5 + 75*T^4 + 100*T^3 + 50*T^2 + 75*T + 75)*x^2, 7)1271"""12721273S = SpecialCubicQuotientRing(Q)1274x, _ = S.gens()1275base_ring = S.base_ring()12761277x_to_p_less_1 = x**(p-1)1278x_to_p = x_to_p_less_1 * x12791280# compute frobQ = Q(x^p)1281x_to_p_squared = x_to_p * x_to_p1282x_to_p_cubed = x_to_p_squared * x_to_p1283frobQ = x_to_p_cubed + Q[1]*x_to_p + Q[0]*S(1)1284# anticipating the day when p = 3 is supported:1285# frobQ = x_to_p_cubed + Q[2]*x_to_p_squared + Q[1]*x_to_p + Q[0]*S(1)12861287E = frobQ - S(1).shift(p) # E = Q(x^p) - Q(x)^p12881289offset = int(((2*M-3)*p-1)/2)1290term = p * x_to_p_less_11291F0 = term.shift((M-2)*p)12921293# todo: Possible speedup idea, perhaps by a factor of 2, but1294# it requires a lot of work:1295# Note that p divides E, so p^k divides E^k. So when we are1296# working with high powers of E, we're doing a lot more work1297# in the multiplications than we need to. To take advantage of1298# this we would need some protocol for "lowering the precision"1299# of a SpecialCubicQuotientRing. This would be quite messy to1300# do properly over an arbitrary base ring. Perhaps it is1301# feasible to do for the most common case (i.e. Z/p^nZ).1302# (but it probably won't save much time unless p^n is very1303# large, because the machine word size is probably pretty1304# big anyway.)13051306for k in range(int(1), int(M-1)):1307term = term * E1308c = base_ring(binomial(-Integer(1)/2, k))1309F0 += (term * c).shift((M-k-2)*p)13101311return F0, F0 * x_to_p, offset131213131314def adjusted_prec(p, prec):1315r"""1316Computes how much precision is required in matrix_of_frobenius to1317get an answer correct to prec `p`-adic digits.13181319The issue is that the algorithm used in1320:func:`matrix_of_frobenius` sometimes performs divisions by `p`,1321so precision is lost during the algorithm.13221323The estimate returned by this function is based on Kedlaya's result1324(Lemmas 2 and 3 of [Ked2001]_),1325which implies that if we start with `M` `p`-adic1326digits, the total precision loss is at most1327`1 + \lfloor \log_p(2M-3) \rfloor` `p`-adic1328digits. (This estimate is somewhat less than the amount you would1329expect by naively counting the number of divisions by1330`p`.)13311332INPUT:13331334- ``p`` -- a prime = 513351336- ``prec`` -- integer, desired output precision, = 113371338OUTPUT: adjusted precision (usually slightly more than prec)1339"""13401341# initial estimate:1342if prec <= 2:1343adjusted = 21344else:1345adjusted = prec + int(log(2*prec - 3, p)) - 113461347# increase it until we have enough1348while adjusted - int(log(2*adjusted - 3, p)) - 1 < prec:1349adjusted += 113501351return adjusted135213531354def matrix_of_frobenius(Q, p, M, trace=None, compute_exact_forms=False):1355"""1356Computes the matrix of Frobenius on Monsky-Washnitzer cohomology,1357with respect to the basis `(dx/y, x dx/y)`.13581359INPUT:13601361- ``Q`` -- cubic polynomial `Q(x) = x^3 + ax + b`1362defining an elliptic curve `E` by1363`y^2 = Q(x)`. The coefficient ring of `Q` should be a1364`\ZZ/(p^M)\ZZ`-algebra in which the matrix of1365frobenius will be constructed.13661367- ``p`` -- prime = 5 for which E has good reduction13681369- ``M`` -- integer = 2; `p` -adic precision of1370the coefficient ring13711372- ``trace`` -- (optional) the trace of the matrix, if1373known in advance. This is easy to compute because it is just the1374`a_p` of the curve. If the trace is supplied,1375matrix_of_frobenius will use it to speed the computation (i.e. we1376know the determinant is `p`, so we have two conditions, so1377really only column of the matrix needs to be computed. it is1378actually a little more complicated than that, but that's the basic1379idea.) If trace=None, then both columns will be computed1380independently, and you can get a strong indication of correctness1381by verifying the trace afterwards.13821383.. warning::13841385THE RESULT WILL NOT NECESSARILY BE CORRECT TO M p-ADIC1386DIGITS. If you want prec digits of precision, you need to use1387the function adjusted_prec(), and then you need to reduce the1388answer mod `p^{\mathrm{prec}}` at the end.13891390OUTPUT:139113922x2 matrix of frobenius on Monsky-Washnitzer cohomology,1393with entries in the coefficient ring of Q.13941395EXAMPLES:13961397A simple example::13981399sage: p = 51400sage: prec = 31401sage: M = monsky_washnitzer.adjusted_prec(p, prec)1402sage: M140351404sage: R.<x> = PolynomialRing(Integers(p**M))1405sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 - x + R(1/4), p, M)1406sage: A1407[3090 187]1408[2945 408]14091410But the result is only accurate to prec digits::14111412sage: B = A.change_ring(Integers(p**prec))1413sage: B1414[90 62]1415[70 33]14161417Check trace (123 = -2 mod 125) and determinant::14181419sage: B.det()142051421sage: B.trace()14221231423sage: EllipticCurve([-1, 1/4]).ap(5)1424-214251426Try using the trace to speed up the calculation::14271428sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 - x + R(1/4),1429... p, M, -2)1430sage: A1431[2715 187]1432[1445 408]14331434Hmmm... it looks different, but that's because the trace of our1435first answer was only -2 modulo `5^3`, not -2 modulo1436`5^5`. So the right answer is::14371438sage: A.change_ring(Integers(p**prec))1439[90 62]1440[70 33]14411442Check it works with only one digit of precision::14431444sage: p = 51445sage: prec = 11446sage: M = monsky_washnitzer.adjusted_prec(p, prec)1447sage: R.<x> = PolynomialRing(Integers(p**M))1448sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 - x + R(1/4), p, M)1449sage: A.change_ring(Integers(p))1450[0 2]1451[0 3]14521453Here is an example that is particularly badly conditioned for1454using the trace trick::14551456sage: p = 111457sage: prec = 31458sage: M = monsky_washnitzer.adjusted_prec(p, prec)1459sage: R.<x> = PolynomialRing(Integers(p**M))1460sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 + 7*x + 8, p, M)1461sage: A.change_ring(Integers(p**prec))1462[1144 176]1463[ 847 185]14641465The problem here is that the top-right entry is divisible by 11,1466and the bottom-left entry is divisible by `11^2`. So when1467you apply the trace trick, neither `F(dx/y)` nor1468`F(x dx/y)` is enough to compute the whole matrix to the1469desired precision, even if you try increasing the target precision1470by one. Nevertheless, ``matrix_of_frobenius`` knows1471how to get the right answer by evaluating `F((x+1) dx/y)`1472instead::14731474sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 + 7*x + 8, p, M, -2)1475sage: A.change_ring(Integers(p**prec))1476[1144 176]1477[ 847 185]14781479The running time is about ``O(p*prec**2)`` (times some logarithmic1480factors), so it is feasible to run on fairly large primes, or1481precision (or both?!?!)::14821483sage: p = 100071484sage: prec = 21485sage: M = monsky_washnitzer.adjusted_prec(p, prec)1486sage: R.<x> = PolynomialRing(Integers(p**M))1487sage: A = monsky_washnitzer.matrix_of_frobenius( # long time1488... x^3 - x + R(1/4), p, M) # long time1489sage: B = A.change_ring(Integers(p**prec)); B # long time1490[74311982 57996908]1491[95877067 25828133]1492sage: B.det() # long time1493100071494sage: B.trace() # long time1495661496sage: EllipticCurve([-1, 1/4]).ap(10007) # long time14976614981499::15001501sage: p = 51502sage: prec = 3001503sage: M = monsky_washnitzer.adjusted_prec(p, prec)1504sage: R.<x> = PolynomialRing(Integers(p**M))1505sage: A = monsky_washnitzer.matrix_of_frobenius( # long time1506... x^3 - x + R(1/4), p, M) # long time1507sage: B = A.change_ring(Integers(p**prec)) # long time1508sage: B.det() # long time150951510sage: -B.trace() # long time151121512sage: EllipticCurve([-1, 1/4]).ap(5) # long time1513-215141515Let us check consistency of the results for a range of precisions::15161517sage: p = 51518sage: max_prec = 601519sage: M = monsky_washnitzer.adjusted_prec(p, max_prec)1520sage: R.<x> = PolynomialRing(Integers(p**M))1521sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 - x + R(1/4), p, M) # long time1522sage: A = A.change_ring(Integers(p**max_prec)) # long time1523sage: result = [] # long time1524sage: for prec in range(1, max_prec): # long time1525... M = monsky_washnitzer.adjusted_prec(p, prec) # long time1526... R.<x> = PolynomialRing(Integers(p^M),'x') # long time1527... B = monsky_washnitzer.matrix_of_frobenius( # long time1528... x^3 - x + R(1/4), p, M) # long time1529... B = B.change_ring(Integers(p**prec)) # long time1530... result.append(B == A.change_ring( # long time1531... Integers(p**prec))) # long time1532sage: result == [True] * (max_prec - 1) # long time1533True15341535The remaining examples discuss what happens when you take the1536coefficient ring to be a power series ring; i.e. in effect you're1537looking at a family of curves.15381539The code does in fact work...15401541::15421543sage: p = 111544sage: prec = 31545sage: M = monsky_washnitzer.adjusted_prec(p, prec)1546sage: S.<t> = PowerSeriesRing(Integers(p**M), default_prec=4)1547sage: a = 7 + t + 3*t^21548sage: b = 8 - 6*t + 17*t^21549sage: R.<x> = PolynomialRing(S)1550sage: Q = x**3 + a*x + b1551sage: A = monsky_washnitzer.matrix_of_frobenius(Q, p, M) # long time1552sage: B = A.change_ring(PowerSeriesRing(Integers(p**prec), 't', default_prec=4)) # long time1553sage: B # long time1554[1144 + 264*t + 841*t^2 + 1025*t^3 + O(t^4) 176 + 1052*t + 216*t^2 + 523*t^3 + O(t^4)]1555[ 847 + 668*t + 81*t^2 + 424*t^3 + O(t^4) 185 + 341*t + 171*t^2 + 642*t^3 + O(t^4)]15561557The trace trick should work for power series rings too, even in the1558badly- conditioned case. Unfortunately I don't know how to compute1559the trace in advance, so I'm not sure exactly how this would help.1560Also, I suspect the running time will be dominated by the1561expansion, so the trace trick won't really speed things up anyway.1562Another problem is that the determinant is not always p::15631564sage: B.det() # long time156511 + 484*t^2 + 451*t^3 + O(t^4)15661567However, it appears that the determinant always has the property1568that if you substitute t - 11t, you do get the constant series p1569(mod p\*\*prec). Similarly for the trace. And since the parameter1570only really makes sense when it is divisible by p anyway, perhaps1571this isn't a problem after all.1572"""15731574M = int(M)1575if M < 2:1576raise ValueError("M (=%s) must be at least 2" % M)15771578base_ring = Q.base_ring()15791580# Expand out frobenius of dx/y and x dx/y.1581# (You can substitute frobenius_expansion_by_series here, that will work1582# as well. See its docstring for some performance notes.)1583F0, F1, offset = frobenius_expansion_by_newton(Q, p, M)1584#F0, F1, offset = frobenius_expansion_by_series(Q, p, M)15851586if compute_exact_forms:1587# we need to do all the work to get the exact expressions f such that F(x^i dx/y) = df + \sum a_i x^i dx/y1588F0_coeffs = transpose_list(F0.coeffs())1589F0_reduced, f_0 = reduce_all(Q, p, F0_coeffs, offset, True)15901591F1_coeffs = transpose_list(F1.coeffs())1592F1_reduced, f_1 = reduce_all(Q, p, F1_coeffs, offset, True)15931594elif M == 2:1595# This implies that only one digit of precision is valid, so we only need1596# to reduce the second column. Also, the trace doesn't help at all.15971598F0_reduced = [base_ring(0), base_ring(0)]15991600F1_coeffs = transpose_list(F1.coeffs())1601F1_reduced = reduce_all(Q, p, F1_coeffs, offset)16021603elif trace is None:1604# No trace provided, just reduce F(dx/y) and F(x dx/y) separately.16051606F0_coeffs = transpose_list(F0.coeffs())1607F0_reduced = reduce_all(Q, p, F0_coeffs, offset)16081609F1_coeffs = transpose_list(F1.coeffs())1610F1_reduced = reduce_all(Q, p, F1_coeffs, offset)16111612else:1613# Trace has been provided.16141615# In most cases this can be used to quickly compute F(dx/y) from1616# F(x dx/y). However, if we're unlucky, the (dx/y)-component of1617# F(x dx/y) (i.e. the top-right corner of the matrix) may be divisible1618# by p, in which case there isn't enough information to get the1619# (x dx/y)-component of F(dx/y) to the desired precision. When this1620# happens, it turns out that F((x+1) dx/y) always *does* give enough1621# information (together with the trace) to get both columns to the1622# desired precision.16231624# First however we need a quick way of telling whether the top-right1625# corner is divisible by p, i.e. we want to compute the second column1626# of the matrix mod p. We could do this by just running the entire1627# algorithm with M = 2 (which assures precision 1). Luckily, we've1628# already done most of the work by computing F1 to high precision; so1629# all we need to do is extract the coefficients that would correspond1630# to the first term of the series, and run the reduction on them.16311632# todo: actually we only need to do this reduction step mod p^2, not1633# mod p^M, which is what the code currently does. If the base ring1634# is Integers(p^M), then it is easy. Otherwise it is tricky to construct1635# the right ring, I don't know how to do it.16361637F1_coeffs = transpose_list(F1.coeffs())1638F1_modp_coeffs = F1_coeffs[int((M-2)*p):]1639# make a copy, because reduce_all will destroy the coefficients:1640F1_modp_coeffs = [[cell for cell in row] for row in F1_modp_coeffs]1641F1_modp_offset = offset - (M-2)*p1642F1_modp_reduced = reduce_all(Q, p, F1_modp_coeffs, F1_modp_offset)16431644if F1_modp_reduced[0].is_unit():1645# If the first entry is invertible mod p, then F(x dx/y) is sufficient1646# to get the whole matrix.16471648F1_reduced = reduce_all(Q, p, F1_coeffs, offset)16491650F0_reduced = [base_ring(trace) - F1_reduced[1], None]1651# using that the determinant is p:1652F0_reduced[1] = (F0_reduced[0] * F1_reduced[1] - base_ring(p)) \1653/ F1_reduced[0]16541655else:1656# If the first entry is zero mod p, then F((x+1) dx/y) will be sufficient1657# to get the whole matrix. (Here we are using the fact that the second1658# entry *cannot* be zero mod p. This is guaranteed by some results in1659# section 3.2 of ``Computation of p-adic Heights and Log Convergence''1660# by Mazur, Stein, Tate. But let's quickly check it anyway :-))1661msg = "The second entry in the second column "1662msg += "should be invertible mod p!"1663assert F1_modp_reduced[1].is_unit(), msg16641665G0_coeffs = transpose_list((F0 + F1).coeffs())1666G0_reduced = reduce_all(Q, p, G0_coeffs, offset)16671668# Now G0_reduced expresses F((x+1) dx/y) in terms of dx/y and x dx/y.1669# Re-express this in terms of (x+1) dx/y and x dx/y.1670H0_reduced = [G0_reduced[0], G0_reduced[1] - G0_reduced[0]]16711672# The thing we're about to divide by better be a unit.1673msg = "The second entry in this column "1674msg += "should be invertible mod p!"1675assert H0_reduced[1].is_unit(), msg16761677# Figure out the second column using the trace...1678H1_reduced = [None, base_ring(trace) - H0_reduced[0]]1679# ... and using that the determinant is p:1680H1_reduced[0] = (H0_reduced[0] * H1_reduced[1] - base_ring(p)) \1681/ H0_reduced[1]16821683# Finally, change back to the usual basis (dx/y, x dx/y)1684F1_reduced = [H1_reduced[0],1685H1_reduced[0] + H1_reduced[1]]1686F0_reduced = [H0_reduced[0] - F1_reduced[0],1687H0_reduced[0] + H0_reduced[1] - F1_reduced[1]]16881689# One more sanity check: our final result should be congruent mod p1690# to the approximation we used earlier.1691msg = "The output matrix is not congruent mod p "1692msg += "to the approximation found earlier!"1693assert not (1694(F1_reduced[0] - F1_modp_reduced[0]).is_unit() or1695(F1_reduced[1] - F1_modp_reduced[1]).is_unit() or1696F0_reduced[0].is_unit() or F0_reduced[1].is_unit()), msg16971698if compute_exact_forms:1699return matrix(base_ring, 2, 2, [F0_reduced[0], F1_reduced[0],1700F0_reduced[1], F1_reduced[1]]), f_0, f_11701else:1702return matrix(base_ring, 2, 2, [F0_reduced[0], F1_reduced[0],1703F0_reduced[1], F1_reduced[1]])170417051706#*****************************************************************************1707# This is a generalization of the above functionality for hyperelliptic curves.1708#1709# THIS IS A WORK IN PROGRESS.1710#1711# I tried to embed must stuff into the rings themselves rather than1712# just extract and manipulate lists of coefficients. Hence the implementations1713# below are much less optimized, so are much slower, but should hopefully be1714# easier to follow. (E.g. one can print/make sense of intermediate results.)1715#1716# AUTHOR:1717# -- Robert Bradshaw (2007-04)1718#1719#*****************************************************************************17201721import weakref17221723from sage.schemes.hyperelliptic_curves.all import is_HyperellipticCurve, HyperellipticCurve1724from sage.rings.padics.all import pAdicField1725from sage.rings.all import QQ17261727from sage.rings.laurent_series_ring import is_LaurentSeriesRing1728from sage.rings.integral_domain import is_IntegralDomain17291730from sage.modules.all import FreeModule, is_FreeModuleElement17311732from sage.misc.profiler import Profiler1733from sage.misc.misc import repr_lincomb173417351736def matrix_of_frobenius_hyperelliptic(Q, p=None, prec=None, M=None):1737"""1738Computes the matrix of Frobenius on Monsky-Washnitzer cohomology,1739with respect to the basis `(dx/2y, x dx/2y, ...x^{d-2} dx/2y)`, where1740`d` is the degree of `Q`.17411742INPUT:17431744- ``Q`` -- monic polynomial `Q(x)`17451746- ``p`` -- prime `\geq 5` for which `E` has good reduction17471748- ``prec`` -- (optional) `p`-adic precision of the coefficient ring17491750- ``M`` -- (optional) adjusted `p`-adic precision of the coefficint ring17511752OUTPUT:17531754`(d-1)` x `(d-1)` matrix `M` of Frobenius on Monsky-Washnitzer cohomology,1755and list of differentials \{f_i \} such that17561757.. math::17581759\phi^* (x^i dx/2y) = df_i + M[i]*vec(dx/2y, ..., x^{d-2} dx/2y)17601761EXAMPLES::17621763sage: p = 51764sage: prec = 31765sage: R.<x> = QQ['x']1766sage: A,f = monsky_washnitzer.matrix_of_frobenius_hyperelliptic(x^5 - 2*x + 3, p, prec)1767sage: A1768[ 4*5 + O(5^3) 5 + 2*5^2 + O(5^3) 2 + 3*5 + 2*5^2 + O(5^3) 2 + 5 + 5^2 + O(5^3)]1769[ 3*5 + 5^2 + O(5^3) 3*5 + O(5^3) 4*5 + O(5^3) 2 + 5^2 + O(5^3)]1770[ 4*5 + 4*5^2 + O(5^3) 3*5 + 2*5^2 + O(5^3) 5 + 3*5^2 + O(5^3) 2*5 + 2*5^2 + O(5^3)]1771[ 5^2 + O(5^3) 5 + 4*5^2 + O(5^3) 4*5 + 3*5^2 + O(5^3) 2*5 + O(5^3)]17721773"""1774prof = Profiler()1775prof("setup")1776if p is None:1777try:1778K = Q.base_ring()1779p = K.prime()1780prec = K.precision_cap()1781except AttributeError:1782raise ValueError("p and prec must be specified if Q is not "1783"defined over a p-adic ring")1784if M is None:1785M = adjusted_prec(p, prec)1786extra_prec_ring = Integers(p**M)1787# extra_prec_ring = pAdicField(p, M) # SLOW!17881789real_prec_ring = pAdicField(p, prec) # pAdicField(p, prec) # To capped absolute?1790S = SpecialHyperellipticQuotientRing(Q, extra_prec_ring, True)1791MW = S.monsky_washnitzer()1792prof("frob basis elements")1793F = MW.frob_basis_elements(M, p)17941795prof("rationalize")1796# do reduction over Q in case we have non-integral entries (and it is so much faster than padics)1797rational_S = S.change_ring(QQ)1798# this is a hack until pAdics are fast1799# (They are in the latest development bundle, but its not standard and I'd need to merge.1800# (it will periodically cast into this ring to reduce coefficient size)1801rational_S._prec_cap = p**M1802rational_S._p = p1803# S._p = p1804# rational_S(F[0]).reduce_fast()1805# prof("reduce others")18061807# rational_S = S.change_ring(pAdicField(p, M))1808F = [rational_S(F_i) for F_i in F]18091810prof("reduce")1811reduced = [F_i.reduce_fast(True) for F_i in F]1812# reduced = [F_i.reduce() for F_i in F]18131814#print reduced[0][0].diff() - F[0]18151816# but the coeffs are WAY more precision than they need to be1817# print reduced[0][1]18181819prof("make matrix")1820# now take care of precision capping1821M = matrix(real_prec_ring, [a for f, a in reduced])1822for i in range(M.ncols()):1823for j in range(M.nrows()):1824M[i, j] = M[i, j].add_bigoh(prec)1825# print prof1826return M.transpose(), [f for f, a in reduced]182718281829# For uniqueness (as many of the non-trivial calculations are cached along the way).18301831_special_ring_cache = {}1832_mw_cache = {}183318341835def SpecialHyperellipticQuotientRing(*args):1836if args in _special_ring_cache:1837R = _special_ring_cache[args]()1838if R is not None:1839return R1840R = SpecialHyperellipticQuotientRing_class(*args)1841_special_ring_cache[args] = weakref.ref(R)1842return R184318441845def MonskyWashnitzerDifferentialRing(base_ring):1846if base_ring in _mw_cache:1847R = _mw_cache[base_ring]()1848if R is not None:1849return R18501851R = MonskyWashnitzerDifferentialRing_class(base_ring)1852_mw_cache[base_ring] = weakref.ref(R)1853return R185418551856class SpecialHyperellipticQuotientRing_class(CommutativeAlgebra):18571858_p = None18591860def __init__(self, Q, R=None, invert_y=True):1861if R is None:1862R = Q.base_ring()18631864# Trac ticket #9138: CommutativeAlgebra.__init__ must not be1865# done so early. It tries to register a coercion, but that1866# requires the hash bein available. But the hash, in its1867# default implementation, relies on the string representation,1868# which is not available at this point.1869#CommutativeAlgebra.__init__(self, R) # moved to below.18701871x = PolynomialRing(R, 'xx').gen(0)1872if is_EllipticCurve(Q):1873E = Q1874if E.a1() != 0 or E.a2() != 0:1875raise NotImplementedError("Curve must be in Weierstrass "1876"normal form.")1877Q = -E.change_ring(R).defining_polynomial()(x, 0, 1)1878self._curve = E18791880elif is_HyperellipticCurve(Q):1881C = Q1882if C.hyperelliptic_polynomials()[1] != 0:1883raise NotImplementedError("Curve must be of form y^2 = Q(x).")1884Q = C.hyperelliptic_polynomials()[0].change_ring(R)1885self._curve = C18861887if is_Polynomial(Q):1888self._Q = Q.change_ring(R)1889self._coeffs = self._Q.coeffs()1890if self._coeffs.pop() != 1:1891raise NotImplementedError("Polynomial must be monic.")1892if not hasattr(self, '_curve'):1893if self._Q.degree() == 3:1894ainvs = [0, self._Q[2], 0, self._Q[1], self._Q[0]]1895self._curve = EllipticCurve(ainvs)1896else:1897self._curve = HyperellipticCurve(self._Q)18981899else:1900raise NotImplementedError("Must be an elliptic curve or polynomial "1901"Q for y^2 = Q(x)\n(Got element of %s)" % Q.parent())19021903self._n = int(Q.degree())1904self._series_ring = (LaurentSeriesRing if invert_y else PolynomialRing)(R, 'y')1905self._series_ring_y = self._series_ring.gen(0)1906self._series_ring_0 = self._series_ring(0)19071908# Trac ticket #9138: Initialise the commutative algebra here!1909# Below, we do self(self._poly_ring.gen(0)), which requires1910# the initialisation being finished.1911CommutativeAlgebra.__init__(self, R)19121913self._poly_ring = PolynomialRing(self._series_ring, 'x')19141915self._x = self(self._poly_ring.gen(0))1916self._y = self(self._series_ring.gen(0))19171918self._Q_coeffs = Q.change_ring(self._series_ring).list()1919self._dQ = Q.derivative().change_ring(self)(self._x)1920self._monsky_washnitzer = MonskyWashnitzerDifferentialRing(self)19211922self._monomial_diffs = {}1923self._monomial_diff_coeffs = {}19241925def _repr_(self):1926"""1927String representation19281929EXAMPLES::19301931sage: R.<x> = QQ['x']1932sage: E = HyperellipticCurve(x^5-3*x+1)1933sage: x,y = E.monsky_washnitzer_gens()1934sage: x.parent() # indirect doctest1935SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 3*x + 1) over Rational Field1936"""1937y_inverse = ",y^-1" if is_LaurentSeriesRing(self._series_ring) else ""1938return "SpecialHyperellipticQuotientRing K[x,y%s] / (y^2 = %s) over %s" % (y_inverse, self._Q, self.base_ring())19391940def base_extend(self, R):1941"""1942Return the base extension of ``self`` to the ring ``R`` if possible.19431944EXAMPLES::19451946sage: R.<x> = QQ['x']1947sage: E = HyperellipticCurve(x^5-3*x+1)1948sage: x,y = E.monsky_washnitzer_gens()1949sage: x.parent().base_extend(UniversalCyclotomicField())1950SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 3*x + 1) over Universal Cyclotomic Field1951sage: x.parent().base_extend(ZZ)1952Traceback (most recent call last):1953...1954TypeError: no such base extension1955"""1956if R.has_coerce_map_from(self.base_ring()):1957return self.change_ring(R)1958else:1959raise TypeError("no such base extension")19601961def change_ring(self, R):1962"""1963Return the analog of ``self`` over the ring ``R``19641965EXAMPLES::19661967sage: R.<x> = QQ['x']1968sage: E = HyperellipticCurve(x^5-3*x+1)1969sage: x,y = E.monsky_washnitzer_gens()1970sage: x.parent().change_ring(ZZ)1971SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 3*x + 1) over Integer Ring1972"""1973return SpecialHyperellipticQuotientRing(self._Q, R, is_LaurentSeriesRing(self._series_ring))19741975def __call__(self, val, offset=0, check=True):1976if isinstance(val, SpecialHyperellipticQuotientElement) and val.parent() is self:1977if offset == 0:1978return val1979else:1980return val << offset1981elif isinstance(val, MonskyWashnitzerDifferential):1982return self._monsky_washnitzer(val)1983return SpecialHyperellipticQuotientElement(self, val, offset, check)19841985def gens(self):1986"""1987Return the generators of ``self``19881989EXAMPLES::19901991sage: R.<x> = QQ['x']1992sage: E = HyperellipticCurve(x^5-3*x+1)1993sage: x,y = E.monsky_washnitzer_gens()1994sage: x.parent().gens()1995(x, y*1)1996"""1997return self._x, self._y19981999def x(self):2000r"""2001Return the generator `x` of ``self``20022003EXAMPLES::20042005sage: R.<x> = QQ['x']2006sage: E = HyperellipticCurve(x^5-3*x+1)2007sage: x,y = E.monsky_washnitzer_gens()2008sage: x.parent().x()2009x2010"""2011return self._x20122013def y(self):2014r"""2015Return the generator `y` of ``self``20162017EXAMPLES::20182019sage: R.<x> = QQ['x']2020sage: E = HyperellipticCurve(x^5-3*x+1)2021sage: x,y = E.monsky_washnitzer_gens()2022sage: x.parent().y()2023y*12024"""2025return self._y20262027def monomial(self, i, j, b=None):2028"""2029Returns `b y^j x^i`, computed quickly.20302031EXAMPLES::20322033sage: R.<x> = QQ['x']2034sage: E = HyperellipticCurve(x^5-3*x+1)2035sage: x,y = E.monsky_washnitzer_gens()2036sage: x.parent().monomial(4,5)2037y^5*x^42038"""2039i = int(i)2040j = int(j)20412042if 0 < i and i < self._n:2043if b is None:2044by_to_j = self._series_ring_y << (j-1)2045else:2046by_to_j = self._series_ring(b) << j2047v = [self._series_ring_0] * self._n2048v[i] = by_to_j2049return self(v)2050else:2051return (self._x ** i) << j if b is None else self.base_ring()(b) * (self._x ** i) << j20522053def monomial_diff_coeffs(self, i, j):2054r"""2055The key here is that the formula for `d(x^iy^j)` is messy2056in terms of `i`, but varies nicely with `j`.20572058.. math::20592060d(x^iy^j) = y^{j-1} (2ix^{i-1}y^2 + j (A_i(x) + B_i(x)y^2)) \frac{dx}{2y}206120622063Where `A,B` have degree at most `n-1` for each2064`i`. Pre-compute `A_i, B_i` for each `i`2065the "hard" way, and the rest are easy.2066"""2067try:2068return self._monomial_diff_coeffs[i, j]2069except KeyError:2070pass2071if i < self._n:2072try:2073A, B, two_i_x_to_i = self._precomputed_diff_coeffs[i]2074except AttributeError:2075self._precomputed_diff_coeffs = self._precompute_monomial_diffs()2076A, B, two_i_x_to_i = self._precomputed_diff_coeffs[i]2077if i == 0:2078return j*A, j*B2079else:2080return j*A, j*B + two_i_x_to_i2081else:2082dg = self.monomial(i, j).diff()2083coeffs = [dg.extract_pow_y(j-1), dg.extract_pow_y(j+1)]2084self._monomial_diff_coeffs[i, j] = coeffs2085return coeffs20862087def monomial_diff_coeffs_matrices(self):2088self.monomial_diff_coeffs(0, 0) # precompute stuff2089R = self.base_ring()2090mat_1 = matrix(R, self._n, self._n)2091mat_2 = matrix(R, self._n, self._n)2092for i in range(self._n):2093mat_1[i] = self._precomputed_diff_coeffs[i][1]2094mat_2[i] = self._precomputed_diff_coeffs[i][2]2095return mat_1.transpose(), mat_2.transpose()20962097def _precompute_monomial_diffs(self):2098x, y = self.gens()2099R = self.base_ring()2100V = FreeModule(R, self.degree())2101As = []2102for i in range(self.degree()):2103dg = self.monomial(i, 1).diff()2104two_i_x_to_i = R(2*i) * x**(i-1) * y*y if i > 0 else self(0)2105A = dg - self._monsky_washnitzer(two_i_x_to_i)2106As.append((V(A.extract_pow_y(0)), V(A.extract_pow_y(2)), V(two_i_x_to_i.extract_pow_y(2))))2107return As21082109def Q(self):2110"""2111Return the defining polynomial of the underlying hyperelliptic curve.21122113EXAMPLES::21142115sage: R.<x> = QQ['x']2116sage: E = HyperellipticCurve(x^5-2*x+1)2117sage: x,y = E.monsky_washnitzer_gens()2118sage: x.parent().Q()2119x^5 - 2*x + 12120"""2121return self._Q21222123def curve(self):2124"""2125Return the underlying hyperelliptic curve.21262127EXAMPLES::21282129sage: R.<x> = QQ['x']2130sage: E = HyperellipticCurve(x^5-3*x+1)2131sage: x,y = E.monsky_washnitzer_gens()2132sage: x.parent().curve()2133Hyperelliptic Curve over Rational Field defined by y^2 = x^5 - 3*x + 12134"""2135return self._curve21362137def degree(self):2138"""2139Return the degree of the underlying hyperelliptic curve.21402141EXAMPLES::21422143sage: R.<x> = QQ['x']2144sage: E = HyperellipticCurve(x^5-3*x+1)2145sage: x,y = E.monsky_washnitzer_gens()2146sage: x.parent().degree()214752148"""2149return self._n21502151def prime(self):2152return self._p21532154def monsky_washnitzer(self):2155return self._monsky_washnitzer21562157def is_field(self, proof=True):2158"""2159Return False as ``self`` is not a field.21602161EXAMPLES::21622163sage: R.<x> = QQ['x']2164sage: E = HyperellipticCurve(x^5-3*x+1)2165sage: x,y = E.monsky_washnitzer_gens()2166sage: x.parent().is_field()2167False2168"""2169return False217021712172class SpecialHyperellipticQuotientElement(CommutativeAlgebraElement):21732174def __init__(self, parent, val=0, offset=0, check=True):2175"""2176Elements in the Hyperelliptic quotient ring21772178EXAMPLES::21792180sage: R.<x> = QQ['x']2181sage: E = HyperellipticCurve(x^5-36*x+1)2182sage: x,y = E.monsky_washnitzer_gens()2183sage: MW = x.parent()2184sage: MW(x+x**2+y-77) # indirect doctest2185-(77-y)*1 + x + x^22186"""2187CommutativeAlgebraElement.__init__(self, parent)2188if not check:2189self._f = parent._poly_ring(val, check=False)2190return2191if isinstance(val, SpecialHyperellipticQuotientElement):2192R = parent.base_ring()2193self._f = parent._poly_ring([a.change_ring(R) for a in val._f])2194return2195if isinstance(val, tuple):2196val, offset = val2197if isinstance(val, list) and len(val) > 0 and is_FreeModuleElement(val[0]):2198val = transpose_list(val)2199self._f = parent._poly_ring(val)2200if offset != 0:2201self._f = self._f.parent()([a << offset for a in self._f], check=False)22022203def __cmp__(self, other):2204"""2205Compares the elements22062207EXAMPLES::22082209sage: R.<x> = QQ['x']2210sage: E = HyperellipticCurve(x^5-36*x+1)2211sage: x,y = E.monsky_washnitzer_gens()2212sage: x.__cmp__(x)221302214sage: x.__cmp__(y)221512216"""2217return cmp(self._f, other._f)22182219def change_ring(self, R):2220"""2221Return the same element after changing the base ring to R22222223EXAMPLES::22242225sage: R.<x> = QQ['x']2226sage: E = HyperellipticCurve(x^5-36*x+1)2227sage: x,y = E.monsky_washnitzer_gens()2228sage: MW = x.parent()2229sage: z = MW(x+x**2+y-77)2230sage: z.change_ring(AA).parent()2231SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 36*x + 1) over Algebraic Real Field2232"""2233return self.parent().change_ring(R)(self)22342235def __call__(self, *x):2236"""2237Evaluate ``self`` at given arguments22382239EXAMPLES::22402241sage: R.<x> = QQ['x']2242sage: E = HyperellipticCurve(x^5-36*x+1)2243sage: x,y = E.monsky_washnitzer_gens()2244sage: MW = x.parent()2245sage: z = MW(x+x**2+y-77); z2246-(77-y)*1 + x + x^22247sage: z(66)22484345 + y2249sage: z(5,4)2250-432251"""2252return self._f(*x)22532254def __invert__(self):2255"""2256Return the inverse of the element22572258The general element in our ring is not invertible, but `y` may2259be. We do not want to pass to the fraction field.22602261EXAMPLES::22622263sage: R.<x> = QQ['x']2264sage: E = HyperellipticCurve(x^5-36*x+1)2265sage: x,y = E.monsky_washnitzer_gens()2266sage: MW = x.parent()2267sage: z = y**(-1) # indirect doctest2268sage: z.parent()2269SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 36*x + 1) over Rational Field22702271sage: z = (x+y)**(-1) # indirect doctest2272Traceback (most recent call last):2273...2274ZeroDivisionError: Element not invertible2275"""2276if self._f.degree() == 0 and self._f[0].is_unit():2277return SpecialHyperellipticQuotientElement(self.parent(), ~self._f[0])2278else:2279raise ZeroDivisionError("Element not invertible")22802281def __nonzero__(self):2282"""2283Return True iff ``self`` is not zero.22842285EXAMPLES::22862287sage: R.<x> = QQ['x']2288sage: E = HyperellipticCurve(x^5-3*x+1)2289sage: x,y = E.monsky_washnitzer_gens()2290sage: x.__nonzero__()2291True2292"""2293return not not self._f22942295def __eq__(self, other):2296"""2297Return True iff ``self`` is equal to ``other``22982299EXAMPLES::23002301sage: R.<x> = QQ['x']2302sage: E = HyperellipticCurve(x^5-3*x+1)2303sage: x,y = E.monsky_washnitzer_gens()2304sage: x == y # indirect doctest2305False2306"""2307if not isinstance(other, SpecialHyperellipticQuotientElement):2308other = self.parent()(other)2309return self._f == other._f23102311def _add_(self, other):2312"""2313Return the sum of two elements23142315EXAMPLES::23162317sage: R.<x> = QQ['x']2318sage: E = HyperellipticCurve(x^5-36*x+1)2319sage: x,y = E.monsky_washnitzer_gens()2320sage: x+y2321y*1 + x2322"""2323return SpecialHyperellipticQuotientElement(self.parent(), self._f + other._f)23242325def _sub_(self, other):2326"""2327Return the difference of two elements23282329EXAMPLES::23302331sage: R.<x> = QQ['x']2332sage: E = HyperellipticCurve(x^5-36*x+1)2333sage: x,y = E.monsky_washnitzer_gens()2334sage: y-x2335y*1 - x2336"""2337return SpecialHyperellipticQuotientElement(self.parent(), self._f - other._f)23382339def _mul_(self, other):2340"""2341Return the product of two elements23422343EXAMPLES::23442345sage: R.<x> = QQ['x']2346sage: E = HyperellipticCurve(x^5-36*x+1)2347sage: x,y = E.monsky_washnitzer_gens()2348sage: y*x2349y*x2350"""2351# over laurent series, addition and subtraction can be2352# expensive, and the degree of this poly is small enough that2353# Karatsuba actually hurts significantly in some cases2354if self._f[0].valuation() + other._f[0].valuation() > -200:2355prod = self._f._mul_generic(other._f)2356else:2357prod = self._f * other._f2358v = prod.list()2359parent = self.parent()2360Q_coeffs = parent._Q_coeffs2361n = len(Q_coeffs) - 12362y2 = self.parent()._series_ring_y << 12363for i in range(len(v)-1, n-1, -1):2364for j in range(n):2365v[i-n+j] -= Q_coeffs[j] * v[i]2366v[i-n] += y2 * v[i]2367return SpecialHyperellipticQuotientElement(parent, v[0:n])23682369def _rmul_(self, c):2370coeffs = self._f.list(copy=False)2371return self.parent()([c*a for a in coeffs], check=False)23722373def _lmul_(self, c):2374coeffs = self._f.list(copy=False)2375return self.parent()([a*c for a in coeffs], check=False)23762377def __lshift__(self, k):2378coeffs = self._f.list(copy=False)2379return self.parent()([a << k for a in coeffs], check=False)23802381def __rshift__(self, k):2382coeffs = self._f.list(copy=False)2383return self.parent()([a >> k for a in coeffs], check=False)23842385def truncate_neg(self, n):2386"""2387Return ``self`` minus its terms of degree less than `n` wrt `y`.23882389EXAMPLES::23902391sage: R.<x> = QQ['x']2392sage: E = HyperellipticCurve(x^5-3*x+1)2393sage: x,y = E.monsky_washnitzer_gens()2394sage: (x+3*y+7*x*2*y**4).truncate_neg(1)23953*y*1 + 14*y^4*x2396"""2397coeffs = self._f.list(copy=False)2398return self.parent()([a.truncate_neg(n) for a in coeffs], check=False)23992400def _repr_(self):2401"""2402Return a string representation of ``self``.24032404EXAMPLES::24052406sage: R.<x> = QQ['x']2407sage: E = HyperellipticCurve(x^5-3*x+1)2408sage: x,y = E.monsky_washnitzer_gens()2409sage: (x+3*y)._repr_()2410'3*y*1 + x'2411"""2412x = PolynomialRing(QQ, 'x').gen(0)2413coeffs = self._f.list()2414return repr_lincomb([(x**i, coeffs[i]) for i in range(len(coeffs))])24152416def _latex_(self):2417"""2418Return a LateX string for ``self``.24192420EXAMPLES::24212422sage: R.<x> = QQ['x']2423sage: E = HyperellipticCurve(x^5-3*x+1)2424sage: x,y = E.monsky_washnitzer_gens()2425sage: (x+3*y)._latex_()2426'3y1 + x'2427"""2428x = PolynomialRing(QQ, 'x').gen(0)2429coeffs = self._f.list()2430return repr_lincomb([(x**i, coeffs[i]) for i in range(len(coeffs))], is_latex=True)24312432def diff(self):2433"""2434Return the differential of ``self``24352436EXAMPLES::24372438sage: R.<x> = QQ['x']2439sage: E = HyperellipticCurve(x^5-3*x+1)2440sage: x,y = E.monsky_washnitzer_gens()2441sage: (x+3*y).diff()2442(-(9-2*y)*1 + 15*x^4) dx/2y2443"""2444# try:2445# return self._diff_x2446# except AttributeError:2447# pass24482449# d(self) = A dx + B dy2450# = (2y A + BQ') dx/2y2451parent = self.parent()2452R = parent.base_ring()2453x, y = parent.gens()2454v = self._f.list()2455n = len(v)2456A = parent([R(i) * v[i] for i in range(1, n)])2457B = parent([a.derivative() for a in v])2458dQ = parent._dQ2459return parent._monsky_washnitzer((R(2) * A << 1) + dQ * B)2460# self._diff = self.parent()._monsky_washnitzer(two_y * A + dQ * B)2461# return self._diff24622463def extract_pow_y(self, k):2464r"""2465Return the coefficients of `y^k` in ``self`` as a list24662467EXAMPLES::24682469sage: R.<x> = QQ['x']2470sage: E = HyperellipticCurve(x^5-3*x+1)2471sage: x,y = E.monsky_washnitzer_gens()2472sage: (x+3*y+9*x*y).extract_pow_y(1)2473[3, 9, 0, 0, 0]2474"""2475v = [a[k] for a in self._f.list()]2476while len(v) < self.parent()._n:2477v.append(0)2478return v24792480def min_pow_y(self):2481"""2482Return the minimal degree of ``self`` w.r.t. y24832484EXAMPLES::24852486sage: R.<x> = QQ['x']2487sage: E = HyperellipticCurve(x^5-3*x+1)2488sage: x,y = E.monsky_washnitzer_gens()2489sage: (x+3*y).min_pow_y()249002491"""2492if self._f.degree() == -1:2493return 02494return min([a.valuation() for a in self._f.list()])24952496def max_pow_y(self):2497"""2498Return the maximal degree of ``self`` w.r.t. y24992500EXAMPLES::25012502sage: R.<x> = QQ['x']2503sage: E = HyperellipticCurve(x^5-3*x+1)2504sage: x,y = E.monsky_washnitzer_gens()2505sage: (x+3*y).max_pow_y()250612507"""2508if self._f.degree() == -1:2509return 02510return max([a.degree() for a in self._f.list()])25112512def coeffs(self, R=None):2513"""2514Returns the raw coefficients of this element.25152516INPUT:25172518- ``R`` -- an (optional) base-ring in which to cast the coefficients25192520OUTPUT:25212522- ``coeffs`` -- a list of coefficients of powers of `x` for each power2523of `y`25242525- ``n`` -- an offset indicating the power of `y` of the first list2526element25272528EXAMPLES::25292530sage: R.<x> = QQ['x']2531sage: E = HyperellipticCurve(x^5-3*x+1)2532sage: x,y = E.monsky_washnitzer_gens()2533sage: x.coeffs()2534([(0, 1, 0, 0, 0)], 0)2535sage: y.coeffs()2536([(0, 0, 0, 0, 0), (1, 0, 0, 0, 0)], 0)25372538sage: a = sum(n*x^n for n in range(5)); a2539x + 2*x^2 + 3*x^3 + 4*x^42540sage: a.coeffs()2541([(0, 1, 2, 3, 4)], 0)2542sage: a.coeffs(Qp(7))2543([(0, 1 + O(7^20), 2 + O(7^20), 3 + O(7^20), 4 + O(7^20))], 0)2544sage: (a*y).coeffs()2545([(0, 0, 0, 0, 0), (0, 1, 2, 3, 4)], 0)2546sage: (a*y^-2).coeffs()2547([(0, 1, 2, 3, 4), (0, 0, 0, 0, 0), (0, 0, 0, 0, 0)], -2)25482549Note that the coefficient list is transposed compared to how they2550are stored and printed::25512552sage: a*y^-22553(y^-2)*x + (2*y^-2)*x^2 + (3*y^-2)*x^3 + (4*y^-2)*x^425542555A more complicated example::25562557sage: a = x^20*y^-3 - x^11*y^2; a2558(y^-3-4*y^-1+6*y-4*y^3+y^5)*1 - (12*y^-3-36*y^-1+36*y+y^2-12*y^3-2*y^4+y^6)*x + (54*y^-3-108*y^-1+54*y+6*y^2-6*y^4)*x^2 - (108*y^-3-108*y^-1+9*y^2)*x^3 + (81*y^-3)*x^42559sage: raw, offset = a.coeffs()2560sage: a.min_pow_y()2561-32562sage: offset2563-32564sage: raw2565[(1, -12, 54, -108, 81),2566(0, 0, 0, 0, 0),2567(-4, 36, -108, 108, 0),2568(0, 0, 0, 0, 0),2569(6, -36, 54, 0, 0),2570(0, -1, 6, -9, 0),2571(-4, 12, 0, 0, 0),2572(0, 2, -6, 0, 0),2573(1, 0, 0, 0, 0),2574(0, -1, 0, 0, 0)]2575sage: sum(c * x^i * y^(j+offset) for j, L in enumerate(raw) for i, c in enumerate(L)) == a2576True25772578Can also be used to construct elements::25792580sage: a.parent()(raw, offset) == a2581True2582"""2583zero = self.base_ring()(0) if R is None else R(0)2584y_offset = min(self.min_pow_y(), 0)2585y_degree = max(self.max_pow_y(), 0)2586coeffs = []2587n = y_degree - y_offset + 12588for a in self._f.list():2589k = a.valuation()2590if k is Infinity: k = 02591k -= y_offset2592z = a.list()2593coeffs.append([zero] * k + z + [zero]*(n - len(z) - k))2594while len(coeffs) < self.parent().degree():2595coeffs.append([zero] * n)2596V = FreeModule(self.base_ring() if R is None else R, self.parent().degree())2597coeffs = transpose_list(coeffs)2598return [V(a) for a in coeffs], y_offset259926002601class MonskyWashnitzerDifferentialRing_class(Module):26022603def __init__(self, base_ring):2604r"""2605Class for the ring of Monsky--Washnitzer differentials over a given2606base ring.2607"""2608Module.__init__(self, base_ring)2609self._cache = {}26102611def invariant_differential(self):2612"""2613Returns `dx/2y` as an element of self.26142615EXAMPLES::26162617sage: R.<x> = QQ['x']2618sage: C = HyperellipticCurve(x^5-4*x+4)2619sage: MW = C.invariant_differential().parent()2620sage: MW.invariant_differential()26211 dx/2y2622"""2623return self(1)26242625def __call__(self, val, offset=0):2626return MonskyWashnitzerDifferential(self, val, offset)26272628def base_extend(self, R):2629"""2630Return a new differential ring which is self base-extended to `R`26312632INPUT:26332634- ``R`` -- ring26352636OUTPUT:26372638Self, base-extended to `R`.26392640EXAMPLES::26412642sage: R.<x> = QQ['x']2643sage: C = HyperellipticCurve(x^5-4*x+4)2644sage: MW = C.invariant_differential().parent()2645sage: MW.base_ring()2646SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 4*x + 4) over Rational Field2647sage: MW.base_extend(Qp(5,5)).base_ring()2648SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = (1 + O(5^5))*x^5 + (1 + 4*5 + 4*5^2 + 4*5^3 + 4*5^4 + O(5^5))*x + (4 + O(5^5))) over 5-adic Field with capped relative precision 52649"""2650return MonskyWashnitzerDifferentialRing(self.base_ring().base_extend(R))26512652def change_ring(self, R):2653"""2654Returns a new differential ring which is self with the coefficient2655ring changed to `R`.26562657INPUT:26582659- ``R`` -- ring of coefficients26602661OUTPUT:26622663Self, with the coefficient ring changed to `R`.26642665EXAMPLES::26662667sage: R.<x> = QQ['x']2668sage: C = HyperellipticCurve(x^5-4*x+4)2669sage: MW = C.invariant_differential().parent()2670sage: MW.base_ring()2671SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 4*x + 4) over Rational Field2672sage: MW.change_ring(Qp(5,5)).base_ring()2673SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = (1 + O(5^5))*x^5 + (1 + 4*5 + 4*5^2 + 4*5^3 + 4*5^4 + O(5^5))*x + (4 + O(5^5))) over 5-adic Field with capped relative precision 52674"""2675return MonskyWashnitzerDifferentialRing(self.base_ring().change_ring(R))26762677def degree(self):2678"""2679Returns the degree of `Q(x)`, where the model of the underlying2680hyperelliptic curve of self is given by `y^2 = Q(x)`.26812682EXAMPLES::26832684sage: R.<x> = QQ['x']2685sage: C = HyperellipticCurve(x^5-4*x+4)2686sage: MW = C.invariant_differential().parent()2687sage: MW.Q()2688x^5 - 4*x + 42689sage: MW.degree()269052691"""2692return self.base_ring().degree()26932694def dimension(self):2695"""2696Returns the dimension of self.26972698EXAMPLES::26992700sage: R.<x> = QQ['x']2701sage: C = HyperellipticCurve(x^5-4*x+4)2702sage: K = Qp(7,5)2703sage: CK = C.change_ring(K)2704sage: MW = CK.invariant_differential().parent()2705sage: MW.dimension()270642707"""2708return self.base_ring().degree()-127092710def Q(self):2711"""2712Returns `Q(x)` where the model of the underlying hyperelliptic curve2713of self is given by `y^2 = Q(x)`.27142715EXAMPLES::27162717sage: R.<x> = QQ['x']2718sage: C = HyperellipticCurve(x^5-4*x+4)2719sage: MW = C.invariant_differential().parent()2720sage: MW.Q()2721x^5 - 4*x + 42722"""2723return self.base_ring().Q()27242725def x_to_p(self, p):2726"""2727Returns and caches `x^p`, reduced via the relations coming from the2728defining polynomial of the hyperelliptic curve.27292730EXAMPLES::27312732sage: R.<x> = QQ['x']2733sage: C = HyperellipticCurve(x^5-4*x+4)2734sage: MW = C.invariant_differential().parent()2735sage: MW.x_to_p(3)2736x^32737sage: MW.x_to_p(5)2738-(4-y^2)*1 + 4*x2739sage: MW.x_to_p(101) is MW.x_to_p(101)2740True2741"""2742try:2743return self._cache["x_to_p", p]2744except KeyError:2745x_to_p = self.base_ring().x() ** p2746self._cache["x_to_p", p] = x_to_p2747return x_to_p27482749def frob_Q(self, p):2750"""2751Returns and caches `Q(x^p)`, which is used in computing the image of2752`y` under a `p`-power lift of Frobenius to `A^{\dagger}`.27532754EXAMPLES::27552756sage: R.<x> = QQ['x']2757sage: C = HyperellipticCurve(x^5-4*x+4)2758sage: MW = C.invariant_differential().parent()2759sage: MW.frob_Q(3)2760-(60-48*y^2+12*y^4-y^6)*1 + (192-96*y^2+12*y^4)*x - (192-48*y^2)*x^2 + 60*x^32761sage: MW.Q()(MW.x_to_p(3))2762-(60-48*y^2+12*y^4-y^6)*1 + (192-96*y^2+12*y^4)*x - (192-48*y^2)*x^2 + 60*x^32763sage: MW.frob_Q(11) is MW.frob_Q(11)2764True2765"""2766try:2767return self._cache["frobQ", p]2768except KeyError:2769x_to_p = self.x_to_p(p)2770frobQ = self.base_ring()._Q.change_ring(self.base_ring())(x_to_p)2771self._cache["frobQ", p] = frobQ2772return frobQ27732774def frob_invariant_differential(self, prec, p):2775r"""2776Kedlaya's algorithm allows us to calculate the action of Frobenius on2777the Monsky-Washnitzer cohomology. First we lift `\phi` to `A^{\dagger}`2778by setting27792780.. math::27812782\phi(x) = x^p27832784\phi(y) = y^p \sqrt{1 + \frac{Q(x^p) - Q(x)^p}{Q(x)^p}}.27852786Pulling back the differential `dx/2y`, we get27872788.. math::27892790\phi^*(dx/2y) = px^{p-1} y(\phi(y))^{-1} dx/2y2791= px^{p-1} y^{1-p} \sqrt{1+ \frac{Q(x^p) - Q(x)^p}{Q(x)^p}} dx/2y27922793Use Newton's method to calculate the square root.27942795EXAMPLES::27962797sage: R.<x> = QQ['x']2798sage: C = HyperellipticCurve(x^5-4*x+4)2799sage: prec = 22800sage: p = 72801sage: MW = C.invariant_differential().parent()2802sage: MW.frob_invariant_differential(prec,p)2803((67894400*y^-20-81198880*y^-18+40140800*y^-16-10035200*y^-14+1254400*y^-12-62720*y^-10)*1 - (119503944*y^-20-116064242*y^-18+43753472*y^-16-7426048*y^-14+514304*y^-12-12544*y^-10+1568*y^-8-70*y^-6-7*y^-4)*x + (78905288*y^-20-61014016*y^-18+16859136*y^-16-2207744*y^-14+250880*y^-12-37632*y^-10+3136*y^-8-70*y^-6)*x^2 - (39452448*y^-20-26148752*y^-18+8085490*y^-16-2007040*y^-14+376320*y^-12-37632*y^-10+1568*y^-8)*x^3 + (21102144*y^-20-18120592*y^-18+8028160*y^-16-2007040*y^-14+250880*y^-12-12544*y^-10)*x^4) dx/2y2804"""2805prof = Profiler()2806prof("setup")2807# TODO, would it be useful to be able to take Frobenius of any element? Less efficient?2808x, y = self.base_ring().gens()2809prof("x_to_p")2810x_to_p_less_1 = x**(p-1)2811x_to_p = x*x_to_p_less_128122813# cache for future use2814self._cache["x_to_p", p] = x_to_p28152816prof("frob_Q")2817a = self.frob_Q(p) >> 2*p # frobQ * y^{-2p}28182819prof("sqrt")28202821# Q = self.base_ring()._Q2822# three_halves = Q.parent().base_ring()(Rational((3,2)))2823# one_half = Q.parent().base_ring()(Rational((1,2)))2824three_halves = self.base_ring()._series_ring.base_ring()(Rational((3, 2)))2825one_half = self.base_ring()._series_ring.base_ring()(Rational((1, 2)))2826half_a = a._rmul_(one_half)28272828# We are solving for t = a^{-1/2} = (F_pQ y^{-p})^{-1/2}2829# Newton's method converges because we know the root is in the same residue class as 1.28302831# t = self.base_ring()(1)2832t = self.base_ring()(three_halves) - half_a2833# first iteration trivial, start with prec 228342835for cur_prec in newton_method_sizes(prec)[2:]:2836# newton_method_sizes = [1, 2, ...]2837y_prec = -(2*cur_prec-1)*p+12838# binomial expansion is $\sum p^{k+1} y^{-(2k+1)p+1} f(x)$2839# so if we are only correct mod p^prec,2840# can ignore y powers less than y_prec2841t_cube = (t*t*t).truncate_neg(y_prec)2842t = t._rmul_(three_halves) - (half_a * t_cube).truncate_neg(y_prec)2843# t = (3/2) t - (1/2) a t^32844# print "a =", a2845# print "t =", t2846# prof("verify")2847# print "a*t^2 =", a * t**228482849prof("compose")2850F_dx_y = (p * x_to_p_less_1 * t) >> (p-1) # px^{p-1} sqrt(a) * y^{-p+1}28512852# print "-----", F_dx_y2853# print "-----", x_to_p * F_dx_y2854prof("done")2855# print prof2856return MonskyWashnitzerDifferential(self, F_dx_y)28572858def frob_basis_elements(self, prec, p):2859"""2860Returns the action of a `p`-power lift of Frobenius on the basis28612862.. math::28632864\{ dx/2y, x dx/2y, ..., x^{d-2} dx/2y \}28652866where `d` is the degree of the underlying hyperelliptic curve.28672868EXAMPLES::28692870sage: R.<x> = QQ['x']2871sage: C = HyperellipticCurve(x^5-4*x+4)2872sage: prec = 12873sage: p = 52874sage: MW = C.invariant_differential().parent()2875sage: MW.frob_basis_elements(prec,p)2876[((92000*y^-14-74200*y^-12+32000*y^-10-8000*y^-8+1000*y^-6-50*y^-4)*1 - (194400*y^-14-153600*y^-12+57600*y^-10-9600*y^-8+600*y^-6)*x + (204800*y^-14-153600*y^-12+38400*y^-10-3200*y^-8)*x^2 - (153600*y^-14-76800*y^-12+9600*y^-10)*x^3 + (63950*y^-14-18550*y^-12+1600*y^-10-400*y^-8+50*y^-6+5*y^-4)*x^4) dx/2y, (-(1391200*y^-14-941400*y^-12+302000*y^-10-76800*y^-8+14400*y^-6-1320*y^-4+30*y^-2)*1 + (2168800*y^-14-1402400*y^-12+537600*y^-10-134400*y^-8+16800*y^-6-720*y^-4)*x - (1596800*y^-14-1433600*y^-12+537600*y^-10-89600*y^-8+5600*y^-6)*x^2 + (1433600*y^-14-1075200*y^-12+268800*y^-10-22400*y^-8)*x^3 - (870200*y^-14-445350*y^-12+63350*y^-10-3200*y^-8+600*y^-6-30*y^-4-5*y^-2)*x^4) dx/2y, ((19488000*y^-14-15763200*y^-12+4944400*y^-10-913800*y^-8+156800*y^-6-22560*y^-4+1480*y^-2-10)*1 - (28163200*y^-14-18669600*y^-12+5774400*y^-10-1433600*y^-8+268800*y^-6-25440*y^-4+760*y^-2)*x + (15062400*y^-14-12940800*y^-12+5734400*y^-10-1433600*y^-8+179200*y^-6-8480*y^-4)*x^2 - (12121600*y^-14-11468800*y^-12+4300800*y^-10-716800*y^-8+44800*y^-6)*x^3 + (9215200*y^-14-6952400*y^-12+1773950*y^-10-165750*y^-8+5600*y^-6-720*y^-4+10*y^-2+5)*x^4) dx/2y, (-(225395200*y^-14-230640000*y^-12+91733600*y^-10-18347400*y^-8+2293600*y^-6-280960*y^-4+31520*y^-2-1480-10*y^2)*1 + (338048000*y^-14-277132800*y^-12+89928000*y^-10-17816000*y^-8+3225600*y^-6-472320*y^-4+34560*y^-2-720)*x - (172902400*y^-14-141504000*y^-12+58976000*y^-10-17203200*y^-8+3225600*y^-6-314880*y^-4+11520*y^-2)*x^2 + (108736000*y^-14-109760000*y^-12+51609600*y^-10-12902400*y^-8+1612800*y^-6-78720*y^-4)*x^3 - (85347200*y^-14-82900000*y^-12+31251400*y^-10-5304150*y^-8+367350*y^-6-8480*y^-4+760*y^-2+10-5*y^2)*x^4) dx/2y]2877"""2878F_i = self.frob_invariant_differential(prec, p)2879x_to_p = self.x_to_p(p)2880F = [F_i]2881for i in range(1, self.degree()-1):2882F_i *= x_to_p2883F.append(F_i)2884return F28852886def helper_matrix(self):2887"""2888We use this to solve for the linear combination of2889`x^i y^j` needed to clear all terms with2890`y^{j-1}`.28912892EXAMPLES::28932894sage: R.<x> = QQ['x']2895sage: C = HyperellipticCurve(x^5-4*x+4)2896sage: MW = C.invariant_differential().parent()2897sage: MW.helper_matrix()2898[ 256/2101 320/2101 400/2101 500/2101 625/2101]2899[-625/8404 -64/2101 -80/2101 -100/2101 -125/2101]2900[-125/2101 -625/8404 -64/2101 -80/2101 -100/2101]2901[-100/2101 -125/2101 -625/8404 -64/2101 -80/2101]2902[ -80/2101 -100/2101 -125/2101 -625/8404 -64/2101]2903"""2904try:2905return self._helper_matrix2906except AttributeError:2907pass29082909# The smallest y term of (1/j) d(x^i y^j) is constant for all j.2910L = []2911x, y = self.base_ring().gens()2912n = self.degree()2913for i in range(n):2914L.append((y*x**i).diff().extract_pow_y(0))2915A = matrix(L).transpose()2916if not is_IntegralDomain(A.base_ring()):2917# must be using integer_mod or something to approximate2918self._helper_matrix = (~A.change_ring(QQ)).change_ring(A.base_ring())2919else:2920self._helper_matrix = ~A2921return self._helper_matrix292229232924class MonskyWashnitzerDifferential(ModuleElement):29252926def __init__(self, parent, val=0, offset=0):2927r"""2928Create an element of the Monsky-Washnitzer ring of differentials, of2929the form `F dx/2y`.29302931INPUT:29322933- ``parent`` -- Monsky-Washnitzer differential ring (instance of class2934:class:`~MonskyWashnitzerDifferentialRing_class`29352936- ``val`` -- element of the base ring, or list of coefficients29372938- ``offset`` -- if non-zero, shift val by `y^\text{offset}` (default 0)29392940EXAMPLES::29412942sage: R.<x> = QQ['x']2943sage: C = HyperellipticCurve(x^5 - 4*x + 4)2944sage: x,y = C.monsky_washnitzer_gens()2945sage: MW = C.invariant_differential().parent()2946sage: sage.schemes.elliptic_curves.monsky_washnitzer.MonskyWashnitzerDifferential(MW, x)2947x dx/2y2948sage: sage.schemes.elliptic_curves.monsky_washnitzer.MonskyWashnitzerDifferential(MW, y)2949y*1 dx/2y2950sage: sage.schemes.elliptic_curves.monsky_washnitzer.MonskyWashnitzerDifferential(MW, x, 10)2951y^10*x dx/2y2952"""2953ModuleElement.__init__(self, parent)2954if isinstance(val, MonskyWashnitzerDifferential):2955val = val._coeff2956self._coeff = self.parent().base_ring()(val, offset)29572958def _add_(left, right):2959"""2960Returns the sum of left and right, both elements of the2961Monsky-Washnitzer ring of differentials.29622963EXAMPLES::29642965sage: R.<x> = QQ['x']2966sage: C = HyperellipticCurve(x^5-4*x + 4)2967sage: x,y = C.monsky_washnitzer_gens()2968sage: w = C.invariant_differential()2969sage: w + w29702*1 dx/2y2971sage: x*w + w2972(1 + x) dx/2y2973sage: x*w + y*w2974(y*1 + x) dx/2y2975"""2976return MonskyWashnitzerDifferential(left.parent(),2977left._coeff + right._coeff)29782979def _sub_(left, right):2980"""2981Returns the difference of left and right, both elements of the2982Monsky-Washnitzer ring of differentials.29832984EXAMPLES::29852986sage: R.<x> = QQ['x']2987sage: C = HyperellipticCurve(x^5-4*x+4)2988sage: x,y = C.monsky_washnitzer_gens()2989sage: w = C.invariant_differential()2990sage: w-w29910 dx/2y2992sage: x*w-w2993(-1 + x) dx/2y2994sage: w - x*w - y*w2995((1-y)*1 - x) dx/2y2996"""2997return MonskyWashnitzerDifferential(left.parent(),2998left._coeff - right._coeff)29993000def __neg__(self):3001"""3002Returns the additive inverse of self.30033004EXAMPLES::30053006sage: R.<x> = QQ['x']3007sage: C = HyperellipticCurve(x^5-4*x+4)3008sage: x,y = C.monsky_washnitzer_gens()3009sage: w = C.invariant_differential()3010sage: -w3011-1 dx/2y3012sage: -((y-x)*w)3013(-y*1 + x) dx/2y3014"""3015return MonskyWashnitzerDifferential(self.parent(), -self._coeff)30163017def _lmul_(self, a):3018"""3019Returns `self * a`.30203021EXAMPLES::30223023sage: R.<x> = QQ['x']3024sage: C = HyperellipticCurve(x^5-4*x+4)3025sage: x,y = C.monsky_washnitzer_gens()3026sage: w = C.invariant_differential()3027sage: w*x3028x dx/2y3029sage: (w*x)*230302*x dx/2y3031sage: w*y3032y*1 dx/2y3033sage: w*(x+y)3034(y*1 + x) dx/2y3035"""3036return MonskyWashnitzerDifferential(self.parent(), self._coeff * a)30373038def _rmul_(self, a):3039"""3040Returns `a * self`.30413042EXAMPLES::30433044sage: R.<x> = QQ['x']3045sage: C = HyperellipticCurve(x^5-4*x+4)3046sage: x,y = C.monsky_washnitzer_gens()3047sage: w = C.invariant_differential()3048sage: x*w3049x dx/2y3050sage: 2*(x*w)30512*x dx/2y3052sage: y*w3053y*1 dx/2y3054sage: (x+y)*w3055(y*1 + x) dx/2y3056"""3057return MonskyWashnitzerDifferential(self.parent(), a * self._coeff)30583059def coeff(self):3060r"""3061Returns `A`, where this element is `A dx/2y`.30623063EXAMPLES::30643065sage: R.<x> = QQ['x']3066sage: C = HyperellipticCurve(x^5-4*x+4)3067sage: x,y = C.monsky_washnitzer_gens()3068sage: w = C.invariant_differential()3069sage: w30701 dx/2y3071sage: w.coeff()307213073sage: (x*y*w).coeff()3074y*x3075"""3076return self._coeff30773078def __nonzero__(self):3079"""3080EXAMPLES::30813082sage: R.<x> = QQ['x']3083sage: C = HyperellipticCurve(x^5-4*x+4)3084sage: x,y = C.monsky_washnitzer_gens()3085sage: w = C.invariant_differential()3086sage: not w3087False3088sage: not 0*w3089True3090sage: not x*y*w3091False3092"""3093return not not self._coeff30943095def _repr_(self):3096"""3097EXAMPLES::30983099sage: R.<x> = QQ['x']3100sage: C = HyperellipticCurve(x^5-4*x+4)3101sage: x,y = C.monsky_washnitzer_gens()3102sage: w = C.invariant_differential()3103sage: w31041 dx/2y3105sage: (2*x+y)*w3106(y*1 + 2*x) dx/2y3107"""3108s = self._coeff._repr_()3109if s.find("+") != -1 or s.find("-") > 0:3110s = "(%s)" % s3111return s + " dx/2y"31123113def _latex_(self):3114"""3115Returns the latex representation of self.31163117EXAMPLES::31183119sage: R.<x> = QQ['x']3120sage: C = HyperellipticCurve(x^5-4*x+4)3121sage: x,y = C.monsky_washnitzer_gens()3122sage: w = C.invariant_differential()3123sage: latex(w)31241 \frac{dx}{2y}3125sage: latex(x*w)3126x \frac{dx}{2y}3127"""3128s = self._coeff._latex_()3129if s.find("+") != -1 or s.find("-") > 0:3130s = "\\left(%s\\right)" % s3131return s + " \\frac{dx}{2y}"31323133def extract_pow_y(self, k):3134"""3135Returns the power of `y` in `A` where self is `A dx/2y`.31363137EXAMPLES::31383139sage: R.<x> = QQ['x']3140sage: C = HyperellipticCurve(x^5-3*x+1)3141sage: x,y = C.monsky_washnitzer_gens()3142sage: A = y^5 - x*y^33143sage: A.extract_pow_y(5)3144[1, 0, 0, 0, 0]3145sage: (A * C.invariant_differential()).extract_pow_y(5)3146[1, 0, 0, 0, 0]3147"""3148return self._coeff.extract_pow_y(k)31493150def min_pow_y(self):3151"""3152Returns the minimum power of `y` in `A` where self is `A dx/2y`.31533154EXAMPLES::31553156sage: R.<x> = QQ['x']3157sage: C = HyperellipticCurve(x^5-3*x+1)3158sage: x,y = C.monsky_washnitzer_gens()3159sage: w = y^5 * C.invariant_differential()3160sage: w.min_pow_y()316153162sage: w = (x^2*y^4 + y^5) * C.invariant_differential()3163sage: w.min_pow_y()316443165"""3166return self._coeff.min_pow_y()31673168def max_pow_y(self):3169"""3170Returns the maximum power of `y` in `A` where self is `A dx/2y`.31713172EXAMPLES::31733174sage: R.<x> = QQ['x']3175sage: C = HyperellipticCurve(x^5-3*x+1)3176sage: x,y = C.monsky_washnitzer_gens()3177sage: w = y^5 * C.invariant_differential()3178sage: w.max_pow_y()317953180sage: w = (x^2*y^4 + y^5) * C.invariant_differential()3181sage: w.max_pow_y()318253183"""3184return self._coeff.max_pow_y()31853186def reduce_neg_y(self):3187"""3188Use homology relations to eliminate negative powers of `y`.31893190EXAMPLES::31913192sage: R.<x> = QQ['x']3193sage: C = HyperellipticCurve(x^5-3*x+1)3194sage: x,y = C.monsky_washnitzer_gens()3195sage: (y^-1).diff().reduce_neg_y()3196((y^-1)*1, 0 dx/2y)3197sage: (y^-5*x^2+y^-1*x).diff().reduce_neg_y()3198((y^-1)*x + (y^-5)*x^2, 0 dx/2y)3199"""3200S = self.parent().base_ring()3201R = S.base_ring()3202M = self.parent().helper_matrix()3203p = S._p3204n = S.degree()3205x, y = S.gens()3206f = S(0)3207reduced = self3208for j in range(self.min_pow_y()+1, 0):3209if p is not None and p.divides(j):3210cs = [a/j for a in reduced.extract_pow_y(j-1)]3211else:3212j_inverse = ~R(j)3213cs = [a*j_inverse for a in reduced.extract_pow_y(j-1)]3214lin_comb = M * vector(M.base_ring(), cs)3215# print "j =", j, "b =", cs, "lin_comb =", lin_comb3216g = self.parent().base_ring()(0)3217if not lin_comb.is_zero():3218for i in range(n):3219if lin_comb[i] != 0:3220g += S.monomial(i, j, lin_comb[i])3221if not g.is_zero():3222f += g3223reduced -= g.diff()3224# print g, g.diff()3225# print "reduced", reduced32263227return f, reduced32283229def reduce_neg_y_fast(self, even_degree_only=False):3230"""3231Use homology relations to eliminate negative powers of `y`.32323233EXAMPLES::32343235sage: R.<x> = QQ['x']3236sage: E = HyperellipticCurve(x^5-3*x+1)3237sage: x, y = E.monsky_washnitzer_gens()3238sage: (y^-1).diff().reduce_neg_y_fast()3239((y^-1)*1, 0 dx/2y)3240sage: (y^-5*x^2+y^-1*x).diff().reduce_neg_y_fast()3241((y^-1)*x + (y^-5)*x^2, 0 dx/2y)32423243It leaves non-negative powers of `y` alone::32443245sage: y.diff()3246(-3*1 + 5*x^4) dx/2y3247sage: y.diff().reduce_neg_y_fast()3248(0, (-3*1 + 5*x^4) dx/2y)3249"""3250# prof = Profiler()3251# prof("reduce setup")3252S = self.parent().base_ring()3253R = S.base_ring()3254M = self.parent().helper_matrix()32553256# prof("extract coeffs")3257coeffs, offset = self.coeffs(R)3258V = coeffs[0].parent()32593260if offset == 0:3261return S(0), self32623263# prof("loop %s"%self.min_pow_y())3264forms = []3265p = S._p3266for j in range(self.min_pow_y()+1, 0):3267if (even_degree_only and j % 2 == 0) or coeffs[j-offset-1].is_zero():3268forms.append(V(0))3269else:3270# this is a total hack to deal with the fact that we're using3271# rational numbers to approximate fixed precision p-adics3272if p is not None and j % 3 == 1:3273try:3274v = coeffs[j-offset-1]3275for kk in range(len(v)):3276a = v[kk]3277ppow = p**max(-a.valuation(S._p), 0)3278v[kk] = ((a * ppow) % S._prec_cap) / ppow3279except AttributeError:3280pass3281lin_comb = ~R(j) * (M * coeffs[j-offset-1])3282forms.append(lin_comb)3283for i in lin_comb.nonzero_positions():3284# g = lin_comb[i] x^i y^j3285# self -= dg3286coeffs[j-offset+1] -= lin_comb[i] * S.monomial_diff_coeffs(i, j)[1]32873288# prof("recreate forms")3289f = S(forms, offset+1)3290reduced = S._monsky_washnitzer(coeffs[-1-offset:], -1)3291# print self - f.diff() - reduced3292# prof("done")3293# print prof3294return f, reduced32953296def reduce_neg_y_faster(self, even_degree_only=False):3297"""3298Use homology relations to eliminate negative powers of `y`.32993300EXAMPLES::33013302sage: R.<x> = QQ['x']3303sage: C = HyperellipticCurve(x^5-3*x+1)3304sage: x,y = C.monsky_washnitzer_gens()3305sage: (y^-1).diff().reduce_neg_y()3306((y^-1)*1, 0 dx/2y)3307sage: (y^-5*x^2+y^-1*x).diff().reduce_neg_y_faster()3308((y^-1)*x + (y^-5)*x^2, 0 dx/2y)3309"""3310# Timings indicate that this is not any faster after all...33113312S = self.parent().base_ring()3313R = S.base_ring()3314M = self.parent().helper_matrix()33153316coeffs, offset = self.coeffs(R)3317V = coeffs[0].parent()3318zeroV = V(0)33193320if offset == 0:3321return S(0), self33223323# See monomial_diff_coeffs3324# this is the B_i and x_to_i contributions respectively for all i3325d_mat_1, d_mat_2 = S.monomial_diff_coeffs_matrices()33263327forms = []3328for j in range(self.min_pow_y()+1, 0):3329if coeffs[j-offset-1].is_zero():3330forms.append(zeroV)3331else:3332# this is a total hack to deal with the fact that we're using3333# rational numbers to approximate fixed precision p-adics3334if j % 3 == 0:3335try:3336v = coeffs[j-offset-1]3337for kk in range(len(v)):3338a = v[kk]3339ppow = S._p**max(-a.valuation(S._p), 0)3340v[kk] = ((a * ppow) % S._prec_cap) / ppow3341except AttributeError:3342pass3343j_inverse = ~R(j)3344lin_comb = (M * coeffs[j-offset-1])3345forms.append(j_inverse * lin_comb)3346coeffs[j-offset+1] -= (d_mat_1 + j_inverse * d_mat_2) * lin_comb33473348f = S(forms, offset+1)3349reduced = S._monsky_washnitzer(coeffs[-1-offset:], -1)3350# reduced = self - f.diff()3351return f, reduced33523353def reduce_pos_y(self):3354"""3355Use homology relations to eliminate positive powers of `y`.33563357EXAMPLES::33583359sage: R.<x> = QQ['x']3360sage: C = HyperellipticCurve(x^3-4*x+4)3361sage: x,y = C.monsky_washnitzer_gens()3362sage: (y^2).diff().reduce_pos_y()3363(y^2*1, 0 dx/2y)3364sage: (y^2*x).diff().reduce_pos_y()3365(y^2*x, 0 dx/2y)3366sage: (y^92*x).diff().reduce_pos_y()3367(y^92*x, 0 dx/2y)3368sage: w = (y^3 + x).diff()3369sage: w += w.parent()(x)3370sage: w.reduce_pos_y_fast()3371(y^3*1 + x, x dx/2y)3372"""3373S = self.parent().base_ring()3374n = S.Q().degree()3375x, y = S.gens()3376f = S(0)3377reduced = self3378for j in range(self.max_pow_y(), 0, -1):3379for i in range(n-1, -1, -1):3380c = reduced.extract_pow_y(j)[i]3381# print "x^%s y^%s"%(i,j), c3382if c != 0:3383g = S.monomial(0, j+1) if i == n-1 else S.monomial(i+1, j-1)3384dg = g.diff()3385# print reduced, " - ", dg3386denom = dg.extract_pow_y(j)[i]3387c /= denom3388c = g.parent()(c)3389f += c * g3390reduced -= c * dg33913392return f, reduced33933394def reduce_pos_y_fast(self, even_degree_only=False):3395"""3396Use homology relations to eliminate positive powers of `y`.33973398EXAMPLES::33993400sage: R.<x> = QQ['x']3401sage: E = HyperellipticCurve(x^3-4*x+4)3402sage: x, y = E.monsky_washnitzer_gens()3403sage: y.diff().reduce_pos_y_fast()3404(y*1, 0 dx/2y)3405sage: (y^2).diff().reduce_pos_y_fast()3406(y^2*1, 0 dx/2y)3407sage: (y^2*x).diff().reduce_pos_y_fast()3408(y^2*x, 0 dx/2y)3409sage: (y^92*x).diff().reduce_pos_y_fast()3410(y^92*x, 0 dx/2y)3411sage: w = (y^3 + x).diff()3412sage: w += w.parent()(x)3413sage: w.reduce_pos_y_fast()3414(y^3*1 + x, x dx/2y)3415"""3416S = self.parent().base_ring()3417R = S.base_ring()3418n = S.Q().degree()34193420coeffs, offset = self.coeffs(R)3421V = coeffs[0].parent()3422zeroV = V(0)3423forms = [V(0), V(0)]34243425for j in range(self.max_pow_y(), -1, -1):34263427if (even_degree_only and j % 2 == 1) or (j > 0 and coeffs[j-offset].is_zero()):3428forms.append(zeroV)3429continue34303431form = V(0)3432i = n-13433c = coeffs[j-offset][i]3434if c != 0:3435dg_coeffs = S.monomial_diff_coeffs(0, j+1)[0]3436c /= dg_coeffs[i]3437forms[len(forms)-2][0] = c3438# self -= c d(y^{j+1})3439coeffs[j-offset] -= c*dg_coeffs34403441if j == 0:3442# the others are basis elements3443break34443445for i in range(n-2, -1, -1):3446c = coeffs[j-offset][i]3447if c != 0:3448dg_coeffs = S.monomial_diff_coeffs(i+1, j-1)3449denom = dg_coeffs[1][i]3450c /= denom3451form[i+1] = c3452# self -= c d(x^{i+1} y^{j-1})3453coeffs[j-offset] -= c*dg_coeffs[1]3454coeffs[j-offset-2] -= c*dg_coeffs[0]3455forms.append(form)34563457forms.reverse()3458f = S(forms)3459reduced = self.parent()(coeffs[:1-offset], offset)3460return f, reduced34613462def reduce(self):3463"""3464Use homology relations to find `a` and `f` such that this element is3465equal to `a + df`, where `a` is given in terms of the `x^i dx/2y`.34663467EXAMPLES::34683469sage: R.<x> = QQ['x']3470sage: C = HyperellipticCurve(x^5-4*x+4)3471sage: x,y = C.monsky_washnitzer_gens()3472sage: w = (y*x).diff()3473sage: w.reduce()3474(y*x, 0 dx/2y)34753476sage: w = x^4 * C.invariant_differential()3477sage: w.reduce()3478(1/5*y*1, 4/5*1 dx/2y)34793480sage: w = sum(QQ.random_element() * x^i * y^j for i in [0..4] for j in [-3..3]) * C.invariant_differential()3481sage: f, a = w.reduce()3482sage: f.diff() + a - w34830 dx/2y3484"""3485# print "max_pow_y = ", self.max_pow_y(), "min_pow_y = ", self.min_pow_y()3486n = self.parent().base_ring().Q().degree()3487f1, a = self.reduce_neg_y()3488f2, a = a.reduce_pos_y()3489f = f1 + f234903491c = a.extract_pow_y(0)[n-1]3492if c != 0:3493x, y = self.parent().base_ring().gens()3494g = y3495dg = g.diff()3496c = g.parent()(c/dg.extract_pow_y(0)[n-1])3497f += c * g3498a -= c * dg3499# print g, dg35003501return f, a35023503def reduce_fast(self, even_degree_only=False):3504"""3505Use homology relations to find `a` and `f` such that this element is3506equal to `a + df`, where `a` is given in terms of the `x^i dx/2y`.35073508EXAMPLES::35093510sage: R.<x> = QQ['x']3511sage: E = HyperellipticCurve(x^3-4*x+4)3512sage: x, y = E.monsky_washnitzer_gens()3513sage: x.diff().reduce_fast()3514(x, (0, 0))3515sage: y.diff().reduce_fast()3516(y*1, (0, 0))3517sage: (y^-1).diff().reduce_fast()3518((y^-1)*1, (0, 0))3519sage: (y^-11).diff().reduce_fast()3520((y^-11)*1, (0, 0))3521sage: (x*y^2).diff().reduce_fast()3522(y^2*x, (0, 0))3523"""3524# print "max_pow_y = ", self.max_pow_y(), "min_pow_y = ", self.min_pow_y()35253526f1, reduced = self.reduce_neg_y_fast(even_degree_only)3527f2, reduced = reduced.reduce_pos_y_fast(even_degree_only)3528# f1, reduced = self.reduce_neg_y()3529# f2, reduced = reduced.reduce_pos_y()3530v = reduced.extract_pow_y(0)3531v.pop()3532V = FreeModule(self.base_ring().base_ring(), len(v))3533return f1+f2, V(v)35343535def coeffs(self, R=None):3536"""3537Used to obtain the raw coefficients of a differential, see3538:meth:`SpecialHyperellipticQuotientElement.coeffs`35393540INPUT:35413542- R -- An (optional) base ring in which to cast the coefficients35433544OUTPUT:35453546The raw coefficients of $A$ where self is $A dx/2y$.35473548EXAMPLES::35493550sage: R.<x> = QQ['x']3551sage: C = HyperellipticCurve(x^5-4*x+4)3552sage: x,y = C.monsky_washnitzer_gens()3553sage: w = C.invariant_differential()3554sage: w.coeffs()3555([(1, 0, 0, 0, 0)], 0)3556sage: (x*w).coeffs()3557([(0, 1, 0, 0, 0)], 0)3558sage: (y*w).coeffs()3559([(0, 0, 0, 0, 0), (1, 0, 0, 0, 0)], 0)3560sage: (y^-2*w).coeffs()3561([(1, 0, 0, 0, 0), (0, 0, 0, 0, 0), (0, 0, 0, 0, 0)], -2)3562"""3563return self._coeff.coeffs(R)35643565def coleman_integral(self, P, Q):3566r"""3567Computes the definite integral of self from $P$ to $Q$.35683569INPUT:35703571- P, Q -- Two points on the underlying curve35723573OUTPUT:35743575`\int_P^Q \text{self}`35763577EXAMPLES::35783579sage: K = pAdicField(5,7)3580sage: E = EllipticCurve(K,[-31/3,-2501/108]) #11a3581sage: P = E(K(14/3), K(11/2))3582sage: w = E.invariant_differential()3583sage: w.coleman_integral(P,2*P)3584O(5^6)35853586sage: Q = E.lift_x(3)3587sage: w.coleman_integral(P,Q)35882*5 + 4*5^2 + 3*5^3 + 4*5^4 + 3*5^5 + O(5^6)3589sage: w.coleman_integral(2*P,Q)35902*5 + 4*5^2 + 3*5^3 + 4*5^4 + 3*5^5 + O(5^6)3591sage: (2*w).coleman_integral(P, Q) == 2*(w.coleman_integral(P, Q))3592True3593"""3594return self.parent().base_ring().curve().coleman_integral(self, P, Q)35953596integrate = coleman_integral359735983599