Frobenius on Monsky-Washnitzer cohomology of a hyperelliptic curve¶
This module provides hypellfrob()
, that is a wrapper for the matrix()
function in hypellfrob.cpp
.
hypellfrob.cpp
is a C++ program for computing the zeta function of a
hyperelliptic curve over a largish prime finite field, based on the method
described in the paper [Harv2007]. More precisely, it computes the matrix of
Frobenius on the Monsky-Washnitzer cohomology of the curve; the zeta function
can be recovered via the characteristic polynomial of the matrix.
AUTHORS:
David Harvey (2007-05): initial version
David Harvey (2007-12): rewrote for
hypellfrob
version 2.0Alex J. Best (2018-02): added wrapper
- sage.schemes.hyperelliptic_curves.hypellfrob.hypellfrob(p, N, Q)[source]¶
Compute the matrix of Frobenius acting on the Monsky-Washnitzer cohomology of a hyperelliptic curve \(y^2 = Q(x)\), with respect to the basis \(x^i dx/y\), \(0 \leq i < 2g\).
INPUT:
p
– a primeQ
– a monic polynomial in \(\ZZ[x]\) of odd degree; must have no multiple roots mod \(p\)N
– precision parameter; the output matrix will be correct modulo \(p^N\)
The prime \(p\) should satisfy \(p > (2g+1)(2N-1)\), where \(g = \left(\deg Q - 1\right) / 2\) is the genus of the curve.
ALGORITHM:
Described in [Harv2007], Section 7. Running time is theoretically \(\widetilde{O}(p^{1/2} N^{5/2} g^3)\).
EXAMPLES:
sage: from sage.schemes.hyperelliptic_curves.hypellfrob import hypellfrob sage: R.<x> = PolynomialRing(ZZ) sage: f = x^5 + 2*x^2 + x + 1; p = 101 sage: M = hypellfrob(p, 4, f); M [ 91844754 + O(101^4) 38295665 + O(101^4) 44498269 + O(101^4) 11854028 + O(101^4)] [ 93514789 + O(101^4) 48987424 + O(101^4) 53287857 + O(101^4) 61431148 + O(101^4)] [ 77916046 + O(101^4) 60656459 + O(101^4) 101244586 + O(101^4) 56237448 + O(101^4)] [ 58643832 + O(101^4) 81727988 + O(101^4) 85294589 + O(101^4) 70104432 + O(101^4)] sage: -M.trace() 7 + O(101^4) sage: sum(legendre_symbol(f(i), p) for i in range(p)) 7 sage: ZZ(M.det()) 10201 sage: M = hypellfrob(p, 1, f); M [ O(101) O(101) 93 + O(101) 62 + O(101)] [ O(101) O(101) 55 + O(101) 19 + O(101)] [ O(101) O(101) 65 + O(101) 42 + O(101)] [ O(101) O(101) 89 + O(101) 29 + O(101)]
>>> from sage.all import * >>> from sage.schemes.hyperelliptic_curves.hypellfrob import hypellfrob >>> R = PolynomialRing(ZZ, names=('x',)); (x,) = R._first_ngens(1) >>> f = x**Integer(5) + Integer(2)*x**Integer(2) + x + Integer(1); p = Integer(101) >>> M = hypellfrob(p, Integer(4), f); M [ 91844754 + O(101^4) 38295665 + O(101^4) 44498269 + O(101^4) 11854028 + O(101^4)] [ 93514789 + O(101^4) 48987424 + O(101^4) 53287857 + O(101^4) 61431148 + O(101^4)] [ 77916046 + O(101^4) 60656459 + O(101^4) 101244586 + O(101^4) 56237448 + O(101^4)] [ 58643832 + O(101^4) 81727988 + O(101^4) 85294589 + O(101^4) 70104432 + O(101^4)] >>> -M.trace() 7 + O(101^4) >>> sum(legendre_symbol(f(i), p) for i in range(p)) 7 >>> ZZ(M.det()) 10201 >>> M = hypellfrob(p, Integer(1), f); M [ O(101) O(101) 93 + O(101) 62 + O(101)] [ O(101) O(101) 55 + O(101) 19 + O(101)] [ O(101) O(101) 65 + O(101) 42 + O(101)] [ O(101) O(101) 89 + O(101) 29 + O(101)]
Todo
Remove the restriction on \(p\). Probably by merging in Robert’s code, which eventually needs a fast C++/NTL implementation.
- sage.schemes.hyperelliptic_curves.hypellfrob.interval_products(M0, M1, target)[source]¶
Given matrices \(M(t)\) with entries linear in \(t\) over \(\ZZ/N\ZZ\) and a list of integers \(a_0 < b_0 \le a_1 < b_1 \le \cdots \le a_n < b_n\), compute the matrices \(\prod_{t = a_i + 1}^{b_i} M(t)\) for \(i = 0\) to \(n\).
INPUT:
M0
,M1
– matrices over \(\ZZ/N\ZZ\), so that \(M(t) = M_0 + M_1t\)target
– list of integers \(a_0, b_0, \dots, a_n, b_n\)
ALGORITHM:
Described in [Harv2007], Theorem 10. Based on the work of Bostan-Gaudry-Schost [BGS2007].
EXAMPLES:
sage: from sage.schemes.hyperelliptic_curves.hypellfrob import interval_products sage: interval_products(Matrix(Integers(9), 2,2, [1,0,1,0]), ....: Matrix(Integers(9), 2, 2, [1, 1, 0, 2]),[0,2,2,4]) [ [7 8] [5 4] [5 1], [2 7] ] sage: [prod(Matrix(Integers(9), 2, 2, [t + 1, t, 1, 2*t]) ....: for t in range(2*i + 1, 2*i + 1 + 2)) for i in range(2)] [ [7 8] [5 4] [5 1], [2 7] ]
>>> from sage.all import * >>> from sage.schemes.hyperelliptic_curves.hypellfrob import interval_products >>> interval_products(Matrix(Integers(Integer(9)), Integer(2),Integer(2), [Integer(1),Integer(0),Integer(1),Integer(0)]), ... Matrix(Integers(Integer(9)), Integer(2), Integer(2), [Integer(1), Integer(1), Integer(0), Integer(2)]),[Integer(0),Integer(2),Integer(2),Integer(4)]) [ [7 8] [5 4] [5 1], [2 7] ] >>> [prod(Matrix(Integers(Integer(9)), Integer(2), Integer(2), [t + Integer(1), t, Integer(1), Integer(2)*t]) ... for t in range(Integer(2)*i + Integer(1), Integer(2)*i + Integer(1) + Integer(2))) for i in range(Integer(2))] [ [7 8] [5 4] [5 1], [2 7] ]
An example with larger modulus:
sage: interval_products(Matrix(Integers(3^8), 1, 1, [1]), ....: Matrix(Integers(3^8), 1, 1, [1]), [2,4]) [[20]] sage: [prod(Matrix(Integers(3^8), 1, 1, [t + 1]) for t in range(3,5))] [[20]]
>>> from sage.all import * >>> interval_products(Matrix(Integers(Integer(3)**Integer(8)), Integer(1), Integer(1), [Integer(1)]), ... Matrix(Integers(Integer(3)**Integer(8)), Integer(1), Integer(1), [Integer(1)]), [Integer(2),Integer(4)]) [[20]] >>> [prod(Matrix(Integers(Integer(3)**Integer(8)), Integer(1), Integer(1), [t + Integer(1)]) for t in range(Integer(3),Integer(5)))] [[20]]
An even larger modulus:
sage: interval_products(Matrix(Integers(3^18), 1, 1, [1]), ....: Matrix(Integers(3^18), 1, 1, [1]), [2,4]) [[20]] sage: [prod(Matrix(Integers(3^18), 1, 1, [t + 1]) for t in range(3,5))] [[20]]
>>> from sage.all import * >>> interval_products(Matrix(Integers(Integer(3)**Integer(18)), Integer(1), Integer(1), [Integer(1)]), ... Matrix(Integers(Integer(3)**Integer(18)), Integer(1), Integer(1), [Integer(1)]), [Integer(2),Integer(4)]) [[20]] >>> [prod(Matrix(Integers(Integer(3)**Integer(18)), Integer(1), Integer(1), [t + Integer(1)]) for t in range(Integer(3),Integer(5)))] [[20]]