Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/schemes/elliptic_curves/monsky_washnitzer.py
8820 views
1
r"""
2
Computation of Frobenius matrix on Monsky-Washnitzer cohomology
3
4
The most interesting functions to be exported here are
5
:func:`matrix_of_frobenius` and :func:`adjusted_prec`.
6
7
Currently this code is limited to the case `p \geq 5` (no
8
`GF(p^n)` for `n > 1`), and only handles the
9
elliptic curve case (not more general hyperelliptic curves).
10
11
REFERENCES:
12
13
.. [Ked2001] Kedlaya, K., "Counting points on hyperelliptic curves using
14
Monsky-Washnitzer cohomology", J. Ramanujan Math. Soc. 16 (2001) no
15
4, 323-338
16
17
.. [Edix] Edixhoven, B., "Point counting after Kedlaya", EIDMA-Stieltjes
18
graduate course, Lieden (lecture notes?).
19
20
AUTHORS:
21
22
- David Harvey and Robert Bradshaw: initial code developed at the 2006
23
MSRI graduate workshop, working with Jennifer Balakrishnan and Liang
24
Xiao
25
26
- David Harvey (2006-08): cleaned up, rewrote some chunks, lots more
27
documentation, added Newton iteration method, added more complete
28
'trace trick', integrated better into Sage.
29
30
- David Harvey (2007-02): added algorithm with sqrt(p) complexity
31
(removed in May 2007 due to better C++ implementation)
32
33
- Robert Bradshaw (2007-03): keep track of exact form in reduction
34
algorithms
35
36
- Robert Bradshaw (2007-04): generalization to hyperelliptic curves
37
"""
38
39
#*****************************************************************************
40
# Copyright (C) 2006 William Stein <[email protected]>
41
# 2006 Robert Bradshaw <[email protected]>
42
# 2006 David Harvey <[email protected]>
43
#
44
# Distributed under the terms of the GNU General Public License (GPL)
45
# http://www.gnu.org/licenses/
46
#*****************************************************************************
47
48
49
from sage.rings.all import Integers, Integer, PolynomialRing, PowerSeriesRing, Rationals, Rational, LaurentSeriesRing
50
51
from sage.rings.polynomial.polynomial_element import is_Polynomial
52
53
from sage.modules.module import Module
54
from sage.structure.element import ModuleElement
55
from sage.matrix.all import matrix
56
from sage.modules.all import vector
57
from sage.rings.ring import CommutativeAlgebra
58
from sage.structure.element import CommutativeAlgebraElement
59
from sage.rings.infinity import Infinity
60
61
from sage.rings.arith import binomial, integer_ceil as ceil
62
from sage.misc.functional import log
63
from sage.misc.misc import newton_method_sizes
64
65
from ell_generic import is_EllipticCurve
66
from constructor import EllipticCurve
67
68
69
class SpecialCubicQuotientRing(CommutativeAlgebra):
70
r"""
71
Specialised class for representing the quotient ring
72
`R[x,T]/(T - x^3 - ax - b)`, where `R` is an
73
arbitrary commutative base ring (in which 2 and 3 are invertible),
74
`a` and `b` are elements of that ring.
75
76
Polynomials are represented internally in the form
77
`p_0 + p_1 x + p_2 x^2` where the `p_i` are
78
polynomials in `T`. Multiplication of polynomials always
79
reduces high powers of `x` (i.e. beyond `x^2`) to
80
powers of `T`.
81
82
Hopefully this ring is faster than a general quotient ring because
83
it uses the special structure of this ring to speed multiplication
84
(which is the dominant operation in the frobenius matrix
85
calculation). I haven't actually tested this theory though...
86
87
.. TODO::
88
89
Eventually we will want to run this in characteristic 3, so we
90
need to: (a) Allow `Q(x)` to contain an `x^2` term, and (b) Remove
91
the requirement that 3 be invertible. Currently this is used in
92
the Toom-Cook algorithm to speed multiplication.
93
94
EXAMPLES::
95
96
sage: B.<t> = PolynomialRing(Integers(125))
97
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
98
sage: R
99
SpecialCubicQuotientRing over Ring of integers modulo 125 with polynomial T = x^3 + 124*x + 94
100
101
Get generators::
102
103
sage: x, T = R.gens()
104
sage: x
105
(0) + (1)*x + (0)*x^2
106
sage: T
107
(T) + (0)*x + (0)*x^2
108
109
Coercions::
110
111
sage: R(7)
112
(7) + (0)*x + (0)*x^2
113
114
Create elements directly from polynomials::
115
116
sage: A, z = R.poly_ring().objgen()
117
sage: A
118
Univariate Polynomial Ring in T over Ring of integers modulo 125
119
sage: R.create_element(z^2, z+1, 3)
120
(T^2) + (T + 1)*x + (3)*x^2
121
122
Some arithmetic::
123
124
sage: x^3
125
(T + 31) + (1)*x + (0)*x^2
126
sage: 3 * x**15 * T**2 + x - T
127
(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^2
128
129
Retrieve coefficients (output is zero-padded)::
130
131
sage: x^10
132
(3*T^2 + 61*T + 8) + (T^3 + 93*T^2 + 12*T + 40)*x + (3*T^2 + 61*T + 9)*x^2
133
sage: (x^10).coeffs()
134
[[8, 61, 3, 0], [40, 12, 93, 1], [9, 61, 3, 0]]
135
136
.. TODO::
137
138
write an example checking multiplication of these polynomials
139
against Sage's ordinary quotient ring arithmetic. I can't seem
140
to get the quotient ring stuff happening right now...
141
"""
142
def __init__(self, Q, laurent_series=False):
143
"""
144
Constructor.
145
146
INPUT:
147
148
- ``Q`` -- a polynomial of the form
149
`Q(x) = x^3 + ax + b`, where `a`, `b` belong to a ring in which
150
2, 3 are invertible.
151
152
- ``laurent_series`` -- whether or not to allow
153
negative powers of `T` (default=False)
154
155
EXAMPLES::
156
157
sage: B.<t> = PolynomialRing(Integers(125))
158
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
159
sage: R
160
SpecialCubicQuotientRing over Ring of integers modulo 125 with polynomial T = x^3 + 124*x + 94
161
162
::
163
164
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 + 2*t^2 - t + B(1/4))
165
Traceback (most recent call last):
166
...
167
ValueError: Q (=t^3 + 2*t^2 + 124*t + 94) must be of the form x^3 + ax + b
168
169
::
170
171
sage: B.<t> = PolynomialRing(Integers(10))
172
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + 1)
173
Traceback (most recent call last):
174
...
175
ArithmeticError: 2 and 3 must be invertible in the coefficient ring (=Ring of integers modulo 10) of Q
176
"""
177
if not is_Polynomial(Q):
178
raise TypeError("Q (=%s) must be a polynomial" % Q)
179
180
if Q.degree() != 3 or not Q[2].is_zero():
181
raise ValueError("Q (=%s) must be of the form x^3 + ax + b" % Q)
182
183
base_ring = Q.parent().base_ring()
184
185
if not base_ring(6).is_unit():
186
raise ArithmeticError("2 and 3 must be invertible in the "
187
"coefficient ring (=%s) of Q" % base_ring)
188
189
# CommutativeAlgebra.__init__ tries to establish a coercion
190
# from the base ring, by trac ticket #9138. The corresponding
191
# hom set is cached. In order to use self as cache key, its
192
# string representation is used. In otder to get the string
193
# representation, we need to know the attributes _a and
194
# _b. Hence, in #9138, we have to move CommutativeAlgebra.__init__
195
# further down:
196
self._a = Q[1]
197
self._b = Q[0]
198
if laurent_series:
199
self._poly_ring = LaurentSeriesRing(base_ring, 'T') # R[T]
200
else:
201
self._poly_ring = PolynomialRing(base_ring, 'T') # R[T]
202
self._poly_generator = self._poly_ring.gen(0) # the generator T
203
CommutativeAlgebra.__init__(self, base_ring)
204
205
# Precompute a matrix that is used in the Toom-Cook multiplication.
206
# This is where we need 2 and 3 invertible.
207
208
# (a good description of Toom-Cook is online at:
209
# http://www.gnu.org/software/gmp/manual/html_node/Toom-Cook-3-Way-Multiplication.html)
210
211
self._speedup_matrix = (matrix(Integers(), 3, 3, [2, 4, 8,
212
1, 1, 1,
213
8, 4, 2])**(-1)).change_ring(base_ring).list()
214
215
def __repr__(self):
216
"""
217
String representation
218
219
EXAMPLES::
220
221
sage: B.<t> = PolynomialRing(Integers(125))
222
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
223
sage: print R
224
SpecialCubicQuotientRing over Ring of integers modulo 125 with polynomial T = x^3 + 124*x + 94
225
"""
226
return "SpecialCubicQuotientRing over %s with polynomial T = %s" % \
227
(self.base_ring(), PolynomialRing(self.base_ring(), 'x')(
228
[self._b, self._a, 0, 1]))
229
230
def poly_ring(self):
231
"""
232
Return the underlying polynomial ring in `T`.
233
234
EXAMPLES::
235
236
sage: B.<t> = PolynomialRing(Integers(125))
237
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
238
sage: R.poly_ring()
239
Univariate Polynomial Ring in T over Ring of integers modulo 125
240
"""
241
return self._poly_ring
242
243
def gens(self):
244
"""
245
Return a list [x, T] where x and T are the generators of the ring
246
(as element *of this ring*).
247
248
.. note::
249
250
I have no idea if this is compatible with the usual Sage
251
'gens' interface.
252
253
EXAMPLES::
254
255
sage: B.<t> = PolynomialRing(Integers(125))
256
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
257
sage: x, T = R.gens()
258
sage: x
259
(0) + (1)*x + (0)*x^2
260
sage: T
261
(T) + (0)*x + (0)*x^2
262
"""
263
return [SpecialCubicQuotientRingElement(self, self._poly_ring(0),
264
self._poly_ring(1),
265
self._poly_ring(0),
266
check=False),
267
SpecialCubicQuotientRingElement(self, self._poly_generator,
268
self._poly_ring(0),
269
self._poly_ring(0),
270
check=False)]
271
272
def create_element(self, p0, p1, p2, check=True):
273
"""
274
Creates the element `p_0 + p_1*x + p_2*x^2`, where the `p_i`
275
are polynomials in `T`.
276
277
INPUT:
278
279
- ``p0, p1, p2`` -- coefficients; must be coercible
280
into poly_ring()
281
282
- ``check`` -- bool (default True): whether to carry
283
out coercion
284
285
EXAMPLES::
286
287
sage: B.<t> = PolynomialRing(Integers(125))
288
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
289
sage: A, z = R.poly_ring().objgen()
290
sage: R.create_element(z^2, z+1, 3)
291
(T^2) + (T + 1)*x + (3)*x^2
292
"""
293
return SpecialCubicQuotientRingElement(self, p0, p1, p2, check)
294
295
def __call__(self, value):
296
"""
297
EXAMPLES::
298
299
sage: B.<t> = PolynomialRing(Integers(125))
300
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
301
sage: R(3)
302
(3) + (0)*x + (0)*x^2
303
"""
304
return self._coerce_(value)
305
306
def _coerce_impl(self, value):
307
"""
308
EXAMPLES::
309
310
sage: B.<t> = PolynomialRing(Integers(125))
311
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
312
sage: R._coerce_impl(3)
313
(3) + (0)*x + (0)*x^2
314
"""
315
# coerce to underlying polynomial ring (possibly via base ring):
316
value = self._poly_ring._coerce_(value)
317
318
return SpecialCubicQuotientRingElement(self, value, self._poly_ring(0),
319
self._poly_ring(0), check=False)
320
321
322
class SpecialCubicQuotientRingElement(CommutativeAlgebraElement):
323
"""
324
An element of a SpecialCubicQuotientRing.
325
"""
326
def __init__(self, parent, p0, p1, p2, check=True):
327
"""
328
Constructs the element `p_0 + p_1*x + p_2*x^2`, where
329
the `p_i` are polynomials in `T`.
330
331
INPUT:
332
333
- ``parent`` -- a SpecialCubicQuotientRing
334
335
- ``p0, p1, p2`` -- coefficients; must be coercible
336
into parent.poly_ring()
337
338
- ``check`` -- bool (default True): whether to carry
339
out coercion
340
341
EXAMPLES::
342
343
sage: B.<t> = PolynomialRing(Integers(125))
344
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
345
sage: from sage.schemes.elliptic_curves.monsky_washnitzer import SpecialCubicQuotientRingElement
346
sage: SpecialCubicQuotientRingElement(R, 2, 3, 4)
347
(2) + (3)*x + (4)*x^2
348
"""
349
if not isinstance(parent, SpecialCubicQuotientRing):
350
raise TypeError("parent (=%s) must be a SpecialCubicQuotientRing" % parent)
351
352
CommutativeAlgebraElement.__init__(self, parent)
353
354
if check:
355
poly_ring = parent.poly_ring()
356
p0 = poly_ring(p0)
357
p1 = poly_ring(p1)
358
p2 = poly_ring(p2)
359
360
self._triple = (p0, p1, p2)
361
362
def coeffs(self):
363
"""
364
Returns list of three lists of coefficients, corresponding to the
365
`x^0`, `x^1`, `x^2` coefficients. The lists
366
are zero padded to the same length. The list entries belong to the
367
base ring.
368
369
EXAMPLES::
370
371
sage: B.<t> = PolynomialRing(Integers(125))
372
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
373
sage: p = R.create_element(t, t^2 - 2, 3)
374
sage: p.coeffs()
375
[[0, 1, 0], [123, 0, 1], [3, 0, 0]]
376
"""
377
coeffs = [column.coeffs() for column in self._triple]
378
degree = max([len(x) for x in coeffs])
379
base_ring = self.parent().base_ring()
380
for column in coeffs:
381
column.extend([base_ring(0)] * (degree - len(column)))
382
return coeffs
383
384
def __nonzero__(self):
385
"""
386
EXAMPLES::
387
388
sage: B.<t> = PolynomialRing(Integers(125))
389
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
390
sage: x, T = R.gens()
391
sage: not x
392
False
393
sage: not T
394
False
395
sage: not R.create_element(0, 0, 0)
396
True
397
"""
398
return not not self._triple[0] or not not self._triple[1] or not not self._triple[2]
399
400
def __cmp__(self, other):
401
"""
402
EXAMPLES::
403
404
sage: B.<t> = PolynomialRing(Integers(125))
405
sage: x, t = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4)).gens()
406
sage: x == t
407
False
408
sage: x == x
409
True
410
sage: x == x + x - x
411
True
412
"""
413
return cmp(self._triple, other._triple)
414
415
def _repr_(self):
416
"""
417
EXAMPLES::
418
419
sage: B.<t> = PolynomialRing(Integers(125))
420
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
421
sage: x, T = R.gens()
422
sage: x + T*x - 2*T^2
423
(123*T^2) + (T + 1)*x + (0)*x^2
424
"""
425
return "(%s) + (%s)*x + (%s)*x^2" % self._triple
426
427
def _latex_(self):
428
"""
429
EXAMPLES::
430
431
sage: B.<t> = PolynomialRing(Integers(125))
432
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
433
sage: x, T = R.gens()
434
sage: f = x + T*x - 2*T^2
435
sage: latex(f)
436
(123 T^{2}) + (T + 1)x + (0)x^2
437
"""
438
return "(%s) + (%s)x + (%s)x^2" % \
439
tuple([column._latex_() for column in self._triple])
440
441
def _add_(self, other):
442
"""
443
EXAMPLES::
444
445
sage: B.<t> = PolynomialRing(Integers(125))
446
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
447
sage: f = R.create_element(2, t, t^2 - 3)
448
sage: g = R.create_element(3 + t, -t, t)
449
sage: f + g
450
(T + 5) + (0)*x + (T^2 + T + 122)*x^2
451
"""
452
return SpecialCubicQuotientRingElement(self.parent(),
453
self._triple[0] + other._triple[0],
454
self._triple[1] + other._triple[1],
455
self._triple[2] + other._triple[2],
456
check=False)
457
458
def _sub_(self, other):
459
"""
460
EXAMPLES::
461
462
sage: B.<t> = PolynomialRing(Integers(125))
463
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
464
sage: f = R.create_element(2, t, t^2 - 3)
465
sage: g = R.create_element(3 + t, -t, t)
466
sage: f - g
467
(124*T + 124) + (2*T)*x + (T^2 + 124*T + 122)*x^2
468
"""
469
return SpecialCubicQuotientRingElement(self.parent(),
470
self._triple[0] - other._triple[0],
471
self._triple[1] - other._triple[1],
472
self._triple[2] - other._triple[2],
473
check=False)
474
475
def shift(self, n):
476
"""
477
Returns this element multiplied by `T^n`.
478
479
EXAMPLES::
480
481
sage: B.<t> = PolynomialRing(Integers(125))
482
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
483
sage: f = R.create_element(2, t, t^2 - 3)
484
sage: f
485
(2) + (T)*x + (T^2 + 122)*x^2
486
sage: f.shift(1)
487
(2*T) + (T^2)*x + (T^3 + 122*T)*x^2
488
sage: f.shift(2)
489
(2*T^2) + (T^3)*x + (T^4 + 122*T^2)*x^2
490
"""
491
return SpecialCubicQuotientRingElement(self.parent(),
492
self._triple[0].shift(n),
493
self._triple[1].shift(n),
494
self._triple[2].shift(n),
495
check=False)
496
497
def scalar_multiply(self, scalar):
498
"""
499
Multiplies this element by a scalar, i.e. just multiply each
500
coefficient of `x^j` by the scalar.
501
502
INPUT:
503
504
- ``scalar`` -- either an element of base_ring, or an
505
element of poly_ring.
506
507
EXAMPLES::
508
509
sage: B.<t> = PolynomialRing(Integers(125))
510
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
511
sage: x, T = R.gens()
512
sage: f = R.create_element(2, t, t^2 - 3)
513
sage: f
514
(2) + (T)*x + (T^2 + 122)*x^2
515
sage: f.scalar_multiply(2)
516
(4) + (2*T)*x + (2*T^2 + 119)*x^2
517
sage: f.scalar_multiply(t)
518
(2*T) + (T^2)*x + (T^3 + 122*T)*x^2
519
"""
520
scalar = self.parent()._poly_ring(scalar)
521
return SpecialCubicQuotientRingElement(self.parent(),
522
scalar * self._triple[0],
523
scalar * self._triple[1],
524
scalar * self._triple[2],
525
check=False)
526
527
def square(self):
528
"""
529
EXAMPLES::
530
531
sage: B.<t> = PolynomialRing(Integers(125))
532
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
533
sage: x, T = R.gens()
534
535
::
536
537
sage: f = R.create_element(1 + 2*t + 3*t^2, 4 + 7*t + 9*t^2, 3 + 5*t + 11*t^2)
538
sage: f.square()
539
(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^2
540
"""
541
return self * self
542
543
def _mul_(self, other):
544
"""
545
EXAMPLES::
546
547
sage: B.<t> = PolynomialRing(Integers(125))
548
sage: R = monsky_washnitzer.SpecialCubicQuotientRing(t^3 - t + B(1/4))
549
sage: x, T = R.gens()
550
551
::
552
553
sage: f = R.create_element(1 + 2*t + 3*t^2, 4 + 7*t + 9*t^2, 3 + 5*t + 11*t^2)
554
sage: g = R.create_element(4 + 3*t + 7*t^2, 2 + 3*t + t^2, 8 + 4*t + 6*t^2)
555
sage: f * g
556
(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^2
557
"""
558
if not isinstance(other, SpecialCubicQuotientRingElement):
559
return self.scalar_multiply(other)
560
561
# Here we do Toom-Cook three-way multiplication, which reduces the
562
# naive 9 polynomial multiplications to only 5 polynomial multiplications.
563
564
a0, a1, a2 = self._triple
565
b0, b1, b2 = other._triple
566
M = self.parent()._speedup_matrix
567
568
if self is other:
569
# faster method if we're squaring
570
p0 = a0 * a0
571
temp = a0 + 2*a1 + 4*a2
572
p1 = temp * temp
573
temp = a0 + a1 + a2
574
p2 = temp * temp
575
temp = 4*a0 + 2*a1 + a2
576
p3 = temp * temp
577
p4 = a2 * a2
578
579
else:
580
p0 = a0 * b0
581
p1 = (a0 + 2*a1 + 4*a2) * (b0 + 2*b1 + 4*b2)
582
p2 = (a0 + a1 + a2) * (b0 + b1 + b2)
583
p3 = (4*a0 + 2*a1 + a2) * (4*b0 + 2*b1 + b2)
584
p4 = a2 * b2
585
586
q1 = p1 - p0 - 16*p4
587
q2 = p2 - p0 - p4
588
q3 = p3 - 16*p0 - p4
589
590
c0 = p0
591
c1 = M[0]*q1 + M[1]*q2 + M[2]*q3
592
c2 = M[3]*q1 + M[4]*q2 + M[5]*q3
593
c3 = M[6]*q1 + M[7]*q2 + M[8]*q3
594
c4 = p4
595
596
# Now the product is c0 + c1 x + c2 x^2 + c3 x^3 + c4 x^4.
597
# We need to reduce mod y = x^3 + ax + b and return result.
598
599
parent = self.parent()
600
T = parent._poly_generator
601
b = parent._b
602
a = parent._a
603
604
# todo: These lines are necessary to get binop stuff working
605
# for certain base rings, e.g. when we compute b*c3 in the
606
# final line. They shouldn't be necessary. Need to fix this
607
# somewhere else in Sage.
608
a = parent._poly_ring(a)
609
b = parent._poly_ring(b)
610
611
return SpecialCubicQuotientRingElement(parent,
612
-b*c3 + c0 + c3*T,
613
-b*c4 - a*c3 + c1 + c4*T,
614
-a*c4 + c2,
615
check=False)
616
617
618
def transpose_list(input):
619
"""
620
INPUT:
621
622
- ``input`` -- a list of lists, each list of the same
623
length
624
625
OUTPUT:
626
627
- ``output`` -- a list of lists such that output[i][j]
628
= input[j][i]
629
630
EXAMPLES::
631
632
sage: from sage.schemes.elliptic_curves.monsky_washnitzer import transpose_list
633
sage: L = [[1, 2], [3, 4], [5, 6]]
634
sage: transpose_list(L)
635
[[1, 3, 5], [2, 4, 6]]
636
"""
637
h = len(input)
638
w = len(input[0])
639
640
output = []
641
for i in range(w):
642
row = []
643
for j in range(h):
644
row.append(input[j][i])
645
output.append(row)
646
return output
647
648
649
def helper_matrix(Q):
650
"""
651
Computes the (constant) matrix used to calculate the linear
652
combinations of the `d(x^i y^j)` needed to eliminate the
653
negative powers of `y` in the cohomology (i.e. in
654
reduce_negative()).
655
656
INPUT:
657
658
- ``Q`` -- cubic polynomial
659
660
EXAMPLES::
661
662
sage: t = polygen(QQ,'t')
663
sage: from sage.schemes.elliptic_curves.monsky_washnitzer import helper_matrix
664
sage: helper_matrix(t**3-4*t-691)
665
[ 64/12891731 -16584/12891731 4297329/12891731]
666
[ 6219/12891731 -32/12891731 8292/12891731]
667
[ -24/12891731 6219/12891731 -32/12891731]
668
"""
669
a = Q[1]
670
b = Q[0]
671
672
# Discriminant (should be invertible for a curve of good reduction)
673
D = 4*a**3 + 27*b**2
674
Dinv = D**(-1) # NB do not use 1/D
675
676
# This is the inverse of the matrix
677
# [ a, -3b, 0 ]
678
# [ 0, -2a, -3b ]
679
# [ 3, 0, -2a ]
680
681
return Dinv * matrix([[4*a**2, -6*b*a, 9*b**2],
682
[-9*b, -2*a**2, 3*b*a],
683
[6*a, -9*b, -2*a**2]])
684
685
686
def lift(x):
687
r"""
688
Tries to call x.lift(), presumably from the `p`-adics to ZZ.
689
690
If this fails, it assumes the input is a power series, and tries to
691
lift it to a power series over QQ.
692
693
This function is just a very kludgy solution to the problem of
694
trying to make the reduction code (below) work over both Zp and
695
Zp[[t]].
696
697
EXAMPLES::
698
699
sage: from sage.schemes.elliptic_curves.monsky_washnitzer import lift
700
sage: l = lift(Qp(13)(131)); l
701
131
702
sage: l.parent()
703
Integer Ring
704
705
sage: x=PowerSeriesRing(Qp(17),'x').gen()
706
sage: l = lift(4+5*x+17*x**6); l
707
4 + 5*t + 17*t^6
708
sage: l.parent()
709
Power Series Ring in t over Rational Field
710
"""
711
try:
712
return x.lift()
713
except AttributeError:
714
return PowerSeriesRing(Rationals(), "t")(x.list(), x.prec())
715
716
717
def reduce_negative(Q, p, coeffs, offset, exact_form=None):
718
"""
719
Applies cohomology relations to incorporate negative powers of
720
`y` into the `y^0` term.
721
722
INPUT:
723
724
- ``p`` -- prime
725
726
- ``Q`` -- cubic polynomial
727
728
- ``coeffs`` -- list of length 3 lists. The
729
`i^{th}` list [a, b, c] represents
730
`y^{2(i - offset)} (a + bx + cx^2) dx/y`.
731
732
- ``offset`` -- nonnegative integer
733
734
OUTPUT: The reduction is performed in-place. The output is placed
735
in coeffs[offset]. Note that coeffs[i] will be meaningless for i
736
offset after this function is finished.
737
738
EXAMPLE::
739
740
sage: R.<x> = Integers(5^3)['x']
741
sage: Q = x^3 - x + R(1/4)
742
sage: coeffs = [[10, 15, 20], [1, 2, 3], [4, 5, 6], [7, 8, 9]]
743
sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]
744
sage: monsky_washnitzer.reduce_negative(Q, 5, coeffs, 3)
745
sage: coeffs[3]
746
[28, 52, 9]
747
748
::
749
750
sage: R.<x> = Integers(7^3)['x']
751
sage: Q = x^3 - x + R(1/4)
752
sage: coeffs = [[7, 14, 21], [1, 2, 3], [4, 5, 6], [7, 8, 9]]
753
sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]
754
sage: monsky_washnitzer.reduce_negative(Q, 7, coeffs, 3)
755
sage: coeffs[3]
756
[245, 332, 9]
757
"""
758
759
m = helper_matrix(Q).list()
760
base_ring = Q.base_ring()
761
next_a = coeffs[0]
762
763
if exact_form is not None:
764
x = exact_form.parent().gen(0)
765
y = exact_form.parent()(exact_form.parent().base_ring().gen(0))
766
767
try:
768
three_j_plus_5 = 5 - base_ring(6*offset)
769
three_j_plus_7 = 7 - base_ring(6*offset)
770
six = base_ring(6)
771
772
for i in range(0, offset):
773
774
j = 2*(i-offset)
775
a = next_a
776
next_a = coeffs[i+1]
777
778
# todo: the following divisions will sometimes involve
779
# a division by (a power of) p. In all cases, we know (from
780
# Kedlaya's estimates) that the answer should be p-integral.
781
# However, since we're working over $Z/p^k Z$, we're not allowed
782
# to "divide by p". So currently we lift to Q, divide, and coerce
783
# back. Eventually, when pAdicInteger is implemented, and plays
784
# nicely with pAdicField, we should reimplement this stuff
785
# using pAdicInteger.
786
787
if (p.divides(j+1)):
788
# need to lift here to perform the division
789
a[0] = base_ring(lift(a[0]) / (j+1))
790
a[1] = base_ring(lift(a[1]) / (j+1))
791
a[2] = base_ring(lift(a[2]) / (j+1))
792
else:
793
j_plus_1_inv = ~base_ring(j+1)
794
a[0] = a[0] * j_plus_1_inv
795
a[1] = a[1] * j_plus_1_inv
796
a[2] = a[2] * j_plus_1_inv
797
798
c1 = m[3]*a[0] + m[4]*a[1] + m[5]*a[2]
799
c2 = m[6]*a[0] + m[7]*a[1] + m[8]*a[2]
800
next_a[0] = next_a[0] - three_j_plus_5 * c1
801
next_a[1] = next_a[1] - three_j_plus_7 * c2
802
803
three_j_plus_7 = three_j_plus_7 + six
804
three_j_plus_5 = three_j_plus_5 + six
805
806
if exact_form is not None:
807
c0 = m[0]*a[0] + m[1]*a[1] + m[2]*a[2]
808
exact_form += (c0 + c1*x + c2 * x**2) * y**(j+1)
809
810
except NotImplementedError:
811
raise NotImplementedError("It looks like you've found a "
812
"non-integral matrix of Frobenius! "
813
"(Q=%s, p=%s)\nTime to write a paper." % (Q, p))
814
815
coeffs[int(offset)] = next_a
816
817
return exact_form
818
819
820
def reduce_positive(Q, p, coeffs, offset, exact_form=None):
821
"""
822
Applies cohomology relations to incorporate positive powers of
823
`y` into the `y^0` term.
824
825
INPUT:
826
827
- ``Q`` -- cubic polynomial
828
829
- ``coeffs`` -- list of length 3 lists. The
830
`i^{th}` list [a, b, c] represents
831
`y^{2(i - offset)} (a + bx + cx^2) dx/y`.
832
833
- ``offset`` -- nonnegative integer
834
835
OUTPUT: The reduction is performed in-place. The output is placed
836
in coeffs[offset]. Note that coeffs[i] will be meaningless for i
837
offset after this function is finished.
838
839
EXAMPLE::
840
841
sage: R.<x> = Integers(5^3)['x']
842
sage: Q = x^3 - x + R(1/4)
843
844
::
845
846
sage: coeffs = [[1, 2, 3], [10, 15, 20]]
847
sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]
848
sage: monsky_washnitzer.reduce_positive(Q, 5, coeffs, 0)
849
sage: coeffs[0]
850
[16, 102, 88]
851
852
::
853
854
sage: coeffs = [[9, 8, 7], [10, 15, 20]]
855
sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]
856
sage: monsky_washnitzer.reduce_positive(Q, 5, coeffs, 0)
857
sage: coeffs[0]
858
[24, 108, 92]
859
"""
860
861
base_ring = Q.base_ring()
862
next_a = coeffs[len(coeffs) - 1]
863
864
Qa = Q[1]
865
Qb = Q[0]
866
867
A = 2*Qa
868
B = 3*Qb
869
870
offset = Integer(offset)
871
872
if exact_form is not None:
873
x = exact_form.parent().gen(0)
874
y = exact_form.parent().base_ring().gen(0)
875
# y = exact_form.parent()(exact_form.parent().base_ring().gen(0))
876
877
for i in range(len(coeffs)-1, offset, -1):
878
j = 2*(i-offset) - 2
879
a = next_a
880
next_a = coeffs[i-1]
881
882
a[0] = a[0] - Qa*a[2]/3 # subtract d(y^j + 3)
883
if exact_form is not None:
884
exact_form += Q.base_ring()(a[2].lift() / (3*j+9)) * y**(j+3)
885
886
# todo: see comments about pAdicInteger in reduceNegative()
887
888
# subtract off c1 of d(x y^j + 1), and
889
if p.divides(3*j + 5):
890
c1 = base_ring(lift(a[0]) / (3*j + 5))
891
else:
892
c1 = a[0] / (3*j + 5)
893
894
# subtract off c2 of d(x^2 y^j + 1)
895
if p.divides(3*j + 7):
896
c2 = base_ring(lift(a[1]) / (3*j + 7))
897
else:
898
c2 = a[1] / (3*j + 7)
899
900
next_a[0] = next_a[0] + B*c1*(j+1)
901
next_a[1] = next_a[1] + A*c1*(j+1) + B*c2*(j+1)
902
next_a[2] = next_a[2] + A*c2*(j+1)
903
904
if exact_form is not None:
905
exact_form += (c1*x + c2 * x**2) * y**(j+1)
906
907
coeffs[int(offset)] = next_a
908
909
return exact_form
910
911
912
def reduce_zero(Q, coeffs, offset, exact_form=None):
913
"""
914
Applies cohomology relation to incorporate `x^2 y^0` term
915
into `x^0 y^0` and `x^1 y^0` terms.
916
917
INPUT:
918
919
- ``Q`` -- cubic polynomial
920
921
- ``coeffs`` -- list of length 3 lists. The
922
`i^{th}` list [a, b, c] represents
923
`y^{2(i - offset)} (a + bx + cx^2) dx/y`.
924
925
- ``offset`` -- nonnegative integer
926
927
OUTPUT: The reduction is performed in-place. The output is placed
928
in coeffs[offset]. This method completely ignores coeffs[i] for i
929
!= offset.
930
931
EXAMPLE::
932
933
sage: R.<x> = Integers(5^3)['x']
934
sage: Q = x^3 - x + R(1/4)
935
sage: coeffs = [[1, 2, 3], [4, 5, 6], [7, 8, 9]]
936
sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]
937
sage: monsky_washnitzer.reduce_zero(Q, coeffs, 1)
938
sage: coeffs[1]
939
[6, 5, 0]
940
"""
941
942
a = coeffs[int(offset)]
943
if a[2] == 0:
944
return exact_form
945
946
Qa = Q[1]
947
948
a[0] = a[0] - a[2]*Qa/3 # $3x^2 dx/y = -a dx/y$
949
950
coeffs[int(offset)] = a
951
952
if exact_form is not None:
953
y = exact_form.parent()(exact_form.parent().base_ring().gen(0))
954
exact_form += Q.base_ring()(a[2] / 3) * y
955
956
a[2] = 0
957
958
coeffs[int(offset)] = a
959
return exact_form
960
961
962
def reduce_all(Q, p, coeffs, offset, compute_exact_form=False):
963
"""
964
Applies cohomology relations to reduce all terms to a linear
965
combination of `dx/y` and `x dx/y`.
966
967
INPUT:
968
969
- ``Q`` -- cubic polynomial
970
971
- ``coeffs`` -- list of length 3 lists. The
972
`i^{th}` list [a, b, c] represents
973
`y^{2(i - offset)} (a + bx + cx^2) dx/y`.
974
975
- ``offset`` -- nonnegative integer
976
977
978
OUTPUT:
979
980
981
- ``A, B`` - pair such that the input differential is
982
cohomologous to (A + Bx) dx/y.
983
984
985
.. note::
986
987
The algorithm operates in-place, so the data in coeffs is
988
destroyed.
989
990
EXAMPLE::
991
992
sage: R.<x> = Integers(5^3)['x']
993
sage: Q = x^3 - x + R(1/4)
994
sage: coeffs = [[1, 2, 3], [4, 5, 6], [7, 8, 9]]
995
sage: coeffs = [[R.base_ring()(a) for a in row] for row in coeffs]
996
sage: monsky_washnitzer.reduce_all(Q, 5, coeffs, 1)
997
(21, 106)
998
"""
999
1000
R = Q.base_ring()
1001
1002
if compute_exact_form:
1003
# exact_form = SpecialCubicQuotientRing(Q, laurent_series=True)(0)
1004
exact_form = PolynomialRing(LaurentSeriesRing(Q.base_ring(), 'y'), 'x')(0)
1005
# t = (Q.base_ring().order().factor())[0]
1006
# from sage.rings.padics.qp import pAdicField
1007
# exact_form = PolynomialRing(LaurentSeriesRing(pAdicField(p, t[1]), 'y'), 'x')(0)
1008
else:
1009
exact_form = None
1010
1011
while len(coeffs) <= offset:
1012
coeffs.append([R(0), R(0), R(0)])
1013
1014
exact_form = reduce_negative(Q, p, coeffs, offset, exact_form)
1015
exact_form = reduce_positive(Q, p, coeffs, offset, exact_form)
1016
exact_form = reduce_zero(Q, coeffs, offset, exact_form)
1017
1018
if exact_form is None:
1019
return coeffs[int(offset)][0], coeffs[int(offset)][1]
1020
else:
1021
return (coeffs[int(offset)][0], coeffs[int(offset)][1]), exact_form
1022
1023
1024
def frobenius_expansion_by_newton(Q, p, M):
1025
r"""
1026
Computes the action of Frobenius on `dx/y` and on
1027
`x dx/y`, using Newton's method (as suggested in Kedlaya's
1028
paper [Ked2001]_).
1029
1030
(This function does *not* yet use the cohomology relations - that
1031
happens afterwards in the "reduction" step.)
1032
1033
More specifically, it finds `F_0` and `F_1` in
1034
the quotient ring `R[x, T]/(T - Q(x))`, such that
1035
1036
.. math::
1037
1038
F( dx/y) = T^{-r} F0 dx/y, \text{\ and\ } F(x dx/y) = T^{-r} F1 dx/y
1039
1040
where
1041
1042
.. math::
1043
1044
r = ( (2M-3)p - 1 )/2.
1045
1046
1047
(Here `T` is `y^2 = z^{-2}`, and `R` is the
1048
coefficient ring of `Q`.)
1049
1050
`F_0` and `F_1` are computed in the
1051
SpecialCubicQuotientRing associated to `Q`, so all powers
1052
of `x^j` for `j \geq 3` are reduced to powers of
1053
`T`.
1054
1055
INPUT:
1056
1057
- ``Q`` -- cubic polynomial of the form
1058
`Q(x) = x^3 + ax + b`, whose coefficient ring is a
1059
`Z/(p^M)Z`-algebra
1060
1061
- ``p`` -- residue characteristic of the p-adic field
1062
1063
- ``M`` -- p-adic precision of the coefficient ring
1064
(this will be used to determine the number of Newton iterations)
1065
1066
OUTPUT:
1067
1068
- ``F0, F1`` - elements of
1069
SpecialCubicQuotientRing(Q), as described above
1070
1071
- ``r`` - non-negative integer, as described above
1072
1073
EXAMPLES::
1074
1075
sage: from sage.schemes.elliptic_curves.monsky_washnitzer import frobenius_expansion_by_newton
1076
sage: R.<x> = Integers(5^3)['x']
1077
sage: Q = x^3 - x + R(1/4)
1078
sage: frobenius_expansion_by_newton(Q,5,3)
1079
((25*T^5 + 75*T^3 + 100*T^2 + 100*T + 100) + (5*T^6 + 80*T^5 + 100*T^3
1080
+ 25*T + 50)*x + (55*T^5 + 50*T^4 + 75*T^3 + 25*T^2 + 25*T + 25)*x^2,
1081
(5*T^8 + 15*T^7 + 95*T^6 + 10*T^5 + 25*T^4 + 25*T^3 + 100*T^2 + 50)
1082
+ (65*T^7 + 55*T^6 + 70*T^5 + 100*T^4 + 25*T^2 + 100*T)*x
1083
+ (15*T^6 + 115*T^5 + 75*T^4 + 100*T^3 + 50*T^2 + 75*T + 75)*x^2, 7)
1084
"""
1085
1086
S = SpecialCubicQuotientRing(Q)
1087
x, _ = S.gens() # T = y^2
1088
base_ring = S.base_ring()
1089
1090
# When we compute Frob(1/y) we actually only need precision M-1, since
1091
# we're going to multiply by p at the end anyway.
1092
M = float(M - 1)
1093
1094
# Kedlaya sets s = Q(x^p)/T^p = 1 + p T^{-p} E, where
1095
# E = (Q(x^p) - Q(x)^p) / p (has integral coefficients).
1096
# Then he computes s^{-1/2} in S, using Newton's method to find
1097
# successive approximations. We follow this plan, but we normalise our
1098
# approximations so that we only ever need positive powers of T.
1099
1100
# Start by setting r = Q(x^p)/2 = 1/2 T^p s.
1101
# (The 1/2 is for convenience later on.)
1102
x_to_p_less_one = x**(p-1)
1103
x_to_p = x_to_p_less_one * x
1104
x_to_p_cubed = x_to_p.square() * x_to_p
1105
r = (base_ring(1) / base_ring(2)) * (x_to_p_cubed + Q[1]*x_to_p + S(Q[0]))
1106
1107
# todo: this next loop would be clearer if it used the newton_method_sizes()
1108
# function
1109
1110
# We will start with a hard-coded initial approximation, which we provide
1111
# up to precision 3. First work out what precision is best to start with.
1112
if M <= 3:
1113
initial_precision = M
1114
elif ceil(log(M/2, 2)) == ceil(log(M/3, 2)):
1115
# In this case there is no advantage to starting with precision three,
1116
# because we'll overshoot at the end. E.g. suppose the final precision
1117
# is 8. If we start with precision 2, we need two iterations to get us
1118
# to 8. If we start at precision 3, we will still need two iterations,
1119
# but we do more work along the way. So may as well start with only 2.
1120
initial_precision = 2
1121
else:
1122
initial_precision = 3
1123
1124
# Now compute the first approximation. In the main loop below, X is the
1125
# normalised approximation, and k is the precision. More specifically,
1126
# X = T^{p(k-1)} x_i, where x_i is an approximation to s^{-1/2}, and the
1127
# approximation is correct mod p^k.
1128
if initial_precision == 1:
1129
k = 1
1130
X = S(1)
1131
elif initial_precision == 2:
1132
# approximation is 3/2 - 1/2 s
1133
k = 2
1134
X = S(base_ring(3) / base_ring(2)).shift(p) - r
1135
elif initial_precision == 3:
1136
# approximation is (15 - 10 s + 3 s^2) / 8
1137
k = 3
1138
X = (base_ring(1) / base_ring(8)) * (S(15).shift(2*p)
1139
- (base_ring(20) * r).shift(p) +
1140
(base_ring(12) * r.square()))
1141
# The key to the following calculation is that the T^{-m} coefficient
1142
# of every x_i is divisible by p^(ceil(m/p)) (for m >= 0). Therefore if
1143
# we are only expecting an answer correct mod p^k, we can truncate
1144
# beyond the T^{-(k-1)p} term without any problems.
1145
1146
# todo: what would be really nice is to be able to work in a lower
1147
# precision *coefficient ring* when we start the iteration, and move up to
1148
# higher precision rings as the iteration proceeds. This would be feasible
1149
# over Integers(p**n), but quite complicated (maybe impossible) over a more
1150
# general base ring. This might give a decent constant factor speedup;
1151
# or it might not, depending on how much the last iteration dominates the
1152
# whole runtime. My guess is that it isn't worth the effort.
1153
1154
three_halves = base_ring(3) / base_ring(2)
1155
1156
# Newton iteration loop
1157
while k < M:
1158
# target_k = k' = precision we want our answer to be after this iteration
1159
target_k = 2*k
1160
1161
# This prevents us overshooting. For example if the current precision
1162
# is 3 and we want to get to 10, we're better off going up to 5
1163
# instead of 6, because it is less work to get from 5 to 10 than it
1164
# is to get from 6 to 10.
1165
if ceil(log(M/target_k, 2)) == ceil(log(M/(target_k-1), 2)):
1166
target_k -= 1
1167
1168
# temp = T^{p(3k-2)} 1/2 s x_i^3
1169
temp = X.square() * (X * r)
1170
1171
# We know that the final result is only going to be correct mod
1172
# p^(target_k), so we might as well truncate the extraneous terms now.
1173
# temp = T^{p(k'-1)} 1/2 s x_i^3
1174
temp = temp.shift(-p*(3*k - target_k - 1))
1175
1176
# X = T^{p(k'-1)} (3/2 x_i - 1/2 s x_i^3)
1177
# = T^{p(k'-1)} x_{i+1}
1178
X = (three_halves * X).shift(p*(target_k - k)) - temp
1179
1180
k = target_k
1181
1182
# Now k should equal M, since we're up to the correct precision
1183
assert k == M, "Oops, something went wrong in the iteration"
1184
1185
# We should have s^{-1/2} correct to precision M.
1186
# The following line can be uncommented to verify this.
1187
# (It is a slow verification though, can double the whole computation time.)
1188
1189
#assert (p * X.square() * r * base_ring(2)).coeffs() == \
1190
# R(p).shift(p*(2*M - 1)).coeffs()
1191
1192
# Finally incorporate frobenius of dx and x dx, and choose offset that
1193
# compensates for our normalisations by powers of T.
1194
F0 = base_ring(p) * x_to_p_less_one * X
1195
F1 = F0 * x_to_p
1196
offset = ((2*k-1)*p - 1)/2
1197
1198
return F0, F1, offset
1199
1200
1201
def frobenius_expansion_by_series(Q, p, M):
1202
r"""
1203
Computes the action of Frobenius on `dx/y` and on `x dx/y`, using a
1204
series expansion.
1205
1206
(This function computes the same thing as
1207
frobenius_expansion_by_newton(), using a different method.
1208
Theoretically the Newton method should be asymptotically faster,
1209
when the precision gets large. However, in practice, this functions
1210
seems to be marginally faster for moderate precision, so I'm
1211
keeping it here until I figure out exactly why it is faster.)
1212
1213
(This function does *not* yet use the cohomology relations - that
1214
happens afterwards in the "reduction" step.)
1215
1216
More specifically, it finds F0 and F1 in the quotient ring
1217
`R[x, T]/(T - Q(x))`, such that
1218
`F( dx/y) = T^{-r} F0 dx/y`, and
1219
`F(x dx/y) = T^{-r} F1 dx/y` where
1220
`r = ( (2M-3)p - 1 )/2`. (Here `T` is `y^2 = z^{-2}`,
1221
and `R` is the coefficient ring of `Q`.)
1222
1223
`F_0` and `F_1` are computed in the
1224
SpecialCubicQuotientRing associated to `Q`, so all powers
1225
of `x^j` for `j \geq 3` are reduced to powers of
1226
`T`.
1227
1228
It uses the sum
1229
1230
.. math::
1231
1232
F0 = \sum_{k=0}^{M-2} {-1/2 \choose k} p x^{p-1} E^k T^{(M-2-k)p}
1233
1234
and
1235
1236
.. math::
1237
1238
F1 = x^p F0,
1239
1240
where `E = Q(x^p) - Q(x)^p`.
1241
1242
INPUT:
1243
1244
- ``Q`` -- cubic polynomial of the form
1245
`Q(x) = x^3 + ax + b`, whose coefficient ring is a
1246
`\ZZ/(p^M)\ZZ` -algebra
1247
1248
- ``p`` -- residue characteristic of the `p`-adic field
1249
1250
- ``M`` -- `p`-adic precision of the coefficient ring
1251
(this will be used to determine the number of terms in the
1252
series)
1253
1254
OUTPUT:
1255
1256
- ``F0, F1`` - elements of
1257
SpecialCubicQuotientRing(Q), as described above
1258
1259
- ``r`` - non-negative integer, as described above
1260
1261
EXAMPLES::
1262
1263
sage: from sage.schemes.elliptic_curves.monsky_washnitzer import frobenius_expansion_by_series
1264
sage: R.<x> = Integers(5^3)['x']
1265
sage: Q = x^3 - x + R(1/4)
1266
sage: frobenius_expansion_by_series(Q,5,3)
1267
((25*T^5 + 75*T^3 + 100*T^2 + 100*T + 100) + (5*T^6 + 80*T^5 + 100*T^3
1268
+ 25*T + 50)*x + (55*T^5 + 50*T^4 + 75*T^3 + 25*T^2 + 25*T + 25)*x^2,
1269
(5*T^8 + 15*T^7 + 95*T^6 + 10*T^5 + 25*T^4 + 25*T^3 + 100*T^2 + 50)
1270
+ (65*T^7 + 55*T^6 + 70*T^5 + 100*T^4 + 25*T^2 + 100*T)*x
1271
+ (15*T^6 + 115*T^5 + 75*T^4 + 100*T^3 + 50*T^2 + 75*T + 75)*x^2, 7)
1272
"""
1273
1274
S = SpecialCubicQuotientRing(Q)
1275
x, _ = S.gens()
1276
base_ring = S.base_ring()
1277
1278
x_to_p_less_1 = x**(p-1)
1279
x_to_p = x_to_p_less_1 * x
1280
1281
# compute frobQ = Q(x^p)
1282
x_to_p_squared = x_to_p * x_to_p
1283
x_to_p_cubed = x_to_p_squared * x_to_p
1284
frobQ = x_to_p_cubed + Q[1]*x_to_p + Q[0]*S(1)
1285
# anticipating the day when p = 3 is supported:
1286
# frobQ = x_to_p_cubed + Q[2]*x_to_p_squared + Q[1]*x_to_p + Q[0]*S(1)
1287
1288
E = frobQ - S(1).shift(p) # E = Q(x^p) - Q(x)^p
1289
1290
offset = int(((2*M-3)*p-1)/2)
1291
term = p * x_to_p_less_1
1292
F0 = term.shift((M-2)*p)
1293
1294
# todo: Possible speedup idea, perhaps by a factor of 2, but
1295
# it requires a lot of work:
1296
# Note that p divides E, so p^k divides E^k. So when we are
1297
# working with high powers of E, we're doing a lot more work
1298
# in the multiplications than we need to. To take advantage of
1299
# this we would need some protocol for "lowering the precision"
1300
# of a SpecialCubicQuotientRing. This would be quite messy to
1301
# do properly over an arbitrary base ring. Perhaps it is
1302
# feasible to do for the most common case (i.e. Z/p^nZ).
1303
# (but it probably won't save much time unless p^n is very
1304
# large, because the machine word size is probably pretty
1305
# big anyway.)
1306
1307
for k in range(int(1), int(M-1)):
1308
term = term * E
1309
c = base_ring(binomial(-Integer(1)/2, k))
1310
F0 += (term * c).shift((M-k-2)*p)
1311
1312
return F0, F0 * x_to_p, offset
1313
1314
1315
def adjusted_prec(p, prec):
1316
r"""
1317
Computes how much precision is required in matrix_of_frobenius to
1318
get an answer correct to prec `p`-adic digits.
1319
1320
The issue is that the algorithm used in
1321
:func:`matrix_of_frobenius` sometimes performs divisions by `p`,
1322
so precision is lost during the algorithm.
1323
1324
The estimate returned by this function is based on Kedlaya's result
1325
(Lemmas 2 and 3 of [Ked2001]_),
1326
which implies that if we start with `M` `p`-adic
1327
digits, the total precision loss is at most
1328
`1 + \lfloor \log_p(2M-3) \rfloor` `p`-adic
1329
digits. (This estimate is somewhat less than the amount you would
1330
expect by naively counting the number of divisions by
1331
`p`.)
1332
1333
INPUT:
1334
1335
- ``p`` -- a prime = 5
1336
1337
- ``prec`` -- integer, desired output precision, = 1
1338
1339
OUTPUT: adjusted precision (usually slightly more than prec)
1340
"""
1341
1342
# initial estimate:
1343
if prec <= 2:
1344
adjusted = 2
1345
else:
1346
adjusted = prec + int(log(2*prec - 3, p)) - 1
1347
1348
# increase it until we have enough
1349
while adjusted - int(log(2*adjusted - 3, p)) - 1 < prec:
1350
adjusted += 1
1351
1352
return adjusted
1353
1354
1355
def matrix_of_frobenius(Q, p, M, trace=None, compute_exact_forms=False):
1356
"""
1357
Computes the matrix of Frobenius on Monsky-Washnitzer cohomology,
1358
with respect to the basis `(dx/y, x dx/y)`.
1359
1360
INPUT:
1361
1362
- ``Q`` -- cubic polynomial `Q(x) = x^3 + ax + b`
1363
defining an elliptic curve `E` by
1364
`y^2 = Q(x)`. The coefficient ring of `Q` should be a
1365
`\ZZ/(p^M)\ZZ`-algebra in which the matrix of
1366
frobenius will be constructed.
1367
1368
- ``p`` -- prime = 5 for which E has good reduction
1369
1370
- ``M`` -- integer = 2; `p` -adic precision of
1371
the coefficient ring
1372
1373
- ``trace`` -- (optional) the trace of the matrix, if
1374
known in advance. This is easy to compute because it is just the
1375
`a_p` of the curve. If the trace is supplied,
1376
matrix_of_frobenius will use it to speed the computation (i.e. we
1377
know the determinant is `p`, so we have two conditions, so
1378
really only column of the matrix needs to be computed. it is
1379
actually a little more complicated than that, but that's the basic
1380
idea.) If trace=None, then both columns will be computed
1381
independently, and you can get a strong indication of correctness
1382
by verifying the trace afterwards.
1383
1384
.. warning::
1385
1386
THE RESULT WILL NOT NECESSARILY BE CORRECT TO M p-ADIC
1387
DIGITS. If you want prec digits of precision, you need to use
1388
the function adjusted_prec(), and then you need to reduce the
1389
answer mod `p^{\mathrm{prec}}` at the end.
1390
1391
OUTPUT:
1392
1393
2x2 matrix of frobenius on Monsky-Washnitzer cohomology,
1394
with entries in the coefficient ring of Q.
1395
1396
EXAMPLES:
1397
1398
A simple example::
1399
1400
sage: p = 5
1401
sage: prec = 3
1402
sage: M = monsky_washnitzer.adjusted_prec(p, prec)
1403
sage: M
1404
5
1405
sage: R.<x> = PolynomialRing(Integers(p**M))
1406
sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 - x + R(1/4), p, M)
1407
sage: A
1408
[3090 187]
1409
[2945 408]
1410
1411
But the result is only accurate to prec digits::
1412
1413
sage: B = A.change_ring(Integers(p**prec))
1414
sage: B
1415
[90 62]
1416
[70 33]
1417
1418
Check trace (123 = -2 mod 125) and determinant::
1419
1420
sage: B.det()
1421
5
1422
sage: B.trace()
1423
123
1424
sage: EllipticCurve([-1, 1/4]).ap(5)
1425
-2
1426
1427
Try using the trace to speed up the calculation::
1428
1429
sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 - x + R(1/4),
1430
... p, M, -2)
1431
sage: A
1432
[2715 187]
1433
[1445 408]
1434
1435
Hmmm... it looks different, but that's because the trace of our
1436
first answer was only -2 modulo `5^3`, not -2 modulo
1437
`5^5`. So the right answer is::
1438
1439
sage: A.change_ring(Integers(p**prec))
1440
[90 62]
1441
[70 33]
1442
1443
Check it works with only one digit of precision::
1444
1445
sage: p = 5
1446
sage: prec = 1
1447
sage: M = monsky_washnitzer.adjusted_prec(p, prec)
1448
sage: R.<x> = PolynomialRing(Integers(p**M))
1449
sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 - x + R(1/4), p, M)
1450
sage: A.change_ring(Integers(p))
1451
[0 2]
1452
[0 3]
1453
1454
Here is an example that is particularly badly conditioned for
1455
using the trace trick::
1456
1457
sage: p = 11
1458
sage: prec = 3
1459
sage: M = monsky_washnitzer.adjusted_prec(p, prec)
1460
sage: R.<x> = PolynomialRing(Integers(p**M))
1461
sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 + 7*x + 8, p, M)
1462
sage: A.change_ring(Integers(p**prec))
1463
[1144 176]
1464
[ 847 185]
1465
1466
The problem here is that the top-right entry is divisible by 11,
1467
and the bottom-left entry is divisible by `11^2`. So when
1468
you apply the trace trick, neither `F(dx/y)` nor
1469
`F(x dx/y)` is enough to compute the whole matrix to the
1470
desired precision, even if you try increasing the target precision
1471
by one. Nevertheless, ``matrix_of_frobenius`` knows
1472
how to get the right answer by evaluating `F((x+1) dx/y)`
1473
instead::
1474
1475
sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 + 7*x + 8, p, M, -2)
1476
sage: A.change_ring(Integers(p**prec))
1477
[1144 176]
1478
[ 847 185]
1479
1480
The running time is about ``O(p*prec**2)`` (times some logarithmic
1481
factors), so it is feasible to run on fairly large primes, or
1482
precision (or both?!?!)::
1483
1484
sage: p = 10007
1485
sage: prec = 2
1486
sage: M = monsky_washnitzer.adjusted_prec(p, prec)
1487
sage: R.<x> = PolynomialRing(Integers(p**M))
1488
sage: A = monsky_washnitzer.matrix_of_frobenius( # long time
1489
... x^3 - x + R(1/4), p, M) # long time
1490
sage: B = A.change_ring(Integers(p**prec)); B # long time
1491
[74311982 57996908]
1492
[95877067 25828133]
1493
sage: B.det() # long time
1494
10007
1495
sage: B.trace() # long time
1496
66
1497
sage: EllipticCurve([-1, 1/4]).ap(10007) # long time
1498
66
1499
1500
::
1501
1502
sage: p = 5
1503
sage: prec = 300
1504
sage: M = monsky_washnitzer.adjusted_prec(p, prec)
1505
sage: R.<x> = PolynomialRing(Integers(p**M))
1506
sage: A = monsky_washnitzer.matrix_of_frobenius( # long time
1507
... x^3 - x + R(1/4), p, M) # long time
1508
sage: B = A.change_ring(Integers(p**prec)) # long time
1509
sage: B.det() # long time
1510
5
1511
sage: -B.trace() # long time
1512
2
1513
sage: EllipticCurve([-1, 1/4]).ap(5) # long time
1514
-2
1515
1516
Let us check consistency of the results for a range of precisions::
1517
1518
sage: p = 5
1519
sage: max_prec = 60
1520
sage: M = monsky_washnitzer.adjusted_prec(p, max_prec)
1521
sage: R.<x> = PolynomialRing(Integers(p**M))
1522
sage: A = monsky_washnitzer.matrix_of_frobenius(x^3 - x + R(1/4), p, M) # long time
1523
sage: A = A.change_ring(Integers(p**max_prec)) # long time
1524
sage: result = [] # long time
1525
sage: for prec in range(1, max_prec): # long time
1526
... M = monsky_washnitzer.adjusted_prec(p, prec) # long time
1527
... R.<x> = PolynomialRing(Integers(p^M),'x') # long time
1528
... B = monsky_washnitzer.matrix_of_frobenius( # long time
1529
... x^3 - x + R(1/4), p, M) # long time
1530
... B = B.change_ring(Integers(p**prec)) # long time
1531
... result.append(B == A.change_ring( # long time
1532
... Integers(p**prec))) # long time
1533
sage: result == [True] * (max_prec - 1) # long time
1534
True
1535
1536
The remaining examples discuss what happens when you take the
1537
coefficient ring to be a power series ring; i.e. in effect you're
1538
looking at a family of curves.
1539
1540
The code does in fact work...
1541
1542
::
1543
1544
sage: p = 11
1545
sage: prec = 3
1546
sage: M = monsky_washnitzer.adjusted_prec(p, prec)
1547
sage: S.<t> = PowerSeriesRing(Integers(p**M), default_prec=4)
1548
sage: a = 7 + t + 3*t^2
1549
sage: b = 8 - 6*t + 17*t^2
1550
sage: R.<x> = PolynomialRing(S)
1551
sage: Q = x**3 + a*x + b
1552
sage: A = monsky_washnitzer.matrix_of_frobenius(Q, p, M) # long time
1553
sage: B = A.change_ring(PowerSeriesRing(Integers(p**prec), 't', default_prec=4)) # long time
1554
sage: B # long time
1555
[1144 + 264*t + 841*t^2 + 1025*t^3 + O(t^4) 176 + 1052*t + 216*t^2 + 523*t^3 + O(t^4)]
1556
[ 847 + 668*t + 81*t^2 + 424*t^3 + O(t^4) 185 + 341*t + 171*t^2 + 642*t^3 + O(t^4)]
1557
1558
The trace trick should work for power series rings too, even in the
1559
badly- conditioned case. Unfortunately I don't know how to compute
1560
the trace in advance, so I'm not sure exactly how this would help.
1561
Also, I suspect the running time will be dominated by the
1562
expansion, so the trace trick won't really speed things up anyway.
1563
Another problem is that the determinant is not always p::
1564
1565
sage: B.det() # long time
1566
11 + 484*t^2 + 451*t^3 + O(t^4)
1567
1568
However, it appears that the determinant always has the property
1569
that if you substitute t - 11t, you do get the constant series p
1570
(mod p\*\*prec). Similarly for the trace. And since the parameter
1571
only really makes sense when it is divisible by p anyway, perhaps
1572
this isn't a problem after all.
1573
"""
1574
1575
M = int(M)
1576
if M < 2:
1577
raise ValueError("M (=%s) must be at least 2" % M)
1578
1579
base_ring = Q.base_ring()
1580
1581
# Expand out frobenius of dx/y and x dx/y.
1582
# (You can substitute frobenius_expansion_by_series here, that will work
1583
# as well. See its docstring for some performance notes.)
1584
F0, F1, offset = frobenius_expansion_by_newton(Q, p, M)
1585
#F0, F1, offset = frobenius_expansion_by_series(Q, p, M)
1586
1587
if compute_exact_forms:
1588
# 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/y
1589
F0_coeffs = transpose_list(F0.coeffs())
1590
F0_reduced, f_0 = reduce_all(Q, p, F0_coeffs, offset, True)
1591
1592
F1_coeffs = transpose_list(F1.coeffs())
1593
F1_reduced, f_1 = reduce_all(Q, p, F1_coeffs, offset, True)
1594
1595
elif M == 2:
1596
# This implies that only one digit of precision is valid, so we only need
1597
# to reduce the second column. Also, the trace doesn't help at all.
1598
1599
F0_reduced = [base_ring(0), base_ring(0)]
1600
1601
F1_coeffs = transpose_list(F1.coeffs())
1602
F1_reduced = reduce_all(Q, p, F1_coeffs, offset)
1603
1604
elif trace is None:
1605
# No trace provided, just reduce F(dx/y) and F(x dx/y) separately.
1606
1607
F0_coeffs = transpose_list(F0.coeffs())
1608
F0_reduced = reduce_all(Q, p, F0_coeffs, offset)
1609
1610
F1_coeffs = transpose_list(F1.coeffs())
1611
F1_reduced = reduce_all(Q, p, F1_coeffs, offset)
1612
1613
else:
1614
# Trace has been provided.
1615
1616
# In most cases this can be used to quickly compute F(dx/y) from
1617
# F(x dx/y). However, if we're unlucky, the (dx/y)-component of
1618
# F(x dx/y) (i.e. the top-right corner of the matrix) may be divisible
1619
# by p, in which case there isn't enough information to get the
1620
# (x dx/y)-component of F(dx/y) to the desired precision. When this
1621
# happens, it turns out that F((x+1) dx/y) always *does* give enough
1622
# information (together with the trace) to get both columns to the
1623
# desired precision.
1624
1625
# First however we need a quick way of telling whether the top-right
1626
# corner is divisible by p, i.e. we want to compute the second column
1627
# of the matrix mod p. We could do this by just running the entire
1628
# algorithm with M = 2 (which assures precision 1). Luckily, we've
1629
# already done most of the work by computing F1 to high precision; so
1630
# all we need to do is extract the coefficients that would correspond
1631
# to the first term of the series, and run the reduction on them.
1632
1633
# todo: actually we only need to do this reduction step mod p^2, not
1634
# mod p^M, which is what the code currently does. If the base ring
1635
# is Integers(p^M), then it is easy. Otherwise it is tricky to construct
1636
# the right ring, I don't know how to do it.
1637
1638
F1_coeffs = transpose_list(F1.coeffs())
1639
F1_modp_coeffs = F1_coeffs[int((M-2)*p):]
1640
# make a copy, because reduce_all will destroy the coefficients:
1641
F1_modp_coeffs = [[cell for cell in row] for row in F1_modp_coeffs]
1642
F1_modp_offset = offset - (M-2)*p
1643
F1_modp_reduced = reduce_all(Q, p, F1_modp_coeffs, F1_modp_offset)
1644
1645
if F1_modp_reduced[0].is_unit():
1646
# If the first entry is invertible mod p, then F(x dx/y) is sufficient
1647
# to get the whole matrix.
1648
1649
F1_reduced = reduce_all(Q, p, F1_coeffs, offset)
1650
1651
F0_reduced = [base_ring(trace) - F1_reduced[1], None]
1652
# using that the determinant is p:
1653
F0_reduced[1] = (F0_reduced[0] * F1_reduced[1] - base_ring(p)) \
1654
/ F1_reduced[0]
1655
1656
else:
1657
# If the first entry is zero mod p, then F((x+1) dx/y) will be sufficient
1658
# to get the whole matrix. (Here we are using the fact that the second
1659
# entry *cannot* be zero mod p. This is guaranteed by some results in
1660
# section 3.2 of ``Computation of p-adic Heights and Log Convergence''
1661
# by Mazur, Stein, Tate. But let's quickly check it anyway :-))
1662
msg = "The second entry in the second column "
1663
msg += "should be invertible mod p!"
1664
assert F1_modp_reduced[1].is_unit(), msg
1665
1666
G0_coeffs = transpose_list((F0 + F1).coeffs())
1667
G0_reduced = reduce_all(Q, p, G0_coeffs, offset)
1668
1669
# Now G0_reduced expresses F((x+1) dx/y) in terms of dx/y and x dx/y.
1670
# Re-express this in terms of (x+1) dx/y and x dx/y.
1671
H0_reduced = [G0_reduced[0], G0_reduced[1] - G0_reduced[0]]
1672
1673
# The thing we're about to divide by better be a unit.
1674
msg = "The second entry in this column "
1675
msg += "should be invertible mod p!"
1676
assert H0_reduced[1].is_unit(), msg
1677
1678
# Figure out the second column using the trace...
1679
H1_reduced = [None, base_ring(trace) - H0_reduced[0]]
1680
# ... and using that the determinant is p:
1681
H1_reduced[0] = (H0_reduced[0] * H1_reduced[1] - base_ring(p)) \
1682
/ H0_reduced[1]
1683
1684
# Finally, change back to the usual basis (dx/y, x dx/y)
1685
F1_reduced = [H1_reduced[0],
1686
H1_reduced[0] + H1_reduced[1]]
1687
F0_reduced = [H0_reduced[0] - F1_reduced[0],
1688
H0_reduced[0] + H0_reduced[1] - F1_reduced[1]]
1689
1690
# One more sanity check: our final result should be congruent mod p
1691
# to the approximation we used earlier.
1692
msg = "The output matrix is not congruent mod p "
1693
msg += "to the approximation found earlier!"
1694
assert not (
1695
(F1_reduced[0] - F1_modp_reduced[0]).is_unit() or
1696
(F1_reduced[1] - F1_modp_reduced[1]).is_unit() or
1697
F0_reduced[0].is_unit() or F0_reduced[1].is_unit()), msg
1698
1699
if compute_exact_forms:
1700
return matrix(base_ring, 2, 2, [F0_reduced[0], F1_reduced[0],
1701
F0_reduced[1], F1_reduced[1]]), f_0, f_1
1702
else:
1703
return matrix(base_ring, 2, 2, [F0_reduced[0], F1_reduced[0],
1704
F0_reduced[1], F1_reduced[1]])
1705
1706
1707
#*****************************************************************************
1708
# This is a generalization of the above functionality for hyperelliptic curves.
1709
#
1710
# THIS IS A WORK IN PROGRESS.
1711
#
1712
# I tried to embed must stuff into the rings themselves rather than
1713
# just extract and manipulate lists of coefficients. Hence the implementations
1714
# below are much less optimized, so are much slower, but should hopefully be
1715
# easier to follow. (E.g. one can print/make sense of intermediate results.)
1716
#
1717
# AUTHOR:
1718
# -- Robert Bradshaw (2007-04)
1719
#
1720
#*****************************************************************************
1721
1722
import weakref
1723
1724
from sage.schemes.hyperelliptic_curves.all import is_HyperellipticCurve, HyperellipticCurve
1725
from sage.rings.padics.all import pAdicField
1726
from sage.rings.all import QQ
1727
1728
from sage.rings.laurent_series_ring import is_LaurentSeriesRing
1729
from sage.rings.integral_domain import is_IntegralDomain
1730
1731
from sage.modules.all import FreeModule, is_FreeModuleElement
1732
1733
from sage.misc.profiler import Profiler
1734
from sage.misc.misc import repr_lincomb
1735
1736
1737
def matrix_of_frobenius_hyperelliptic(Q, p=None, prec=None, M=None):
1738
"""
1739
Computes the matrix of Frobenius on Monsky-Washnitzer cohomology,
1740
with respect to the basis `(dx/2y, x dx/2y, ...x^{d-2} dx/2y)`, where
1741
`d` is the degree of `Q`.
1742
1743
INPUT:
1744
1745
- ``Q`` -- monic polynomial `Q(x)`
1746
1747
- ``p`` -- prime `\geq 5` for which `E` has good reduction
1748
1749
- ``prec`` -- (optional) `p`-adic precision of the coefficient ring
1750
1751
- ``M`` -- (optional) adjusted `p`-adic precision of the coefficint ring
1752
1753
OUTPUT:
1754
1755
`(d-1)` x `(d-1)` matrix `M` of Frobenius on Monsky-Washnitzer cohomology,
1756
and list of differentials \{f_i \} such that
1757
1758
.. math::
1759
1760
\phi^* (x^i dx/2y) = df_i + M[i]*vec(dx/2y, ..., x^{d-2} dx/2y)
1761
1762
EXAMPLES::
1763
1764
sage: p = 5
1765
sage: prec = 3
1766
sage: R.<x> = QQ['x']
1767
sage: A,f = monsky_washnitzer.matrix_of_frobenius_hyperelliptic(x^5 - 2*x + 3, p, prec)
1768
sage: A
1769
[ 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)]
1770
[ 3*5 + 5^2 + O(5^3) 3*5 + O(5^3) 4*5 + O(5^3) 2 + 5^2 + O(5^3)]
1771
[ 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)]
1772
[ 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)]
1773
1774
"""
1775
prof = Profiler()
1776
prof("setup")
1777
if p is None:
1778
try:
1779
K = Q.base_ring()
1780
p = K.prime()
1781
prec = K.precision_cap()
1782
except AttributeError:
1783
raise ValueError("p and prec must be specified if Q is not "
1784
"defined over a p-adic ring")
1785
if M is None:
1786
M = adjusted_prec(p, prec)
1787
extra_prec_ring = Integers(p**M)
1788
# extra_prec_ring = pAdicField(p, M) # SLOW!
1789
1790
real_prec_ring = pAdicField(p, prec) # pAdicField(p, prec) # To capped absolute?
1791
S = SpecialHyperellipticQuotientRing(Q, extra_prec_ring, True)
1792
MW = S.monsky_washnitzer()
1793
prof("frob basis elements")
1794
F = MW.frob_basis_elements(M, p)
1795
1796
prof("rationalize")
1797
# do reduction over Q in case we have non-integral entries (and it is so much faster than padics)
1798
rational_S = S.change_ring(QQ)
1799
# this is a hack until pAdics are fast
1800
# (They are in the latest development bundle, but its not standard and I'd need to merge.
1801
# (it will periodically cast into this ring to reduce coefficient size)
1802
rational_S._prec_cap = p**M
1803
rational_S._p = p
1804
# S._p = p
1805
# rational_S(F[0]).reduce_fast()
1806
# prof("reduce others")
1807
1808
# rational_S = S.change_ring(pAdicField(p, M))
1809
F = [rational_S(F_i) for F_i in F]
1810
1811
prof("reduce")
1812
reduced = [F_i.reduce_fast(True) for F_i in F]
1813
# reduced = [F_i.reduce() for F_i in F]
1814
1815
#print reduced[0][0].diff() - F[0]
1816
1817
# but the coeffs are WAY more precision than they need to be
1818
# print reduced[0][1]
1819
1820
prof("make matrix")
1821
# now take care of precision capping
1822
M = matrix(real_prec_ring, [a for f, a in reduced])
1823
for i in range(M.ncols()):
1824
for j in range(M.nrows()):
1825
M[i, j] = M[i, j].add_bigoh(prec)
1826
# print prof
1827
return M.transpose(), [f for f, a in reduced]
1828
1829
1830
# For uniqueness (as many of the non-trivial calculations are cached along the way).
1831
1832
_special_ring_cache = {}
1833
_mw_cache = {}
1834
1835
1836
def SpecialHyperellipticQuotientRing(*args):
1837
if args in _special_ring_cache:
1838
R = _special_ring_cache[args]()
1839
if R is not None:
1840
return R
1841
R = SpecialHyperellipticQuotientRing_class(*args)
1842
_special_ring_cache[args] = weakref.ref(R)
1843
return R
1844
1845
1846
def MonskyWashnitzerDifferentialRing(base_ring):
1847
if base_ring in _mw_cache:
1848
R = _mw_cache[base_ring]()
1849
if R is not None:
1850
return R
1851
1852
R = MonskyWashnitzerDifferentialRing_class(base_ring)
1853
_mw_cache[base_ring] = weakref.ref(R)
1854
return R
1855
1856
1857
class SpecialHyperellipticQuotientRing_class(CommutativeAlgebra):
1858
1859
_p = None
1860
1861
def __init__(self, Q, R=None, invert_y=True):
1862
if R is None:
1863
R = Q.base_ring()
1864
1865
# Trac ticket #9138: CommutativeAlgebra.__init__ must not be
1866
# done so early. It tries to register a coercion, but that
1867
# requires the hash bein available. But the hash, in its
1868
# default implementation, relies on the string representation,
1869
# which is not available at this point.
1870
#CommutativeAlgebra.__init__(self, R) # moved to below.
1871
1872
x = PolynomialRing(R, 'xx').gen(0)
1873
if is_EllipticCurve(Q):
1874
E = Q
1875
if E.a1() != 0 or E.a2() != 0:
1876
raise NotImplementedError("Curve must be in Weierstrass "
1877
"normal form.")
1878
Q = -E.change_ring(R).defining_polynomial()(x, 0, 1)
1879
self._curve = E
1880
1881
elif is_HyperellipticCurve(Q):
1882
C = Q
1883
if C.hyperelliptic_polynomials()[1] != 0:
1884
raise NotImplementedError("Curve must be of form y^2 = Q(x).")
1885
Q = C.hyperelliptic_polynomials()[0].change_ring(R)
1886
self._curve = C
1887
1888
if is_Polynomial(Q):
1889
self._Q = Q.change_ring(R)
1890
self._coeffs = self._Q.coeffs()
1891
if self._coeffs.pop() != 1:
1892
raise NotImplementedError("Polynomial must be monic.")
1893
if not hasattr(self, '_curve'):
1894
if self._Q.degree() == 3:
1895
ainvs = [0, self._Q[2], 0, self._Q[1], self._Q[0]]
1896
self._curve = EllipticCurve(ainvs)
1897
else:
1898
self._curve = HyperellipticCurve(self._Q)
1899
1900
else:
1901
raise NotImplementedError("Must be an elliptic curve or polynomial "
1902
"Q for y^2 = Q(x)\n(Got element of %s)" % Q.parent())
1903
1904
self._n = int(Q.degree())
1905
self._series_ring = (LaurentSeriesRing if invert_y else PolynomialRing)(R, 'y')
1906
self._series_ring_y = self._series_ring.gen(0)
1907
self._series_ring_0 = self._series_ring(0)
1908
1909
# Trac ticket #9138: Initialise the commutative algebra here!
1910
# Below, we do self(self._poly_ring.gen(0)), which requires
1911
# the initialisation being finished.
1912
CommutativeAlgebra.__init__(self, R)
1913
1914
self._poly_ring = PolynomialRing(self._series_ring, 'x')
1915
1916
self._x = self(self._poly_ring.gen(0))
1917
self._y = self(self._series_ring.gen(0))
1918
1919
self._Q_coeffs = Q.change_ring(self._series_ring).list()
1920
self._dQ = Q.derivative().change_ring(self)(self._x)
1921
self._monsky_washnitzer = MonskyWashnitzerDifferentialRing(self)
1922
1923
self._monomial_diffs = {}
1924
self._monomial_diff_coeffs = {}
1925
1926
def _repr_(self):
1927
"""
1928
String representation
1929
1930
EXAMPLES::
1931
1932
sage: R.<x> = QQ['x']
1933
sage: E = HyperellipticCurve(x^5-3*x+1)
1934
sage: x,y = E.monsky_washnitzer_gens()
1935
sage: x.parent() # indirect doctest
1936
SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 3*x + 1) over Rational Field
1937
"""
1938
y_inverse = ",y^-1" if is_LaurentSeriesRing(self._series_ring) else ""
1939
return "SpecialHyperellipticQuotientRing K[x,y%s] / (y^2 = %s) over %s" % (y_inverse, self._Q, self.base_ring())
1940
1941
def base_extend(self, R):
1942
"""
1943
Return the base extension of ``self`` to the ring ``R`` if possible.
1944
1945
EXAMPLES::
1946
1947
sage: R.<x> = QQ['x']
1948
sage: E = HyperellipticCurve(x^5-3*x+1)
1949
sage: x,y = E.monsky_washnitzer_gens()
1950
sage: x.parent().base_extend(UniversalCyclotomicField())
1951
SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 3*x + 1) over Universal Cyclotomic Field
1952
sage: x.parent().base_extend(ZZ)
1953
Traceback (most recent call last):
1954
...
1955
TypeError: no such base extension
1956
"""
1957
if R.has_coerce_map_from(self.base_ring()):
1958
return self.change_ring(R)
1959
else:
1960
raise TypeError("no such base extension")
1961
1962
def change_ring(self, R):
1963
"""
1964
Return the analog of ``self`` over the ring ``R``
1965
1966
EXAMPLES::
1967
1968
sage: R.<x> = QQ['x']
1969
sage: E = HyperellipticCurve(x^5-3*x+1)
1970
sage: x,y = E.monsky_washnitzer_gens()
1971
sage: x.parent().change_ring(ZZ)
1972
SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 3*x + 1) over Integer Ring
1973
"""
1974
return SpecialHyperellipticQuotientRing(self._Q, R, is_LaurentSeriesRing(self._series_ring))
1975
1976
def __call__(self, val, offset=0, check=True):
1977
if isinstance(val, SpecialHyperellipticQuotientElement) and val.parent() is self:
1978
if offset == 0:
1979
return val
1980
else:
1981
return val << offset
1982
elif isinstance(val, MonskyWashnitzerDifferential):
1983
return self._monsky_washnitzer(val)
1984
return SpecialHyperellipticQuotientElement(self, val, offset, check)
1985
1986
def gens(self):
1987
"""
1988
Return the generators of ``self``
1989
1990
EXAMPLES::
1991
1992
sage: R.<x> = QQ['x']
1993
sage: E = HyperellipticCurve(x^5-3*x+1)
1994
sage: x,y = E.monsky_washnitzer_gens()
1995
sage: x.parent().gens()
1996
(x, y*1)
1997
"""
1998
return self._x, self._y
1999
2000
def x(self):
2001
r"""
2002
Return the generator `x` of ``self``
2003
2004
EXAMPLES::
2005
2006
sage: R.<x> = QQ['x']
2007
sage: E = HyperellipticCurve(x^5-3*x+1)
2008
sage: x,y = E.monsky_washnitzer_gens()
2009
sage: x.parent().x()
2010
x
2011
"""
2012
return self._x
2013
2014
def y(self):
2015
r"""
2016
Return the generator `y` of ``self``
2017
2018
EXAMPLES::
2019
2020
sage: R.<x> = QQ['x']
2021
sage: E = HyperellipticCurve(x^5-3*x+1)
2022
sage: x,y = E.monsky_washnitzer_gens()
2023
sage: x.parent().y()
2024
y*1
2025
"""
2026
return self._y
2027
2028
def monomial(self, i, j, b=None):
2029
"""
2030
Returns `b y^j x^i`, computed quickly.
2031
2032
EXAMPLES::
2033
2034
sage: R.<x> = QQ['x']
2035
sage: E = HyperellipticCurve(x^5-3*x+1)
2036
sage: x,y = E.monsky_washnitzer_gens()
2037
sage: x.parent().monomial(4,5)
2038
y^5*x^4
2039
"""
2040
i = int(i)
2041
j = int(j)
2042
2043
if 0 < i and i < self._n:
2044
if b is None:
2045
by_to_j = self._series_ring_y << (j-1)
2046
else:
2047
by_to_j = self._series_ring(b) << j
2048
v = [self._series_ring_0] * self._n
2049
v[i] = by_to_j
2050
return self(v)
2051
else:
2052
return (self._x ** i) << j if b is None else self.base_ring()(b) * (self._x ** i) << j
2053
2054
def monomial_diff_coeffs(self, i, j):
2055
r"""
2056
The key here is that the formula for `d(x^iy^j)` is messy
2057
in terms of `i`, but varies nicely with `j`.
2058
2059
.. math::
2060
2061
d(x^iy^j) = y^{j-1} (2ix^{i-1}y^2 + j (A_i(x) + B_i(x)y^2)) \frac{dx}{2y}
2062
2063
2064
Where `A,B` have degree at most `n-1` for each
2065
`i`. Pre-compute `A_i, B_i` for each `i`
2066
the "hard" way, and the rest are easy.
2067
"""
2068
try:
2069
return self._monomial_diff_coeffs[i, j]
2070
except KeyError:
2071
pass
2072
if i < self._n:
2073
try:
2074
A, B, two_i_x_to_i = self._precomputed_diff_coeffs[i]
2075
except AttributeError:
2076
self._precomputed_diff_coeffs = self._precompute_monomial_diffs()
2077
A, B, two_i_x_to_i = self._precomputed_diff_coeffs[i]
2078
if i == 0:
2079
return j*A, j*B
2080
else:
2081
return j*A, j*B + two_i_x_to_i
2082
else:
2083
dg = self.monomial(i, j).diff()
2084
coeffs = [dg.extract_pow_y(j-1), dg.extract_pow_y(j+1)]
2085
self._monomial_diff_coeffs[i, j] = coeffs
2086
return coeffs
2087
2088
def monomial_diff_coeffs_matrices(self):
2089
self.monomial_diff_coeffs(0, 0) # precompute stuff
2090
R = self.base_ring()
2091
mat_1 = matrix(R, self._n, self._n)
2092
mat_2 = matrix(R, self._n, self._n)
2093
for i in range(self._n):
2094
mat_1[i] = self._precomputed_diff_coeffs[i][1]
2095
mat_2[i] = self._precomputed_diff_coeffs[i][2]
2096
return mat_1.transpose(), mat_2.transpose()
2097
2098
def _precompute_monomial_diffs(self):
2099
x, y = self.gens()
2100
R = self.base_ring()
2101
V = FreeModule(R, self.degree())
2102
As = []
2103
for i in range(self.degree()):
2104
dg = self.monomial(i, 1).diff()
2105
two_i_x_to_i = R(2*i) * x**(i-1) * y*y if i > 0 else self(0)
2106
A = dg - self._monsky_washnitzer(two_i_x_to_i)
2107
As.append((V(A.extract_pow_y(0)), V(A.extract_pow_y(2)), V(two_i_x_to_i.extract_pow_y(2))))
2108
return As
2109
2110
def Q(self):
2111
"""
2112
Return the defining polynomial of the underlying hyperelliptic curve.
2113
2114
EXAMPLES::
2115
2116
sage: R.<x> = QQ['x']
2117
sage: E = HyperellipticCurve(x^5-2*x+1)
2118
sage: x,y = E.monsky_washnitzer_gens()
2119
sage: x.parent().Q()
2120
x^5 - 2*x + 1
2121
"""
2122
return self._Q
2123
2124
def curve(self):
2125
"""
2126
Return the underlying hyperelliptic curve.
2127
2128
EXAMPLES::
2129
2130
sage: R.<x> = QQ['x']
2131
sage: E = HyperellipticCurve(x^5-3*x+1)
2132
sage: x,y = E.monsky_washnitzer_gens()
2133
sage: x.parent().curve()
2134
Hyperelliptic Curve over Rational Field defined by y^2 = x^5 - 3*x + 1
2135
"""
2136
return self._curve
2137
2138
def degree(self):
2139
"""
2140
Return the degree of the underlying hyperelliptic curve.
2141
2142
EXAMPLES::
2143
2144
sage: R.<x> = QQ['x']
2145
sage: E = HyperellipticCurve(x^5-3*x+1)
2146
sage: x,y = E.monsky_washnitzer_gens()
2147
sage: x.parent().degree()
2148
5
2149
"""
2150
return self._n
2151
2152
def prime(self):
2153
return self._p
2154
2155
def monsky_washnitzer(self):
2156
return self._monsky_washnitzer
2157
2158
def is_field(self, proof=True):
2159
"""
2160
Return False as ``self`` is not a field.
2161
2162
EXAMPLES::
2163
2164
sage: R.<x> = QQ['x']
2165
sage: E = HyperellipticCurve(x^5-3*x+1)
2166
sage: x,y = E.monsky_washnitzer_gens()
2167
sage: x.parent().is_field()
2168
False
2169
"""
2170
return False
2171
2172
2173
class SpecialHyperellipticQuotientElement(CommutativeAlgebraElement):
2174
2175
def __init__(self, parent, val=0, offset=0, check=True):
2176
"""
2177
Elements in the Hyperelliptic quotient ring
2178
2179
EXAMPLES::
2180
2181
sage: R.<x> = QQ['x']
2182
sage: E = HyperellipticCurve(x^5-36*x+1)
2183
sage: x,y = E.monsky_washnitzer_gens()
2184
sage: MW = x.parent()
2185
sage: MW(x+x**2+y-77) # indirect doctest
2186
-(77-y)*1 + x + x^2
2187
"""
2188
CommutativeAlgebraElement.__init__(self, parent)
2189
if not check:
2190
self._f = parent._poly_ring(val, check=False)
2191
return
2192
if isinstance(val, SpecialHyperellipticQuotientElement):
2193
R = parent.base_ring()
2194
self._f = parent._poly_ring([a.change_ring(R) for a in val._f])
2195
return
2196
if isinstance(val, tuple):
2197
val, offset = val
2198
if isinstance(val, list) and len(val) > 0 and is_FreeModuleElement(val[0]):
2199
val = transpose_list(val)
2200
self._f = parent._poly_ring(val)
2201
if offset != 0:
2202
self._f = self._f.parent()([a << offset for a in self._f], check=False)
2203
2204
def __cmp__(self, other):
2205
"""
2206
Compares the elements
2207
2208
EXAMPLES::
2209
2210
sage: R.<x> = QQ['x']
2211
sage: E = HyperellipticCurve(x^5-36*x+1)
2212
sage: x,y = E.monsky_washnitzer_gens()
2213
sage: x.__cmp__(x)
2214
0
2215
sage: x.__cmp__(y)
2216
1
2217
"""
2218
return cmp(self._f, other._f)
2219
2220
def change_ring(self, R):
2221
"""
2222
Return the same element after changing the base ring to R
2223
2224
EXAMPLES::
2225
2226
sage: R.<x> = QQ['x']
2227
sage: E = HyperellipticCurve(x^5-36*x+1)
2228
sage: x,y = E.monsky_washnitzer_gens()
2229
sage: MW = x.parent()
2230
sage: z = MW(x+x**2+y-77)
2231
sage: z.change_ring(AA).parent()
2232
SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 36*x + 1) over Algebraic Real Field
2233
"""
2234
return self.parent().change_ring(R)(self)
2235
2236
def __call__(self, *x):
2237
"""
2238
Evaluate ``self`` at given arguments
2239
2240
EXAMPLES::
2241
2242
sage: R.<x> = QQ['x']
2243
sage: E = HyperellipticCurve(x^5-36*x+1)
2244
sage: x,y = E.monsky_washnitzer_gens()
2245
sage: MW = x.parent()
2246
sage: z = MW(x+x**2+y-77); z
2247
-(77-y)*1 + x + x^2
2248
sage: z(66)
2249
4345 + y
2250
sage: z(5,4)
2251
-43
2252
"""
2253
return self._f(*x)
2254
2255
def __invert__(self):
2256
"""
2257
Return the inverse of the element
2258
2259
The general element in our ring is not invertible, but `y` may
2260
be. We do not want to pass to the fraction field.
2261
2262
EXAMPLES::
2263
2264
sage: R.<x> = QQ['x']
2265
sage: E = HyperellipticCurve(x^5-36*x+1)
2266
sage: x,y = E.monsky_washnitzer_gens()
2267
sage: MW = x.parent()
2268
sage: z = y**(-1) # indirect doctest
2269
sage: z.parent()
2270
SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 36*x + 1) over Rational Field
2271
2272
sage: z = (x+y)**(-1) # indirect doctest
2273
Traceback (most recent call last):
2274
...
2275
ZeroDivisionError: Element not invertible
2276
"""
2277
if self._f.degree() == 0 and self._f[0].is_unit():
2278
return SpecialHyperellipticQuotientElement(self.parent(), ~self._f[0])
2279
else:
2280
raise ZeroDivisionError("Element not invertible")
2281
2282
def __nonzero__(self):
2283
"""
2284
Return True iff ``self`` is not zero.
2285
2286
EXAMPLES::
2287
2288
sage: R.<x> = QQ['x']
2289
sage: E = HyperellipticCurve(x^5-3*x+1)
2290
sage: x,y = E.monsky_washnitzer_gens()
2291
sage: x.__nonzero__()
2292
True
2293
"""
2294
return not not self._f
2295
2296
def __eq__(self, other):
2297
"""
2298
Return True iff ``self`` is equal to ``other``
2299
2300
EXAMPLES::
2301
2302
sage: R.<x> = QQ['x']
2303
sage: E = HyperellipticCurve(x^5-3*x+1)
2304
sage: x,y = E.monsky_washnitzer_gens()
2305
sage: x == y # indirect doctest
2306
False
2307
"""
2308
if not isinstance(other, SpecialHyperellipticQuotientElement):
2309
other = self.parent()(other)
2310
return self._f == other._f
2311
2312
def _add_(self, other):
2313
"""
2314
Return the sum of two elements
2315
2316
EXAMPLES::
2317
2318
sage: R.<x> = QQ['x']
2319
sage: E = HyperellipticCurve(x^5-36*x+1)
2320
sage: x,y = E.monsky_washnitzer_gens()
2321
sage: x+y
2322
y*1 + x
2323
"""
2324
return SpecialHyperellipticQuotientElement(self.parent(), self._f + other._f)
2325
2326
def _sub_(self, other):
2327
"""
2328
Return the difference of two elements
2329
2330
EXAMPLES::
2331
2332
sage: R.<x> = QQ['x']
2333
sage: E = HyperellipticCurve(x^5-36*x+1)
2334
sage: x,y = E.monsky_washnitzer_gens()
2335
sage: y-x
2336
y*1 - x
2337
"""
2338
return SpecialHyperellipticQuotientElement(self.parent(), self._f - other._f)
2339
2340
def _mul_(self, other):
2341
"""
2342
Return the product of two elements
2343
2344
EXAMPLES::
2345
2346
sage: R.<x> = QQ['x']
2347
sage: E = HyperellipticCurve(x^5-36*x+1)
2348
sage: x,y = E.monsky_washnitzer_gens()
2349
sage: y*x
2350
y*x
2351
"""
2352
# over laurent series, addition and subtraction can be
2353
# expensive, and the degree of this poly is small enough that
2354
# Karatsuba actually hurts significantly in some cases
2355
if self._f[0].valuation() + other._f[0].valuation() > -200:
2356
prod = self._f._mul_generic(other._f)
2357
else:
2358
prod = self._f * other._f
2359
v = prod.list()
2360
parent = self.parent()
2361
Q_coeffs = parent._Q_coeffs
2362
n = len(Q_coeffs) - 1
2363
y2 = self.parent()._series_ring_y << 1
2364
for i in range(len(v)-1, n-1, -1):
2365
for j in range(n):
2366
v[i-n+j] -= Q_coeffs[j] * v[i]
2367
v[i-n] += y2 * v[i]
2368
return SpecialHyperellipticQuotientElement(parent, v[0:n])
2369
2370
def _rmul_(self, c):
2371
coeffs = self._f.list(copy=False)
2372
return self.parent()([c*a for a in coeffs], check=False)
2373
2374
def _lmul_(self, c):
2375
coeffs = self._f.list(copy=False)
2376
return self.parent()([a*c for a in coeffs], check=False)
2377
2378
def __lshift__(self, k):
2379
coeffs = self._f.list(copy=False)
2380
return self.parent()([a << k for a in coeffs], check=False)
2381
2382
def __rshift__(self, k):
2383
coeffs = self._f.list(copy=False)
2384
return self.parent()([a >> k for a in coeffs], check=False)
2385
2386
def truncate_neg(self, n):
2387
"""
2388
Return ``self`` minus its terms of degree less than `n` wrt `y`.
2389
2390
EXAMPLES::
2391
2392
sage: R.<x> = QQ['x']
2393
sage: E = HyperellipticCurve(x^5-3*x+1)
2394
sage: x,y = E.monsky_washnitzer_gens()
2395
sage: (x+3*y+7*x*2*y**4).truncate_neg(1)
2396
3*y*1 + 14*y^4*x
2397
"""
2398
coeffs = self._f.list(copy=False)
2399
return self.parent()([a.truncate_neg(n) for a in coeffs], check=False)
2400
2401
def _repr_(self):
2402
"""
2403
Return a string representation of ``self``.
2404
2405
EXAMPLES::
2406
2407
sage: R.<x> = QQ['x']
2408
sage: E = HyperellipticCurve(x^5-3*x+1)
2409
sage: x,y = E.monsky_washnitzer_gens()
2410
sage: (x+3*y)._repr_()
2411
'3*y*1 + x'
2412
"""
2413
x = PolynomialRing(QQ, 'x').gen(0)
2414
coeffs = self._f.list()
2415
return repr_lincomb([(x**i, coeffs[i]) for i in range(len(coeffs))])
2416
2417
def _latex_(self):
2418
"""
2419
Return a LateX string for ``self``.
2420
2421
EXAMPLES::
2422
2423
sage: R.<x> = QQ['x']
2424
sage: E = HyperellipticCurve(x^5-3*x+1)
2425
sage: x,y = E.monsky_washnitzer_gens()
2426
sage: (x+3*y)._latex_()
2427
'3y1 + x'
2428
"""
2429
x = PolynomialRing(QQ, 'x').gen(0)
2430
coeffs = self._f.list()
2431
return repr_lincomb([(x**i, coeffs[i]) for i in range(len(coeffs))], is_latex=True)
2432
2433
def diff(self):
2434
"""
2435
Return the differential of ``self``
2436
2437
EXAMPLES::
2438
2439
sage: R.<x> = QQ['x']
2440
sage: E = HyperellipticCurve(x^5-3*x+1)
2441
sage: x,y = E.monsky_washnitzer_gens()
2442
sage: (x+3*y).diff()
2443
(-(9-2*y)*1 + 15*x^4) dx/2y
2444
"""
2445
# try:
2446
# return self._diff_x
2447
# except AttributeError:
2448
# pass
2449
2450
# d(self) = A dx + B dy
2451
# = (2y A + BQ') dx/2y
2452
parent = self.parent()
2453
R = parent.base_ring()
2454
x, y = parent.gens()
2455
v = self._f.list()
2456
n = len(v)
2457
A = parent([R(i) * v[i] for i in range(1, n)])
2458
B = parent([a.derivative() for a in v])
2459
dQ = parent._dQ
2460
return parent._monsky_washnitzer((R(2) * A << 1) + dQ * B)
2461
# self._diff = self.parent()._monsky_washnitzer(two_y * A + dQ * B)
2462
# return self._diff
2463
2464
def extract_pow_y(self, k):
2465
r"""
2466
Return the coefficients of `y^k` in ``self`` as a list
2467
2468
EXAMPLES::
2469
2470
sage: R.<x> = QQ['x']
2471
sage: E = HyperellipticCurve(x^5-3*x+1)
2472
sage: x,y = E.monsky_washnitzer_gens()
2473
sage: (x+3*y+9*x*y).extract_pow_y(1)
2474
[3, 9, 0, 0, 0]
2475
"""
2476
v = [a[k] for a in self._f.list()]
2477
while len(v) < self.parent()._n:
2478
v.append(0)
2479
return v
2480
2481
def min_pow_y(self):
2482
"""
2483
Return the minimal degree of ``self`` w.r.t. y
2484
2485
EXAMPLES::
2486
2487
sage: R.<x> = QQ['x']
2488
sage: E = HyperellipticCurve(x^5-3*x+1)
2489
sage: x,y = E.monsky_washnitzer_gens()
2490
sage: (x+3*y).min_pow_y()
2491
0
2492
"""
2493
if self._f.degree() == -1:
2494
return 0
2495
return min([a.valuation() for a in self._f.list()])
2496
2497
def max_pow_y(self):
2498
"""
2499
Return the maximal degree of ``self`` w.r.t. y
2500
2501
EXAMPLES::
2502
2503
sage: R.<x> = QQ['x']
2504
sage: E = HyperellipticCurve(x^5-3*x+1)
2505
sage: x,y = E.monsky_washnitzer_gens()
2506
sage: (x+3*y).max_pow_y()
2507
1
2508
"""
2509
if self._f.degree() == -1:
2510
return 0
2511
return max([a.degree() for a in self._f.list()])
2512
2513
def coeffs(self, R=None):
2514
"""
2515
Returns the raw coefficients of this element.
2516
2517
INPUT:
2518
2519
- ``R`` -- an (optional) base-ring in which to cast the coefficients
2520
2521
OUTPUT:
2522
2523
- ``coeffs`` -- a list of coefficients of powers of `x` for each power
2524
of `y`
2525
2526
- ``n`` -- an offset indicating the power of `y` of the first list
2527
element
2528
2529
EXAMPLES::
2530
2531
sage: R.<x> = QQ['x']
2532
sage: E = HyperellipticCurve(x^5-3*x+1)
2533
sage: x,y = E.monsky_washnitzer_gens()
2534
sage: x.coeffs()
2535
([(0, 1, 0, 0, 0)], 0)
2536
sage: y.coeffs()
2537
([(0, 0, 0, 0, 0), (1, 0, 0, 0, 0)], 0)
2538
2539
sage: a = sum(n*x^n for n in range(5)); a
2540
x + 2*x^2 + 3*x^3 + 4*x^4
2541
sage: a.coeffs()
2542
([(0, 1, 2, 3, 4)], 0)
2543
sage: a.coeffs(Qp(7))
2544
([(0, 1 + O(7^20), 2 + O(7^20), 3 + O(7^20), 4 + O(7^20))], 0)
2545
sage: (a*y).coeffs()
2546
([(0, 0, 0, 0, 0), (0, 1, 2, 3, 4)], 0)
2547
sage: (a*y^-2).coeffs()
2548
([(0, 1, 2, 3, 4), (0, 0, 0, 0, 0), (0, 0, 0, 0, 0)], -2)
2549
2550
Note that the coefficient list is transposed compared to how they
2551
are stored and printed::
2552
2553
sage: a*y^-2
2554
(y^-2)*x + (2*y^-2)*x^2 + (3*y^-2)*x^3 + (4*y^-2)*x^4
2555
2556
A more complicated example::
2557
2558
sage: a = x^20*y^-3 - x^11*y^2; a
2559
(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^4
2560
sage: raw, offset = a.coeffs()
2561
sage: a.min_pow_y()
2562
-3
2563
sage: offset
2564
-3
2565
sage: raw
2566
[(1, -12, 54, -108, 81),
2567
(0, 0, 0, 0, 0),
2568
(-4, 36, -108, 108, 0),
2569
(0, 0, 0, 0, 0),
2570
(6, -36, 54, 0, 0),
2571
(0, -1, 6, -9, 0),
2572
(-4, 12, 0, 0, 0),
2573
(0, 2, -6, 0, 0),
2574
(1, 0, 0, 0, 0),
2575
(0, -1, 0, 0, 0)]
2576
sage: sum(c * x^i * y^(j+offset) for j, L in enumerate(raw) for i, c in enumerate(L)) == a
2577
True
2578
2579
Can also be used to construct elements::
2580
2581
sage: a.parent()(raw, offset) == a
2582
True
2583
"""
2584
zero = self.base_ring()(0) if R is None else R(0)
2585
y_offset = min(self.min_pow_y(), 0)
2586
y_degree = max(self.max_pow_y(), 0)
2587
coeffs = []
2588
n = y_degree - y_offset + 1
2589
for a in self._f.list():
2590
k = a.valuation()
2591
if k is Infinity: k = 0
2592
k -= y_offset
2593
z = a.list()
2594
coeffs.append([zero] * k + z + [zero]*(n - len(z) - k))
2595
while len(coeffs) < self.parent().degree():
2596
coeffs.append([zero] * n)
2597
V = FreeModule(self.base_ring() if R is None else R, self.parent().degree())
2598
coeffs = transpose_list(coeffs)
2599
return [V(a) for a in coeffs], y_offset
2600
2601
2602
class MonskyWashnitzerDifferentialRing_class(Module):
2603
2604
def __init__(self, base_ring):
2605
r"""
2606
Class for the ring of Monsky--Washnitzer differentials over a given
2607
base ring.
2608
"""
2609
Module.__init__(self, base_ring)
2610
self._cache = {}
2611
2612
def invariant_differential(self):
2613
"""
2614
Returns `dx/2y` as an element of self.
2615
2616
EXAMPLES::
2617
2618
sage: R.<x> = QQ['x']
2619
sage: C = HyperellipticCurve(x^5-4*x+4)
2620
sage: MW = C.invariant_differential().parent()
2621
sage: MW.invariant_differential()
2622
1 dx/2y
2623
"""
2624
return self(1)
2625
2626
def __call__(self, val, offset=0):
2627
return MonskyWashnitzerDifferential(self, val, offset)
2628
2629
def base_extend(self, R):
2630
"""
2631
Return a new differential ring which is self base-extended to `R`
2632
2633
INPUT:
2634
2635
- ``R`` -- ring
2636
2637
OUTPUT:
2638
2639
Self, base-extended to `R`.
2640
2641
EXAMPLES::
2642
2643
sage: R.<x> = QQ['x']
2644
sage: C = HyperellipticCurve(x^5-4*x+4)
2645
sage: MW = C.invariant_differential().parent()
2646
sage: MW.base_ring()
2647
SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 4*x + 4) over Rational Field
2648
sage: MW.base_extend(Qp(5,5)).base_ring()
2649
SpecialHyperellipticQuotientRing 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 5
2650
"""
2651
return MonskyWashnitzerDifferentialRing(self.base_ring().base_extend(R))
2652
2653
def change_ring(self, R):
2654
"""
2655
Returns a new differential ring which is self with the coefficient
2656
ring changed to `R`.
2657
2658
INPUT:
2659
2660
- ``R`` -- ring of coefficients
2661
2662
OUTPUT:
2663
2664
Self, with the coefficient ring changed to `R`.
2665
2666
EXAMPLES::
2667
2668
sage: R.<x> = QQ['x']
2669
sage: C = HyperellipticCurve(x^5-4*x+4)
2670
sage: MW = C.invariant_differential().parent()
2671
sage: MW.base_ring()
2672
SpecialHyperellipticQuotientRing K[x,y,y^-1] / (y^2 = x^5 - 4*x + 4) over Rational Field
2673
sage: MW.change_ring(Qp(5,5)).base_ring()
2674
SpecialHyperellipticQuotientRing 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 5
2675
"""
2676
return MonskyWashnitzerDifferentialRing(self.base_ring().change_ring(R))
2677
2678
def degree(self):
2679
"""
2680
Returns the degree of `Q(x)`, where the model of the underlying
2681
hyperelliptic curve of self is given by `y^2 = Q(x)`.
2682
2683
EXAMPLES::
2684
2685
sage: R.<x> = QQ['x']
2686
sage: C = HyperellipticCurve(x^5-4*x+4)
2687
sage: MW = C.invariant_differential().parent()
2688
sage: MW.Q()
2689
x^5 - 4*x + 4
2690
sage: MW.degree()
2691
5
2692
"""
2693
return self.base_ring().degree()
2694
2695
def dimension(self):
2696
"""
2697
Returns the dimension of self.
2698
2699
EXAMPLES::
2700
2701
sage: R.<x> = QQ['x']
2702
sage: C = HyperellipticCurve(x^5-4*x+4)
2703
sage: K = Qp(7,5)
2704
sage: CK = C.change_ring(K)
2705
sage: MW = CK.invariant_differential().parent()
2706
sage: MW.dimension()
2707
4
2708
"""
2709
return self.base_ring().degree()-1
2710
2711
def Q(self):
2712
"""
2713
Returns `Q(x)` where the model of the underlying hyperelliptic curve
2714
of self is given by `y^2 = Q(x)`.
2715
2716
EXAMPLES::
2717
2718
sage: R.<x> = QQ['x']
2719
sage: C = HyperellipticCurve(x^5-4*x+4)
2720
sage: MW = C.invariant_differential().parent()
2721
sage: MW.Q()
2722
x^5 - 4*x + 4
2723
"""
2724
return self.base_ring().Q()
2725
2726
def x_to_p(self, p):
2727
"""
2728
Returns and caches `x^p`, reduced via the relations coming from the
2729
defining polynomial of the hyperelliptic curve.
2730
2731
EXAMPLES::
2732
2733
sage: R.<x> = QQ['x']
2734
sage: C = HyperellipticCurve(x^5-4*x+4)
2735
sage: MW = C.invariant_differential().parent()
2736
sage: MW.x_to_p(3)
2737
x^3
2738
sage: MW.x_to_p(5)
2739
-(4-y^2)*1 + 4*x
2740
sage: MW.x_to_p(101) is MW.x_to_p(101)
2741
True
2742
"""
2743
try:
2744
return self._cache["x_to_p", p]
2745
except KeyError:
2746
x_to_p = self.base_ring().x() ** p
2747
self._cache["x_to_p", p] = x_to_p
2748
return x_to_p
2749
2750
def frob_Q(self, p):
2751
"""
2752
Returns and caches `Q(x^p)`, which is used in computing the image of
2753
`y` under a `p`-power lift of Frobenius to `A^{\dagger}`.
2754
2755
EXAMPLES::
2756
2757
sage: R.<x> = QQ['x']
2758
sage: C = HyperellipticCurve(x^5-4*x+4)
2759
sage: MW = C.invariant_differential().parent()
2760
sage: MW.frob_Q(3)
2761
-(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^3
2762
sage: MW.Q()(MW.x_to_p(3))
2763
-(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^3
2764
sage: MW.frob_Q(11) is MW.frob_Q(11)
2765
True
2766
"""
2767
try:
2768
return self._cache["frobQ", p]
2769
except KeyError:
2770
x_to_p = self.x_to_p(p)
2771
frobQ = self.base_ring()._Q.change_ring(self.base_ring())(x_to_p)
2772
self._cache["frobQ", p] = frobQ
2773
return frobQ
2774
2775
def frob_invariant_differential(self, prec, p):
2776
r"""
2777
Kedlaya's algorithm allows us to calculate the action of Frobenius on
2778
the Monsky-Washnitzer cohomology. First we lift `\phi` to `A^{\dagger}`
2779
by setting
2780
2781
.. math::
2782
2783
\phi(x) = x^p
2784
2785
\phi(y) = y^p \sqrt{1 + \frac{Q(x^p) - Q(x)^p}{Q(x)^p}}.
2786
2787
Pulling back the differential `dx/2y`, we get
2788
2789
.. math::
2790
2791
\phi^*(dx/2y) = px^{p-1} y(\phi(y))^{-1} dx/2y
2792
= px^{p-1} y^{1-p} \sqrt{1+ \frac{Q(x^p) - Q(x)^p}{Q(x)^p}} dx/2y
2793
2794
Use Newton's method to calculate the square root.
2795
2796
EXAMPLES::
2797
2798
sage: R.<x> = QQ['x']
2799
sage: C = HyperellipticCurve(x^5-4*x+4)
2800
sage: prec = 2
2801
sage: p = 7
2802
sage: MW = C.invariant_differential().parent()
2803
sage: MW.frob_invariant_differential(prec,p)
2804
((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/2y
2805
"""
2806
prof = Profiler()
2807
prof("setup")
2808
# TODO, would it be useful to be able to take Frobenius of any element? Less efficient?
2809
x, y = self.base_ring().gens()
2810
prof("x_to_p")
2811
x_to_p_less_1 = x**(p-1)
2812
x_to_p = x*x_to_p_less_1
2813
2814
# cache for future use
2815
self._cache["x_to_p", p] = x_to_p
2816
2817
prof("frob_Q")
2818
a = self.frob_Q(p) >> 2*p # frobQ * y^{-2p}
2819
2820
prof("sqrt")
2821
2822
# Q = self.base_ring()._Q
2823
# three_halves = Q.parent().base_ring()(Rational((3,2)))
2824
# one_half = Q.parent().base_ring()(Rational((1,2)))
2825
three_halves = self.base_ring()._series_ring.base_ring()(Rational((3, 2)))
2826
one_half = self.base_ring()._series_ring.base_ring()(Rational((1, 2)))
2827
half_a = a._rmul_(one_half)
2828
2829
# We are solving for t = a^{-1/2} = (F_pQ y^{-p})^{-1/2}
2830
# Newton's method converges because we know the root is in the same residue class as 1.
2831
2832
# t = self.base_ring()(1)
2833
t = self.base_ring()(three_halves) - half_a
2834
# first iteration trivial, start with prec 2
2835
2836
for cur_prec in newton_method_sizes(prec)[2:]:
2837
# newton_method_sizes = [1, 2, ...]
2838
y_prec = -(2*cur_prec-1)*p+1
2839
# binomial expansion is $\sum p^{k+1} y^{-(2k+1)p+1} f(x)$
2840
# so if we are only correct mod p^prec,
2841
# can ignore y powers less than y_prec
2842
t_cube = (t*t*t).truncate_neg(y_prec)
2843
t = t._rmul_(three_halves) - (half_a * t_cube).truncate_neg(y_prec)
2844
# t = (3/2) t - (1/2) a t^3
2845
# print "a =", a
2846
# print "t =", t
2847
# prof("verify")
2848
# print "a*t^2 =", a * t**2
2849
2850
prof("compose")
2851
F_dx_y = (p * x_to_p_less_1 * t) >> (p-1) # px^{p-1} sqrt(a) * y^{-p+1}
2852
2853
# print "-----", F_dx_y
2854
# print "-----", x_to_p * F_dx_y
2855
prof("done")
2856
# print prof
2857
return MonskyWashnitzerDifferential(self, F_dx_y)
2858
2859
def frob_basis_elements(self, prec, p):
2860
"""
2861
Returns the action of a `p`-power lift of Frobenius on the basis
2862
2863
.. math::
2864
2865
\{ dx/2y, x dx/2y, ..., x^{d-2} dx/2y \}
2866
2867
where `d` is the degree of the underlying hyperelliptic curve.
2868
2869
EXAMPLES::
2870
2871
sage: R.<x> = QQ['x']
2872
sage: C = HyperellipticCurve(x^5-4*x+4)
2873
sage: prec = 1
2874
sage: p = 5
2875
sage: MW = C.invariant_differential().parent()
2876
sage: MW.frob_basis_elements(prec,p)
2877
[((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]
2878
"""
2879
F_i = self.frob_invariant_differential(prec, p)
2880
x_to_p = self.x_to_p(p)
2881
F = [F_i]
2882
for i in range(1, self.degree()-1):
2883
F_i *= x_to_p
2884
F.append(F_i)
2885
return F
2886
2887
def helper_matrix(self):
2888
"""
2889
We use this to solve for the linear combination of
2890
`x^i y^j` needed to clear all terms with
2891
`y^{j-1}`.
2892
2893
EXAMPLES::
2894
2895
sage: R.<x> = QQ['x']
2896
sage: C = HyperellipticCurve(x^5-4*x+4)
2897
sage: MW = C.invariant_differential().parent()
2898
sage: MW.helper_matrix()
2899
[ 256/2101 320/2101 400/2101 500/2101 625/2101]
2900
[-625/8404 -64/2101 -80/2101 -100/2101 -125/2101]
2901
[-125/2101 -625/8404 -64/2101 -80/2101 -100/2101]
2902
[-100/2101 -125/2101 -625/8404 -64/2101 -80/2101]
2903
[ -80/2101 -100/2101 -125/2101 -625/8404 -64/2101]
2904
"""
2905
try:
2906
return self._helper_matrix
2907
except AttributeError:
2908
pass
2909
2910
# The smallest y term of (1/j) d(x^i y^j) is constant for all j.
2911
L = []
2912
x, y = self.base_ring().gens()
2913
n = self.degree()
2914
for i in range(n):
2915
L.append((y*x**i).diff().extract_pow_y(0))
2916
A = matrix(L).transpose()
2917
if not is_IntegralDomain(A.base_ring()):
2918
# must be using integer_mod or something to approximate
2919
self._helper_matrix = (~A.change_ring(QQ)).change_ring(A.base_ring())
2920
else:
2921
self._helper_matrix = ~A
2922
return self._helper_matrix
2923
2924
2925
class MonskyWashnitzerDifferential(ModuleElement):
2926
2927
def __init__(self, parent, val=0, offset=0):
2928
r"""
2929
Create an element of the Monsky-Washnitzer ring of differentials, of
2930
the form `F dx/2y`.
2931
2932
INPUT:
2933
2934
- ``parent`` -- Monsky-Washnitzer differential ring (instance of class
2935
:class:`~MonskyWashnitzerDifferentialRing_class`
2936
2937
- ``val`` -- element of the base ring, or list of coefficients
2938
2939
- ``offset`` -- if non-zero, shift val by `y^\text{offset}` (default 0)
2940
2941
EXAMPLES::
2942
2943
sage: R.<x> = QQ['x']
2944
sage: C = HyperellipticCurve(x^5 - 4*x + 4)
2945
sage: x,y = C.monsky_washnitzer_gens()
2946
sage: MW = C.invariant_differential().parent()
2947
sage: sage.schemes.elliptic_curves.monsky_washnitzer.MonskyWashnitzerDifferential(MW, x)
2948
x dx/2y
2949
sage: sage.schemes.elliptic_curves.monsky_washnitzer.MonskyWashnitzerDifferential(MW, y)
2950
y*1 dx/2y
2951
sage: sage.schemes.elliptic_curves.monsky_washnitzer.MonskyWashnitzerDifferential(MW, x, 10)
2952
y^10*x dx/2y
2953
"""
2954
ModuleElement.__init__(self, parent)
2955
if isinstance(val, MonskyWashnitzerDifferential):
2956
val = val._coeff
2957
self._coeff = self.parent().base_ring()(val, offset)
2958
2959
def _add_(left, right):
2960
"""
2961
Returns the sum of left and right, both elements of the
2962
Monsky-Washnitzer ring of differentials.
2963
2964
EXAMPLES::
2965
2966
sage: R.<x> = QQ['x']
2967
sage: C = HyperellipticCurve(x^5-4*x + 4)
2968
sage: x,y = C.monsky_washnitzer_gens()
2969
sage: w = C.invariant_differential()
2970
sage: w + w
2971
2*1 dx/2y
2972
sage: x*w + w
2973
(1 + x) dx/2y
2974
sage: x*w + y*w
2975
(y*1 + x) dx/2y
2976
"""
2977
return MonskyWashnitzerDifferential(left.parent(),
2978
left._coeff + right._coeff)
2979
2980
def _sub_(left, right):
2981
"""
2982
Returns the difference of left and right, both elements of the
2983
Monsky-Washnitzer ring of differentials.
2984
2985
EXAMPLES::
2986
2987
sage: R.<x> = QQ['x']
2988
sage: C = HyperellipticCurve(x^5-4*x+4)
2989
sage: x,y = C.monsky_washnitzer_gens()
2990
sage: w = C.invariant_differential()
2991
sage: w-w
2992
0 dx/2y
2993
sage: x*w-w
2994
(-1 + x) dx/2y
2995
sage: w - x*w - y*w
2996
((1-y)*1 - x) dx/2y
2997
"""
2998
return MonskyWashnitzerDifferential(left.parent(),
2999
left._coeff - right._coeff)
3000
3001
def __neg__(self):
3002
"""
3003
Returns the additive inverse of self.
3004
3005
EXAMPLES::
3006
3007
sage: R.<x> = QQ['x']
3008
sage: C = HyperellipticCurve(x^5-4*x+4)
3009
sage: x,y = C.monsky_washnitzer_gens()
3010
sage: w = C.invariant_differential()
3011
sage: -w
3012
-1 dx/2y
3013
sage: -((y-x)*w)
3014
(-y*1 + x) dx/2y
3015
"""
3016
return MonskyWashnitzerDifferential(self.parent(), -self._coeff)
3017
3018
def _lmul_(self, a):
3019
"""
3020
Returns `self * a`.
3021
3022
EXAMPLES::
3023
3024
sage: R.<x> = QQ['x']
3025
sage: C = HyperellipticCurve(x^5-4*x+4)
3026
sage: x,y = C.monsky_washnitzer_gens()
3027
sage: w = C.invariant_differential()
3028
sage: w*x
3029
x dx/2y
3030
sage: (w*x)*2
3031
2*x dx/2y
3032
sage: w*y
3033
y*1 dx/2y
3034
sage: w*(x+y)
3035
(y*1 + x) dx/2y
3036
"""
3037
return MonskyWashnitzerDifferential(self.parent(), self._coeff * a)
3038
3039
def _rmul_(self, a):
3040
"""
3041
Returns `a * self`.
3042
3043
EXAMPLES::
3044
3045
sage: R.<x> = QQ['x']
3046
sage: C = HyperellipticCurve(x^5-4*x+4)
3047
sage: x,y = C.monsky_washnitzer_gens()
3048
sage: w = C.invariant_differential()
3049
sage: x*w
3050
x dx/2y
3051
sage: 2*(x*w)
3052
2*x dx/2y
3053
sage: y*w
3054
y*1 dx/2y
3055
sage: (x+y)*w
3056
(y*1 + x) dx/2y
3057
"""
3058
return MonskyWashnitzerDifferential(self.parent(), a * self._coeff)
3059
3060
def coeff(self):
3061
r"""
3062
Returns `A`, where this element is `A dx/2y`.
3063
3064
EXAMPLES::
3065
3066
sage: R.<x> = QQ['x']
3067
sage: C = HyperellipticCurve(x^5-4*x+4)
3068
sage: x,y = C.monsky_washnitzer_gens()
3069
sage: w = C.invariant_differential()
3070
sage: w
3071
1 dx/2y
3072
sage: w.coeff()
3073
1
3074
sage: (x*y*w).coeff()
3075
y*x
3076
"""
3077
return self._coeff
3078
3079
def __nonzero__(self):
3080
"""
3081
EXAMPLES::
3082
3083
sage: R.<x> = QQ['x']
3084
sage: C = HyperellipticCurve(x^5-4*x+4)
3085
sage: x,y = C.monsky_washnitzer_gens()
3086
sage: w = C.invariant_differential()
3087
sage: not w
3088
False
3089
sage: not 0*w
3090
True
3091
sage: not x*y*w
3092
False
3093
"""
3094
return not not self._coeff
3095
3096
def _repr_(self):
3097
"""
3098
EXAMPLES::
3099
3100
sage: R.<x> = QQ['x']
3101
sage: C = HyperellipticCurve(x^5-4*x+4)
3102
sage: x,y = C.monsky_washnitzer_gens()
3103
sage: w = C.invariant_differential()
3104
sage: w
3105
1 dx/2y
3106
sage: (2*x+y)*w
3107
(y*1 + 2*x) dx/2y
3108
"""
3109
s = self._coeff._repr_()
3110
if s.find("+") != -1 or s.find("-") > 0:
3111
s = "(%s)" % s
3112
return s + " dx/2y"
3113
3114
def _latex_(self):
3115
"""
3116
Returns the latex representation of self.
3117
3118
EXAMPLES::
3119
3120
sage: R.<x> = QQ['x']
3121
sage: C = HyperellipticCurve(x^5-4*x+4)
3122
sage: x,y = C.monsky_washnitzer_gens()
3123
sage: w = C.invariant_differential()
3124
sage: latex(w)
3125
1 \frac{dx}{2y}
3126
sage: latex(x*w)
3127
x \frac{dx}{2y}
3128
"""
3129
s = self._coeff._latex_()
3130
if s.find("+") != -1 or s.find("-") > 0:
3131
s = "\\left(%s\\right)" % s
3132
return s + " \\frac{dx}{2y}"
3133
3134
def extract_pow_y(self, k):
3135
"""
3136
Returns the power of `y` in `A` where self is `A dx/2y`.
3137
3138
EXAMPLES::
3139
3140
sage: R.<x> = QQ['x']
3141
sage: C = HyperellipticCurve(x^5-3*x+1)
3142
sage: x,y = C.monsky_washnitzer_gens()
3143
sage: A = y^5 - x*y^3
3144
sage: A.extract_pow_y(5)
3145
[1, 0, 0, 0, 0]
3146
sage: (A * C.invariant_differential()).extract_pow_y(5)
3147
[1, 0, 0, 0, 0]
3148
"""
3149
return self._coeff.extract_pow_y(k)
3150
3151
def min_pow_y(self):
3152
"""
3153
Returns the minimum power of `y` in `A` where self is `A dx/2y`.
3154
3155
EXAMPLES::
3156
3157
sage: R.<x> = QQ['x']
3158
sage: C = HyperellipticCurve(x^5-3*x+1)
3159
sage: x,y = C.monsky_washnitzer_gens()
3160
sage: w = y^5 * C.invariant_differential()
3161
sage: w.min_pow_y()
3162
5
3163
sage: w = (x^2*y^4 + y^5) * C.invariant_differential()
3164
sage: w.min_pow_y()
3165
4
3166
"""
3167
return self._coeff.min_pow_y()
3168
3169
def max_pow_y(self):
3170
"""
3171
Returns the maximum power of `y` in `A` where self is `A dx/2y`.
3172
3173
EXAMPLES::
3174
3175
sage: R.<x> = QQ['x']
3176
sage: C = HyperellipticCurve(x^5-3*x+1)
3177
sage: x,y = C.monsky_washnitzer_gens()
3178
sage: w = y^5 * C.invariant_differential()
3179
sage: w.max_pow_y()
3180
5
3181
sage: w = (x^2*y^4 + y^5) * C.invariant_differential()
3182
sage: w.max_pow_y()
3183
5
3184
"""
3185
return self._coeff.max_pow_y()
3186
3187
def reduce_neg_y(self):
3188
"""
3189
Use homology relations to eliminate negative powers of `y`.
3190
3191
EXAMPLES::
3192
3193
sage: R.<x> = QQ['x']
3194
sage: C = HyperellipticCurve(x^5-3*x+1)
3195
sage: x,y = C.monsky_washnitzer_gens()
3196
sage: (y^-1).diff().reduce_neg_y()
3197
((y^-1)*1, 0 dx/2y)
3198
sage: (y^-5*x^2+y^-1*x).diff().reduce_neg_y()
3199
((y^-1)*x + (y^-5)*x^2, 0 dx/2y)
3200
"""
3201
S = self.parent().base_ring()
3202
R = S.base_ring()
3203
M = self.parent().helper_matrix()
3204
p = S._p
3205
n = S.degree()
3206
x, y = S.gens()
3207
f = S(0)
3208
reduced = self
3209
for j in range(self.min_pow_y()+1, 0):
3210
if p is not None and p.divides(j):
3211
cs = [a/j for a in reduced.extract_pow_y(j-1)]
3212
else:
3213
j_inverse = ~R(j)
3214
cs = [a*j_inverse for a in reduced.extract_pow_y(j-1)]
3215
lin_comb = M * vector(M.base_ring(), cs)
3216
# print "j =", j, "b =", cs, "lin_comb =", lin_comb
3217
g = self.parent().base_ring()(0)
3218
if not lin_comb.is_zero():
3219
for i in range(n):
3220
if lin_comb[i] != 0:
3221
g += S.monomial(i, j, lin_comb[i])
3222
if not g.is_zero():
3223
f += g
3224
reduced -= g.diff()
3225
# print g, g.diff()
3226
# print "reduced", reduced
3227
3228
return f, reduced
3229
3230
def reduce_neg_y_fast(self, even_degree_only=False):
3231
"""
3232
Use homology relations to eliminate negative powers of `y`.
3233
3234
EXAMPLES::
3235
3236
sage: R.<x> = QQ['x']
3237
sage: E = HyperellipticCurve(x^5-3*x+1)
3238
sage: x, y = E.monsky_washnitzer_gens()
3239
sage: (y^-1).diff().reduce_neg_y_fast()
3240
((y^-1)*1, 0 dx/2y)
3241
sage: (y^-5*x^2+y^-1*x).diff().reduce_neg_y_fast()
3242
((y^-1)*x + (y^-5)*x^2, 0 dx/2y)
3243
3244
It leaves non-negative powers of `y` alone::
3245
3246
sage: y.diff()
3247
(-3*1 + 5*x^4) dx/2y
3248
sage: y.diff().reduce_neg_y_fast()
3249
(0, (-3*1 + 5*x^4) dx/2y)
3250
"""
3251
# prof = Profiler()
3252
# prof("reduce setup")
3253
S = self.parent().base_ring()
3254
R = S.base_ring()
3255
M = self.parent().helper_matrix()
3256
3257
# prof("extract coeffs")
3258
coeffs, offset = self.coeffs(R)
3259
V = coeffs[0].parent()
3260
3261
if offset == 0:
3262
return S(0), self
3263
3264
# prof("loop %s"%self.min_pow_y())
3265
forms = []
3266
p = S._p
3267
for j in range(self.min_pow_y()+1, 0):
3268
if (even_degree_only and j % 2 == 0) or coeffs[j-offset-1].is_zero():
3269
forms.append(V(0))
3270
else:
3271
# this is a total hack to deal with the fact that we're using
3272
# rational numbers to approximate fixed precision p-adics
3273
if p is not None and j % 3 == 1:
3274
try:
3275
v = coeffs[j-offset-1]
3276
for kk in range(len(v)):
3277
a = v[kk]
3278
ppow = p**max(-a.valuation(S._p), 0)
3279
v[kk] = ((a * ppow) % S._prec_cap) / ppow
3280
except AttributeError:
3281
pass
3282
lin_comb = ~R(j) * (M * coeffs[j-offset-1])
3283
forms.append(lin_comb)
3284
for i in lin_comb.nonzero_positions():
3285
# g = lin_comb[i] x^i y^j
3286
# self -= dg
3287
coeffs[j-offset+1] -= lin_comb[i] * S.monomial_diff_coeffs(i, j)[1]
3288
3289
# prof("recreate forms")
3290
f = S(forms, offset+1)
3291
reduced = S._monsky_washnitzer(coeffs[-1-offset:], -1)
3292
# print self - f.diff() - reduced
3293
# prof("done")
3294
# print prof
3295
return f, reduced
3296
3297
def reduce_neg_y_faster(self, even_degree_only=False):
3298
"""
3299
Use homology relations to eliminate negative powers of `y`.
3300
3301
EXAMPLES::
3302
3303
sage: R.<x> = QQ['x']
3304
sage: C = HyperellipticCurve(x^5-3*x+1)
3305
sage: x,y = C.monsky_washnitzer_gens()
3306
sage: (y^-1).diff().reduce_neg_y()
3307
((y^-1)*1, 0 dx/2y)
3308
sage: (y^-5*x^2+y^-1*x).diff().reduce_neg_y_faster()
3309
((y^-1)*x + (y^-5)*x^2, 0 dx/2y)
3310
"""
3311
# Timings indicate that this is not any faster after all...
3312
3313
S = self.parent().base_ring()
3314
R = S.base_ring()
3315
M = self.parent().helper_matrix()
3316
3317
coeffs, offset = self.coeffs(R)
3318
V = coeffs[0].parent()
3319
zeroV = V(0)
3320
3321
if offset == 0:
3322
return S(0), self
3323
3324
# See monomial_diff_coeffs
3325
# this is the B_i and x_to_i contributions respectively for all i
3326
d_mat_1, d_mat_2 = S.monomial_diff_coeffs_matrices()
3327
3328
forms = []
3329
for j in range(self.min_pow_y()+1, 0):
3330
if coeffs[j-offset-1].is_zero():
3331
forms.append(zeroV)
3332
else:
3333
# this is a total hack to deal with the fact that we're using
3334
# rational numbers to approximate fixed precision p-adics
3335
if j % 3 == 0:
3336
try:
3337
v = coeffs[j-offset-1]
3338
for kk in range(len(v)):
3339
a = v[kk]
3340
ppow = S._p**max(-a.valuation(S._p), 0)
3341
v[kk] = ((a * ppow) % S._prec_cap) / ppow
3342
except AttributeError:
3343
pass
3344
j_inverse = ~R(j)
3345
lin_comb = (M * coeffs[j-offset-1])
3346
forms.append(j_inverse * lin_comb)
3347
coeffs[j-offset+1] -= (d_mat_1 + j_inverse * d_mat_2) * lin_comb
3348
3349
f = S(forms, offset+1)
3350
reduced = S._monsky_washnitzer(coeffs[-1-offset:], -1)
3351
# reduced = self - f.diff()
3352
return f, reduced
3353
3354
def reduce_pos_y(self):
3355
"""
3356
Use homology relations to eliminate positive powers of `y`.
3357
3358
EXAMPLES::
3359
3360
sage: R.<x> = QQ['x']
3361
sage: C = HyperellipticCurve(x^3-4*x+4)
3362
sage: x,y = C.monsky_washnitzer_gens()
3363
sage: (y^2).diff().reduce_pos_y()
3364
(y^2*1, 0 dx/2y)
3365
sage: (y^2*x).diff().reduce_pos_y()
3366
(y^2*x, 0 dx/2y)
3367
sage: (y^92*x).diff().reduce_pos_y()
3368
(y^92*x, 0 dx/2y)
3369
sage: w = (y^3 + x).diff()
3370
sage: w += w.parent()(x)
3371
sage: w.reduce_pos_y_fast()
3372
(y^3*1 + x, x dx/2y)
3373
"""
3374
S = self.parent().base_ring()
3375
n = S.Q().degree()
3376
x, y = S.gens()
3377
f = S(0)
3378
reduced = self
3379
for j in range(self.max_pow_y(), 0, -1):
3380
for i in range(n-1, -1, -1):
3381
c = reduced.extract_pow_y(j)[i]
3382
# print "x^%s y^%s"%(i,j), c
3383
if c != 0:
3384
g = S.monomial(0, j+1) if i == n-1 else S.monomial(i+1, j-1)
3385
dg = g.diff()
3386
# print reduced, " - ", dg
3387
denom = dg.extract_pow_y(j)[i]
3388
c /= denom
3389
c = g.parent()(c)
3390
f += c * g
3391
reduced -= c * dg
3392
3393
return f, reduced
3394
3395
def reduce_pos_y_fast(self, even_degree_only=False):
3396
"""
3397
Use homology relations to eliminate positive powers of `y`.
3398
3399
EXAMPLES::
3400
3401
sage: R.<x> = QQ['x']
3402
sage: E = HyperellipticCurve(x^3-4*x+4)
3403
sage: x, y = E.monsky_washnitzer_gens()
3404
sage: y.diff().reduce_pos_y_fast()
3405
(y*1, 0 dx/2y)
3406
sage: (y^2).diff().reduce_pos_y_fast()
3407
(y^2*1, 0 dx/2y)
3408
sage: (y^2*x).diff().reduce_pos_y_fast()
3409
(y^2*x, 0 dx/2y)
3410
sage: (y^92*x).diff().reduce_pos_y_fast()
3411
(y^92*x, 0 dx/2y)
3412
sage: w = (y^3 + x).diff()
3413
sage: w += w.parent()(x)
3414
sage: w.reduce_pos_y_fast()
3415
(y^3*1 + x, x dx/2y)
3416
"""
3417
S = self.parent().base_ring()
3418
R = S.base_ring()
3419
n = S.Q().degree()
3420
3421
coeffs, offset = self.coeffs(R)
3422
V = coeffs[0].parent()
3423
zeroV = V(0)
3424
forms = [V(0), V(0)]
3425
3426
for j in range(self.max_pow_y(), -1, -1):
3427
3428
if (even_degree_only and j % 2 == 1) or (j > 0 and coeffs[j-offset].is_zero()):
3429
forms.append(zeroV)
3430
continue
3431
3432
form = V(0)
3433
i = n-1
3434
c = coeffs[j-offset][i]
3435
if c != 0:
3436
dg_coeffs = S.monomial_diff_coeffs(0, j+1)[0]
3437
c /= dg_coeffs[i]
3438
forms[len(forms)-2][0] = c
3439
# self -= c d(y^{j+1})
3440
coeffs[j-offset] -= c*dg_coeffs
3441
3442
if j == 0:
3443
# the others are basis elements
3444
break
3445
3446
for i in range(n-2, -1, -1):
3447
c = coeffs[j-offset][i]
3448
if c != 0:
3449
dg_coeffs = S.monomial_diff_coeffs(i+1, j-1)
3450
denom = dg_coeffs[1][i]
3451
c /= denom
3452
form[i+1] = c
3453
# self -= c d(x^{i+1} y^{j-1})
3454
coeffs[j-offset] -= c*dg_coeffs[1]
3455
coeffs[j-offset-2] -= c*dg_coeffs[0]
3456
forms.append(form)
3457
3458
forms.reverse()
3459
f = S(forms)
3460
reduced = self.parent()(coeffs[:1-offset], offset)
3461
return f, reduced
3462
3463
def reduce(self):
3464
"""
3465
Use homology relations to find `a` and `f` such that this element is
3466
equal to `a + df`, where `a` is given in terms of the `x^i dx/2y`.
3467
3468
EXAMPLES::
3469
3470
sage: R.<x> = QQ['x']
3471
sage: C = HyperellipticCurve(x^5-4*x+4)
3472
sage: x,y = C.monsky_washnitzer_gens()
3473
sage: w = (y*x).diff()
3474
sage: w.reduce()
3475
(y*x, 0 dx/2y)
3476
3477
sage: w = x^4 * C.invariant_differential()
3478
sage: w.reduce()
3479
(1/5*y*1, 4/5*1 dx/2y)
3480
3481
sage: w = sum(QQ.random_element() * x^i * y^j for i in [0..4] for j in [-3..3]) * C.invariant_differential()
3482
sage: f, a = w.reduce()
3483
sage: f.diff() + a - w
3484
0 dx/2y
3485
"""
3486
# print "max_pow_y = ", self.max_pow_y(), "min_pow_y = ", self.min_pow_y()
3487
n = self.parent().base_ring().Q().degree()
3488
f1, a = self.reduce_neg_y()
3489
f2, a = a.reduce_pos_y()
3490
f = f1 + f2
3491
3492
c = a.extract_pow_y(0)[n-1]
3493
if c != 0:
3494
x, y = self.parent().base_ring().gens()
3495
g = y
3496
dg = g.diff()
3497
c = g.parent()(c/dg.extract_pow_y(0)[n-1])
3498
f += c * g
3499
a -= c * dg
3500
# print g, dg
3501
3502
return f, a
3503
3504
def reduce_fast(self, even_degree_only=False):
3505
"""
3506
Use homology relations to find `a` and `f` such that this element is
3507
equal to `a + df`, where `a` is given in terms of the `x^i dx/2y`.
3508
3509
EXAMPLES::
3510
3511
sage: R.<x> = QQ['x']
3512
sage: E = HyperellipticCurve(x^3-4*x+4)
3513
sage: x, y = E.monsky_washnitzer_gens()
3514
sage: x.diff().reduce_fast()
3515
(x, (0, 0))
3516
sage: y.diff().reduce_fast()
3517
(y*1, (0, 0))
3518
sage: (y^-1).diff().reduce_fast()
3519
((y^-1)*1, (0, 0))
3520
sage: (y^-11).diff().reduce_fast()
3521
((y^-11)*1, (0, 0))
3522
sage: (x*y^2).diff().reduce_fast()
3523
(y^2*x, (0, 0))
3524
"""
3525
# print "max_pow_y = ", self.max_pow_y(), "min_pow_y = ", self.min_pow_y()
3526
3527
f1, reduced = self.reduce_neg_y_fast(even_degree_only)
3528
f2, reduced = reduced.reduce_pos_y_fast(even_degree_only)
3529
# f1, reduced = self.reduce_neg_y()
3530
# f2, reduced = reduced.reduce_pos_y()
3531
v = reduced.extract_pow_y(0)
3532
v.pop()
3533
V = FreeModule(self.base_ring().base_ring(), len(v))
3534
return f1+f2, V(v)
3535
3536
def coeffs(self, R=None):
3537
"""
3538
Used to obtain the raw coefficients of a differential, see
3539
:meth:`SpecialHyperellipticQuotientElement.coeffs`
3540
3541
INPUT:
3542
3543
- R -- An (optional) base ring in which to cast the coefficients
3544
3545
OUTPUT:
3546
3547
The raw coefficients of $A$ where self is $A dx/2y$.
3548
3549
EXAMPLES::
3550
3551
sage: R.<x> = QQ['x']
3552
sage: C = HyperellipticCurve(x^5-4*x+4)
3553
sage: x,y = C.monsky_washnitzer_gens()
3554
sage: w = C.invariant_differential()
3555
sage: w.coeffs()
3556
([(1, 0, 0, 0, 0)], 0)
3557
sage: (x*w).coeffs()
3558
([(0, 1, 0, 0, 0)], 0)
3559
sage: (y*w).coeffs()
3560
([(0, 0, 0, 0, 0), (1, 0, 0, 0, 0)], 0)
3561
sage: (y^-2*w).coeffs()
3562
([(1, 0, 0, 0, 0), (0, 0, 0, 0, 0), (0, 0, 0, 0, 0)], -2)
3563
"""
3564
return self._coeff.coeffs(R)
3565
3566
def coleman_integral(self, P, Q):
3567
r"""
3568
Computes the definite integral of self from $P$ to $Q$.
3569
3570
INPUT:
3571
3572
- P, Q -- Two points on the underlying curve
3573
3574
OUTPUT:
3575
3576
`\int_P^Q \text{self}`
3577
3578
EXAMPLES::
3579
3580
sage: K = pAdicField(5,7)
3581
sage: E = EllipticCurve(K,[-31/3,-2501/108]) #11a
3582
sage: P = E(K(14/3), K(11/2))
3583
sage: w = E.invariant_differential()
3584
sage: w.coleman_integral(P,2*P)
3585
O(5^6)
3586
3587
sage: Q = E.lift_x(3)
3588
sage: w.coleman_integral(P,Q)
3589
2*5 + 4*5^2 + 3*5^3 + 4*5^4 + 3*5^5 + O(5^6)
3590
sage: w.coleman_integral(2*P,Q)
3591
2*5 + 4*5^2 + 3*5^3 + 4*5^4 + 3*5^5 + O(5^6)
3592
sage: (2*w).coleman_integral(P, Q) == 2*(w.coleman_integral(P, Q))
3593
True
3594
"""
3595
return self.parent().base_ring().curve().coleman_integral(self, P, Q)
3596
3597
integrate = coleman_integral
3598
3599