Path: blob/master/src/sage/modular/modform/vm_basis.py
8820 views
"""1The Victor Miller Basis23This module contains functions for quick calculation of a basis of4`q`-expansions for the space of modular forms of level 1 and any weight. The5basis returned is the Victor Miller basis, which is the unique basis of6elliptic modular forms `f_1, \dots, f_d` for which `a_i(f_j) = \delta_{ij}`7for `1 \le i, j \le d` (where `d` is the dimension of the space).89This basis is calculated using a standard set of generators for the ring of10modular forms, using the fast multiplication algorithms for polynomials and11power series provided by the FLINT library. (This is far quicker than using12modular symbols).1314TESTS::1516sage: ModularSymbols(1, 36, 1).cuspidal_submodule().q_expansion_basis(30) == victor_miller_basis(36, 30, cusp_only=True)17True18"""1920#*****************************************************************************21# Copyright (C) 2006 William Stein <[email protected]>22#23# Distributed under the terms of the GNU General Public License (GPL)24# as published by the Free Software Foundation; either version 2 of25# the License, or (at your option) any later version.26# http://www.gnu.org/licenses/27#*****************************************************************************2829import math3031from sage.rings.all import QQ, ZZ, Integer, \32PolynomialRing, PowerSeriesRing, O as bigO33from sage.structure.all import Sequence34from sage.libs.flint.fmpz_poly import Fmpz_poly35from sage.misc.all import verbose3637from eis_series_cython import eisenstein_series_poly3839def victor_miller_basis(k, prec=10, cusp_only=False, var='q'):40r"""41Compute and return the Victor Miller basis for modular forms of42weight `k` and level 1 to precision `O(q^{prec})`. If43``cusp_only`` is True, return only a basis for the cuspidal44subspace.4546INPUT:4748- ``k`` -- an integer4950- ``prec`` -- (default: 10) a positive integer5152- ``cusp_only`` -- bool (default: False)5354- ``var`` -- string (default: 'q')5556OUTPUT:5758A sequence whose entries are power series in ``ZZ[[var]]``.5960EXAMPLES::6162sage: victor_miller_basis(1, 6)63[]64sage: victor_miller_basis(0, 6)65[661 + O(q^6)67]68sage: victor_miller_basis(2, 6)69[]70sage: victor_miller_basis(4, 6)71[721 + 240*q + 2160*q^2 + 6720*q^3 + 17520*q^4 + 30240*q^5 + O(q^6)73]7475sage: victor_miller_basis(6, 6, var='w')76[771 - 504*w - 16632*w^2 - 122976*w^3 - 532728*w^4 - 1575504*w^5 + O(w^6)78]7980sage: victor_miller_basis(6, 6)81[821 - 504*q - 16632*q^2 - 122976*q^3 - 532728*q^4 - 1575504*q^5 + O(q^6)83]84sage: victor_miller_basis(12, 6)85[861 + 196560*q^2 + 16773120*q^3 + 398034000*q^4 + 4629381120*q^5 + O(q^6),87q - 24*q^2 + 252*q^3 - 1472*q^4 + 4830*q^5 + O(q^6)88]8990sage: victor_miller_basis(12, 6, cusp_only=True)91[92q - 24*q^2 + 252*q^3 - 1472*q^4 + 4830*q^5 + O(q^6)93]94sage: victor_miller_basis(24, 6, cusp_only=True)95[96q + 195660*q^3 + 12080128*q^4 + 44656110*q^5 + O(q^6),97q^2 - 48*q^3 + 1080*q^4 - 15040*q^5 + O(q^6)98]99sage: victor_miller_basis(24, 6)100[1011 + 52416000*q^3 + 39007332000*q^4 + 6609020221440*q^5 + O(q^6),102q + 195660*q^3 + 12080128*q^4 + 44656110*q^5 + O(q^6),103q^2 - 48*q^3 + 1080*q^4 - 15040*q^5 + O(q^6)104]105sage: victor_miller_basis(32, 6)106[1071 + 2611200*q^3 + 19524758400*q^4 + 19715347537920*q^5 + O(q^6),108q + 50220*q^3 + 87866368*q^4 + 18647219790*q^5 + O(q^6),109q^2 + 432*q^3 + 39960*q^4 - 1418560*q^5 + O(q^6)110]111112sage: victor_miller_basis(40,200)[1:] == victor_miller_basis(40,200,cusp_only=True)113True114sage: victor_miller_basis(200,40)[1:] == victor_miller_basis(200,40,cusp_only=True)115True116117AUTHORS:118119- William Stein, Craig Citro: original code120121- Martin Raum (2009-08-02): use FLINT for polynomial arithmetic (instead of NTL)122"""123k = Integer(k)124if k%2 == 1 or k==2:125return Sequence([])126elif k < 0:127raise ValueError, "k must be non-negative"128elif k == 0:129return Sequence([PowerSeriesRing(ZZ,var)(1).add_bigoh(prec)], cr=True)130e = k.mod(12)131if e == 2: e += 12132n = (k-e) // 12133134if n == 0 and cusp_only:135return Sequence([])136137# If prec is less than or equal to the dimension of the space of138# cusp forms, which is just n, then we know the answer, and we139# simply return it.140if prec <= n:141q = PowerSeriesRing(ZZ,var).gen(0)142err = bigO(q**prec)143ls = [0] * (n+1)144if not cusp_only:145ls[0] = 1 + err146for i in range(1,prec):147ls[i] = q**i + err148for i in range(prec,n+1):149ls[i] = err150return Sequence(ls, cr=True)151152F6 = eisenstein_series_poly(6,prec)153154if e == 0:155A = Fmpz_poly(1)156elif e == 4:157A = eisenstein_series_poly(4,prec)158elif e == 6:159A = F6160elif e == 8:161A = eisenstein_series_poly(8,prec)162elif e == 10:163A = eisenstein_series_poly(10,prec)164else: # e == 14165A = eisenstein_series_poly(14,prec)166167if A[0] == -1 :168A = -A169170if n == 0:171return Sequence([PowerSeriesRing(ZZ,var)(A.list()).add_bigoh(prec)],cr=True)172173F6_squared = F6**2174F6_squared._unsafe_mutate_truncate(prec)175D = _delta_poly(prec)176Fprod = F6_squared177Dprod = D178179if cusp_only:180ls = [Fmpz_poly(0)] + [A] * n181else:182ls = [A] * (n+1)183184for i in xrange(1,n+1):185ls[n-i] *= Fprod186ls[i] *= Dprod187ls[n-i]._unsafe_mutate_truncate(prec)188ls[i]._unsafe_mutate_truncate(prec)189190Fprod *= F6_squared191Dprod *= D192Fprod._unsafe_mutate_truncate(prec)193Dprod._unsafe_mutate_truncate(prec)194195196P = PowerSeriesRing(ZZ,var)197if cusp_only :198for i in xrange(1,n+1) :199for j in xrange(1, i) :200ls[j] = ls[j] - ls[j][i]*ls[i]201202return Sequence(map(lambda l: P(l.list()).add_bigoh(prec), ls[1:]),cr=True)203else :204for i in xrange(1,n+1) :205for j in xrange(i) :206ls[j] = ls[j] - ls[j][i]*ls[i]207208return Sequence(map(lambda l: P(l.list()).add_bigoh(prec), ls), cr=True)209210def _delta_poly(prec=10):211"""212Return the q-expansion of Delta as a FLINT polynomial. Used internally by213the :func:`~delta_qexp` function. See the docstring of :func:`~delta_qexp`214for more information.215216INPUT:217218- ``prec`` -- integer; the absolute precision of the output219220OUTPUT:221222the q-expansion of Delta to precision ``prec``, as a FLINT223:class:`~sage.libs.flint.fmpz_poly.Fmpz_poly` object.224225EXAMPLES::226227sage: from sage.modular.modform.vm_basis import _delta_poly228sage: _delta_poly(7)2297 0 1 -24 252 -1472 4830 -6048230"""231if prec <= 0:232raise ValueError, "prec must be positive"233v = [0] * prec234235# Let F = \sum_{n >= 0} (-1)^n (2n+1) q^(floor(n(n+1)/2)).236# Then delta is F^8.237238# First compute F^2 directly by naive polynomial multiplication,239# since F is very sparse.240241stop = int((-1+math.sqrt(1+8*prec))/2.0)242# make list of index/value pairs for the sparse poly243values = [(n*(n+1)//2, ((-2*n-1) if (n & 1) else (2*n+1))) \244for n in xrange(stop+1)]245246for (i1, v1) in values:247for (i2, v2) in values:248try:249v[i1 + i2] += v1 * v2250except IndexError:251break252253f = Fmpz_poly(v)254t = verbose('made series')255f = f*f256f._unsafe_mutate_truncate(prec)257t = verbose('squared (2 of 3)', t)258f = f*f259f._unsafe_mutate_truncate(prec - 1)260t = verbose('squared (3 of 3)', t)261f = f.left_shift(1)262t = verbose('shifted', t)263264return f265266def _delta_poly_modulo(N, prec=10):267"""268Return the q-expansion of `\Delta` modulo `N`. Used internally by269the :func:`~delta_qexp` function. See the docstring of :func:`~delta_qexp`270for more information.271272INPUT:273274- `N` -- positive integer modulo which we want to compute `\Delta`275276- ``prec`` -- integer; the absolute precision of the output277278OUTPUT:279280the polynomial of degree ``prec``-1 which is the truncation281of `\Delta` modulo `N`, as an element of the polynomial282ring in `q` over the integers modulo `N`.283284EXAMPLES::285286sage: from sage.modular.modform.vm_basis import _delta_poly_modulo287sage: _delta_poly_modulo(5, 7)2882*q^6 + 3*q^4 + 2*q^3 + q^2 + q289sage: _delta_poly_modulo(10, 12)2902*q^11 + 7*q^9 + 6*q^7 + 2*q^6 + 8*q^4 + 2*q^3 + 6*q^2 + q291"""292if prec <= 0:293raise ValueError( "prec must be positive" )294v = [0] * prec295296# Let F = \sum_{n >= 0} (-1)^n (2n+1) q^(floor(n(n+1)/2)).297# Then delta is F^8.298299stop = int((-1+math.sqrt(8*prec))/2.0)300301for n in xrange(stop+1):302v[n*(n+1)//2] = ((N-1)*(2*n+1) if (n & 1) else (2*n+1))303304from sage.rings.all import Integers305306P = PolynomialRing(Integers(N), 'q')307f = P(v)308t = verbose('made series')309# fast way of computing f*f truncated at prec310f = f._mul_trunc(f, prec)311t = verbose('squared (1 of 3)', t)312f = f._mul_trunc(f, prec)313t = verbose('squared (2 of 3)', t)314f = f._mul_trunc(f, prec - 1)315t = verbose('squared (3 of 3)', t)316f = f.shift(1)317t = verbose('shifted', t)318319return f320321322def delta_qexp(prec=10, var='q', K=ZZ) :323"""324Return the `q`-expansion of the weight 12 cusp form `\Delta` as a power325series with coefficients in the ring K (`= \ZZ` by default).326327INPUT:328329- ``prec`` -- integer (default 10), the absolute precision of the output330(must be positive)331332- ``var`` -- string (default: 'q'), variable name333334- ``K`` -- ring (default: `\ZZ`), base ring of answer335336OUTPUT:337338a power series over K in the variable ``var``339340ALGORITHM:341342Compute the theta series343344.. math::345346\sum_{n \ge 0} (-1)^n (2n+1) q^{n(n+1)/2},347348a very simple explicit modular form whose 8th power is `\Delta`. Then349compute the 8th power. All computations are done over `\ZZ` or `\ZZ`350modulo `N` depending on the characteristic of the given coefficient351ring `K`, and coerced into `K` afterwards.352353EXAMPLES::354355sage: delta_qexp(7)356q - 24*q^2 + 252*q^3 - 1472*q^4 + 4830*q^5 - 6048*q^6 + O(q^7)357sage: delta_qexp(7,'z')358z - 24*z^2 + 252*z^3 - 1472*z^4 + 4830*z^5 - 6048*z^6 + O(z^7)359sage: delta_qexp(-3)360Traceback (most recent call last):361...362ValueError: prec must be positive363sage: delta_qexp(20, K = GF(3))364q + q^4 + 2*q^7 + 2*q^13 + q^16 + 2*q^19 + O(q^20)365sage: delta_qexp(20, K = GF(3^5, 'a'))366q + q^4 + 2*q^7 + 2*q^13 + q^16 + 2*q^19 + O(q^20)367sage: delta_qexp(10, K = IntegerModRing(60))368q + 36*q^2 + 12*q^3 + 28*q^4 + 30*q^5 + 12*q^6 + 56*q^7 + 57*q^9 + O(q^10)369370TESTS:371372Test algorithm with modular arithmetic (see also #11804)::373374sage: delta_qexp(10^4).change_ring(GF(13)) == delta_qexp(10^4, K=GF(13))375True376sage: delta_qexp(1000).change_ring(IntegerModRing(5^100)) == delta_qexp(1000, K=IntegerModRing(5^100))377True378379AUTHORS:380381- William Stein: original code382383- David Harvey (2007-05): sped up first squaring step384385- Martin Raum (2009-08-02): use FLINT for polynomial arithmetic (instead of NTL)386"""387R = PowerSeriesRing(K, var)388if K in (ZZ, QQ):389return R(_delta_poly(prec).list(), prec, check=False)390ch = K.characteristic()391if ch > 0 and prec > 150:392return R(_delta_poly_modulo(ch, prec), prec, check=False)393else:394# compute over ZZ and coerce395return R(_delta_poly(prec).list(), prec, check=True)396397398