Hot-keys on this page
r m x p toggle line displays
j k next/prev highlighted chunk
0 (zero) top of page
1 (one) first highlighted chunk
# -*- coding: utf-8 -*- ######################################################################### # Copyright (C) 2011 Cameron Franc and Marc Masdeu # # Distributed under the terms of the GNU General Public License (GPL) # # http://www.gnu.org/licenses/ ######################################################################### r""" Quotients of the Bruhat-Tits tree
This package contains all the functionality described and developed in [FM]_. It allows for computations with fundamental domains of the Bruhat-Tits tree, under the action of arithmetic groups arising from units in definite quaternion algebras.
EXAMPLES:
Create the quotient attached to a maximal order of the quaternion algebra of discriminant `13`, at the prime `p = 5`. ::
sage: Y = BruhatTitsQuotient(5, 13)
We can query for its genus, as well as get it back as a graph::
sage: Y.genus() 5 sage: Y.get_graph() Multi-graph on 2 vertices
The rest of functionality can be found in the docstrings below.
REFERENCES:
.. [FM] Computing fundamental domains for the Bruhat-Tits tree for `\textrm{GL}_2(\QQ_p)`, `p`-adic automorphic forms, and the canonical embedding of Shimura curves Cameron Franc, Marc Masdeu LMS Journal of Computation and Mathematics (2014), volume 17, issue 01, pp. 1-23. """
from __future__ import print_function from __future__ import absolute_import from sage.rings.integer import Integer from sage.matrix.constructor import Matrix from sage.matrix.matrix_space import MatrixSpace from sage.structure.sage_object import SageObject from sage.rings.all import ZZ, Zmod, QQ from sage.misc.latex import latex from sage.rings.padics.precision_error import PrecisionError import collections from sage.misc.misc_c import prod from sage.structure.unique_representation import UniqueRepresentation from sage.misc.cachefunc import cached_method from sage.arith.all import gcd, xgcd, kronecker_symbol, fundamental_discriminant from sage.rings.padics.all import Qp, Zp from sage.rings.finite_rings.finite_field_constructor import GF from sage.algebras.quatalg.all import QuaternionAlgebra from sage.quadratic_forms.all import QuadraticForm from sage.graphs.all import Graph from sage.libs.all import pari from sage.interfaces.all import magma from copy import copy from sage.plot.colors import rainbow from sage.rings.number_field.all import NumberField from sage.modular.arithgroup.all import Gamma0 from sage.misc.lazy_attribute import lazy_attribute from sage.modular.dirichlet import DirichletGroup from sage.modular.arithgroup.congroup_gammaH import GammaH_constructor from sage.misc.misc import verbose
class DoubleCosetReduction(SageObject): r""" Edges in the Bruhat-Tits tree are represented by cosets of matrices in `GL_2`. Given a matrix `x` in `GL_2`, this class computes and stores the data corresponding to the double coset representation of `x` in terms of a fundamental domain of edges for the action of the arithmetic group `\Gamma`.
More precisely:
Initialized with an element `x` of `GL_2(\ZZ)`, finds elements `\gamma` in `\Gamma`, `t` and an edge `e` such that `get=x`. It stores these values as members ``gamma``, ``label`` and functions ``self.sign()``, ``self.t()`` and ``self.igamma()``, satisfying:
- if ``self.sign() == +1``: ``igamma() * edge_list[label].rep * t() == x``
- if ``self.sign() == -1``: ``igamma() * edge_list[label].opposite.rep * t() == x``
It also stores a member called power so that:
``p**(2*power) = gamma.reduced_norm()``
The usual decomposition `get=x` would be:
- g = gamma / (p ** power)
- e = edge_list[label]
- t' = t * p ** power
Here usual denotes that we have rescaled gamma to have unit determinant, and so that the result is honestly an element of the arithmetic quaternion group under consideration. In practice we store integral multiples and keep track of the powers of `p`.
INPUT:
- ``Y`` - BruhatTitsQuotient object in which to work - ``x`` - Something coercible into a matrix in `GL_2(\ZZ)`. In principle we should allow elements in `GL_2(\QQ_p)`, but it is enough to work with integral entries - ``extrapow`` - gets added to the power attribute, and it is used for the Hecke action.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import DoubleCosetReduction sage: Y = BruhatTitsQuotient(5, 13) sage: x = Matrix(ZZ,2,2,[123,153,1231,1231]) sage: d = DoubleCosetReduction(Y,x) sage: d.sign() -1 sage: d.igamma()*Y._edge_list[d.label - len(Y.get_edge_list())].opposite.rep*d.t() == x True sage: x = Matrix(ZZ,2,2,[1423,113553,11231,12313]) sage: d = DoubleCosetReduction(Y,x) sage: d.sign() 1 sage: d.igamma()*Y._edge_list[d.label].rep*d.t() == x True
AUTHORS:
- Cameron Franc (2012-02-20) - Marc Masdeu """
def __init__(self, Y, x, extrapow=0): r""" Initialize and compute the reduction as a double coset.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import DoubleCosetReduction sage: Y = BruhatTitsQuotient(5, 13) sage: x = Matrix(ZZ,2,2,[123,153,1231,1231]) sage: d = DoubleCosetReduction(Y,x) sage: TestSuite(d).run() """ valuation=valuation)
# The label will encode whether it is an edge or its opposite !
def _repr_(self): r""" Return the representation of self as a string.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import DoubleCosetReduction sage: Y = BruhatTitsQuotient(5, 13) sage: x = Matrix(ZZ,2,2,[123,153,1231,1231]) sage: DoubleCosetReduction(Y,x) Double coset data (-1, [(4), (5), (-4), (-4)], 8) """
def __eq__(self, other): """ Return self == other
TESTS::
sage: from sage.modular.btquotients.btquotient import DoubleCosetReduction sage: Y = BruhatTitsQuotient(5, 13) sage: x = Matrix(ZZ,2,2,[123,153,1231,1231]) sage: d1 = DoubleCosetReduction(Y,x) sage: d1 == d1 True """ return False return False return False return False return False return False return False return False return False
def __ne__(self, other): """ Return self != other
TESTS::
sage: from sage.modular.btquotients.btquotient import DoubleCosetReduction sage: Y = BruhatTitsQuotient(5, 13) sage: x = Matrix(ZZ,2,2,[123,153,1231,1231]) sage: d1 = DoubleCosetReduction(Y,x) sage: d1 != d1 False """
def sign(self): r""" The direction of the edge.
The Bruhat Tits quotients are directed graphs but we only store half the edges (we treat them more like unordered graphs). The sign tells whether the matrix self.x is equivalent to the representative in the quotient (sign = +1), or to the opposite of one of the representatives (sign = -1).
OUTPUT :
an int that is +1 or -1 according to the sign of self
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import DoubleCosetReduction sage: Y = BruhatTitsQuotient(3, 11) sage: x = Matrix(ZZ,2,2,[123,153,1231,1231]) sage: d = DoubleCosetReduction(Y,x) sage: d.sign() -1 sage: d.igamma()*Y._edge_list[d.label - len(Y.get_edge_list())].opposite.rep*d.t() == x True sage: x = Matrix(ZZ,2,2,[1423,113553,11231,12313]) sage: d = DoubleCosetReduction(Y,x) sage: d.sign() 1 sage: d.igamma()*Y._edge_list[d.label].rep*d.t() == x True """ else:
def igamma(self, embedding=None, scale=1): r""" Image under gamma.
Elements of the arithmetic group can be regarded as elements of the global quaternion order, and hence may be represented exactly. This function computes the image of such an element under the local splitting and returns the corresponding `p`-adic approximation.
INPUT:
- ``embedding`` - an integer, or a function (default: none). If ``embedding`` is None, then the image of ``self.gamma`` under the local splitting associated to ``self.Y`` is used. If ``embedding`` is an integer, then the precision of the local splitting of self.Y is raised (if necessary) to be larger than this integer, and this new local splitting is used. If a function is passed, then map ``self.gamma`` under ``embedding``. - ``scale`` -- (default: 1) scaling factor applied to the output
OUTPUT:
a 2x2 matrix with `p`-adic entries encoding the image of ``self`` under the local splitting
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import DoubleCosetReduction sage: Y = BruhatTitsQuotient(7, 11) sage: d = DoubleCosetReduction(Y,Matrix(ZZ,2,2,[123,45,88,1])) sage: d.igamma() [6 + 6*7 + 6*7^2 + 6*7^3 + 6*7^4 + O(7^5) O(7^5)] [ O(7^5) 6 + 6*7 + 6*7^2 + 6*7^3 + 6*7^4 + O(7^5)] sage: d.igamma(embedding = 7) [6 + 6*7 + 6*7^2 + 6*7^3 + 6*7^4 + 6*7^5 + 6*7^6 + O(7^7) O(7^7)] [ O(7^7) 6 + 6*7 + 6*7^2 + 6*7^3 + 6*7^4 + 6*7^5 + 6*7^6 + O(7^7)] """ else: # The user wants higher precision # The user knows what she is doing, so let it go prec=prec)
def t(self, prec=None): r""" Return the 't part' of the decomposition using the rest of the data.
INPUT:
- ``prec`` - a `p`-adic precision that t will be computed to. Defaults to the default working precision of self.
OUTPUT:
a 2x2 `p`-adic matrix with entries of precision ``prec`` that is the 't-part' of the decomposition of self
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import DoubleCosetReduction sage: Y = BruhatTitsQuotient(5, 13) sage: x = Matrix(ZZ,2,2,[123,153,1231,1232]) sage: d = DoubleCosetReduction(Y,x) sage: t = d.t(20) sage: t[1,0].valuation() > 0 True """ # assert self._cached_t[1, 0].valuation()>self._cached_t[1,1].valuation() else: # assert self._cached_t[1, 0].valuation()>self._cached_t[1,1].valuation() for xx in self._cached_t.list()])
class BruhatTitsTree(SageObject, UniqueRepresentation): r""" An implementation of the Bruhat-Tits tree for `GL_2(\QQ_p)`.
INPUT:
- ``p`` - a prime number. The corresponding tree is then `p+1` regular
EXAMPLES:
We create the tree for `GL_2(\QQ_5)`::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: p = 5 sage: T = BruhatTitsTree(p) sage: m = Matrix(ZZ,2,2,[p**5,p**2,p**3,1+p+p*3]) sage: e = T.edge(m); e [ 0 25] [625 21] sage: v0 = T.origin(e); v0 [ 25 0] [ 21 125] sage: v1 = T.target(e); v1 [ 25 0] [ 21 625] sage: T.origin(T.opposite(e)) == v1 True sage: T.target(T.opposite(e)) == v0 True
A value error is raised if a prime is not passed::
sage: T = BruhatTitsTree(4) Traceback (most recent call last): ... ValueError: Input (4) must be prime
AUTHORS:
- Marc Masdeu (2012-02-20) """ def __init__(self, p): """ Initialize a BruhatTitsTree object for a given prime `p`
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: T = BruhatTitsTree(17) sage: TestSuite(T).run() """
def target(self, e, normalized=False): r""" Return the target vertex of the edge represented by the input matrix e.
INPUT:
- ``e`` - a 2x2 matrix with integer entries
- ``normalized`` - boolean (default: false). If True then the input matrix is assumed to be normalized.
OUTPUT:
- ``e`` - 2x2 integer matrix representing the target of the input edge
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: T = BruhatTitsTree(7) sage: T.target(Matrix(ZZ,2,2,[1,5,8,9])) [1 0] [0 1] """ #then the normalized target vertex is also M and we save some #row reductions with a simple return return e else: #must normalize the target vertex representative
def origin(self, e, normalized=False): r""" Return the origin vertex of the edge represented by the input matrix e.
INPUT:
- ``e`` - a 2x2 matrix with integer entries
- ``normalized`` - boolean (default: false). If True then the input matrix M is assumed to be normalized
OUTPUT:
- ``e`` - A 2x2 integer matrix
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: T = BruhatTitsTree(7) sage: T.origin(Matrix(ZZ,2,2,[1,5,8,9])) [1 0] [1 7] """ #then normalize else: x = copy(e)
def edge(self, M): r""" Normalize a matrix to the correct normalized edge representative.
INPUT:
- ``M`` - a 2x2 integer matrix
OUTPUT:
- ``newM`` - a 2x2 integer matrix
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: T = BruhatTitsTree(3) sage: T.edge( Matrix(ZZ,2,2,[0,-1,3,0]) ) [0 1] [3 0] """ # M_orig = M
""" Naively approximate a p-adic integer by a positive integer.
INPUT:
- ``a`` - a p-adic integer.
OUTPUT:
An integer.
EXAMPLES::
sage: x = Zp(3)(-17) sage: lift(x) 3486784384 """ except AttributeError: return ZZ(a)
else: # assert self.is_in_group(M_orig.inverse()*newM, as_edge = True)
def vertex(self, M): r""" Normalize a matrix to the corresponding normalized vertex representative
INPUT:
- ``M`` - 2x2 integer matrix
OUTPUT:
- a 2x2 integer matrix
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: p = 5 sage: T = BruhatTitsTree(p) sage: m = Matrix(ZZ,2,2,[p**5,p**2,p**3,1+p+p*3]) sage: e = T.edge(m) sage: t = m.inverse()*e sage: scaling = Qp(p,20)(t.determinant()).sqrt() sage: t = 1/scaling * t sage: min([t[ii,jj].valuation(p) for ii in range(2) for jj in range(2)]) >= 0 True sage: t[1,0].valuation(p) > 0 True """ # M_orig = M
except AttributeError: return ZZ(a)
# r = ZZ(r) % bigpower # assert self.is_in_group(M_orig.inverse()*newM, as_edge=False)
def edges_leaving_origin(self): r""" Find normalized representatives for the `p+1` edges leaving the origin vertex corresponding to the homothety class of `\ZZ_p^2`. These are cached.
OUTPUT:
- A list of size `p+1` of 2x2 integer matrices
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: T = BruhatTitsTree(3) sage: T.edges_leaving_origin() [ [0 1] [3 0] [0 1] [0 1] [3 0], [0 1], [3 1], [3 2] ] """
def edge_between_vertices(self, v1, v2, normalized=False): r""" Compute the normalized matrix rep. for the edge passing between two vertices.
INPUT:
- ``v1`` - 2x2 integer matrix
- ``v2`` - 2x2 integer matrix
- ``normalized`` - boolean (default: False), whether the vertices are normalized.
OUTPUT:
- 2x2 integer matrix, representing the edge from ``v1`` to ``v2``. If ``v1`` and ``v2`` are not at distance `1`, raise a ``ValueError``.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: p = 7 sage: T = BruhatTitsTree(p) sage: v1 = T.vertex(Matrix(ZZ,2,2,[p,0,0,1])); v1 [7 0] [0 1] sage: v2 = T.vertex(Matrix(ZZ,2,2,[p,1,0,1])); v2 [1 0] [1 7] sage: T.edge_between_vertices(v1,v2) Traceback (most recent call last): ... ValueError: Vertices are not adjacent.
sage: v3 = T.vertex(Matrix(ZZ,2,2,[1,0,0,1])); v3 [1 0] [0 1] sage: T.edge_between_vertices(v1,v3) [0 1] [1 0] """ v22 = v2 else:
def leaving_edges(self, M): r""" Return edges leaving a vertex
INPUT:
- ``M`` - 2x2 integer matrix
OUTPUT:
List of size `p+1` of 2x2 integer matrices
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: p = 7 sage: T = BruhatTitsTree(p) sage: T.leaving_edges(Matrix(ZZ,2,2,[1,0,0,1])) [ [0 1] [7 0] [0 1] [0 1] [0 1] [0 1] [0 1] [0 1] [7 0], [0 1], [7 1], [7 4], [7 5], [7 2], [7 3], [7 6] ] """
def opposite(self, e): r""" This function returns the edge oriented oppositely to a given edge.
INPUT:
- ``e`` - 2x2 integer matrix
OUTPUT:
2x2 integer matrix
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: p = 7 sage: T = BruhatTitsTree(p) sage: e = Matrix(ZZ,2,2,[1,0,0,1]) sage: T.opposite(e) [0 1] [7 0] sage: T.opposite(T.opposite(e)) == e True """
def entering_edges(self, v): r""" This function returns the edges entering a given vertex.
INPUT:
- ``v`` - 2x2 integer matrix
OUTPUT:
A list of size `p+1` of 2x2 integer matrices
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: p = 7 sage: T = BruhatTitsTree(p) sage: T.entering_edges(Matrix(ZZ,2,2,[1,0,0,1])) [ [1 0] [0 1] [1 0] [1 0] [1 0] [1 0] [1 0] [1 0] [0 1], [1 0], [1 1], [4 1], [5 1], [2 1], [3 1], [6 1] ] """
def subdivide(self, edgelist, level): r""" (Ordered) edges of self may be regarded as open balls in `P^1(\QQ_p)`. Given a list of edges, this function return a list of edges corresponding to the level-th subdivision of the corresponding opens. That is, each open ball of the input is broken up into `p^{\mbox{level}}` subballs of equal radius.
INPUT:
- ``edgelist`` - a list of edges
- ``level`` - an integer
OUTPUT:
A list of 2x2 integer matrices
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: p = 3 sage: T = BruhatTitsTree(p) sage: T.subdivide([Matrix(ZZ,2,2,[p,0,0,1])],2) [ [27 0] [0 9] [0 9] [0 3] [0 3] [0 3] [0 3] [0 3] [0 3] [ 0 1], [3 1], [3 2], [9 1], [9 4], [9 7], [9 2], [9 5], [9 8] ] """ return [] else:
def get_balls(self, center=1, level=1): r""" Return a decomposition of `P^1(\QQ_p)` into compact open balls.
Each vertex in the Bruhat-Tits tree gives a decomposition of `P^1(\QQ_p)` into `p+1` open balls. Each of these balls may be further subdivided, to get a finer decomposition.
This function returns the decomposition of `P^1(\QQ_p)` corresponding to ``center`` into `(p+1)p^{\mbox{level}}` balls.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: p = 2 sage: T = BruhatTitsTree(p) sage: T.get_balls(Matrix(ZZ,2,2,[p,0,0,1]),1) [ [0 1] [0 1] [8 0] [0 4] [0 2] [0 2] [2 0], [2 1], [0 1], [2 1], [4 1], [4 3] ] """
def find_path(self, v, boundary=None): r""" Compute a path from a vertex to a given set of so-called boundary vertices, whose interior must contain the origin vertex. In the case that the boundary is not specified, it computes the geodesic between the given vertex and the origin. In the case that the boundary contains more than one vertex, it computes the geodesic to some point of the boundary.
INPUT:
- ``v`` - a 2x2 matrix representing a vertex ``boundary``
- a list of matrices (default: None). If ommitted, finds the geodesic from ``v`` to the central vertex.
OUTPUT:
An ordered list of vertices describing the geodesic from ``v`` to ``boundary``, followed by the vertex in the boundary that is closest to ``v``.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: p = 3 sage: T = BruhatTitsTree(p) sage: T.find_path( Matrix(ZZ,2,2,[p^4,0,0,1]) ) ( [[81 0] [ 0 1], [27 0] [ 0 1], [9 0] [0 1], [3 0] [1 0] [0 1]] , [0 1] ) sage: T.find_path( Matrix(ZZ,2,2,[p^3,0,134,p^2]) ) ( [[27 0] [ 8 9], [27 0] [ 2 3], [27 0] [ 0 1], [9 0] [0 1], [3 0] [1 0] [0 1]] , [0 1] ) """
raise RuntimeError
def find_containing_affinoid(self, z): r""" Return the vertex corresponding to the affinoid in the `p`-adic upper half plane that a given (unramified!) point reduces to.
INPUT:
- ``z`` - an element of an unramified extension of `\QQ_p` that is not contained in `\QQ_p`.
OUTPUT:
A 2x2 integer matrix representing a vertex of ``self``.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: T = BruhatTitsTree(5) sage: K.<a> = Qq(5^2,20) sage: T.find_containing_affinoid(a) [1 0] [0 1] sage: z = 5*a+3 sage: v = T.find_containing_affinoid(z).inverse(); v [ 1 0] [-2/5 1/5]
Note that the translate of ``z`` belongs to the standard affinoid. That is, it is a `p`-adic unit and its reduction modulo `p` is not in `\mathbb{F}_p`::
sage: gz = (v[0,0]*z+v[0,1])/(v[1,0]*z+v[1,1]); gz (a + 1) + O(5^19) sage: gz.valuation() == 0 True """ #Assume z belongs to some extension of QQp. return self.vertex(self._Mat_22([0, 1, p, 0]) * self.find_containing_affinoid(1 / (p * z))) L.append(0)
def find_geodesic(self, v1, v2, normalized=True): r""" This function computes the geodesic between two vertices
INPUT:
- ``v1`` - 2x2 integer matrix representing a vertex
- ``v2`` - 2x2 integer matrix representing a vertex
- ``normalized`` - boolean (default: True)
OUTPUT:
An ordered list of 2x2 integer matrices representing the vertices of the paths joining ``v1`` and ``v2``.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: p = 3 sage: T = BruhatTitsTree(p) sage: v1 = T.vertex( Matrix(ZZ,2,2,[p^3, 0, 1, p^1]) ); v1 [27 0] [ 1 3] sage: v2 = T.vertex( Matrix(ZZ,2,2,[p,2,0,p]) ); v2 [1 0] [6 9] sage: T.find_geodesic(v1,v2) [ [27 0] [27 0] [9 0] [3 0] [1 0] [1 0] [1 0] [ 1 3], [ 0 1], [0 1], [0 1], [0 1], [0 3], [6 9] ] """ v1, v2 = self.vertex(v1), self.vertex(v2)
def find_covering(self, z1, z2, level=0): r""" Compute a covering of `P^1(\QQ_p)` adapted to a certain geodesic in self.
More precisely, the `p`-adic upper half plane points ``z1`` and ``z2`` reduce to vertices `v_1`, `v_2`. The returned covering consists of all the edges leaving the geodesic from `v_1` to `v_2`.
INPUT:
- ``z1``, ``z2`` - unramified algebraic points of h_p
OUTPUT:
a list of 2x2 integer matrices representing edges of self
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import BruhatTitsTree sage: p = 3 sage: K.<a> = Qq(p^2) sage: T = BruhatTitsTree(p) sage: z1 = a + a*p sage: z2 = 1 + a*p + a*p^2 - p^6 sage: T.find_covering(z1,z2) [ [0 1] [3 0] [0 1] [0 1] [0 1] [0 1] [3 0], [0 1], [3 2], [9 1], [9 4], [9 7] ]
.. NOTE::
This function is used to compute certain Coleman integrals on `P^1`. That's why the input consists of two points of the `p`-adic upper half plane, but decomposes `P^1(\QQ_p)`. This decomposition is what allows us to represent the relevant integrand as a locally analytic function. The ``z1`` and ``z2`` appear in the integrand. """ # m = vv.determinant().valuation(self._p)
class Vertex(SageObject): r""" This is a structure to represent vertices of quotients of the Bruhat-Tits tree. It is useful to enrich the representation of the vertex as a matrix with extra data.
INPUT:
- ``p`` - a prime integer.
- ``label`` - An integer which uniquely identifies this vertex.
- ``rep`` - A 2x2 matrix in reduced form representing this vertex.
- ``leaving_edges`` - (default: empty list) A list of edges leaving this vertex.
- ``entering_edges`` - (default: empty list) A list of edges entering this vertex.
- ``determinant`` - (default: None) The determinant of ``rep``, if known.
- ``valuation`` - (default: None) The valuation of the determinant of ``rep``, if known.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import Vertex sage: v1 = Vertex(5,0,Matrix(ZZ,2,2,[1,2,3,18])) sage: v1.rep [ 1 2] [ 3 18] sage: v1.entering_edges []
AUTHORS:
- Marc Masdeu (2012-02-20) """ def __init__(self, p, label, rep, leaving_edges=None, entering_edges=None, determinant=None, valuation=None): """ This initializes a structure to represent vertices of quotients of the Bruhat-Tits tree. It is useful to enrich the representation of the vertex as a matrix with extra data.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import Vertex sage: Y = BruhatTitsQuotient(5,13) sage: v1 = Vertex(5,0,Matrix(ZZ,2,2,[1,2,3,18])) sage: TestSuite(v1).run() """
def _repr_(self): r""" Return the representation of self as a string.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,5) sage: X.get_vertex_list()[0] Vertex of Bruhat-Tits tree for p = 3 """
def __eq__(self, other): """ Return self == other
TESTS::
sage: from sage.modular.btquotients.btquotient import Vertex sage: v1 = Vertex(7,0,Matrix(ZZ,2,2,[1,2,3,18])) sage: v1 == v1 True """ return False return False return False return False return False return False
def __ne__(self, other): """ Return self != other
TESTS::
sage: from sage.modular.btquotients.btquotient import Vertex sage: v1 = Vertex(7,0,Matrix(ZZ,2,2,[1,2,3,18])) sage: v1 != v1 False """
class Edge(SageObject): r""" This is a structure to represent edges of quotients of the Bruhat-Tits tree. It is useful to enrich the representation of an edge as a matrix with extra data.
INPUT:
- ``p`` - a prime integer.
- ``label`` - An integer which uniquely identifies this edge.
- ``rep`` - A 2x2 matrix in reduced form representing this edge.
- ``origin`` - The origin vertex of ``self``.
- ``target`` - The target vertex of ``self``.
- ``links`` - (Default: empty list) A list of elements of `\Gamma` which identify different edges in the Bruhat-Tits tree which are equivalent to ``self``.
- ``opposite`` - (Default: None) The edge opposite to ``self``
- ``determinant`` - (Default: None) The determinant of ``rep``, if known.
- ``valuation`` - (Default: None) The valuation of the determinant of ``rep``, if known.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import Edge, Vertex sage: v1 = Vertex(7,0,Matrix(ZZ,2,2,[1,2,3,18])) sage: v2 = Vertex(7,0,Matrix(ZZ,2,2,[3,2,1,18])) sage: e1 = Edge(7,0,Matrix(ZZ,2,2,[1,2,3,18]),v1,v2) sage: e1.rep [ 1 2] [ 3 18]
AUTHORS:
- Marc Masdeu (2012-02-20) """ def __init__(self, p, label, rep, origin, target, links=None, opposite=None, determinant=None, valuation=None): """ Representation for edges of quotients of the Bruhat-Tits tree. It is useful to enrich the representation of an edge as a matrix with extra data.
EXAMPLES::
sage: from sage.modular.btquotients.btquotient import Edge sage: Y = BruhatTitsQuotient(5,11) sage: el = Y.get_edge_list() sage: e1 = el.pop() sage: e2 = Edge(5,e1.label,e1.rep,e1.origin,e1.target) sage: TestSuite(e2).run() """
def _repr_(self): r""" Return the representation of self as a string.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,5) sage: X.get_edge_list()[0] Edge of Bruhat-Tits tree for p = 3 """
def __eq__(self, other): """ Return self == other
TESTS::
sage: from sage.modular.btquotients.btquotient import Edge,Vertex sage: v1 = Vertex(7,0,Matrix(ZZ,2,2,[1,2,3,18])) sage: v2 = Vertex(7,0,Matrix(ZZ,2,2,[3,2,1,18])) sage: e1 = Edge(7,0,Matrix(ZZ,2,2,[1,2,3,18]),v1,v2) sage: e1 == e1 True """ return False return False return False return False return False return False return False return False return False
def __ne__(self, other): """ Return self != other
TESTS::
sage: from sage.modular.btquotients.btquotient import Edge,Vertex sage: v1 = Vertex(7,0,Matrix(ZZ,2,2,[1,2,3,18])) sage: v2 = Vertex(7,0,Matrix(ZZ,2,2,[3,2,1,18])) sage: e1 = Edge(7,0,Matrix(ZZ,2,2,[1,2,3,18]),v1,v2) sage: e1 != e1 False """
class BruhatTitsQuotient(SageObject, UniqueRepresentation): r""" This function computes the quotient of the Bruhat-Tits tree by an arithmetic quaternionic group. The group in question is the group of norm 1 elements in an Eichler `\ZZ[1/p]`-order of some (tame) level inside of a definite quaternion algebra that is unramified at the prime `p`. Note that this routine relies in Magma in the case `p = 2` or when `N^{+} > 1`.
INPUT:
- ``p`` - a prime number
- ``Nminus`` - squarefree integer divisible by an odd number of distinct primes and relatively prime to p. This is the discriminant of the definite quaternion algebra that one is quotienting by.
- ``Nplus`` - an integer coprime to pNminus (Default: 1). This is the tame level. It need not be squarefree! If Nplus is not 1 then the user currently needs magma installed due to sage's inability to compute well with nonmaximal Eichler orders in rational (definite) quaternion algebras.
- ``character`` - a Dirichlet character (Default: None) of modulus `pN^-N^+`.
- ``use_magma`` - boolean (default: False). If True, uses Magma for quaternion arithmetic.
- ``magma_session`` -- (default: None). If specified, the Magma session to use.
EXAMPLES:
Here is an example without a Dirichlet character::
sage: X = BruhatTitsQuotient(13, 19) sage: X.genus() 19 sage: G = X.get_graph(); G Multi-graph on 4 vertices
And an example with a Dirichlet character::
sage: f = DirichletGroup(6)[1] sage: X = BruhatTitsQuotient(3,2*5*7,character = f) sage: X.genus() 5
.. NOTE::
A sage implementation of Eichler orders in rational quaternions algebras would remove the dependency on magma.
AUTHORS:
- Marc Masdeu (2012-02-20) """ @staticmethod def __classcall__(cls, p, Nminus, Nplus=1, character=None, use_magma=False, seed=None, magma_session=None): """ Ensure that a canonical BruhatTitsQuotient is created.
EXAMPLES:
sage: BruhatTitsQuotient(3,17) is BruhatTitsQuotient(3,17,1) True """ character, use_magma, seed, magma_session)
def __init__(self, p, Nminus, Nplus=1, character=None, use_magma=False, seed=None, magma_session=None): """ Compute the quotient of the Bruhat-Tits tree by an arithmetic quaternionic group.
EXAMPLES::
sage: Y = BruhatTitsQuotient(19,11) sage: TestSuite(Y).run() """ raise ValueError("character must be of squarefree conductor") else:
raise ValueError("p must be a prime") raise ValueError("level must be squarefree") raise ValueError("level and conductor must be coprime")
# if len(Nminus.factor()) % 2 != 1: # raise ValueError("Nminus should be divisible by an odd number of primes")
try: if magma_session is None: self._magma = magma else: self._magma = magma_session magmap = self._magma(p) # print("Warning: this input needs magma to work...") except RuntimeError: raise NotImplementedError('Sage does not know yet how to work with the kind of orders that you are trying to use. Try installing Magma first and set it up so that Sage can use it.')
## This is added for debugging, in order to have reproducible results if seed is not None: self._magma.function_call('SetSeed', seed, nvals=0) self._use_magma = True else:
else: self._Mat_22([0, 1, 0, 0]), self._Mat_22([0, 0, 1, 0]), self._Mat_22([0, 0, 0, 1])] self._Mat_22([0, 1, 0, 0]), self._Mat_22([0, 0, self._p, 0]), self._Mat_22([0, 0, 0, 1])]
def _cache_key(self): r""" Return a hash of self, for using in caching.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,13) sage: X._cache_key() -406423199 # 32-bit 1375458358400022881 # 64-bit
sage: Y = BruhatTitsQuotient(5,13,use_magma = True) # optional - magma sage: Y._cache_key() == X._cache_key() # optional - magma False """
def _repr_(self): r""" Return the representation of self as a string.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,13); X Quotient of the Bruhat Tits tree of GL_2(QQ_5) with discriminant 13 and level 1 """
def __eq__(self, other): r""" Compare self with other.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,13) sage: Y = BruhatTitsQuotient(p = 5, Nminus = 13, Nplus=1,seed = 1231) sage: X == Y True """ return False return False return False return False
def __ne__(self, other): r""" Compare self with other.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,13) sage: Y = BruhatTitsQuotient(p = 5, Nminus = 13, Nplus=1,seed = 1231) sage: X != Y False """
def _latex_(self): r""" Return the LaTeX representation of ``self``.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,13); latex(X) X(5 \cdot 13,1)\otimes_{\mathbb{Z}} \mathbb{F}_{5} """
def get_vertex_dict(self): r""" This function returns the vertices of the quotient viewed as a dict.
OUTPUT:
A python dict with the vertices of the quotient.
EXAMPLES::
sage: X = BruhatTitsQuotient(37,3) sage: X.get_vertex_dict() {[1 0] [0 1]: Vertex of Bruhat-Tits tree for p = 37, [ 1 0] [ 0 37]: Vertex of Bruhat-Tits tree for p = 37} """
def get_vertex_list(self): r""" Return a list of the vertices of the quotient.
OUTPUT:
- A list with the vertices of the quotient.
EXAMPLES::
sage: X = BruhatTitsQuotient(37,3) sage: X.get_vertex_list() [Vertex of Bruhat-Tits tree for p = 37, Vertex of Bruhat-Tits tree for p = 37] """
def get_edge_list(self): r""" Return a list of ``Edge`` which represent a fundamental domain inside the Bruhat-Tits tree for the quotient.
OUTPUT:
A list of ``Edge``.
EXAMPLES::
sage: X = BruhatTitsQuotient(37,3) sage: len(X.get_edge_list()) 8
"""
def get_list(self): r""" Return a list of ``Edge`` which represent a fundamental domain inside the Bruhat-Tits tree for the quotient, together with a list of the opposite edges. This is used to work with automorphic forms.
OUTPUT:
A list of ``Edge``.
EXAMPLES::
sage: X = BruhatTitsQuotient(37,3) sage: len(X.get_list()) 16 """
def get_nontorsion_generators(self): r""" Use a fundamental domain in the Bruhat-Tits tree, and certain gluing data for boundary vertices, in order to compute a collection of generators for the nontorsion part of the arithmetic quaternionic group that one is quotienting by. This is analogous to using a polygonal rep. of a compact real surface to present its fundamental domain.
OUTPUT:
- A generating list of elements of an arithmetic quaternionic group.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,13) sage: len(X.get_nontorsion_generators()) 3 """
@cached_method def get_generators(self): r""" Use a fundamental domain in the Bruhat-Tits tree, and certain gluing data for boundary vertices, in order to compute a collection of generators for the arithmetic quaternionic group that one is quotienting by. This is analogous to using a polygonal rep. of a compact real surface to present its fundamental domain.
OUTPUT:
- A generating list of elements of an arithmetic quaternionic group.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,2) sage: len(X.get_generators()) 2 """
def _compute_invariants(self): """ Compute certain invariants from the level data of the quotient which allow one to compute the genus of the curve.
Details to be found in Theorem 9 of [FM]_.
EXAMPLES::
sage: X = BruhatTitsQuotient(23,11) sage: X._compute_invariants() """ if (f[1] == 1): e4 *= (1 + kronecker_symbol(-4, Integer(f[0]))) e3 *= (1 + kronecker_symbol(-3, Integer(f[0]))) else: if kronecker_symbol(-4, Integer(f[0])) == 1: e4 *= 2 else: e4 = 0 if kronecker_symbol(-3, Integer(f[0])) == 1: e3 *= 2 else: e3 = 0 mu *= 1 + 1 / Integer(f[0])
@lazy_attribute def e3(self): r""" Compute the `e_3` invariant defined by the formula
.. MATH::
e_k =\prod_{\ell\mid pN^-}\left(1-\left(\frac{-3}{\ell}\right)\right)\prod_{\ell \| N^+}\left(1+\left(\frac{-3}{\ell}\right)\right)\prod_{\ell^2\mid N^+} \nu_\ell(3)
OUTPUT:
an integer
EXAMPLES::
sage: X = BruhatTitsQuotient(31,3) sage: X.e3 1 """
@lazy_attribute def e4(self): r""" Compute the `e_4` invariant defined by the formula
.. MATH::
e_k =\prod_{\ell\mid pN^-}\left(1-\left(\frac{-k}{\ell}\right)\right)\prod_{\ell \| N^+}\left(1+\left(\frac{-k}{\ell}\right)\right)\prod_{\ell^2\mid N^+} \nu_\ell(k)
OUTPUT:
an integer
EXAMPLES::
sage: X = BruhatTitsQuotient(31,3) sage: X.e4 2 """
@lazy_attribute def mu(self): """ Compute the mu invariant of self.
OUTPUT:
An integer.
EXAMPLES::
sage: X = BruhatTitsQuotient(29,3) sage: X.mu 2 """
@cached_method def get_num_verts(self): """ Return the number of vertices in the quotient using the formula `V = 2(\mu/12 + e_3/3 + e_4/4)`.
OUTPUT:
- An integer (the number of vertices)
EXAMPLES::
sage: X = BruhatTitsQuotient(29,11) sage: X.get_num_verts() 4 """
@cached_method def get_num_ordered_edges(self): """ Return the number of ordered edges `E` in the quotient using the formula relating the genus `g` with the number of vertices `V` and that of unordered edges `E/2`: `E = 2(g + V - 1)`.
OUTPUT:
- An integer
EXAMPLES::
sage: X = BruhatTitsQuotient(3,2) sage: X.get_num_ordered_edges() 2 """
def genus_no_formula(self): """ Compute the genus of the quotient from the data of the quotient graph. This should agree with self.genus().
OUTPUT:
An integer
EXAMPLES::
sage: X = BruhatTitsQuotient(5,2*3*29) sage: X.genus_no_formula() 17 sage: X.genus_no_formula() == X.genus() True """
@cached_method def genus(self): r""" Compute the genus of the quotient graph using a formula This should agree with self.genus_no_formula().
Compute the genus of the Shimura curve corresponding to this quotient via Cerednik-Drinfeld. It is computed via a formula and not in terms of the quotient graph.
INPUT:
- level: Integer (default: None) a level. By default, use that of ``self``.
- Nplus: Integer (default: None) a conductor. By default, use that of ``self``.
OUTPUT:
An integer equal to the genus
EXAMPLES::
sage: X = BruhatTitsQuotient(3,2*5*31) sage: X.genus() 21 sage: X.genus() == X.genus_no_formula() True """
@cached_method def dimension_harmonic_cocycles(self, k, lev=None, Nplus=None, character=None): r""" Compute the dimension of the space of harmonic cocycles of weight `k` on ``self``.
OUTPUT:
An integer equal to the dimension
EXAMPLES::
sage: X = BruhatTitsQuotient(3,7) sage: [X.dimension_harmonic_cocycles(k) for k in range(2,20,2)] [1, 4, 4, 8, 8, 12, 12, 16, 16]
sage: X = BruhatTitsQuotient(2,5) # optional - magma sage: [X.dimension_harmonic_cocycles(k) for k in range(2,40,2)] # optional - magma [0, 1, 3, 1, 3, 5, 3, 5, 7, 5, 7, 9, 7, 9, 11, 9, 11, 13, 11]
sage: X = BruhatTitsQuotient(7, 2 * 3 * 5) sage: X.dimension_harmonic_cocycles(4) 12 sage: X = BruhatTitsQuotient(7, 2 * 3 * 5 * 11 * 13) sage: X.dimension_harmonic_cocycles(2) 481 sage: X.dimension_harmonic_cocycles(4) 1440 """ else: lev = ZZ(lev) else: Nplus = ZZ(Nplus)
range(lev * Nplus)) else:
return 0
return Gamma0(Nplus).dimension_cusp_forms(k=k)
raise NotImplementedError('The level should be squarefree for ' 'this function to work... Sorry!')
return ZZ(0)
def Nplus(self): r""" Return the tame level `N^+`.
OUTPUT:
An integer equal to `N^+`.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,7,1) sage: X.Nplus() 1 """
def Nminus(self): r""" Return the discriminant of the relevant definite quaternion algebra.
OUTPUT:
An integer equal to `N^-`.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,7) sage: X.Nminus() 7 """
@cached_method def level(self): r""" Return `p N^-`, which is the discriminant of the indefinite quaternion algebra that is uniformed by Cerednik-Drinfeld.
OUTPUT:
An integer equal to `p N^-`.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,7) sage: X.level() 35 """
def prime(self): r""" Return the prime one is working with.
OUTPUT:
An integer equal to the fixed prime `p`
EXAMPLES::
sage: X = BruhatTitsQuotient(5,7) sage: X.prime() 5 """
def get_graph(self): r""" Return the quotient graph (and compute it if needed).
OUTPUT:
A graph representing the quotient of the Bruhat-Tits tree.
EXAMPLES::
sage: X = BruhatTitsQuotient(11,5) sage: X.get_graph() Multi-graph on 2 vertices """
def get_fundom_graph(self): r""" Return the fundamental domain (and computes it if needed).
OUTPUT:
A fundamental domain for the action of `\Gamma`.
EXAMPLES::
sage: X = BruhatTitsQuotient(11,5) sage: X.get_fundom_graph() Graph on 24 vertices """ except AttributeError: self._compute_quotient() return self._Sfun
def plot(self, *args, **kwargs): r""" Plot the quotient graph.
OUTPUT:
A plot of the quotient graph
EXAMPLES::
sage: X = BruhatTitsQuotient(7,23) sage: X.plot() Graphics object consisting of 17 graphics primitives """ vertex_colors[key].append(v) else:
def plot_fundom(self, *args, **kwargs): r""" Plot a fundamental domain.
OUTPUT:
A plot of the fundamental domain.
EXAMPLES::
sage: X = BruhatTitsQuotient(7,23) sage: X.plot_fundom() Graphics object consisting of 88 graphics primitives """ else:
def is_admissible(self, D): r""" Test whether the imaginary quadratic field of discriminant `D` embeds in the quaternion algebra. It furthermore tests the Heegner hypothesis in this setting (e.g., is `p` inert in the field, etc).
INPUT:
- ``D`` - an integer whose squarefree part will define the quadratic field
OUTPUT:
A boolean describing whether the quadratic field is admissible
EXAMPLES::
sage: X = BruhatTitsQuotient(5,7) sage: [X.is_admissible(D) for D in range(-1,-20,-1)] [False, True, False, False, False, False, False, True, False, False, False, False, False, False, False, False, False, True, False] """ if kronecker_symbol(disc, f[0]) != 1: return False
def _local_splitting_map(self, prec): r""" Return an embedding of the definite quaternion algebra into the algebra of 2x2 matrices with coefficients in `\QQ_p`.
INPUT:
- ``prec`` -- Integer. The precision of the splitting.
OUTPUT:
A function giving the splitting.
EXAMPLES::
sage: X = BruhatTitsQuotient(11,3) sage: phi = X._local_splitting_map(10) sage: B.<i,j,k> = QuaternionAlgebra(3) sage: phi(i)**2 == QQ(i**2)*phi(B(1)) True """
def _local_splitting(self, prec): r""" Find an embedding of the definite quaternion algebra into the algebra of 2x2 matrices with coefficients in `\QQ_p`.
INPUT:
- ``prec`` - Integer. The precision of the splitting.
OUTPUT:
- Matrices `I`, `J`, `K` giving the splitting.
EXAMPLES::
sage: X = BruhatTitsQuotient(11,3) sage: phi = X._local_splitting_map(10) sage: B.<i,j,k> = QuaternionAlgebra(3) sage: phi(i)**2 == QQ(i**2)*phi(B(1)) True """ return self._II, self._JJ, self._KK
raise ValueError("must be rational quaternion algebra") raise ValueError("p (=%s) must be an unramified prime" % self._p)
else: else:
def _compute_embedding_matrix(self, prec, force_computation=False): r""" Return a matrix representing the embedding with the given precision.
INPUT:
- ``prec`` - Integer. The precision of the embedding matrix.
EXAMPLES:
Note that the entries of the matrix are elements of Zmod::
sage: X = BruhatTitsQuotient(3,7) sage: A = X.get_embedding_matrix(10) # indirect doctest sage: R = A.base_ring() sage: B = X.get_eichler_order_basis() sage: R(B[0].reduced_trace()) == A[0,0]+A[3,0] True """ if not force_computation: try: return Matrix(Zmod(self._pN), 4, 4, self._cached_Iota0_matrix) except AttributeError: pass
Ord = self.get_eichler_order(magma=True) # force_computation = force_computation) OrdMax = self.get_maximal_order(magma=True)
OBasis = Ord.Basis() verbose('Calling magma: pMatrixRing, args = %s' % [OrdMax, self._p]) M, f, rho = self._magma.function_call('pMatrixRing', args=[OrdMax, self._p], params={'Precision': 2000}, nvals=3) v = [f.Image(OBasis[i]) for i in [1, 2, 3, 4]]
self._cached_Iota0_matrix = [v[kk][ii, jj].sage() for ii in range(1, 3) for jj in range(1, 3) for kk in range(4)] return Matrix(Zmod(self._pN), 4, 4, self._cached_Iota0_matrix) else: [phi(B[kk])[ii, jj] for ii in range(2) for jj in range(2) for kk in range(4)])
@cached_method def get_extra_embedding_matrices(self): r""" Return a list of matrices representing the different embeddings.
.. NOTE::
The precision is very low (currently set to 5 digits), since these embeddings are only used to apply a character.
EXAMPLES:
This portion of the code is only relevant when working with a nontrivial Dirichlet character. If there is no such character then the code returns an empty list. Even if the character is not trivial it might return an empty list::
sage: f = DirichletGroup(6)[1] sage: X = BruhatTitsQuotient(3,2*5*7,character = f) sage: X.get_extra_embedding_matrices() []
::
sage: f = DirichletGroup(6)[1] sage: X = BruhatTitsQuotient(5,2,3, character = f, use_magma=True) # optional - magma sage: X.get_extra_embedding_matrices() # optional - magma [ [1 0 2 0] [0 0 2 0] [0 0 0 0] [1 0 2 2] ] """ n_iters = 0 Ord = self.get_eichler_order(magma=True) OrdMax = self.get_maximal_order(magma=True) OBasis = Ord.Basis() extra_embeddings = [] success = False while not success: success = True for l in self._extra_level: success = False found = False while not found: verbose('Calling magma: pMatrixRing, args = %s' % [OrdMax, l]) M, f, rho = self._magma.function_call('pMatrixRing', args=[OrdMax, l], params={'Precision': 20}, nvals=3) v = [f.Image(OBasis[i]) for i in [1, 2, 3, 4]] if all([Qp(l, 5)(v[kk][2, 1].sage()).valuation() >= 1 for kk in range(4)]) and not all([Qp(l, 5)(v[kk][2, 1].sage()).valuation() >= 2 for kk in range(4)]): found = True success = True else: n_iters += 1 verbose('Restarting magma...') self._magma.quit() self._magma = magma self._magma.function_call('SetSeed', n_iters, nvals=0) self._order_is_initialized = False self._init_order() self._compute_embedding_matrix(self._prec, force_computation=True) Ord = self.get_eichler_order(magma=True) OrdMax = self.get_maximal_order(magma=True) OBasis = Ord.Basis() extra_embeddings = [] success = False break if not success: break mat = Matrix(GF(l), 4, 4, [v[kk][ii, jj].sage() for ii in range(1, 3) for jj in range(1, 3) for kk in range(4)]) extra_embeddings.append(mat) return extra_embeddings
def _increase_precision(self, amount=1): r""" Increase the working precision.
INPUT:
- ``amount`` Integer (default: 1). The amount by which to increase the precision.
EXAMPLES:
sage: X = BruhatTitsQuotient(3,101) sage: X.get_embedding_matrix() [ O(3) 1 + O(3) 1 + O(3) 1 + O(3)] [2 + O(3) O(3) 2 + O(3) 2 + O(3)] [1 + O(3) 1 + O(3) O(3) 2 + O(3)] [1 + O(3) 2 + O(3) 2 + O(3) 2 + O(3)] sage: X._increase_precision(5) sage: X.get_embedding_matrix()[0,0] 2*3^3 + 2*3^5 + O(3^6) """
def get_embedding_matrix(self, prec=None, exact=False): r""" Return the matrix of the embedding.
INPUT:
- ``exact`` boolean (Default: ``False``). If ``True``, return an embedding into a matrix algebra with coefficients in a number field. Otherwise, embed into matrices over `p`-adic numbers.
- ``prec`` Integer (Default: ``None``). If specified, return the matrix with precision ``prec``. Otherwise, return the cached matrix (with the current working precision).
OUTPUT:
- A 4x4 matrix representing the embedding.
EXAMPLES::
sage: X = BruhatTitsQuotient(7,2*3*5) sage: X.get_embedding_matrix(4) [ 1 + O(7^4) 5 + 2*7 + 3*7^3 + O(7^4) 4 + 5*7 + 6*7^2 + 6*7^3 + O(7^4) 6 + 3*7^2 + 4*7^3 + O(7^4)] [ O(7^4) O(7^4) 3 + 7 + O(7^4) 1 + 6*7 + 3*7^2 + 2*7^3 + O(7^4)] [ O(7^4) 2 + 5*7 + 6*7^3 + O(7^4) 3 + 5*7 + 6*7^2 + 6*7^3 + O(7^4) 3 + 3*7 + 3*7^2 + O(7^4)] [ 1 + O(7^4) 3 + 4*7 + 6*7^2 + 3*7^3 + O(7^4) 3 + 7 + O(7^4) 1 + 6*7 + 3*7^2 + 2*7^3 + O(7^4)] sage: X.get_embedding_matrix(3) [ 1 + O(7^4) 5 + 2*7 + 3*7^3 + O(7^4) 4 + 5*7 + 6*7^2 + 6*7^3 + O(7^4) 6 + 3*7^2 + 4*7^3 + O(7^4)] [ O(7^4) O(7^4) 3 + 7 + O(7^4) 1 + 6*7 + 3*7^2 + 2*7^3 + O(7^4)] [ O(7^4) 2 + 5*7 + 6*7^3 + O(7^4) 3 + 5*7 + 6*7^2 + 6*7^3 + O(7^4) 3 + 3*7 + 3*7^2 + O(7^4)] [ 1 + O(7^4) 3 + 4*7 + 6*7^2 + 3*7^3 + O(7^4) 3 + 7 + O(7^4) 1 + 6*7 + 3*7^2 + 2*7^3 + O(7^4)] sage: X.get_embedding_matrix(5) [ 1 + O(7^5) 5 + 2*7 + 3*7^3 + 6*7^4 + O(7^5) 4 + 5*7 + 6*7^2 + 6*7^3 + 6*7^4 + O(7^5) 6 + 3*7^2 + 4*7^3 + 5*7^4 + O(7^5)] [ O(7^5) O(7^5) 3 + 7 + O(7^5) 1 + 6*7 + 3*7^2 + 2*7^3 + 7^4 + O(7^5)] [ O(7^5) 2 + 5*7 + 6*7^3 + 5*7^4 + O(7^5) 3 + 5*7 + 6*7^2 + 6*7^3 + 6*7^4 + O(7^5) 3 + 3*7 + 3*7^2 + 5*7^4 + O(7^5)] [ 1 + O(7^5) 3 + 4*7 + 6*7^2 + 3*7^3 + O(7^5) 3 + 7 + O(7^5) 1 + 6*7 + 3*7^2 + 2*7^3 + 7^4 + O(7^5)] """ try: return self._Iota_exact except AttributeError: raise RuntimeError('Exact splitting not available.') else:
except AttributeError: pass
for ii in range(4) for jj in range(4)])
def embed_quaternion(self, g, exact=False, prec=None): r""" Embed the quaternion element ``g`` into a matrix algebra.
INPUT:
- ``g`` a row vector of size `4` whose entries represent a quaternion in our basis.
- ``exact`` boolean (default: False) - If True, tries to embed ``g`` into a matrix algebra over a number field. If False, the target is the matrix algebra over `\QQ_p`.
OUTPUT:
A 2x2 matrix with coefficients in `\QQ_p` if ``exact`` is False, or a number field if ``exact`` is True.
EXAMPLES::
sage: X = BruhatTitsQuotient(7,2) sage: l = X.get_units_of_order() sage: len(l) 12 sage: l[3] # random [-1] [ 0] [ 1] [ 1] sage: u = X.embed_quaternion(l[3]); u # random [ O(7) 3 + O(7)] [2 + O(7) 6 + O(7)] sage: X._increase_precision(5) sage: v = X.embed_quaternion(l[3]); v # random [ 7 + 3*7^2 + 7^3 + 4*7^4 + O(7^6) 3 + 7 + 3*7^2 + 7^3 + 4*7^4 + O(7^6)] [ 2 + 7 + 3*7^2 + 7^3 + 4*7^4 + O(7^6) 6 + 5*7 + 3*7^2 + 5*7^3 + 2*7^4 + 6*7^5 + O(7^6)] sage: u == v True """ return Matrix(self.get_splitting_field(), 2, 2, (self.get_embedding_matrix(exact=True) * g).list()) else:
embed = embed_quaternion
def get_embedding(self, prec=None): r""" Return a function which embeds quaternions into a matrix algebra.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,3) sage: f = X.get_embedding(prec = 4) sage: b = Matrix(ZZ,4,1,[1,2,3,4]) sage: f(b) [2 + 3*5 + 2*5^2 + 4*5^3 + O(5^4) 3 + 2*5^2 + 4*5^3 + O(5^4)] [ 5 + 5^2 + 3*5^3 + O(5^4) 4 + 5 + 2*5^2 + O(5^4)] """
def get_edge_stabilizers(self): r""" Compute the stabilizers in the arithmetic group of all edges in the Bruhat-Tits tree within a fundamental domain for the quotient graph. The stabilizers of an edge and its opposite are equal, and so we only store half the data.
OUTPUT:
A list of lists encoding edge stabilizers. It contains one entry for each edge. Each entry is a list of data corresponding to the group elements in the stabilizer of the edge. The data consists of: (0) a column matrix representing a quaternion, (1) the power of `p` that one needs to divide by in order to obtain a quaternion of norm 1, and hence an element of the arithmetic group `\Gamma`, (2) a boolean that is only used to compute spaces of modular forms.
EXAMPLES::
sage: X=BruhatTitsQuotient(3,2) sage: s = X.get_edge_stabilizers() sage: len(s) == X.get_num_ordered_edges()/2 True sage: len(s[0]) 3 """ for e in self.get_edge_list()]
def get_stabilizers(self): r""" Compute the stabilizers in the arithmetic group of all edges in the Bruhat-Tits tree within a fundamental domain for the quotient graph. This is similar to get_edge_stabilizers, except that here we also store the stabilizers of the opposites.
OUTPUT:
A list of lists encoding edge stabilizers. It contains one entry for each edge. Each entry is a list of data corresponding to the group elements in the stabilizer of the edge. The data consists of: (0) a column matrix representing a quaternion, (1) the power of `p` that one needs to divide by in order to obtain a quaternion of norm 1, and hence an element of the arithmetic group `\Gamma`, (2) a boolean that is only used to compute spaces of modular forms.
EXAMPLES::
sage: X=BruhatTitsQuotient(3,5) sage: s = X.get_stabilizers() sage: len(s) == X.get_num_ordered_edges() True sage: gamma = X.embed_quaternion(s[1][0][0][0],prec = 20) sage: v = X.get_edge_list()[0].rep sage: X._BT.edge(gamma*v) == v True """
def get_vertex_stabs(self): r""" This function computes the stabilizers in the arithmetic group of all vertices in the Bruhat-Tits tree within a fundamental domain for the quotient graph.
OUTPUT:
A list of vertex stabilizers. Each vertex stabilizer is a finite cyclic subgroup, so we return generators for these subgroups.
EXAMPLES::
sage: X = BruhatTitsQuotient(13,2) sage: S = X.get_vertex_stabs() sage: gamma = X.embed_quaternion(S[0][0][0],prec = 20) sage: v = X.get_vertex_list()[0].rep sage: X._BT.vertex(gamma*v) == v True """ for v in self.get_vertex_list()]
def get_quaternion_algebra(self): r""" Return the underlying quaternion algebra.
OUTPUT:
The underlying definite quaternion algebra
EXAMPLES::
sage: X = BruhatTitsQuotient(5,7) sage: X.get_quaternion_algebra() Quaternion Algebra (-1, -7) with base ring Rational Field """
def get_eichler_order(self, magma=False, force_computation=False): r""" Return the underlying Eichler order of level `N^+`.
OUTPUT:
An Eichler order.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,7) sage: X.get_eichler_order() Order of Quaternion Algebra (-1, -7) with base ring Rational Field with basis (1/2 + 1/2*j, 1/2*i + 1/2*k, j, k) """ if not force_computation: try: return self._Omagma except AttributeError: pass self._init_order() return self._Omagma else:
def get_maximal_order(self, magma=False, force_computation=False): r""" Return the underlying maximal order containing the Eichler order.
OUTPUT:
A maximal order.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,7) sage: X.get_maximal_order() Order of Quaternion Algebra (-1, -7) with base ring Rational Field with basis (1/2 + 1/2*j, 1/2*i + 1/2*k, j, k) """ if not force_computation: try: return self._OMaxmagma except AttributeError: pass self._init_order() return self._OMaxmagma else:
def get_splitting_field(self): r""" Return a quadratic field that splits the quaternion algebra attached to ``self``. Currently requires Magma.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,11) sage: X.get_splitting_field() Traceback (most recent call last): ... NotImplementedError: Sage does not know yet how to work with the kind of orders that you are trying to use. Try installing Magma first and set it up so that Sage can use it.
If we do have Magma installed, then it works::
sage: X = BruhatTitsQuotient(5,11,use_magma=True) # optional - magma sage: X.get_splitting_field() # optional - magma Number Field in a with defining polynomial X1^2 + 11 """ try: return self._FF except AttributeError: pass self._compute_exact_splitting() return self._FF
def get_eichler_order_basis(self): r""" Return a basis for the global Eichler order.
OUTPUT:
Basis for the underlying Eichler order of level Nplus.
EXAMPLES::
sage: X = BruhatTitsQuotient(7,11) sage: X.get_eichler_order_basis() [1/2 + 1/2*j, 1/2*i + 1/2*k, j, k] """
def get_eichler_order_quadform(self): r""" This function return the norm form for the underlying Eichler order of level ``Nplus``. Required for finding elements in the arithmetic subgroup Gamma.
OUTPUT:
The norm form of the underlying Eichler order
EXAMPLES::
sage: X = BruhatTitsQuotient(7,11) sage: X.get_eichler_order_quadform() Quadratic form in 4 variables over Integer Ring with coefficients: [ 3 0 11 0 ] [ * 3 0 11 ] [ * * 11 0 ] [ * * * 11 ] """
def get_eichler_order_quadmatrix(self): r""" This function returns the matrix of the quadratic form of the underlying Eichler order in the fixed basis.
OUTPUT:
A 4x4 integral matrix describing the norm form.
EXAMPLES::
sage: X = BruhatTitsQuotient(7,11) sage: X.get_eichler_order_quadmatrix() [ 6 0 11 0] [ 0 6 0 11] [11 0 22 0] [ 0 11 0 22] """
@cached_method def get_units_of_order(self): r""" Return the units of the underlying Eichler `\ZZ`-order. This is a finite group since the order lives in a definite quaternion algebra over `\QQ`.
OUTPUT:
A list of elements of the global Eichler `\ZZ`-order of level `N^+`.
EXAMPLES::
sage: X = BruhatTitsQuotient(7,11) sage: X.get_units_of_order() [ [ 0] [-2] [-2] [ 0] [ 0] [ 1] [ 1], [ 0] ] """
@cached_method def _get_Up_data(self): r""" Return (compute if necessary) Up data.
The Up data is a vector of length `p`, and each entry consists of the corresponding data for the matrix `[p,a,0,1]` where a varies from 0 to `p-1`. The data is a tuple (acter,edge_images), with edge images being of type ``DoubleCosetReduction``.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,7) sage: [o[0] for o in X._get_Up_data()] [ [1/3 0] [-1/3 1/3] [-2/3 1/3] [ 0 1], [ 1 0], [ 1 0] ] """ [DoubleCosetReduction(self, e.rep * alpha) for e in E] + [DoubleCosetReduction(self, e.opposite.rep * alpha) for e in E]] for alpha in vec_a]
@cached_method def _get_atkin_lehner_data(self, q): r""" Return (and compute if necessary) data to compute the Atkin-Lehner involution.
INPUT:
- ``q`` - integer dividing p*Nminus*Nplus
EXAMPLES::
sage: X = BruhatTitsQuotient(3,5) sage: X._get_atkin_lehner_data(3)[0] [ 2] [ 4] [-3] [-2] """ # self._increase_precision(20)
#print 'Searching for norm', q*self._p**nninc
[DoubleCosetReduction(self, x.adjoint() * e.rep, extrapow=nn) for e in E]] except (PrecisionError, NotImplementedError): self._increase_precision(10)
@cached_method def _get_hecke_data(self, l): r""" Return (and compute if necessary) data to compute the Hecke operator at a prime.
INPUT:
- ``l`` - a prime l.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,17) sage: len(X._get_hecke_data(5)) 2 """ Sset = [] else:
else:
exact=False) * alphamat extrapow=nn) for e in E]
def _find_equivalent_vertex(self, v0, V=None, valuation=None): r""" Find a vertex in ``V`` equivalent to ``v0``.
INPUT:
- ``v0`` -- a 2x2 matrix in `\ZZ_p` representing a vertex in the Bruhat-Tits tree.
- ``V`` -- list (Default: None) If a list of Vertex is given, restrict the search to the vertices in ``V``. Otherwise use all the vertices in a fundamental domain.
- ``valuation`` -- an integer (Default: None): The valuation of the determinant of ``v0``, if known (otherwise it is calculated).
OUTPUT:
A pair ``g``, ``v``, where ``v`` is a Vertex in ``V`` equivalent to ``v0``, and ``g`` is such that `g\cdot v_0= v`.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,7) sage: M = Matrix(ZZ,2,2,[1,3,2,7]) sage: M.set_immutable() sage: X._find_equivalent_vertex(M)[-1] in X.get_vertex_list() True """
def _find_equivalent_edge(self, e0, E=None, valuation=None): r""" Find an edge in ``E`` equivalent to ``e0``.
INPUT:
- ``e0`` -- a 2x2 matrix in `\ZZ_p` representing an edge in the Bruhat-Tits tree.
- ``E`` -- list (Default: None) If a list of Edge is given, restrict the search to the vertices in ``E``. Otherwise use all the edges in a fundamental domain.
- ``valuation`` -- an integer (Default: None): The valuation of the determinant of ``e0``, if known (otherwise it is calculated).
OUTPUT:
A pair ``g``, ``e``, where ``e`` is an Edge in ``E`` equivalent to ``e0``, and ``g`` is such that `g\cdot e_0= e`.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,7) sage: M = Matrix(ZZ,2,2,[1,3,2,7]) sage: M.set_immutable() sage: X._find_equivalent_edge(M)[-1] in X.get_edge_list() True """ else: E = [e.opposite for e in self._edge_list]
def fundom_rep(self, v1): r""" Find an equivalent vertex in the fundamental domain.
INPUT:
- ``v1`` - a 2x2 matrix representing a normalized vertex.
OUTPUT:
A ``Vertex`` equivalent to ``v1``, in the fundamental domain.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,7) sage: M = Matrix(ZZ,2,2,[1,3,2,7]) sage: M.set_immutable() sage: X.fundom_rep(M) Vertex of Bruhat-Tits tree for p = 3 """ # print('v1=',v1) # print('chain =', chain) print('Given vertex:', v0) print('Not equivalent to any existing vertex in the list:') if V is not None: print([ve.label for ve in V]) assert 0 # what the hell is that ?
def _find_lattice(self, v1, v2, as_edges, m): r""" Find the lattice attached to the pair ``v1``,``v2``.
INPUT:
- ``v1``, ``v2`` - 2x2 matrices. They represent either a pair of normalized vertices or a pair of normalized edges.
- ``as_edges`` - boolean. If True, the inputs will be considered as edges instead of vertices.
- ``m`` - integer - The valuation of the determinant of ``v1``*``v2``.
OUTPUT:
A 4x4 integer matrix whose columns encode a lattice and a 4x4 integer matrix encoding a quadratic form.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,17) sage: X._find_lattice(Matrix(ZZ,2,2,[1,2,3,4]),Matrix(ZZ,2,2,[3,2,1,5]), True,0) ( [1 0 0 0] [138 204 -35 102] [2 3 0 0] [204 306 -51 153] [0 0 1 0] [-35 -51 12 -34] [0 0 0 1], [102 153 -34 102] ) """ else: for jj in range(2) for kk in range(2)])
def _stabilizer(self, e, as_edge=True): r""" Find the stabilizer of an edge or vertex.
INPUT:
- ``e`` - A 2x2 matrix representing an edge or vertex
- ``as_edge`` - Boolean (Default = True). Determines whether ``e`` is treated as an edge or vertex
OUTPUT:
A list of data describing the (finite) stabilizing subgroup of ``e``.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,7) sage: X._stabilizer(Matrix(ZZ,2,2,[3,8,2,9]))[0][2] False """ ## Using PARI to get the shortest vector in the lattice (via LLL) ## We used to pass qfminim flag = 2 else:
def _nebentype_check(self, vec, twom, E, A, flag = 2): """ Check if a quaternion maps into a subgroup of matrices determined by a nontrivial Dirichlet character (associated to self). If `N^+ = 1` then the condition is trivially satisfied.
INPUT:
- ``vec`` - 4x1 integer matrix. It encodes the quaternion to test in the basis defined by the columns of E.
- ``twom`` - An integer.
- ``E`` - 4x4 integer matrix. Its columns should form a basis for an order in the quaternion algebra.
- ``A`` - 4x4 integer matrix. It encodes the quadratic form on the order defined by the columns of E.
- ``flag`` - integer (Default = 0). Passed to Pari for finding minimal elements in a positive definite lattice.
OUTPUT:
A pair consisting of a quaternion (represented by a 4x1 column matrix) and a boolean saying whether the quaternion is in the subgroup of `M_2(\Qp)` determined by the Dirichlet character. Note that if `N^+` is trivial then this function aways outputs true.
EXAMPLES::
sage: f = DirichletGroup(6)[1] sage: X = BruhatTitsQuotient(3,2,1,f) sage: e = Matrix(ZZ,2,2,[1,2,5,7]) sage: m = e.determinant().valuation(3) sage: twom = 2*m sage: E,A = X._find_lattice(e,e,True,twom) sage: X._nebentype_check(E**(-1)*Matrix(ZZ,4,1,[1,0,0,0]),twom,E,A) ( [1] [0] [0] [0], True ) """ m = ZZ(twom / 2) mat = pari('qfminim(%s,,%s,flag = %s)' % (A.__pari__(), 1000, flag))[2].sage().transpose() n_vecs = mat.nrows() p = self._p pinv = Zmod(self._character.modulus())(p) ** -1 for jj in range(n_vecs): vect = mat.row(jj).row() vec = vect.transpose() nrd = Integer((vect * A * vec)[0, 0] / 2) if nrd == p ** twom: g = E * vec if prod([self._character(ZZ(pinv ** m * (v * g)[0, 0])) for v in self.get_extra_embedding_matrices()]) == 1: return g, True return None, False
def _are_equivalent(self, v1, v2, as_edges=False, twom=None, check_parity=False): r""" Determine whether two vertices (or edges) of the Bruhat-Tits tree are equivalent under the arithmetic group in question. The computation boils down to an application of the LLL short-vector algorithm to a particular lattice; for details see [FM2014]_.
INPUT:
- ``v1``, ``v2`` - two 2x2 integral matrices representing either vertices or edges
- ``as_edges`` - boolean (Default: False). Tells whether the matrices should be interpreted as edges (if true), or as vertices (if false)
- ``twom`` - integer (Default: None) If specified, indicates the valuation of the determinant of ``v1`` `\times` ``v2``.
OUTPUT:
If the objects are equivalent, returns an element of the arithmetic group Gamma that takes ``v1`` to ``v2``. Otherwise returns False.
EXAMPLES::
sage: X = BruhatTitsQuotient(7,5) sage: M1 = Matrix(ZZ,2,2,[88,3,1,1]) sage: M1.set_immutable() sage: X._are_equivalent(M1,M1) == False False sage: M2 = Matrix(ZZ,2,2,[1,2,8,1]); M2.set_immutable() sage: X._are_equivalent(M1,M2, as_edges=True) sage: X._are_equivalent(M1,M2) == False False """ if twom % 2 != 0: self._cached_equivalent[(v1, v2, as_edges)] = None return None ## Using PARI to get the shortest vector in the lattice (via LLL)
def _compute_exact_splitting(self): r""" Use Magma to calculate a splitting of the order into the Matrix algebra with coefficients in an appropriate number field.
TESTS::
sage: X = BruhatTitsQuotient(3,23,use_magma=True) # optional - magma sage: X._compute_exact_splitting() # optional - magma """ # A = self.get_quaternion_algebra() R = self.get_maximal_order(magma=True) f = R.MatrixRepresentation() self._FF = NumberField(f.Codomain().BaseRing().DefiningPolynomial().sage(), 'a') allmats = [] verbose('Calling magma, compute exact splitting') for kk in range(4): xseq = self._magma('%s(%s)' % (f.name(), R.gen(kk + 1).name())).ElementToSequence() allmats.append(Matrix(self._FF, 2, 2, [self._FF([QQ(xseq[ii + 1][jj + 1]) for jj in range(2)]) for ii in range(4)])) self._Iota_exact = Matrix(self._FF, 4, 4, [self._FF(allmats[kk][ii, jj]) for ii in range(2) for jj in range(2) for kk in range(4)])
def _init_order(self): r""" Initialize the order of the quaternion algebra. Here we possibly use Magma to split it.
EXAMPLES::
sage: X = BruhatTitsQuotient(3,23) sage: X._init_order() """ return verbose('Calling magma, init_order') A = self._magma.QuaternionAlgebra(self._Nminus) g = A.gens() # We store the order because we need to split it OMaxmagma = A.QuaternionOrder(1) Omagma = OMaxmagma.Order(self._Nplus) OBasis = Omagma.Basis() self._A = QuaternionAlgebra((g[0] ** 2).sage(), (g[1] ** 2).sage()) i, j, k = self._A.gens() v = [1] + self._A.gens() self._B = [self._A(sum([OBasis[tt + 1][rr + 1].sage() * v[rr] for rr in range(4)])) for tt in range(4)] self._O = self._A.quaternion_order(self._B) self._Omagma = Omagma self._OMaxmagma = OMaxmagma else: # Note that we can't work with non-maximal orders in sage
for jj in range(4)]).inverse()
def B_one(self): r""" Return the coordinates of `1` in the basis for the quaternion order.
EXAMPLES::
sage: X = BruhatTitsQuotient(7,11) sage: v,pow = X.B_one() sage: X._conv(v) == 1 True """
def _conv(self, v): r""" Return a quaternion having coordinates in the fixed basis for the order given by ``v``.
OUTPUT:
A quaternion.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,7) sage: A = X.get_quaternion_algebra() sage: i,j,k = A.gens() sage: B = X.get_eichler_order_basis() sage: X._conv([1,2,3,4]) == B[0]+2*B[1]+3*B[2]+4*B[3] True """
@cached_method def _find_elements_in_order(self, norm, trace=None, primitive=False): r""" Return elements in the order of the quaternion algebra of specified reduced norm. One may optionally choose to specify the reduced trace.
INPUT:
- ``norm`` - integer. The required reduced norm.
- ``trace`` - integer (Default: None). If specified, returns elements only reduced trace ``trace``.
- ``primitive`` boolean (Default: False). If True, return only elements that cannot be divided by `p`.
EXAMPLES::
sage: X = BruhatTitsQuotient(5,7) sage: X._find_elements_in_order(23) [[2, 9, -1, -5], [0, 8, 0, -5], [-2, 9, 1, -5], [6, 7, -3, -4], [2, 5, -1, -4], [0, 6, -1, -4], [0, 8, -1, -4], [2, 9, -1, -4], [-2, 5, 1, -4], [0, 6, 1, -4], [0, 8, 1, -4], [-2, 9, 1, -4], [-6, 7, 3, -4], [7, 6, -4, -3], [7, 6, -3, -3], [6, 7, -3, -3], [0, 8, 0, -3], [-7, 6, 3, -3], [-6, 7, 3, -3], [-7, 6, 4, -3], [0, 1, -1, -2], [0, 6, -1, -2], [0, 1, 1, -2], [0, 6, 1, -2], [9, 2, -5, -1], [6, 0, -4, -1], [8, 0, -4, -1], [5, 2, -4, -1], [9, 2, -4, -1], [1, 0, -2, -1], [6, 0, -2, -1], [0, -1, -1, -1], [-1, 0, -1, -1], [5, 2, -1, -1], [2, 5, -1, -1], [0, -1, 1, -1], [1, 0, 1, -1], [-5, 2, 1, -1], [-2, 5, 1, -1], [-6, 0, 2, -1], [-1, 0, 2, -1], [-8, 0, 4, -1], [-6, 0, 4, -1], [-9, 2, 4, -1], [-5, 2, 4, -1], [-9, 2, 5, -1], [8, 0, -5, 0], [8, 0, -3, 0]] sage: X._find_elements_in_order(23,1) [[1, 0, -2, -1], [1, 0, 1, -1]] """ verbose('Warning: norm (= %s) is quite large, this may take some time!' % norm)
def _compute_quotient(self, check=True): r""" Compute the quotient graph.
INPUT:
- ``check`` - Boolean (Default = True).
EXAMPLES::
sage: X = BruhatTitsQuotient(11,2) sage: X.get_graph() # indirect doctest Multi-graph on 2 vertices
sage: X = BruhatTitsQuotient(17,19) sage: X.get_graph() # indirect doctest Multi-graph on 4 vertices
The following examples require magma::
sage: X = BruhatTitsQuotient(5,7,12) # optional - magma sage: X.get_graph() # optional - magma Multi-graph on 24 vertices sage: len(X._edge_list) # optional - magma 72
sage: X = BruhatTitsQuotient(2,3,5) # optional - magma sage: X.get_graph() # optional - magma Multi-graph on 4 vertices
sage: X = BruhatTitsQuotient(2,3,35) # optional - magma sage: X.get_graph() # optional - magma Multi-graph on 16 vertices
sage: X = BruhatTitsQuotient(53,11,2) # optional - magma sage: X.get_graph() # optional - magma Multi-graph on 6 vertices
sage: X = BruhatTitsQuotient(2,13,9) # optional - magma sage: X.get_graph() # optional - magma Multi-graph on 24 vertices
AUTHORS:
- Cameron Franc (2012-02-20) - Marc Masdeu """ determinant=1, valuation=0) # total_verts = self.get_num_verts() # total_edges = genus + total_verts -1
valuation=edge_valuation)
else:
else: # The edge is new. # new_parity = new_valuation % 2 #The vertex is also new valuation=new_valuation) #Add the vertex to the list of pending vertices else:
# Add the edge to the list valuation=edge_valuation)
# Add the edge to the graph
# Find the opposite edge # opp_det = opp.determinant()
else:
print('You found a bug! Please report!') print('Computed genus =', computed_genus) print('Theoretical genus =', genus) raise RuntimeError raise RuntimeError('Number of vertices different ' 'from expected.')
def harmonic_cocycle_from_elliptic_curve(self, E, prec=None): r""" Return a harmonic cocycle with the same Hecke eigenvalues as ``E``.
Given an elliptic curve `E` having a conductor `N` of the form `pN^-N^+`, return the harmonic cocycle over ``self`` which is attached to ``E`` via modularity. The result is only well-defined up to scaling.
INPUT:
- ``E`` -- an elliptic curve over the rational numbers
- ``prec`` -- (default: None) If specified, the harmonic cocycle will take values in `\QQ_p` with precision ``prec``. Otherwise it will take values in `\ZZ`.
OUTPUT:
A harmonic cocycle attached via modularity to the given elliptic curve.
EXAMPLES::
sage: E = EllipticCurve('21a1') sage: X = BruhatTitsQuotient(7,3) sage: f = X.harmonic_cocycle_from_elliptic_curve(E,10) sage: T29 = f.parent().hecke_operator(29) sage: T29(f) == E.ap(29) * f True sage: E = EllipticCurve('51a1') sage: X = BruhatTitsQuotient(3,17) sage: f = X.harmonic_cocycle_from_elliptic_curve(E,20) sage: T31 = f.parent().hecke_operator(31) sage: T31(f) == E.ap(31) * f True """ except TypeError: try: N = E.conductor().norm() except ValueError: N = E.conductor().norm(QQ) continue else: Q = F(q).factor()[0][0] Eap = ZZ(Q.norm() + 1 - E.reduction(Q).count_points())
def harmonic_cocycles(self, k, prec=None, basis_matrix=None, base_field=None): r""" Compute the space of harmonic cocycles of a given even weight ``k``.
INPUT:
- ``k`` - integer - The weight. It must be even.
- ``prec`` - integer (default: None). If specified, the precision for the coefficient module
- ``basis_matrix`` - a matrix (default: None).
- ``base_field`` - a ring (default: None)
OUTPUT: A space of harmonic cocycles
EXAMPLES::
sage: X = BruhatTitsQuotient(31,7) sage: H = X.harmonic_cocycles(2,prec=10) sage: H Space of harmonic cocycles of weight 2 on Quotient of the Bruhat Tits tree of GL_2(QQ_31) with discriminant 7 and level 1 sage: H.basis()[0] Harmonic cocycle with values in Sym^0 Q_31^2 """
def padic_automorphic_forms(self, U, prec=None, t=None, R=None, overconvergent=False): r""" The module of (quaternionic) `p`-adic automorphic forms over ``self``.
INPUT:
- ``U`` -- A distributions module or an integer. If ``U`` is a distributions module then this creates the relevant space of automorphic forms. If ``U`` is an integer then the coefficients are the (`U-2`)nd power of the symmetric representation of `GL_2(\QQ_p)`.
- ``prec`` -- A precision (default : None). If not None should be a positive integer.
- ``t`` -- (default : None). The number of additional moments to store. If None, determine it automatically from ``prec``, ``U`` and the ``overconvergent`` flag.
- ``R`` -- (default : None). If specified, coefficient field of the automorphic forms. If not specified it defaults to the base ring of the distributions ``U``, or to `\QQ_p` with the working precision ``prec``.
- ``overconvergent`` -- Boolean (default = False). If True, will construct overconvergent `p`-adic automorphic forms. Otherwise it constructs the finite dimensional space of `p`-adic automorphic forms which is isomorphic to the space of harmonic cocycles.
EXAMPLES::
sage: X = BruhatTitsQuotient(11,5) sage: X.padic_automorphic_forms(2,prec=10) Space of automorphic forms on Quotient of the Bruhat Tits tree of GL_2(QQ_11) with discriminant 5 and level 1 with values in Sym^0 Q_11^2 """ |