-functions from PARI¶
This is a wrapper around the general PARI
REFERENCES:
AUTHORS:
Frédéric Chapoton (2018) interface
- class sage.lfunctions.pari.LFunction(lfun, prec=None, max_im=1)[source]¶
Bases:
SageObjectBuild the
-function from a PARI -function.Rank 1 elliptic curve
We compute with the
-series of a rank curve.sage: E = EllipticCurve('37a') sage: L = E.lseries().dokchitser(algorithm='pari'); L PARI L-function associated to Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field sage: L(1) 0.000000000000000 sage: L.derivative(1) 0.305999773834052 sage: L.derivative(1, 2) 0.373095594536324 sage: L.cost() 50 sage: L.taylor_series(1, 4) 0.000000000000000 + 0.305999773834052*z + 0.186547797268162*z^2 - 0.136791463097188*z^3 + O(z^4) sage: L.check_functional_equation() # abs tol 4e-19 1.08420217248550e-19
Rank 2 elliptic curve
We compute the leading coefficient and Taylor expansion of the
-series of a rank elliptic curve:sage: E = EllipticCurve('389a') sage: L = E.lseries().dokchitser(algorithm='pari') sage: L.cost() 163 sage: L.derivative(1, E.rank()) 1.51863300057685 sage: L.taylor_series(1, 4) ...e-19 + (...e-19)*z + 0.759316500288427*z^2 - 0.430302337583362*z^3 + O(z^4)
Number field
We compute with the Dedekind zeta function of a number field:
sage: x = var('x') sage: K = NumberField(x**4 - x**2 - 1,'a') sage: L = K.zeta_function(algorithm='pari') sage: L.conductor 400 sage: L.cost() 313 sage: L(2) 1.10398438736918 sage: L.taylor_series(2, 3) 1.10398438736918 - 0.215822638498759*z + 0.279836437522536*z^2 + O(z^3)
Ramanujan
L-functionThe coefficients are given by Ramanujan’s tau function:
sage: from sage.lfunctions.pari import lfun_generic, LFunction sage: lf = lfun_generic(conductor=1, gammaV=[0, 1], weight=12, eps=1) sage: tau = pari('k->vector(k,n,(5*sigma(n,3)+7*sigma(n,5))*n/12 - 35*sum(k=1,n-1,(6*k-4*(n-k))*sigma(k,3)*sigma(n-k,5)))') sage: lf.init_coeffs(tau) sage: L = LFunction(lf)
Now we are ready to evaluate, etc.
sage: L(1) 0.0374412812685155 sage: L.taylor_series(1, 3) 0.0374412812685155 + 0.0709221123619322*z + 0.0380744761270520*z^2 + O(z^3)
- Lambda(s)[source]¶
Evaluate the completed
-function at s.EXAMPLES:
sage: from sage.lfunctions.pari import lfun_number_field, LFunction sage: L = LFunction(lfun_number_field(QQ)) sage: L.Lambda(2) 0.523598775598299 sage: L.Lambda(1 - 2) 0.523598775598299
- check_functional_equation()[source]¶
Verify how well numerically the functional equation is satisfied.
If what this function returns does not look like 0 at all, probably the functional equation is wrong, i.e. some of the parameters gammaV, conductor, etc., or the coefficients are wrong.
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_generic, LFunction sage: lf = lfun_generic(conductor=1, gammaV=[0], weight=1, eps=1, ....: poles=[1], residues=[1], v=pari('k->vector(k,n,1)')) sage: L = LFunction(lf) sage: L.check_functional_equation() 4.33680868994202e-19
If we choose the sign in functional equation for the
function incorrectly, the functional equation does not check out:sage: lf = lfun_generic(conductor=1, gammaV=[0], weight=1, ....: eps=-1, poles=[1], residues=[1]) sage: lf.init_coeffs([1]*2000) sage: L = LFunction(lf) sage: L.check_functional_equation() 16.0000000000000
- property conductor¶
Return the conductor.
EXAMPLES:
sage: from sage.lfunctions.pari import * sage: L = LFunction(lfun_number_field(QQ)); L.conductor 1 sage: E = EllipticCurve('11a') sage: L = LFunction(lfun_number_field(E)); L.conductor 11
- cost(domain=None)[source]¶
Return number of coefficients
that are needed in order to perform most relevant -function computations to the desired precision.INPUT:
domain– optional triple (center, width, height)
The domain is then a rectangle around the real point
centerwith size2*widthand2*height.For computation with real arguments, one should set
heightto zero.EXAMPLES:
sage: E = EllipticCurve('11a') sage: L = E.lseries().dokchitser(algorithm='pari') sage: L.cost() 27 sage: E = EllipticCurve('5077a') sage: L = E.lseries().dokchitser(algorithm='pari') sage: L.cost() 591 sage: from sage.lfunctions.pari import lfun_generic, LFunction sage: lf = lfun_generic(conductor=1, gammaV=[0], weight=1, eps=1, ....: poles=[1], residues=[-1], v=pari('k->vector(k,n,1)')) sage: L = LFunction(lf) sage: L.cost() 4
- derivative(s, D=1)[source]¶
Return the derivative of the
-function at point s and order D.INPUT:
s– complex numberD– integer (default: 1)
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_number_field, LFunction sage: L = LFunction(lfun_number_field(QQ)) sage: L.derivative(2) -0.937548254315844
- hardy(t)[source]¶
Evaluate the associated Hardy function at t.
This only works for real t.
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_number_field, LFunction sage: L = LFunction(lfun_number_field(QQ)) sage: L.hardy(2) -0.962008487244041
- taylor_series(s, k=6, var='z')[source]¶
Return the first
terms of the Taylor series expansion of the -series about .This is returned as a formal power series in
var.INPUT:
s– complex number; point about which to expandk– integer (default: 6); series isvar– string (default:'z'); variable of power series
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_number_field, LFunction sage: lf = lfun_number_field(QQ) sage: L = LFunction(lf) sage: L.taylor_series(2, 3) 1.64493406684823 - 0.937548254315844*z + 0.994640117149451*z^2 + O(z^3) sage: E = EllipticCurve('37a') sage: L = E.lseries().dokchitser(algorithm='pari') sage: L.taylor_series(1) 0.000000000000000 + 0.305999773834052*z + 0.186547797268162*z^2 - 0.136791463097188*z^3 + 0.0161066468496401*z^4 + 0.0185955175398802*z^5 + O(z^6)
We compute a Taylor series where each coefficient is to high precision:
sage: E = EllipticCurve('389a') sage: L = E.lseries().dokchitser(200,algorithm='pari') sage: L.taylor_series(1, 3) 2...e-63 + (...e-63)*z + 0.75931650028842677023019260789472201907809751649492435158581*z^2 + O(z^3)
Check that Issue #25402 is fixed:
sage: L = EllipticCurve("24a1").modular_form().lseries() sage: L.taylor_series(-1, 3) 0.000000000000000 - 0.702565506265199*z + 0.638929001045535*z^2 + O(z^3)
- sage.lfunctions.pari.lfun_character(chi)[source]¶
Create the
-function of a primitive Dirichlet character.If the given character is not primitive, it is replaced by its associated primitive character.
OUTPUT: one pari:lfun object
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_character, LFunction sage: chi = DirichletGroup(6).gen().primitive_character() sage: L = LFunction(lfun_character(chi)) sage: L(3) 0.884023811750080
- sage.lfunctions.pari.lfun_delta()[source]¶
Return the
-function of Ramanujan’s Delta modular form.EXAMPLES:
sage: from sage.lfunctions.pari import lfun_delta, LFunction sage: L = LFunction(lfun_delta()) sage: L(1) 0.0374412812685155
- sage.lfunctions.pari.lfun_eisenstein(j, algorithm='mf')[source]¶
Return the
-function of the Eisenstein form .INPUT:
j– an even integer, at least
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_eisenstein, LFunction sage: L = LFunction(lfun_eisenstein(16)) sage: L(1) -0.291657724743874 sage: L = LFunction(lfun_eisenstein(20)) sage: L(2) -5.02355351645998
- sage.lfunctions.pari.lfun_elliptic_curve(E)[source]¶
Create the
-function of an elliptic curve.OUTPUT: one pari:lfun object
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_elliptic_curve, LFunction sage: E = EllipticCurve('11a1') sage: L = LFunction(lfun_elliptic_curve(E)) sage: L(3) 0.752723147389051 sage: L(1) 0.253841860855911
Over number fields:
sage: K.<a> = QuadraticField(2) sage: E = EllipticCurve([1, a]) sage: L = LFunction(lfun_elliptic_curve(E)) sage: L(3) 1.00412346717019
- sage.lfunctions.pari.lfun_eta_quotient(scalings, exponents)[source]¶
Return the
-function of an eta-quotient.This uses pari:lfunetaquo.
INPUT:
scalings– list of integers; the scaling factorsexponents– list of integers; the exponents
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_eta_quotient, LFunction sage: L = LFunction(lfun_eta_quotient([1], [24])) sage: L(1) 0.0374412812685155 sage: lfun_eta_quotient([6], [4]) [[Vecsmall([7]), [Vecsmall([6]), Vecsmall([4]), 0]], 0, [0, 1], 2, 36, 1] sage: lfun_eta_quotient([2, 1, 4], [5, -2, -2]) Traceback (most recent call last): ... PariError: sorry, noncuspidal eta quotient is not yet implemented
- class sage.lfunctions.pari.lfun_generic(conductor, gammaV, weight, eps, poles=[], residues='automatic', *args, **kwds)[source]¶
Bases:
objectCreate a PARI
-function (pari:lfun instance).The arguments are:
lfun_generic(conductor, gammaV, weight, eps, poles, residues, init)where
conductor– integer; the conductorgammaV– list of Gamma-factor parameters, e.g. [0] for Riemann zeta, [0,1] for elliptic curves, (see examples)weight– positive real number, usually an integer e.g. 1 for Riemann zeta, 2 for of curves overeps– complex number; sign in functional equationpoles– (default:[]) list of points where has (simple) poles; only poles with should be includedresidues– vector of residues of in those poles or setresidues='automatic'(default)init– list of coefficients (optional)
RIEMANN ZETA FUNCTION:
We compute with the Riemann Zeta function:
sage: from sage.lfunctions.pari import lfun_generic, LFunction sage: lf = lfun_generic(conductor=1, gammaV=[0], weight=1, ....: eps=1, poles=[1], residues=[1]) sage: lf.init_coeffs([1]*2000)
Now we can wrap this PARI
-function into one Sage -function:sage: L = LFunction(lf); L L-series of conductor 1 and weight 1 sage: L(1) Traceback (most recent call last): ... ArithmeticError: pole here sage: L(2) 1.64493406684823 sage: L.derivative(2) -0.937548254315844 sage: h = RR('0.0000000000001') sage: (zeta(2+h) - zeta(2.))/h -0.937028232783632 sage: L.taylor_series(2, k=5) 1.64493406684823 - 0.937548254315844*z + 0.994640117149451*z^2 - 1.00002430047384*z^3 + 1.00006193307...*z^4 + O(z^5)
- init_coeffs(v, w=1)[source]¶
Set the coefficients
of the -series.If
is not equal to its dual, pass the coefficients of the dual as the second argument.INPUT:
v– list of complex numbers or unary functionw– list of complex numbers or unary function
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_generic, LFunction sage: lf = lfun_generic(conductor=1, gammaV=[0, 1], weight=12, eps=1) sage: pari_coeffs = pari('k->vector(k,n,(5*sigma(n,3)+7*sigma(n,5))*n/12 - 35*sum(k=1,n-1,(6*k-4*(n-k))*sigma(k,3)*sigma(n-k,5)))') sage: lf.init_coeffs(pari_coeffs)
Evaluate the resulting
-function at a point, and compare with the answer that one gets “by definition” (of -function attached to a modular form):sage: L = LFunction(lf) sage: L(14) 0.998583063162746 sage: a = delta_qexp(1000) sage: sum(a[n]/float(n)^14 for n in reversed(range(1,1000))) 0.9985830631627461
Illustrate that one can give a list of complex numbers for v (see Issue #10937):
sage: l2 = lfun_generic(conductor=1, gammaV=[0, 1], weight=12, eps=1) sage: l2.init_coeffs(list(delta_qexp(1000))[1:]) sage: L2 = LFunction(l2) sage: L2(14) 0.998583063162746
- init_empty()[source]¶
Create the Pari object with the Dokchitser parameters only.
This is useful to ask Pari for the number of terms.
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_generic, LFunction sage: lf = lfun_generic(conductor=1, gammaV=[0], weight=1, ....: eps=1, poles=[1], residues=[1]) sage: lf.init_empty()
- sage.lfunctions.pari.lfun_genus2(C)[source]¶
Return the
-function of a curve of genus 2.INPUT:
C– hyperelliptic curve of genus 2
Currently, the model needs to be minimal at 2.
This uses pari:lfungenus2.
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_genus2, LFunction sage: x = polygen(QQ, 'x') sage: C = HyperellipticCurve(x^5 + x + 1) sage: L = LFunction(lfun_genus2(C)) # this one is broken ... sage: L(3) 0.965946926261520 sage: C = HyperellipticCurve(x^2 + x, x^3 + x^2 + 1) sage: L = LFunction(lfun_genus2(C)) sage: L(2) 0.364286342944359
- sage.lfunctions.pari.lfun_hgm(motif, t)[source]¶
Create the
-function of an hypergeometric motive.OUTPUT: one pari:lfun object
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_hgm, LFunction sage: from sage.modular.hypergeometric_motive import HypergeometricData as Hyp sage: H = Hyp(gamma_list=([3,-1,-1,-1])) sage: L = LFunction(lfun_hgm(H, 1/5)) sage: L(3) 0.901925346034773
- sage.lfunctions.pari.lfun_modular_form(f)[source]¶
Return the
-function of the modular form .INPUT:
– a modular form
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_modular_form, LFunction sage: h = Newforms(37)[1] sage: L = LFunction(lfun_modular_form(h)) sage: L(1) 0.725681061936153
- sage.lfunctions.pari.lfun_number_field(K)[source]¶
Create the Dedekind zeta function of a number field.
OUTPUT: one pari:lfun object
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_number_field, LFunction sage: L = LFunction(lfun_number_field(QQ)) sage: L(3) 1.20205690315959 sage: K = NumberField(x**2 - 2, 'a') sage: L = LFunction(lfun_number_field(K)) sage: L(3) 1.15202784126080 sage: L(0) ...0.000000000000000
- sage.lfunctions.pari.lfun_quadratic_form(qf)[source]¶
Return the
-function of a positive definite quadratic form.This uses pari:lfunqf.
EXAMPLES:
sage: from sage.lfunctions.pari import lfun_quadratic_form, LFunction sage: Q = QuadraticForm(ZZ, 2, [2, 3, 4]) sage: L = LFunction(lfun_quadratic_form(Q)) sage: L(3) 0.377597233183583