Numerical computation of newforms

class sage.modular.modform.numerical.NumericalEigenforms(group, weight=2, eps=1e-20, delta=0.01, tp=[2, 3, 5])[source]

Bases: SageObject

numerical_eigenforms(group, weight=2, eps=1e-20, delta=1e-2, tp=[2,3,5]).

INPUT:

  • group – a congruence subgroup of a Dirichlet character of order 1 or 2

  • weight – integer >= 2

  • eps – a small float; abs( ) < eps is what “equal to zero” is interpreted as for floating point numbers

  • delta – a small-ish float; eigenvalues are considered distinct if their difference has absolute value at least delta

  • tp – use the Hecke operators T_p for p in tp when searching for a random Hecke operator with distinct Hecke eigenvalues

OUTPUT: a numerical eigenforms object, with the following useful methods:

  • ap() – return all eigenvalues of \(T_p\)

  • eigenvalues() – list of eigenvalues corresponding to the given list of primes, e.g.,:

    [[eigenvalues of T_2],
     [eigenvalues of T_3],
     [eigenvalues of T_5], ...]
    
  • systems_of_eigenvalues() – list of the systems of eigenvalues of eigenforms such that the chosen random linear combination of Hecke operators has multiplicity 1 eigenvalues.

EXAMPLES:

sage: n = numerical_eigenforms(23)
sage: n == loads(dumps(n))
True
sage: n.ap(2)  # abs tol 1e-12
[3.0, -1.6180339887498947, 0.6180339887498968]
sage: n.systems_of_eigenvalues(7)  # abs tol 2e-12
[
[-1.6180339887498947, 2.2360679774997894, -3.2360679774997894],
[0.6180339887498968, -2.236067977499788, 1.2360679774997936],
[3.0, 4.0, 6.0]
]
sage: n.systems_of_abs(7)  # abs tol 2e-12
[
[0.6180339887498943, 2.2360679774997894, 1.2360679774997887],
[1.6180339887498947, 2.23606797749979, 3.2360679774997894],
[3.0, 4.0, 6.0]
]
sage: n.eigenvalues([2,3,5])  # rel tol 2e-12
[[3.0, -1.6180339887498947, 0.6180339887498968],
 [4.0, 2.2360679774997894, -2.236067977499788],
 [6.0, -3.2360679774997894, 1.2360679774997936]]
>>> from sage.all import *
>>> n = numerical_eigenforms(Integer(23))
>>> n == loads(dumps(n))
True
>>> n.ap(Integer(2))  # abs tol 1e-12
[3.0, -1.6180339887498947, 0.6180339887498968]
>>> n.systems_of_eigenvalues(Integer(7))  # abs tol 2e-12
[
[-1.6180339887498947, 2.2360679774997894, -3.2360679774997894],
[0.6180339887498968, -2.236067977499788, 1.2360679774997936],
[3.0, 4.0, 6.0]
]
>>> n.systems_of_abs(Integer(7))  # abs tol 2e-12
[
[0.6180339887498943, 2.2360679774997894, 1.2360679774997887],
[1.6180339887498947, 2.23606797749979, 3.2360679774997894],
[3.0, 4.0, 6.0]
]
>>> n.eigenvalues([Integer(2),Integer(3),Integer(5)])  # rel tol 2e-12
[[3.0, -1.6180339887498947, 0.6180339887498968],
 [4.0, 2.2360679774997894, -2.236067977499788],
 [6.0, -3.2360679774997894, 1.2360679774997936]]
n = numerical_eigenforms(23)
n == loads(dumps(n))
n.ap(2)  # abs tol 1e-12
n.systems_of_eigenvalues(7)  # abs tol 2e-12
n.systems_of_abs(7)  # abs tol 2e-12
n.eigenvalues([2,3,5])  # rel tol 2e-12
ap(p)[source]

Return a list of the eigenvalues of the Hecke operator \(T_p\) on all the computed eigenforms. The eigenvalues match up between one prime and the next.

INPUT:

  • p – integer; a prime number

OUTPUT: list of double precision complex numbers

EXAMPLES:

sage: n = numerical_eigenforms(11,4)
sage: n.ap(2) # random order
[9.0, 9.0, 2.73205080757, -0.732050807569]
sage: n.ap(3) # random order
[28.0, 28.0, -7.92820323028, 5.92820323028]
sage: m = n.modular_symbols()
sage: x = polygen(QQ, 'x')
sage: m.T(2).charpoly('x').factor()
(x - 9)^2 * (x^2 - 2*x - 2)
sage: m.T(3).charpoly('x').factor()
(x - 28)^2 * (x^2 + 2*x - 47)
>>> from sage.all import *
>>> n = numerical_eigenforms(Integer(11),Integer(4))
>>> n.ap(Integer(2)) # random order
[9.0, 9.0, 2.73205080757, -0.732050807569]
>>> n.ap(Integer(3)) # random order
[28.0, 28.0, -7.92820323028, 5.92820323028]
>>> m = n.modular_symbols()
>>> x = polygen(QQ, 'x')
>>> m.T(Integer(2)).charpoly('x').factor()
(x - 9)^2 * (x^2 - 2*x - 2)
>>> m.T(Integer(3)).charpoly('x').factor()
(x - 28)^2 * (x^2 + 2*x - 47)
n = numerical_eigenforms(11,4)
n.ap(2) # random order
n.ap(3) # random order
m = n.modular_symbols()
x = polygen(QQ, 'x')
m.T(2).charpoly('x').factor()
m.T(3).charpoly('x').factor()
eigenvalues(primes)[source]

Return the eigenvalues of the Hecke operators corresponding to the primes in the input list of primes. The eigenvalues match up between one prime and the next.

INPUT:

  • primes – list of primes

OUTPUT: list of lists of eigenvalues

EXAMPLES:

sage: n = numerical_eigenforms(1,12)
sage: n.eigenvalues([3,5,13])  # rel tol 2.4e-10
[[177148.0, 252.00000000001896], [48828126.0, 4830.000000001376], [1792160394038.0, -577737.9999898539]]
>>> from sage.all import *
>>> n = numerical_eigenforms(Integer(1),Integer(12))
>>> n.eigenvalues([Integer(3),Integer(5),Integer(13)])  # rel tol 2.4e-10
[[177148.0, 252.00000000001896], [48828126.0, 4830.000000001376], [1792160394038.0, -577737.9999898539]]
n = numerical_eigenforms(1,12)
n.eigenvalues([3,5,13])  # rel tol 2.4e-10
level()[source]

Return the level of this set of modular eigenforms.

EXAMPLES:

sage: n = numerical_eigenforms(61) ; n.level()
61
>>> from sage.all import *
>>> n = numerical_eigenforms(Integer(61)) ; n.level()
61
n = numerical_eigenforms(61) ; n.level()
modular_symbols()[source]

Return the space of modular symbols used for computing this set of modular eigenforms.

EXAMPLES:

sage: n = numerical_eigenforms(61) ; n.modular_symbols()
Modular Symbols space of dimension 5 for Gamma_0(61) of weight 2 with sign 1 over Rational Field
>>> from sage.all import *
>>> n = numerical_eigenforms(Integer(61)) ; n.modular_symbols()
Modular Symbols space of dimension 5 for Gamma_0(61) of weight 2 with sign 1 over Rational Field
n = numerical_eigenforms(61) ; n.modular_symbols()
systems_of_abs(bound)[source]

Return the absolute values of all systems of eigenvalues for self for primes up to bound.

EXAMPLES:

sage: numerical_eigenforms(61).systems_of_abs(10)  # rel tol 1e-9
[
[0.3111078174659775, 2.903211925911551, 2.525427560843529, 3.214319743377552],
[1.0, 2.0000000000000027, 3.000000000000003, 1.0000000000000044],
[1.4811943040920152, 0.8060634335253695, 3.1563251746586642, 0.6751308705666477],
[2.170086486626034, 1.7092753594369208, 1.63089761381512, 0.46081112718908984],
[3.0, 4.0, 6.0, 8.0]
]
>>> from sage.all import *
>>> numerical_eigenforms(Integer(61)).systems_of_abs(Integer(10))  # rel tol 1e-9
[
[0.3111078174659775, 2.903211925911551, 2.525427560843529, 3.214319743377552],
[1.0, 2.0000000000000027, 3.000000000000003, 1.0000000000000044],
[1.4811943040920152, 0.8060634335253695, 3.1563251746586642, 0.6751308705666477],
[2.170086486626034, 1.7092753594369208, 1.63089761381512, 0.46081112718908984],
[3.0, 4.0, 6.0, 8.0]
]
numerical_eigenforms(61).systems_of_abs(10)  # rel tol 1e-9
systems_of_eigenvalues(bound)[source]

Return all systems of eigenvalues for self for primes up to bound.

EXAMPLES:

sage: numerical_eigenforms(61).systems_of_eigenvalues(10)  # rel tol 1e-9
[
[-1.4811943040920152, 0.8060634335253695, 3.1563251746586642, 0.6751308705666477],
[-1.0, -2.0000000000000027, -3.000000000000003, 1.0000000000000044],
[0.3111078174659775, 2.903211925911551, -2.525427560843529, -3.214319743377552],
[2.170086486626034, -1.7092753594369208, -1.63089761381512, -0.46081112718908984],
[3.0, 4.0, 6.0, 8.0]
]
>>> from sage.all import *
>>> numerical_eigenforms(Integer(61)).systems_of_eigenvalues(Integer(10))  # rel tol 1e-9
[
[-1.4811943040920152, 0.8060634335253695, 3.1563251746586642, 0.6751308705666477],
[-1.0, -2.0000000000000027, -3.000000000000003, 1.0000000000000044],
[0.3111078174659775, 2.903211925911551, -2.525427560843529, -3.214319743377552],
[2.170086486626034, -1.7092753594369208, -1.63089761381512, -0.46081112718908984],
[3.0, 4.0, 6.0, 8.0]
]
numerical_eigenforms(61).systems_of_eigenvalues(10)  # rel tol 1e-9
weight()[source]

Return the weight of this set of modular eigenforms.

EXAMPLES:

sage: n = numerical_eigenforms(61) ; n.weight()
2
>>> from sage.all import *
>>> n = numerical_eigenforms(Integer(61)) ; n.weight()
2
n = numerical_eigenforms(61) ; n.weight()
sage.modular.modform.numerical.support(v, eps)[source]

Given a vector \(v\) and a threshold eps, return all indices where \(|v|\) is larger than eps.

EXAMPLES:

sage: sage.modular.modform.numerical.support( numerical_eigenforms(61)._easy_vector(), 1.0 )
[]

sage: sage.modular.modform.numerical.support( numerical_eigenforms(61)._easy_vector(), 0.5 )
[0, 4]
>>> from sage.all import *
>>> sage.modular.modform.numerical.support( numerical_eigenforms(Integer(61))._easy_vector(), RealNumber('1.0') )
[]

>>> sage.modular.modform.numerical.support( numerical_eigenforms(Integer(61))._easy_vector(), RealNumber('0.5') )
[0, 4]
sage.modular.modform.numerical.support( numerical_eigenforms(61)._easy_vector(), 1.0 )
sage.modular.modform.numerical.support( numerical_eigenforms(61)._easy_vector(), 0.5 )