Path: blob/master/src/sage/algebras/quatalg/quaternion_algebra.py
8822 views
"""1Quaternion Algebras23AUTHORS:45- Jon Bobber -- 2009 rewrite67- William Stein -- 2009 rewrite89This code is partly based on Sage code by David Kohel from 2005.1011TESTS:1213Pickling test::1415sage: Q.<i,j,k> = QuaternionAlgebra(QQ,-5,-2)16sage: Q == loads(dumps(Q))17True18"""1920########################################################################21# Copyright (C) 2009 William Stein <[email protected]>22# Copyright (C) 2009 Jonathon Bober <[email protected]>23#24# Distributed under the terms of the GNU General Public License (GPL)25#26# This code is distributed in the hope that it will be useful,27# but WITHOUT ANY WARRANTY; without even the implied warranty of28# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU29# General Public License for more details.30#31# The full text of the GPL is available at:32#33# http://www.gnu.org/licenses/34########################################################################353637from sage.rings.arith import (GCD,38hilbert_conductor_inverse, hilbert_conductor,39factor, gcd, lcm, kronecker_symbol, valuation)40from sage.rings.all import RR, Integer41from sage.rings.integer_ring import ZZ42from sage.rings.rational import Rational43from sage.rings.finite_rings.constructor import GF4445from sage.rings.ring import Algebra46from sage.rings.ideal import Ideal_fractional47from sage.rings.rational_field import is_RationalField, QQ48from sage.rings.infinity import infinity49from sage.rings.number_field.number_field import is_NumberField50from sage.structure.parent_gens import ParentWithGens, normalize_names51from sage.matrix.matrix_space import MatrixSpace52from sage.matrix.constructor import diagonal_matrix, matrix53from sage.structure.sequence import Sequence54from sage.structure.element import is_RingElement55from sage.modules.free_module import VectorSpace, FreeModule56from sage.modules.free_module_element import vector5758from operator import itemgetter5960import quaternion_algebra_element61import quaternion_algebra_cython6263from sage.modular.modsym.p1list import P1List6465from sage.misc.cachefunc import cached_method6667from sage.categories.fields import Fields68_Fields = Fields()6970########################################################71# Constructor72########################################################7374_cache = {}7576def QuaternionAlgebra(arg0, arg1=None, arg2=None, names='i,j,k'):77"""78There are three input formats:7980- ``QuaternionAlgebra(a, b)``: quaternion algebra generated by ``i``, ``j``81subject to `i^2 = a`, `j^2 = b`, `j * i = -i * j`.8283- ``QuaternionAlgebra(K, a, b)``: same as above but over a field ``K``.84Here, ``a`` and ``b`` are nonzero elements of a field (``K``) of85characteristic not 2, and we set `k = i * j`.8687- ``QuaternionAlgebra(D)``: a rational quaternion algebra with88discriminant ``D``, where `D > 1` is a squarefree integer.8990EXAMPLES:9192``QuaternionAlgebra(a, b)`` - return quaternion algebra over the93*smallest* field containing the nonzero elements ``a`` and ``b`` with94generators ``i``, ``j``, ``k`` with `i^2=a`, `j^2=b` and `j*i=-i*j`::9596sage: QuaternionAlgebra(-2,-3)97Quaternion Algebra (-2, -3) with base ring Rational Field98sage: QuaternionAlgebra(GF(5)(2), GF(5)(3))99Quaternion Algebra (2, 3) with base ring Finite Field of size 5100sage: QuaternionAlgebra(2, GF(5)(3))101Quaternion Algebra (2, 3) with base ring Finite Field of size 5102sage: QuaternionAlgebra(QQ[sqrt(2)](-1), -5)103Quaternion Algebra (-1, -5) with base ring Number Field in sqrt2 with defining polynomial x^2 - 2104sage: QuaternionAlgebra(sqrt(-1), sqrt(-3))105Quaternion Algebra (I, sqrt(-3)) with base ring Symbolic Ring106sage: QuaternionAlgebra(1r,1)107Quaternion Algebra (1, 1) with base ring Rational Field108109Python ints, longs and floats may be passed to the ``QuaternionAlgebra(a, b)`` constructor, as may110all pairs of nonzero elements of a ring not of characteristic 2. The following tests address111the issues raised in trac ticket 10601::112113sage: QuaternionAlgebra(1r,1)114Quaternion Algebra (1, 1) with base ring Rational Field115sage: QuaternionAlgebra(1,1.0r)116Quaternion Algebra (1.00000000000000, 1.00000000000000) with base ring Real Field with 53 bits of precision117sage: QuaternionAlgebra(0,0)118Traceback (most recent call last):119...120ValueError: a and b must be nonzero121sage: QuaternionAlgebra(GF(2)(1),1)122Traceback (most recent call last):123...124ValueError: a and b must be elements of a ring with characteristic not 2125sage: a = PermutationGroupElement([1,2,3])126sage: QuaternionAlgebra(a, a)127Traceback (most recent call last):128...129ValueError: a and b must be elements of a ring with characteristic not 2130131``QuaternionAlgebra(K, a, b)`` - return quaternion algebra over the132field ``K`` with generators ``i``, ``j``, ``k`` with `i^2=a`, `j^2=b`133and `i*j=-j*i`::134135sage: QuaternionAlgebra(QQ, -7, -21)136Quaternion Algebra (-7, -21) with base ring Rational Field137sage: QuaternionAlgebra(QQ[sqrt(2)], -2,-3)138Quaternion Algebra (-2, -3) with base ring Number Field in sqrt2 with defining polynomial x^2 - 2139140``QuaternionAlgebra(D)`` - ``D`` is a squarefree integer; returns a141rational quaternion algebra of discriminant ``D``::142143sage: QuaternionAlgebra(1)144Quaternion Algebra (-1, 1) with base ring Rational Field145sage: QuaternionAlgebra(2)146Quaternion Algebra (-1, -1) with base ring Rational Field147sage: QuaternionAlgebra(7)148Quaternion Algebra (-1, -7) with base ring Rational Field149sage: QuaternionAlgebra(2*3*5*7)150Quaternion Algebra (-22, 210) with base ring Rational Field151152If the coefficients `a` and `b` in the definition of the quaternion153algebra are not integral, then a slower generic type is used for154arithmetic::155156sage: type(QuaternionAlgebra(-1,-3).0)157<type 'sage.algebras.quatalg.quaternion_algebra_element.QuaternionAlgebraElement_rational_field'>158sage: type(QuaternionAlgebra(-1,-3/2).0)159<type 'sage.algebras.quatalg.quaternion_algebra_element.QuaternionAlgebraElement_generic'>160161Make sure caching is sane::162163sage: A = QuaternionAlgebra(2,3); A164Quaternion Algebra (2, 3) with base ring Rational Field165sage: B = QuaternionAlgebra(GF(5)(2),GF(5)(3)); B166Quaternion Algebra (2, 3) with base ring Finite Field of size 5167sage: A is QuaternionAlgebra(2,3)168True169sage: B is QuaternionAlgebra(GF(5)(2),GF(5)(3))170True171sage: Q = QuaternionAlgebra(2); Q172Quaternion Algebra (-1, -1) with base ring Rational Field173sage: Q is QuaternionAlgebra(QQ,-1,-1)174True175sage: Q is QuaternionAlgebra(-1,-1)176True177sage: Q.<ii,jj,kk> = QuaternionAlgebra(15); Q.variable_names()178('ii', 'jj', 'kk')179sage: QuaternionAlgebra(15).variable_names()180('i', 'j', 'k')181182TESTS:183184Verify that bug found when working on trac 12006 involving coercing invariants185into the base field is fixed::186187sage: Q = QuaternionAlgebra(-1,-1); Q188Quaternion Algebra (-1, -1) with base ring Rational Field189sage: parent(Q._a)190Rational Field191sage: parent(Q._b)192Rational Field193"""194195# QuaternionAlgebra(D)196if arg1 is None and arg2 is None:197K = QQ198D = Integer(arg0)199a, b = hilbert_conductor_inverse(D)200a = Rational(a); b = Rational(b)201202elif arg2 is None:203# If arg0 or arg1 are Python data types, coerce them204# to the relevant Sage types. This is a bit inelegant.205L = []206for a in [arg0,arg1]:207if is_RingElement(a):208L.append(a)209elif isinstance(a, int) or isinstance(a, long):210L.append(Integer(a))211elif isinstance(a, float):212L.append(RR(a))213else:214raise ValueError("a and b must be elements of a ring with characteristic not 2")215216# QuaternionAlgebra(a, b)217v = Sequence(L)218K = v.universe().fraction_field()219a = K(v[0])220b = K(v[1])221222# QuaternionAlgebra(K, a, b)223else:224K = arg0225if K not in _Fields:226raise TypeError("base ring of quaternion algebra must be a field")227a = K(arg1)228b = K(arg2)229230if K.characteristic() == 2:231# Lameness!232raise ValueError("a and b must be elements of a ring with characteristic not 2")233if a == 0 or b == 0:234raise ValueError("a and b must be nonzero")235236global _cache237names = normalize_names(3, names)238key = (K, a, b, names)239if key in _cache:240return _cache[key]241A = QuaternionAlgebra_ab(K, a, b, names=names)242A._key = key243_cache[key] = A244return A245246247248249########################################################250# Classes251########################################################252253def is_QuaternionAlgebra(A):254"""255Return ``True`` if ``A`` is of the QuaternionAlgebra data type.256257EXAMPLES::258259sage: sage.algebras.quatalg.quaternion_algebra.is_QuaternionAlgebra(QuaternionAlgebra(QQ,-1,-1))260True261sage: sage.algebras.quatalg.quaternion_algebra.is_QuaternionAlgebra(ZZ)262False263"""264return isinstance(A, QuaternionAlgebra_abstract)265266class QuaternionAlgebra_abstract(Algebra):267def _repr_(self):268"""269EXAMPLES::270271sage: sage.algebras.quatalg.quaternion_algebra.QuaternionAlgebra_abstract(QQ)._repr_()272'Quaternion Algebra with base ring Rational Field'273"""274return "Quaternion Algebra with base ring %s"%self.base_ring()275276def ngens(self):277"""278Return the number of generators of the quaternion algebra as a K-vector279space, not including 1. This value is always 3: the algebra is spanned280by the standard basis `1`, `i`, `j`, `k`.281282EXAMPLES::283284sage: Q.<i,j,k> = QuaternionAlgebra(QQ,-5,-2)285sage: Q.ngens()2863287sage: Q.gens()288[i, j, k]289"""290return 3291292def basis(self):293"""294Return the fixed basis of ``self``, which is `1`, `i`, `j`, `k`, where295`i`, `j`, `k` are the generators of ``self``.296297EXAMPLES::298299sage: Q.<i,j,k> = QuaternionAlgebra(QQ,-5,-2)300sage: Q.basis()301(1, i, j, k)302303sage: Q.<xyz,abc,theta> = QuaternionAlgebra(GF(9,'a'),-5,-2)304sage: Q.basis()305(1, xyz, abc, theta)306307The basis is cached::308309sage: Q.basis() is Q.basis()310True311"""312try:313return self.__basis314except AttributeError:315self.__basis = tuple([self(1)] + list(self.gens()))316return self.__basis317318def inner_product_matrix(self):319"""320Return the inner product matrix associated to ``self``, i.e. the321Gram matrix of the reduced norm as a quadratic form on ``self``.322The standard basis `1`, `i`, `j`, `k` is orthogonal, so this matrix323is just the diagonal matrix with diagonal entries `2`, `2a`, `2b`,324`2ab`.325326EXAMPLES::327328sage: Q.<i,j,k> = QuaternionAlgebra(-5,-19)329sage: Q.inner_product_matrix()330[ 2 0 0 0]331[ 0 10 0 0]332[ 0 0 38 0]333[ 0 0 0 190]334"""335try: return self.__inner_product_matrix336except AttributeError: pass337338a, b = self._a, self._b339M = diagonal_matrix(self.base_ring(), [2, -2*a, -2*b, 2*a*b])340M.set_immutable()341self.__inner_product_matrix = M342return M343344def is_commutative(self):345"""346Return ``False`` always, since all quaternion algebras are347noncommutative.348349EXAMPLES::350351sage: Q.<i,j,k> = QuaternionAlgebra(QQ, -3,-7)352sage: Q.is_commutative()353False354"""355return False356357def is_division_algebra(self):358"""359Return ``True`` if the quaternion algebra is a division algebra (i.e.360every nonzero element in ``self`` is invertible), and ``False`` if the361quaternion algebra is isomorphic to the 2x2 matrix algebra.362363EXAMPLES::364365sage: QuaternionAlgebra(QQ,-5,-2).is_division_algebra()366True367sage: QuaternionAlgebra(1).is_division_algebra()368False369sage: QuaternionAlgebra(2,9).is_division_algebra()370False371sage: QuaternionAlgebra(RR(2.),1).is_division_algebra()372Traceback (most recent call last):373...374NotImplementedError: base field must be rational numbers375"""376if not is_RationalField(self.base_ring()):377raise NotImplementedError("base field must be rational numbers")378return self.discriminant() != 1379380def is_matrix_ring(self):381"""382Return ``True`` if the quaternion algebra is isomorphic to the 2x2383matrix ring, and ``False`` if ``self`` is a division algebra (i.e.384every nonzero element in ``self`` is invertible).385386EXAMPLES::387388sage: QuaternionAlgebra(QQ,-5,-2).is_matrix_ring()389False390sage: QuaternionAlgebra(1).is_matrix_ring()391True392sage: QuaternionAlgebra(2,9).is_matrix_ring()393True394sage: QuaternionAlgebra(RR(2.),1).is_matrix_ring()395Traceback (most recent call last):396...397NotImplementedError: base field must be rational numbers398399"""400if not is_RationalField(self.base_ring()):401raise NotImplementedError("base field must be rational numbers")402return self.discriminant() == 1403404def is_exact(self):405"""406Return ``True`` if elements of this quaternion algebra are represented407exactly, i.e. there is no precision loss when doing arithmetic. A408quaternion algebra is exact if and only if its base field is409exact.410411EXAMPLES::412413sage: Q.<i,j,k> = QuaternionAlgebra(QQ, -3, -7)414sage: Q.is_exact()415True416sage: Q.<i,j,k> = QuaternionAlgebra(Qp(7), -3, -7)417sage: Q.is_exact()418False419"""420return self.base_ring().is_exact()421422def is_field(self, proof = True):423"""424Return ``False`` always, since all quaternion algebras are425noncommutative and all fields are commutative.426427EXAMPLES::428429sage: Q.<i,j,k> = QuaternionAlgebra(QQ, -3, -7)430sage: Q.is_field()431False432"""433return False434435def is_finite(self):436"""437Return ``True`` if the quaternion algebra is finite as a set.438439Algorithm: A quaternion algebra is finite if and only if the440base field is finite.441442EXAMPLES::443444sage: Q.<i,j,k> = QuaternionAlgebra(QQ, -3, -7)445sage: Q.is_finite()446False447sage: Q.<i,j,k> = QuaternionAlgebra(GF(5), -3, -7)448sage: Q.is_finite()449True450"""451return self.base_ring().is_finite()452453def is_integral_domain(self, proof = True):454"""455Return ``False`` always, since all quaternion algebras are456noncommutative and integral domains are commutative (in Sage).457458EXAMPLES::459460sage: Q.<i,j,k> = QuaternionAlgebra(QQ, -3, -7)461sage: Q.is_integral_domain()462False463"""464return False465466def is_noetherian(self):467"""468Return ``True`` always, since any quaternion algebra is a noetherian469ring (because it is a finitely generated module over a field).470471EXAMPLES::472473sage: Q.<i,j,k> = QuaternionAlgebra(QQ, -3, -7)474sage: Q.is_noetherian()475True476"""477return True478479def order(self):480"""481Return the number of elements of the quaternion algebra, or482``+Infinity`` if the algebra is not finite.483484EXAMPLES::485486sage: Q.<i,j,k> = QuaternionAlgebra(QQ, -3, -7)487sage: Q.order()488+Infinity489sage: Q.<i,j,k> = QuaternionAlgebra(GF(5), -3, -7)490sage: Q.order()491625492"""493return (self.base_ring().order())**4494495def random_element(self, *args, **kwds):496"""497Return a random element of this quaternion algebra.498499The ``args`` and ``kwds`` are passed to the ``random_element`` method500of the base ring.501502EXAMPLES::503504sage: QuaternionAlgebra(QQ[sqrt(2)],-3,7).random_element()505(sqrt2 + 2)*i + (-12*sqrt2 - 2)*j + (-sqrt2 + 1)*k506sage: QuaternionAlgebra(-3,19).random_element()507-1 + 2*i - j - 6/5*k508sage: QuaternionAlgebra(GF(17)(2),3).random_element()50914 + 10*i + 4*j + 7*k510511Specify the numerator and denominator bounds::512513sage: QuaternionAlgebra(-3,19).random_element(10^6,10^6)514-979933/553629 + 255525/657688*i - 3511/6929*j - 700105/258683*k515"""516K = self.base_ring()517return self([ K.random_element(*args, **kwds) for _ in range(4) ])518519def vector_space(self):520"""521Return the vector space associated to ``self`` with inner product given522by the reduced norm.523524EXAMPLES::525526sage: QuaternionAlgebra(-3,19).vector_space()527Ambient quadratic space of dimension 4 over Rational Field528Inner product matrix:529[ 2 0 0 0]530[ 0 6 0 0]531[ 0 0 -38 0]532[ 0 0 0 -114]533"""534try:535return self.__vector_space536except AttributeError:537V = VectorSpace(self.base_ring(), 4, inner_product_matrix = self.inner_product_matrix())538self.__vector_space = V539return V540541542class QuaternionAlgebra_ab(QuaternionAlgebra_abstract):543"""544The quaternion algebra of the form `(a, b/K)`, where `i^2=a`, `j^2 = b`,545and `j*i = -i*j`. ``K`` is a field not of characteristic 2 and ``a``,546``b`` are nonzero elements of ``K``.547548See ``QuaternionAlgebra`` for many more examples.549550EXAMPLES::551552sage: QuaternionAlgebra(QQ, -7, -21) # indirect doctest553Quaternion Algebra (-7, -21) with base ring Rational Field554"""555def __init__(self, base_ring, a, b, names='i,j,k'):556"""557Create the quaternion algebra with `i^2 = a`, `j^2 = b`, and558`i*j = -j*i = k`.559560INPUT:561562- ``base_ring`` - commutative ring563- ``a, b`` - elements of ``base_ring``564- ``names`` - string (optional, default 'i, j, k') names of the generators565566TESTS:567568Test making quaternion elements (using the element constructor)::569570sage: Q.<i,j,k> = QuaternionAlgebra(QQ,-1,-2)571sage: a = Q(2/3); a5722/3573sage: type(a)574<type 'sage.algebras.quatalg.quaternion_algebra_element.QuaternionAlgebraElement_rational_field'>575sage: Q(a)5762/3577sage: Q([1,2,3,4])5781 + 2*i + 3*j + 4*k579sage: Q((1,2,3,4))5801 + 2*i + 3*j + 4*k581sage: Q(-3/5)582-3/5583584The base ring must be a field::585586sage: Q.<ii,jj,kk> = QuaternionAlgebra(ZZ,-5,-19)587Traceback (most recent call last):588...589TypeError: base ring of quaternion algebra must be a field590"""591ParentWithGens.__init__(self, base_ring, names=names)592self._a = a593self._b = b594if base_ring not in _Fields:595raise TypeError("base ring of quaternion algebra must be a field")596if is_RationalField(base_ring) and a.denominator() == 1 and b.denominator() == 1:597element_constructor = quaternion_algebra_element.QuaternionAlgebraElement_rational_field598elif is_NumberField(base_ring) and base_ring.degree() > 2 and base_ring.is_absolute() and \599a.denominator() == 1 and b.denominator() == 1 and base_ring.defining_polynomial().is_monic():600# This QuaternionAlgebraElement_number_field class is not601# designed to work with elements of a quadratic field. To602# do that, the main thing would be to implement603# __getitem__, etc. This would maybe give a factor of 2604# (or more?) speedup. Much care must be taken because the605# underlying representation of quadratic fields is a bit606# tricky.607element_constructor = quaternion_algebra_element.QuaternionAlgebraElement_number_field608else:609element_constructor = quaternion_algebra_element.QuaternionAlgebraElement_generic610self._populate_coercion_lists_(coerce_list=[base_ring], element_constructor=element_constructor)611self._gens = [self([0,1,0,0]), self([0,0,1,0]), self([0,0,0,1])]612613def maximal_order(self, take_shortcuts = True):614"""615Return a maximal order in this quaternion algebra.616617The algorithm used is from [Voi2012].618619INPUT:620621- ``take_shortcuts`` - (default: True) If the discriminant is prime and the invariants of the algebra are of a nice form, use Proposition 5.2 of [Pizer, 1980].622623OUTPUT: a maximal order in this quaternion algebra624625EXAMPLES:626627::628629sage: QuaternionAlgebra(-1,-7).maximal_order()630Order 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)631632sage: QuaternionAlgebra(-1,-1).maximal_order().basis()633(1/2 + 1/2*i + 1/2*j + 1/2*k, i, j, k)634635sage: QuaternionAlgebra(-1,-11).maximal_order().basis()636(1/2 + 1/2*j, 1/2*i + 1/2*k, j, k)637638sage: QuaternionAlgebra(-1,-3).maximal_order().basis()639(1/2 + 1/2*j, 1/2*i + 1/2*k, j, k)640641sage: QuaternionAlgebra(-3,-1).maximal_order().basis()642(1/2 + 1/2*i, 1/2*j - 1/2*k, i, -k)643644sage: QuaternionAlgebra(-2,-5).maximal_order().basis()645(1/2 + 1/2*j + 1/2*k, 1/4*i + 1/2*j + 1/4*k, j, k)646647sage: QuaternionAlgebra(-5,-2).maximal_order().basis()648(1/2 + 1/2*i - 1/2*k, 1/2*i + 1/4*j - 1/4*k, i, -k)649650sage: QuaternionAlgebra(-17,-3).maximal_order().basis()651(1/2 + 1/2*j, 1/2*i + 1/2*k, -1/3*j - 1/3*k, k)652653sage: QuaternionAlgebra(-3,-17).maximal_order().basis()654(1/2 + 1/2*i, 1/2*j - 1/2*k, -1/3*i + 1/3*k, -k)655656sage: QuaternionAlgebra(-17*9,-3).maximal_order().basis()657(1, 1/3*i, 1/6*i + 1/2*j, 1/2 + 1/3*j + 1/18*k)658659sage: QuaternionAlgebra(-2, -389).maximal_order().basis()660(1/2 + 1/2*j + 1/2*k, 1/4*i + 1/2*j + 1/4*k, j, k)661662If you want bases containing 1, switch off ``take_shortcuts``:663664::665666sage: QuaternionAlgebra(-3,-89).maximal_order(take_shortcuts=False)667Order of Quaternion Algebra (-3, -89) with base ring Rational Field with basis (1, 1/2 + 1/2*i, j, 1/2 + 1/6*i + 1/2*j + 1/6*k)668669sage: QuaternionAlgebra(1,1).maximal_order(take_shortcuts=False) # Matrix ring670Order of Quaternion Algebra (1, 1) with base ring Rational Field with basis (1, 1/2 + 1/2*i, j, 1/2*j + 1/2*k)671672sage: QuaternionAlgebra(-22,210).maximal_order(take_shortcuts=False)673Order of Quaternion Algebra (-22, 210) with base ring Rational Field with basis (1, i, 1/2*i + 1/2*j, 1/2 + 17/22*i + 1/44*k)674675sage: for d in ( m for m in range(1, 750) if is_squarefree(m) ): # long time (3s)676... A = QuaternionAlgebra(d)677... R = A.maximal_order(take_shortcuts=False)678... assert A.discriminant() == R.discriminant()679680We don't support number fields other than the rationals yet:681682::683684sage: K = QuadraticField(5)685sage: QuaternionAlgebra(K,-1,-1).maximal_order()686Traceback (most recent call last):687...688NotImplementedError: maximal order only implemented for rational quaternion algebras689690691REFERENCES:692693.. [Piz1980] A. Pizer. An Algorithm for Computing Modular Forms on `\\Gamma_0(N)`,694J. Algebra 64 (1980), 340-390.695696.. [Voi2012] J. Voight. Identifying the matrix ring: algorithms697for quaternion algebras and quadratic forms, to appear.698"""699try: return self.__maximal_order700except AttributeError: pass701702if self.base_ring() != QQ:703raise NotImplementedError("maximal order only implemented for rational quaternion algebras")704705d_A = self.discriminant()706707# The following only works over QQ if the discriminant is prime708# and if the invariants are of the special form709# (every quaternion algebra of prime discriminant has a representation710# of such a form though)711a, b = self.invariants()712if take_shortcuts and d_A.is_prime() and a in ZZ and b in ZZ:713a = ZZ(a)714b = ZZ(b)715i,j,k = self.gens()716717# if necessary, try to swap invariants to match Pizer's paper718if (a != -1 and b == -1) or (b == -2) \719or (a != -1 and a != -2 and (-a) % 8 != 1):720a, b = b, a721i, j = j, i722k = i*j723724basis = []725if (a,b) == (-1,-1):726basis = [(1+i+j+k)/2, i, j, k]727elif a == -1 and (-b).is_prime() and ((-b) % 4 == 3):728basis = [(1+j)/2, (i+k)/2, j, k]729elif a == -2 and (-b).is_prime() and ((-b) % 8 == 5):730basis = [(1+j+k)/2, (i+2*j+k)/4, j, k]731elif (-a).is_prime() and (-b).is_prime():732q = -b733p = -a734735if q % 4 == 3 and kronecker_symbol(p,q) == -1:736a = 0737while (a*a*p + 1)%q != 0:738a += 1739basis = [(1+j)/2, (i+k)/2, -(j+a*k)/q, k]740741if basis:742self.__maximal_order = self.quaternion_order(basis)743return self.__maximal_order744745# The following code should always work (over QQ)746# Start with <1,i,j,k>747R = self.quaternion_order([1] + self.gens())748d_R = R.discriminant()749750e_new_gens = []751752# For each prime at which R is not yet maximal, make it bigger753for (p,p_val) in d_R.factor():754e = R.basis()755while self.quaternion_order(e).discriminant().valuation(p) > d_A.valuation(p):756# Compute a normalized basis at p757f = normalize_basis_at_p(list(e), p)758759# Ensure the basis lies in R by clearing denominators760# (this may make the order smaller at q != p)761# Also saturate the basis (divide out p as far as possible)762V = self.base_ring()**4763A = matrix(self.base_ring(), 4, 4, [ list(g) for g in e ]);764765e_n = []766x_rows = A.solve_left(matrix([ V(vec.coefficient_tuple()) for (vec,val) in f ]), check=False).rows()767denoms = [ x.denominator() for x in x_rows ]768for i in range(4):769vec = f[i][0]770val = f[i][1]771772v = (val/2).floor()773e_n.append(denoms[i] / p**(v) * vec)774775# for e_n to become p-saturated we still need to sort by776# ascending valuation of the quadratic form777lst = sorted(zip(e_n, [f[m][1].mod(2) for m in range(4)]),778key = itemgetter(1))779e_n = list(zip(*lst)[0])780781# Final step: Enlarge the basis at p782if p != 2:783# ensure that v_p(e_n[1]**2) = 0 by swapping basis elements784if ZZ(e_n[1]**2).valuation(p) != 0:785if ZZ(e_n[2]**2).valuation(p) == 0:786e_n[1], e_n[2] = e_n[2], e_n[1]787else:788e_n[1], e_n[3] = e_n[3], e_n[1]789790a = ZZ(e_n[1]**2)791b = ZZ(e_n[2]**2)792793if b.valuation(p) > 0: # if v_p(b) = 0, then already p-maximal794F = ZZ.quo(p)795if F(a).is_square():796x = F(a).sqrt().lift()797if (x**2 - a).mod(p**2) == 0: # make sure v_p(x**2 - a) = 1798x = x + p799g = 1/p*(x - e_n[1])*e_n[2]800e_n[2] = g801e_n[3] = e_n[1]*g802803else: # p == 2804t = e_n[1].reduced_trace()805a = -e_n[1].reduced_norm()806b = ZZ(e_n[2]**2)807808if t.valuation(p) == 0:809if b.valuation(p) > 0:810x = a811if (x**2 - t*x + a).mod(p**2) == 0: # make sure v_p(...) = 1812x = x + p813g = 1/p*(x - e_n[1])*e_n[2]814e_n[2] = g815e_n[3] = e_n[1]*g816817else: # t.valuation(p) > 0818(y,z,w) = maxord_solve_aux_eq(a, b, p)819g = 1/p*(1 + y*e_n[1] + z*e_n[2] + w*e_n[1]*e_n[2])820h = (z*b)*e_n[1] - (y*a)*e_n[2]821e_n[1:4] = [g,h,g*h]822if (1 - a*y**2 - b*z**2 + a*b*w**2).valuation(2) > 2:823e_n = basis_for_quaternion_lattice(list(e) + e_n[1:], reverse=True)824825# e_n now contains elements that locally at p give a bigger order,826# but the basis may be messed up at other primes (it might not even827# be an order). We will join them all together at the end828e = e_n829830e_new_gens.extend(e[1:])831832e_new = basis_for_quaternion_lattice(list(R.basis()) + e_new_gens, reverse=True)833self.__maximal_order = self.quaternion_order(e_new)834return self.__maximal_order835836837def invariants(self):838"""839Return the structural invariants `a`, `b` of this quaternion840algebra: ``self`` is generated by `i`, `j` subject to841`i^2 = a`, `j^2 = b` and `j*i = -i*j`.842843EXAMPLES::844845sage: Q.<i,j,k> = QuaternionAlgebra(15)846sage: Q.invariants()847(-3, 5)848sage: i^2849-3850sage: j^28515852"""853return self._a, self._b854855def __cmp__(self, other):856"""857Compare self and other.858859EXAMPLES::860861sage: cmp(QuaternionAlgebra(-1,-7), QuaternionAlgebra(-1,-7))8620863sage: cmp(QuaternionAlgebra(-1,-7), QuaternionAlgebra(-1,-5))864-1865sage: cmp(QuaternionAlgebra(-1,-7), QuaternionAlgebra(-1,-10))8661867"""868if not isinstance(other, QuaternionAlgebra_abstract):869return cmp(type(self), type(other))870c = cmp(self.base_ring(), other.base_ring())871if c: return c872return cmp((self._a, self._b), (other._a, other._b))873874def __reduce__(self):875"""876Internal method used for pickling.877878TESTS::879880sage: QuaternionAlgebra(QQ,-1,-2).__reduce__()881(<function unpickle_QuaternionAlgebra_v0 at ...>, (Rational Field, -1, -2, ('i', 'j', 'k')))882883Test uniqueness of parent::884885sage: Q = QuaternionAlgebra(QQ,-1,-2)886sage: loads(dumps(Q)) is Q887True888"""889return unpickle_QuaternionAlgebra_v0, self._key890891def gen(self, i=0):892"""893Return the `i^{th}` generator of ``self``.894895INPUT:896897- ``i`` - integer (optional, default 0)898899EXAMPLES::900901sage: Q.<ii,jj,kk> = QuaternionAlgebra(QQ,-1,-2); Q902Quaternion Algebra (-1, -2) with base ring Rational Field903sage: Q.gen(0)904ii905sage: Q.gen(1)906jj907sage: Q.gen(2)908kk909sage: Q.gens()910[ii, jj, kk]911"""912return self._gens[i]913914def _repr_(self):915"""916Print representation.917918TESTS::919920sage: Q.<i,j,k> = QuaternionAlgebra(QQ,-5,-2)921sage: type(Q)922<class 'sage.algebras.quatalg.quaternion_algebra.QuaternionAlgebra_ab'>923sage: Q._repr_()924'Quaternion Algebra (-5, -2) with base ring Rational Field'925sage: Q926Quaternion Algebra (-5, -2) with base ring Rational Field927sage: print Q928Quaternion Algebra (-5, -2) with base ring Rational Field929sage: str(Q)930'Quaternion Algebra (-5, -2) with base ring Rational Field'931"""932return "Quaternion Algebra (%r, %r) with base ring %s"%(self._a, self._b, self.base_ring())933934def inner_product_matrix(self):935"""936Return the inner product matrix associated to ``self``, i.e. the937Gram matrix of the reduced norm as a quadratic form on ``self``.938The standard basis `1`, `i`, `j`, `k` is orthogonal, so this matrix939is just the diagonal matrix with diagonal entries `1`, `a`, `b`, `ab`.940941EXAMPLES::942943sage: Q.<i,j,k> = QuaternionAlgebra(-5,-19)944sage: Q.inner_product_matrix()945[ 2 0 0 0]946[ 0 10 0 0]947[ 0 0 38 0]948[ 0 0 0 190]949950sage: R.<a,b> = QQ[]; Q.<i,j,k> = QuaternionAlgebra(Frac(R),a,b)951sage: Q.inner_product_matrix()952[ 2 0 0 0]953[ 0 -2*a 0 0]954[ 0 0 -2*b 0]955[ 0 0 0 2*a*b]956"""957a, b = self._a, self._b958return diagonal_matrix(self.base_ring(), [2, -2*a, -2*b, 2*a*b])959960def discriminant(self):961"""962Given a quaternion algebra `A` defined over a number field,963return the discriminant of `A`, i.e. the964product of the ramified primes of `A`.965966EXAMPLES::967968sage: QuaternionAlgebra(210,-22).discriminant()969210970sage: QuaternionAlgebra(19).discriminant()97119972973sage: F.<a> = NumberField(x^2-x-1)974sage: B.<i,j,k> = QuaternionAlgebra(F, 2*a,F(-1))975sage: B.discriminant()976Fractional ideal (2)977978sage: QuaternionAlgebra(QQ[sqrt(2)],3,19).discriminant()979Fractional ideal (1)980"""981try: return self.__discriminant982except AttributeError: pass983if not is_RationalField(self.base_ring()):984try:985F = self.base_ring()986self.__discriminant = F.hilbert_conductor(self._a,self._b)987except NotImplementedError: raise "base field must be rational numbers or number field"988else:989self.__discriminant = hilbert_conductor(self._a, self._b)990return self.__discriminant991992def ramified_primes(self):993"""994Return the primes that ramify in this quaternion algebra. Currently995only implemented over the rational numbers.996997EXAMPLES::998999sage: QuaternionAlgebra(QQ, -1, -1).ramified_primes()1000[2]1001"""1002#TODO: more examples10031004return [f[0] for f in factor(self.discriminant())]10051006def _magma_init_(self, magma):1007"""1008Return Magma version of this quaternion algebra.10091010EXAMPLES::10111012sage: Q = QuaternionAlgebra(-1,-1); Q1013Quaternion Algebra (-1, -1) with base ring Rational Field1014sage: Q._magma_init_(magma) # optional - magma1015'QuaternionAlgebra(_sage_[...],-1/1,-1/1)'1016sage: A = magma(Q); A # optional - magma1017Quaternion Algebra with base ring Rational Field, defined by i^2 = -1, j^2 = -11018sage: A.RamifiedPlaces() # optional - magma1019[1020Ideal of Integer Ring generated by 21021]10221023A more complicated example involving a quaternion algebra over a number field::10241025sage: K.<a> = QQ[sqrt(2)]; Q = QuaternionAlgebra(K,-1,a); Q1026Quaternion Algebra (-1, sqrt2) with base ring Number Field in sqrt2 with defining polynomial x^2 - 21027sage: magma(Q) # optional - magma1028Quaternion Algebra with base ring Number Field with defining polynomial x^2 - 2 over the Rational Field, defined by i^2 = -1, j^2 = sqrt21029sage: Q._magma_init_(magma) # optional - magma1030'QuaternionAlgebra(_sage_[...],(_sage_[...]![-1, 0]),(_sage_[...]![0, 1]))'1031"""1032R = magma(self.base_ring())1033return 'QuaternionAlgebra(%s,%s,%s)'%(R.name(),1034self._a._magma_init_(magma),1035self._b._magma_init_(magma))10361037def quaternion_order(self, basis, check=True):1038"""1039Return the order of this quaternion order with given basis.10401041INPUT:10421043- ``basis`` - list of 4 elements of ``self``1044- ``check`` - bool (default: ``True``)10451046EXAMPLES::10471048sage: Q.<i,j,k> = QuaternionAlgebra(-11,-1)1049sage: Q.quaternion_order([1,i,j,k])1050Order of Quaternion Algebra (-11, -1) with base ring Rational Field with basis (1, i, j, k)10511052We test out ``check=False``::10531054sage: Q.quaternion_order([1,i,j,k], check=False)1055Order of Quaternion Algebra (-11, -1) with base ring Rational Field with basis [1, i, j, k]1056sage: Q.quaternion_order([i,j,k], check=False)1057Order of Quaternion Algebra (-11, -1) with base ring Rational Field with basis [i, j, k]1058"""1059return QuaternionOrder(self, basis, check=check)10601061def ideal(self, gens, left_order=None, right_order=None, check=True, **kwds):1062r"""1063Return the quaternion ideal with given gens over `\ZZ`.1064Neither a left or right order structure need be specified.10651066INPUT:10671068- ``gens`` -- a list of elements of this quaternion order10691070- ``check`` -- bool (default: ``True``); if ``False``, then ``gens`` must10714-tuple that forms a Hermite basis for an ideal10721073- ``left_order`` -- a quaternion order or ``None``10741075- ``right_order`` -- a quaternion order or ``None``10761077EXAMPLES::10781079sage: R = QuaternionAlgebra(-11,-1)1080sage: R.ideal([2*a for a in R.basis()])1081Fractional ideal (2, 2*i, 2*j, 2*k)1082"""1083if self.base_ring() == QQ:1084return QuaternionFractionalIdeal_rational(gens, left_order=left_order, right_order=right_order, check=check)1085else:1086raise NotImplementedError("ideal only implemented for quaternion algebras over QQ")10871088@cached_method1089def modp_splitting_data(self, p):1090r"""1091Return mod `p` splitting data for this quaternion algebra at1092the unramified prime `p`. This is `2\times 2`1093matrices `I`, `J`, `K` over the finite field `\GF{p}` such that if1094the quaternion algebra has generators `i, j, k`, then `I^2 =1095i^2`, `J^2 = j^2`, `IJ=K` and `IJ=-JI`.10961097.. note::10981099Currently only implemented when `p` is odd and the base1100ring is `\QQ`.11011102INPUT:11031104- `p` -- unramified odd prime11051106OUTPUT:11071108- 2-tuple of matrices over finite field11091110EXAMPLES::11111112sage: Q = QuaternionAlgebra(-15, -19)1113sage: Q.modp_splitting_data(7)1114(1115[0 6] [6 1] [6 6]1116[1 0], [1 1], [6 1]1117)1118sage: Q.modp_splitting_data(next_prime(10^5))1119(1120[ 0 99988] [97311 4] [99999 59623]1121[ 1 0], [13334 2692], [97311 4]1122)1123sage: I,J,K = Q.modp_splitting_data(23)1124sage: I1125[0 8]1126[1 0]1127sage: I^21128[8 0]1129[0 8]1130sage: J1131[19 2]1132[17 4]1133sage: J^21134[4 0]1135[0 4]1136sage: I*J == -J*I1137True1138sage: I*J == K1139True11401141The following is a good test because of the asserts in the code::11421143sage: v = [Q.modp_splitting_data(p) for p in primes(20,1000)]114411451146Proper error handling::11471148sage: Q.modp_splitting_data(5)1149Traceback (most recent call last):1150...1151NotImplementedError: algorithm for computing local splittings not implemented in general (currently require the first invariant to be coprime to p)11521153sage: Q.modp_splitting_data(2)1154Traceback (most recent call last):1155...1156NotImplementedError: p must be odd1157"""1158if self.base_ring() != QQ:1159raise NotImplementedError("must be rational quaternion algebra")1160p = ZZ(p)1161if not p.is_prime():1162raise ValueError("p (=%s) must be prime"%p)1163if p == 2:1164raise NotImplementedError("p must be odd")1165if self.discriminant() % p == 0:1166raise ValueError("p (=%s) must be an unramified prime"%p)11671168i, j, k = self.gens()1169F = GF(p)1170i2 = F(i*i)1171j2 = F(j*j)11721173M = MatrixSpace(F, 2)1174I = M([0,i2,1,0])1175if i2 == 0:1176raise NotImplementedError("algorithm for computing local splittings not implemented in general (currently require the first invariant to be coprime to p)")1177i2inv = 1/i21178a = None1179for b in list(F):1180if not b: continue1181c = j2 + i2inv * b*b1182if c.is_square():1183a = -c.sqrt()1184break11851186if a is None:1187# do a fallback search, maybe needed in char 3 sometimes.1188for J in M:1189K = I*J1190if J*J == j2 and K == -J*I:1191return I, J, K11921193J = M([a,b,(j2-a*a)/b, -a])1194K = I*J1195assert K == -J*I, "bug in that I,J don't skew commute"1196return I, J, K11971198def modp_splitting_map(self, p):1199r"""1200Return Python map from the (`p`-integral) quaternion algebra to1201the set of `2\times 2` matrices over `\GF{p}`.12021203INPUT:12041205- `p` -- prime number12061207EXAMPLES::12081209sage: Q.<i,j,k> = QuaternionAlgebra(-1, -7)1210sage: f = Q.modp_splitting_map(13)1211sage: a = 2+i-j+3*k; b = 7+2*i-4*j+k1212sage: f(a*b)1213[12 3]1214[10 5]1215sage: f(a)*f(b)1216[12 3]1217[10 5]1218"""1219I, J, K = self.modp_splitting_data(p)1220F = I.base_ring()1221def phi(q):1222v = [F(a) for a in q.coefficient_tuple()]1223return v[0] + I*v[1] + J*v[2] + K*v[3]1224return phi122512261227############################################################1228# Unpickling1229############################################################1230def unpickle_QuaternionAlgebra_v0(*key):1231"""1232The 0th version of pickling for quaternion algebras.12331234EXAMPLES::12351236sage: Q = QuaternionAlgebra(-5,-19)1237sage: f, t = Q.__reduce__()1238sage: sage.algebras.quatalg.quaternion_algebra.unpickle_QuaternionAlgebra_v0(*t)1239Quaternion Algebra (-5, -19) with base ring Rational Field1240sage: loads(dumps(Q)) == Q1241True1242sage: loads(dumps(Q)) is Q1243True1244"""1245return QuaternionAlgebra(*key)124612471248class QuaternionOrder(Algebra):1249"""1250An order in a quaternion algebra.12511252EXAMPLES::12531254sage: QuaternionAlgebra(-1,-7).maximal_order()1255Order 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)1256sage: type(QuaternionAlgebra(-1,-7).maximal_order())1257<class 'sage.algebras.quatalg.quaternion_algebra.QuaternionOrder'>1258"""1259def __init__(self, A, basis, check=True):1260"""1261INPUT:12621263- ``A`` - a quaternion algebra1264- ``basis`` - list of 4 integral quaternions in ``A``1265- ``check`` - whether to do type and other consistency checks12661267.. WARNING::12681269Currently most methods silently assume that the A.base_ring() is QQ12701271EXAMPLES::12721273sage: A.<i,j,k> = QuaternionAlgebra(-3,-5)1274sage: sage.algebras.quatalg.quaternion_algebra.QuaternionOrder(A, [1,i,j,k])1275Order of Quaternion Algebra (-3, -5) with base ring Rational Field with basis (1, i, j, k)1276sage: R = sage.algebras.quatalg.quaternion_algebra.QuaternionOrder(A, [1,2*i,2*j,2*k]); R1277Order of Quaternion Algebra (-3, -5) with base ring Rational Field with basis (1, 2*i, 2*j, 2*k)1278sage: type(R)1279<class 'sage.algebras.quatalg.quaternion_algebra.QuaternionOrder'>12801281Over QQ and number fields it is checked whether the given1282basis actually gives a an order (as a module over the maximal order):12831284sage: A.<i,j,k> = QuaternionAlgebra(-1,-1)1285sage: A.quaternion_order([1,i,j,i-j])1286Traceback (most recent call last):1287...1288ValueError: basis must have rank 41289sage: A.quaternion_order([2,i,j,k])1290Traceback (most recent call last):1291...1292ValueError: lattice must contain 11293sage: A.quaternion_order([1,i/2,j/2,k/2])1294Traceback (most recent call last):1295...1296ValueError: given lattice must be a ring12971298sage: K = QuadraticField(10)1299sage: A.<i,j,k> = QuaternionAlgebra(K,-1,-1)1300sage: A.quaternion_order([1,i,j,k])1301Order of Quaternion Algebra (-1, -1) with base ring Number Field in a with defining polynomial x^2 - 10 with basis (1, i, j, k)1302sage: A.quaternion_order([1,i/2,j,k])1303Traceback (most recent call last):1304...1305ValueError: given lattice must be a ring1306"""1307if check:1308# right data type1309if not isinstance(basis, (list, tuple)):1310raise TypeError("basis must be a list or tuple")1311# right length1312if len(basis) != 4:1313raise ValueError("basis must have length 4")1314# coerce to common parent1315basis = tuple([A(x) for x in basis])13161317# has rank 41318V = A.base_ring()**41319if V.span([ V(x.coefficient_tuple()) for x in basis]).dimension() != 4:1320raise ValueError, "basis must have rank 4"13211322# The additional checks will work over QQ and over number fields,1323# but we can't actually do much with an order defined over a number1324# field13251326if A.base_ring() == QQ: # fast code over QQ1327M = matrix(QQ, 4, 4, [ x.coefficient_tuple() for x in basis])1328v = M.solve_left(V([1,0,0,0]))13291330if v.denominator() != 1:1331raise ValueError, "lattice must contain 1"13321333# check if multiplicatively closed1334M1 = basis_for_quaternion_lattice(basis)1335M2 = basis_for_quaternion_lattice(list(basis) + [ x*y for x in basis for y in basis])1336if M1 != M2:1337raise ValueError, "given lattice must be a ring"13381339if A.base_ring() != QQ: # slow code over number fields (should eventually use PARI's nfhnf)1340O = None1341try:1342O = A.base_ring().maximal_order()1343except AttributeError:1344pass13451346if O:1347M = matrix(A.base_ring(), 4, 4, [ x.coefficient_tuple() for x in basis])1348v = M.solve_left(V([1,0,0,0]))13491350if any([ not a in O for a in v]):1351raise ValueError, "lattice must contain 1"13521353# check if multiplicatively closed1354Y = matrix(QQ, 16, 4, [ (x*y).coefficient_tuple() for x in basis for y in basis])1355X = M.solve_left(Y)1356if any([ not a in O for x in X for a in x ]):1357raise ValueError, "given lattice must be a ring"13581359self.__basis = basis1360self.__quaternion_algebra = A1361ParentWithGens.__init__(self, ZZ, names=None)13621363def gens(self):1364"""1365Return generators for self.13661367EXAMPLES::13681369sage: QuaternionAlgebra(-1,-7).maximal_order().gens()1370(1/2 + 1/2*j, 1/2*i + 1/2*k, j, k)1371"""1372return self.__basis13731374def ngens(self):1375"""1376Return the number of generators (which is 4).13771378EXAMPLES::13791380sage: QuaternionAlgebra(-1,-7).maximal_order().ngens()138141382"""1383return 413841385def gen(self, n):1386"""1387Return the n-th generator.13881389INPUT:13901391- ``n`` - an integer between 0 and 3, inclusive.13921393EXAMPLES::13941395sage: R = QuaternionAlgebra(-11,-1).maximal_order(); R1396Order of Quaternion Algebra (-11, -1) with base ring Rational Field with basis (1/2 + 1/2*i, 1/2*j - 1/2*k, i, -k)1397sage: R.gen(0)13981/2 + 1/2*i1399sage: R.gen(1)14001/2*j - 1/2*k1401sage: R.gen(2)1402i1403sage: R.gen(3)1404-k1405"""1406return self.__basis[n]14071408def __cmp__(self, R):1409"""1410Compare orders self and other. Two orders are equal if they1411have the same basis and are in the same quaternion algebra.14121413EXAMPLES::14141415sage: R = QuaternionAlgebra(-11,-1).maximal_order()1416sage: R == R # indirect doctest1417True1418sage: R == QuaternionAlgebra(-1,-1).maximal_order()1419False1420sage: R==51421False1422"""1423if not isinstance(R, QuaternionOrder):1424return cmp(type(self), type(R))1425c = cmp(self.__quaternion_algebra, R.__quaternion_algebra)1426if c: return c1427return cmp(self.__basis, R.__basis)142814291430def basis(self):1431"""1432Return fix choice of basis for this quaternion order.14331434EXAMPLES::14351436sage: QuaternionAlgebra(-11,-1).maximal_order().basis()1437(1/2 + 1/2*i, 1/2*j - 1/2*k, i, -k)1438"""1439return self.__basis14401441def quaternion_algebra(self):1442"""1443Return ambient quaternion algebra that contains this quaternion order.14441445EXAMPLES::14461447sage: QuaternionAlgebra(-11,-1).maximal_order().quaternion_algebra()1448Quaternion Algebra (-11, -1) with base ring Rational Field1449"""1450return self.__quaternion_algebra14511452def _repr_(self):1453"""1454Return string representation of this order.14551456EXAMPLES::14571458sage: QuaternionAlgebra(-11,-1).maximal_order()._repr_()1459'Order of Quaternion Algebra (-11, -1) with base ring Rational Field with basis (1/2 + 1/2*i, 1/2*j - 1/2*k, i, -k)'1460sage: QuaternionAlgebra(-11,-1).maximal_order()1461Order of Quaternion Algebra (-11, -1) with base ring Rational Field with basis (1/2 + 1/2*i, 1/2*j - 1/2*k, i, -k)1462"""1463return 'Order of %s with basis %s'%(self.quaternion_algebra(), self.basis())14641465def random_element(self, *args, **kwds):1466"""1467Return a random element of this order.14681469The args and kwds are passed to the random_element method of1470the integer ring, and we return an element of the form14711472.. math::14731474ae_1 + be_2 + ce_3 + de_414751476where `e_1`, ..., `e_4` are the basis of this order and `a`,1477`b`, `c`, `d` are random integers.14781479EXAMPLES::14801481sage: QuaternionAlgebra(-11,-1).maximal_order().random_element()1482-4 - 4*i + j - k1483sage: QuaternionAlgebra(-11,-1).maximal_order().random_element(-10,10)1484-9/2 - 7/2*i - 7/2*j - 3/2*k1485"""1486return sum( (ZZ.random_element(*args, **kwds) * b for b in self.basis()) )14871488def intersection(self, other):1489"""1490Return the intersection of this order with other.14911492INPUT:14931494- ``other`` - a quaternion order in the same ambient quaternion algebra14951496OUTPUT: a quaternion order14971498EXAMPLES::14991500sage: R = QuaternionAlgebra(-11,-1).maximal_order()1501sage: R.intersection(R)1502Order of Quaternion Algebra (-11, -1) with base ring Rational Field with basis (1/2 + 1/2*i, i, 1/2*j + 1/2*k, k)15031504We intersect various orders in the quaternion algebra ramified at 11::15051506sage: B = BrandtModule(11,3)1507sage: R = B.maximal_order(); S = B.order_of_level_N()1508sage: R.intersection(S)1509Order of Quaternion Algebra (-1, -11) with base ring Rational Field with basis (1/2 + 1/2*j, 1/2*i + 5/2*k, j, 3*k)1510sage: R.intersection(S) == S1511True1512sage: B = BrandtModule(11,5)1513sage: T = B.order_of_level_N()1514sage: S.intersection(T)1515Order of Quaternion Algebra (-1, -11) with base ring Rational Field with basis (1/2 + 1/2*j, 1/2*i + 23/2*k, j, 15*k)1516"""1517if not isinstance(other, QuaternionOrder):1518raise TypeError("other must be a QuaternionOrder")15191520A = self.quaternion_algebra()1521if other.quaternion_algebra() != A:1522raise ValueError("self and other must be in the same ambient quaternion algebra")15231524V = A.base_ring()**415251526B = V.span([V(list(g)) for g in self.basis()], ZZ)1527C = V.span([V(list(g)) for g in other.basis()], ZZ)15281529# todo -- A(list(e)) could be A(e)1530return QuaternionOrder(A, [A(list(e)) for e in B.intersection(C).basis()])15311532def free_module(self):1533"""1534Return the free `\\ZZ`-module that corresponds to this order1535inside the vector space corresponding to the ambient1536quaternion algebra.15371538OUTPUT: a free `\\ZZ`-module of rank 415391540EXAMPLES::15411542sage: R = QuaternionAlgebra(-11,-1).maximal_order()1543sage: R.basis()1544(1/2 + 1/2*i, 1/2*j - 1/2*k, i, -k)1545sage: R.free_module()1546Free module of degree 4 and rank 4 over Integer Ring1547Echelon basis matrix:1548[1/2 1/2 0 0]1549[ 0 1 0 0]1550[ 0 0 1/2 1/2]1551[ 0 0 0 1]1552"""1553try: return self.__free_module1554except AttributeError: pass1555V = self.quaternion_algebra().base_ring()**41556M = V.span([V(list(g)) for g in self.basis()], ZZ)1557self.__free_module = M1558return M15591560def discriminant(self):1561r"""1562Return the discriminant of this order, which we define as1563`\sqrt{ det ( Tr(e_i \bar{e_j} ) ) }`, where `\{e_i\}` is the1564basis of the order.15651566OUTPUT: rational number15671568EXAMPLES::15691570sage: QuaternionAlgebra(-11,-1).maximal_order().discriminant()1571111572sage: S = BrandtModule(11,5).order_of_level_N()1573sage: S.discriminant()1574551575sage: type(S.discriminant())1576<type 'sage.rings.rational.Rational'>1577"""1578L = []1579for d in self.basis():1580MM = []1581for e in self.basis():1582MM.append( (d * e.conjugate()).reduced_trace() )1583L.append(MM)15841585return (MatrixSpace(QQ, 4, 4)(L)).determinant().sqrt()15861587def left_ideal(self, gens, check=True):1588r"""1589Return the ideal with given gens over `\ZZ`.15901591INPUT:15921593- ``gens`` -- a list of elements of this quaternion order15941595- ``check`` -- bool (default: ``True``); if ``False``, then ``gens`` must15964-tuple that forms a Hermite basis for an ideal15971598EXAMPLES::15991600sage: R = QuaternionAlgebra(-11,-1).maximal_order()1601sage: R.left_ideal([2*a for a in R.basis()])1602Fractional ideal (1 + i, 2*i, j + k, 2*k)1603"""1604if self.base_ring() == ZZ:1605return QuaternionFractionalIdeal_rational(gens, left_order=self, check=check)1606else:1607raise NotImplementedError("ideal only implemented for quaternion algebras over QQ")16081609def right_ideal(self, gens, check=True):1610r"""1611Return the ideal with given gens over `\ZZ`.16121613INPUT:16141615- ``gens`` -- a list of elements of this quaternion order16161617- ``check`` -- bool (default: ``True``); if ``False``, then ``gens`` must16184-tuple that forms a Hermite basis for an ideal16191620EXAMPLES::16211622sage: R = QuaternionAlgebra(-11,-1).maximal_order()1623sage: R.right_ideal([2*a for a in R.basis()])1624Fractional ideal (1 + i, 2*i, j + k, 2*k)1625"""1626if self.base_ring() == ZZ:1627return QuaternionFractionalIdeal_rational(gens, right_order=self, check=check)1628else:1629raise NotImplementedError("ideal only implemented for quaternion algebras over QQ")16301631def unit_ideal(self):1632"""1633Return the unit ideal in this quaternion order.16341635EXAMPLES::16361637sage: R = QuaternionAlgebra(-11,-1).maximal_order()1638sage: I = R.unit_ideal(); I1639Fractional ideal (1/2 + 1/2*i, 1/2*j - 1/2*k, i, -k)1640"""1641if self.base_ring() == ZZ:1642return QuaternionFractionalIdeal_rational(self.basis(), left_order=self, right_order=self, check=False)1643else:1644raise NotImplementedError("ideal only implemented for quaternion algebras over QQ")16451646def quadratic_form(self):1647"""1648Return the normalized quadratic form associated to this quaternion order.16491650OUTPUT: quadratic form16511652EXAMPLES::16531654sage: R = BrandtModule(11,13).order_of_level_N()1655sage: Q = R.quadratic_form(); Q1656Quadratic form in 4 variables over Rational Field with coefficients:1657[ 14 253 55 286 ]1658[ * 1455 506 3289 ]1659[ * * 55 572 ]1660[ * * * 1859 ]1661sage: Q.theta_series(10)16621 + 2*q + 2*q^4 + 4*q^6 + 4*q^8 + 2*q^9 + O(q^10)1663"""1664return self.unit_ideal().quadratic_form()16651666def ternary_quadratic_form(self, include_basis=False):1667"""1668Return the ternary quadratic form associated to this order.16691670INPUT:16711672- ``include_basis`` -- bool (default: False), if True also1673return a basis for the dimension 3 subspace `G`16741675OUTPUT:16761677- QuadraticForm16781679- optional basis for dimension 3 subspace16801681This function computes the positive definition quadratic form1682obtained by letting G be the trace zero subspace of `\ZZ` +16832* ``self``, which has rank 3, and restricting the pairing16841685(x,y) = (x.conjugate()*y).reduced_trace()16861687to `G`.16881689APPLICATIONS: Ternary quadratic forms associated to an order1690in a rational quaternion algebra are useful in computing with1691Gross points, in decided whether quaternion orders have1692embeddings from orders in quadratic imaginary fields, and in1693computing elements of the Kohnen plus subspace of modular1694forms of weight 3/2.16951696EXAMPLES::16971698sage: R = BrandtModule(11,13).order_of_level_N()1699sage: Q = R.ternary_quadratic_form(); Q1700Quadratic form in 3 variables over Rational Field with coefficients:1701[ 5820 1012 13156 ]1702[ * 55 1144 ]1703[ * * 7436 ]1704sage: factor(Q.disc())17052^4 * 11^2 * 13^217061707The following theta series is a modular form of weight 3/2 and level 4*11*13::17081709sage: Q.theta_series(100)17101 + 2*q^23 + 2*q^55 + 2*q^56 + 2*q^75 + 4*q^92 + O(q^100)1711"""1712if self.base_ring() != ZZ:1713raise NotImplementedError("ternary quadratic form of order only implemented for quaternion algebras over QQ")17141715Q = self.quaternion_algebra()1716# 2*R + ZZ1717twoR = self.free_module().scale(2)1718A = twoR.ambient_module()1719Z = twoR.span( [Q(1).coefficient_tuple()], ZZ)1720S = twoR + Z1721# Now we intersect with the trace 0 submodule1722v = [b.reduced_trace() for b in Q.basis()]1723M = matrix(QQ,4,1,v)1724tr0 = M.kernel()1725G = tr0.intersection(S)1726B = [Q(a) for a in G.basis()]1727m = matrix(QQ,[[x.pair(y) for x in B] for y in B])1728from sage.quadratic_forms.quadratic_form import QuadraticForm1729Q = QuadraticForm(m)1730if include_basis:1731return Q, B1732else:1733return Q17341735class QuaternionFractionalIdeal(Ideal_fractional):1736pass17371738class QuaternionFractionalIdeal_rational(QuaternionFractionalIdeal):1739"""1740A fractional ideal in a rational quaternion algebra.1741"""1742def __init__(self, basis, left_order=None, right_order=None, check=True):1743"""1744INPUT:17451746- ``left_order`` -- a quaternion order or ``None``17471748- ``right_order`` -- a quaternion order or ``None``17491750- ``basis`` -- tuple of length 4 of elements in of ambient1751quaternion algebra whose `\\ZZ`-span is an ideal17521753- ``check`` -- bool (default: ``True``); if ``False``, do no type1754checking, and the input basis *must* be in Hermite form.17551756EXAMPLES::17571758sage: R = QuaternionAlgebra(-11,-1).maximal_order()1759sage: R.right_ideal(R.basis())1760Fractional ideal (1/2 + 1/2*i, i, 1/2*j + 1/2*k, k)1761sage: R.right_ideal(tuple(R.basis()), check=False)1762Fractional ideal (1/2 + 1/2*i, 1/2*j - 1/2*k, i, -k)1763"""1764if check:1765if left_order is not None and not isinstance(left_order, QuaternionOrder):1766raise TypeError("left_order must be a quaternion order or None")1767if right_order is not None and not isinstance(right_order, QuaternionOrder):1768raise TypeError("right_order must be a quaternion order or None")1769if not isinstance(basis, (list, tuple)):1770raise TypeError("basis must be a list or tuple")17711772self.__left_order = left_order1773self.__right_order = right_order17741775if check:1776try:1777Q = self.quaternion_order().quaternion_algebra()1778except RuntimeError:1779Q = basis[0].parent()1780basis = tuple([Q(v) for v in1781(QQ**4).span([Q(v).coefficient_tuple() for v in basis], ZZ).basis()])1782self.__basis = basis17831784def scale(self, alpha, left=False):1785r"""1786Scale the fractional ideal self by multiplying the basis by alpha.17871788INPUT:17891790- `\alpha` -- element of quaternion algebra17911792- ``left`` -- bool (default: False); if true multiply1793`\alpha` on the left, otherwise multiply `\alpha` on the1794right.17951796OUTPUT:17971798- a new fractional ideal179918001801EXAMPLES::18021803sage: B = BrandtModule(5,37); I = B.right_ideals()[0]; i,j,k = B.quaternion_algebra().gens(); I1804Fractional ideal (2 + 2*j + 106*k, i + 2*j + 105*k, 4*j + 64*k, 148*k)1805sage: I.scale(i)1806Fractional ideal [2*i + 212*j - 2*k, -2 + 210*j - 2*k, 128*j - 4*k, 296*j]1807sage: I.scale(i, left=True)1808Fractional ideal [2*i - 212*j + 2*k, -2 - 210*j + 2*k, -128*j + 4*k, -296*j]1809sage: I.scale(i, left=False)1810Fractional ideal [2*i + 212*j - 2*k, -2 + 210*j - 2*k, 128*j - 4*k, 296*j]1811sage: i * I.gens()[0]18122*i - 212*j + 2*k1813sage: I.gens()[0] * i18142*i + 212*j - 2*k1815"""18161817Q = self.quaternion_algebra()1818alpha = Q(alpha)1819if left:1820gens = [alpha*b for b in self.basis()]1821else:1822gens = [b*alpha for b in self.basis()]1823return Q.ideal(gens, left_order = self.__left_order,1824right_order = self.__right_order, check=False)18251826def quaternion_algebra(self):1827"""1828Return the ambient quaternion algebra that contains this fractional ideal.18291830OUTPUT: a quaternion algebra18311832EXAMPLES::18331834sage: I = BrandtModule(3,5).right_ideals()[1]; I1835Fractional ideal (2 + 6*j + 4*k, 2*i + 4*j + 34*k, 8*j + 32*k, 40*k)1836sage: I.quaternion_algebra()1837Quaternion Algebra (-1, -3) with base ring Rational Field1838"""1839try: return self.__quaternion_algebra1840except AttributeError: pass1841A = self.__basis[0].parent()1842self.__quaternion_algebra = A1843return A18441845def _compute_order(self, side='left'):1846r"""1847Used internally to compute either the left or right order1848associated to an ideal in a quaternion algebra. If1849action='right', compute the left order, and if action='left'1850compute the right order.18511852INPUT:18531854- ``side`` -- 'left' or 'right'18551856EXAMPLES::18571858sage: R.<i,j,k> = QuaternionAlgebra(-1,-11)1859sage: I = R.ideal([2 + 2*j + 140*k, 2*i + 4*j + 150*k, 8*j + 104*k, 152*k])1860sage: Ol = I._compute_order('left'); Ol1861Order of Quaternion Algebra (-1, -11) with base ring Rational Field with basis (1/2 + 1/2*j + 35*k, 1/4*i + 1/2*j + 75/4*k, j + 32*k, 38*k)1862sage: Or = I._compute_order('right'); Or1863Order of Quaternion Algebra (-1, -11) with base ring Rational Field with basis (1/2 + 1/2*j + 16*k, 1/2*i + 11/2*k, j + 13*k, 19*k)1864sage: Ol.discriminant()18652091866sage: Or.discriminant()18672091868sage: I.left_order() == Ol1869True1870sage: I.right_order() == Or1871True18721873ALGORITHM: Let `b_1, b_2, b_3, b_3` be a basis for this1874fractional ideal `I`, and assume we want to compute the left1875order of `I` in the quaternion algebra `Q`. Then1876multiplication by `b_i` on the right defines a map `B_i:Q \to1877Q`. We have1878`R = B_1^{-1}(I) \cap B_2^{-1}(I) \cap B_3^{-1}(I)\cap B_4^{-1}(I).`1879This is because1880`B_n^{-1}(I) = \{\alpha \in Q : \alpha b_n \in I \},`1881and1882`R = \{\alpha \in Q : \alpha b_n \in I, n=1,2,3,4\}.`1883"""1884if side == 'left': action = 'right'1885elif side == 'right': action = 'left'1886else: ValueError, "side must be 'left' or 'right'"1887Q = self.quaternion_algebra()1888if Q.base_ring() != QQ:1889raise NotImplementedError("computation of left and right orders only implemented over QQ")1890M = [(~b).matrix(action=action) for b in self.basis()]1891B = self.basis_matrix()1892invs = [B*m for m in M]1893# Now intersect the row spans of each matrix in invs1894ISB = [Q(v) for v in intersection_of_row_modules_over_ZZ(invs).row_module(ZZ).basis()]1895return Q.quaternion_order(ISB)18961897def left_order(self):1898"""1899Return the left order associated to this fractional ideal.19001901OUTPUT: an order in a quaternion algebra19021903EXAMPLES::19041905sage: B = BrandtModule(11)1906sage: R = B.maximal_order()1907sage: I = R.unit_ideal()1908sage: I.left_order()1909Order of Quaternion Algebra (-1, -11) with base ring Rational Field with basis (1/2 + 1/2*j, 1/2*i + 1/2*k, j, k)19101911We do a consistency check::19121913sage: B = BrandtModule(11,19); R = B.right_ideals()1914sage: print [r.left_order().discriminant() for r in R]1915[209, 209, 209, 209, 209, 209, 209, 209, 209, 209, 209, 209, 209, 209, 209, 209, 209, 209]1916"""1917if self.__left_order is None:1918self.__left_order = self._compute_order(side='left')1919return self.__left_order19201921def right_order(self):1922"""1923Return the right order associated to this fractional ideal.19241925OUTPUT: an order in a quaternion algebra19261927EXAMPLES::19281929sage: I = BrandtModule(389).right_ideals()[1]; I1930Fractional ideal (2 + 6*j + 2*k, i + 2*j + k, 8*j, 8*k)1931sage: I.right_order()1932Order of Quaternion Algebra (-2, -389) with base ring Rational Field with basis (1/2 + 1/2*j + 1/2*k, 1/4*i + 1/2*j + 1/4*k, j, k)1933sage: I.left_order()1934Order of Quaternion Algebra (-2, -389) with base ring Rational Field with basis (1/2 + 1/2*j + 3/2*k, 1/8*i + 1/4*j + 9/8*k, j + k, 2*k)19351936The following is a big consistency check. We take reps for1937all the right ideal classes of a certain order, take the1938corresponding left orders, then take ideals in the left orders1939and from those compute the right order again::19401941sage: B = BrandtModule(11,19); R = B.right_ideals()1942sage: O = [r.left_order() for r in R]1943sage: J = [O[i].left_ideal(R[i].basis()) for i in range(len(R))]1944sage: len(set(J))1945181946sage: len(set([I.right_order() for I in J]))194711948sage: J[0].right_order() == B.order_of_level_N()1949True1950"""1951if self.__right_order is None:1952self.__right_order = self._compute_order(side='right')1953return self.__right_order19541955def __repr__(self):1956"""1957Return string representation of this quaternion fractional ideal.19581959EXAMPLES::19601961sage: I = BrandtModule(11).right_ideals()[1]1962sage: type(I)1963<class 'sage.algebras.quatalg.quaternion_algebra.QuaternionFractionalIdeal_rational'>1964sage: I.__repr__()1965'Fractional ideal (2 + 6*j + 4*k, 2*i + 4*j + 2*k, 8*j, 8*k)'1966"""1967return 'Fractional ideal %s'%(self.gens(),)19681969def quaternion_order(self):1970"""1971Return the order for which this ideal is a left or right1972fractional ideal. If this ideal has both a left and right1973ideal structure, then the left order is returned. If it has1974neither structure, then an error is raised.19751976OUTPUT: QuaternionOrder19771978EXAMPLES::19791980sage: R = QuaternionAlgebra(-11,-1).maximal_order()1981sage: R.unit_ideal().quaternion_order() is R1982True1983"""1984try: return self.__quaternion_order1985except AttributeError: pass1986if self.__left_order is not None:1987A = self.__left_order1988elif self.__right_order is not None:1989A = self.__right_order1990else:1991raise RuntimeError("unable to determine quaternion order of ideal without known order")1992self.__quaternion_order = A1993return A19941995def ring(self):1996"""1997Return ring that this is a fractional ideal for.19981999EXAMPLES::20002001sage: R = QuaternionAlgebra(-11,-1).maximal_order()2002sage: R.unit_ideal().ring() is R2003True2004"""2005return self.quaternion_order()20062007def basis(self):2008"""2009Return basis for this fractional ideal. The basis is in Hermite form.20102011OUTPUT: tuple20122013EXAMPLES::20142015sage: QuaternionAlgebra(-11,-1).maximal_order().unit_ideal().basis()2016(1/2 + 1/2*i, 1/2*j - 1/2*k, i, -k)2017"""2018return self.__basis20192020def gens(self):2021"""2022Return the generators for this ideal, which are the same as2023the `\\ZZ`-basis for this ideal.20242025EXAMPLES::20262027sage: QuaternionAlgebra(-11,-1).maximal_order().unit_ideal().gens()2028(1/2 + 1/2*i, 1/2*j - 1/2*k, i, -k)2029"""2030return self.__basis20312032def __cmp__(self, right):2033"""2034Compare this fractional quaternion ideal to ``right``. If2035``right`` is not a fractional quaternion ideal a TypeError is2036raised. If the fractional ideals are in different ambient2037quaternion algebras, then the quaternion algebras themselves2038are compared.20392040INPUT:20412042- ``right`` - another fractional quaternion ideal20432044EXAMPLES::20452046sage: I = QuaternionAlgebra(-11,-1).maximal_order().unit_ideal()2047sage: I == I # indirect doctest2048True2049sage: I == 52050False2051"""2052if not isinstance(right, QuaternionFractionalIdeal_rational):2053raise TypeError2054return cmp(self.__basis, right.__basis)20552056def basis_matrix(self):2057"""2058Return basis matrix `M` in Hermite normal form for self as a2059matrix with rational entries.20602061If `Q` is the ambient quaternion algebra, then the `\\ZZ`-span of2062the rows of `M` viewed as linear combinations of Q.basis() =2063`[1,i,j,k]` is the fractional ideal self. Also, M *2064M.denominator() is an integer matrix in Hermite normal form.20652066OUTPUT: matrix over `\\QQ`20672068EXAMPLES::20692070sage: QuaternionAlgebra(-11,-1).maximal_order().unit_ideal().basis_matrix()2071[ 1/2 1/2 0 0]2072[ 0 0 1/2 -1/2]2073[ 0 1 0 0]2074[ 0 0 0 -1]2075"""2076try: return self.__hermite_basis_matrix2077except AttributeError: pass2078B = quaternion_algebra_cython.rational_matrix_from_rational_quaternions(self.__basis)2079self.__hermite_basis_matrix = B2080return B20812082def free_module(self):2083"""2084Return the free module associated to this quaternionic2085fractional ideal, viewed as a submodule of2086``Q.free_module()``, where ``Q`` is the ambient quaternion2087algebra.20882089OUTPUT: free `\\ZZ`-module of rank 4 embedded in an ambient `\\QQ^4`.20902091EXAMPLES::20922093sage: QuaternionAlgebra(-11,-1).maximal_order().unit_ideal().basis_matrix()2094[ 1/2 1/2 0 0]2095[ 0 0 1/2 -1/2]2096[ 0 1 0 0]2097[ 0 0 0 -1]20982099This shows that the issue at trac ticket 6760 is fixed::21002101sage: R.<i,j,k> = QuaternionAlgebra(-1, -13)2102sage: I = R.ideal([2+i, 3*i, 5*j, j+k]); I2103Fractional ideal (2 + i, 3*i, j + k, 5*k)2104sage: I.free_module()2105Free module of degree 4 and rank 4 over Integer Ring2106Echelon basis matrix:2107[2 1 0 0]2108[0 3 0 0]2109[0 0 1 1]2110[0 0 0 5]2111"""2112try: return self.__free_module2113except AttributeError:2114M = self.basis_matrix().row_module(ZZ)2115self.__free_module = M2116return M21172118def theta_series_vector(self, B):2119"""2120Return theta series coefficients of self, as a vector of `B` integers.21212122INPUT:21232124- ``B`` -- positive integer21252126OUTPUT: vector over `\\ZZ` with ``B`` entries21272128EXAMPLES::21292130sage: I = BrandtModule(37).right_ideals()[1]; I2131Fractional ideal (2 + 6*j + 2*k, i + 2*j + k, 8*j, 8*k)2132sage: I.theta_series_vector(5)2133(1, 0, 2, 2, 6)2134sage: I.theta_series_vector(10)2135(1, 0, 2, 2, 6, 4, 8, 6, 10, 10)2136sage: I.theta_series_vector(5)2137(1, 0, 2, 2, 6)2138"""2139B = Integer(B)2140try:2141if len(self.__theta_series_vector)>= B: return self.__theta_series_vector[:B]2142except AttributeError: pass2143V = FreeModule(ZZ, B)2144Q = self.quadratic_form()2145v = V(Q.representation_number_list(B))2146self.__theta_series_vector = v2147return v21482149def quadratic_form(self):2150"""2151Return the normalized quadratic form associated to this quaternion ideal.21522153OUTPUT: quadratic form21542155EXAMPLES::21562157sage: I = BrandtModule(11).right_ideals()[1]2158sage: Q = I.quadratic_form(); Q2159Quadratic form in 4 variables over Rational Field with coefficients:2160[ 18 22 33 22 ]2161[ * 7 22 11 ]2162[ * * 22 0 ]2163[ * * * 22 ]2164sage: Q.theta_series(10)21651 + 12*q^2 + 12*q^3 + 12*q^4 + 12*q^5 + 24*q^6 + 24*q^7 + 36*q^8 + 36*q^9 + O(q^10)2166sage: I.theta_series(10)21671 + 12*q^2 + 12*q^3 + 12*q^4 + 12*q^5 + 24*q^6 + 24*q^7 + 36*q^8 + 36*q^9 + O(q^10)2168"""2169try: return self.__quadratic_form2170except AttributeError: pass2171from sage.quadratic_forms.quadratic_form import QuadraticForm2172# first get the gram matrix2173gram_matrix = self.gram_matrix()2174# rescale so that there are no denominators2175gram_matrix, _ = gram_matrix._clear_denom()2176# Make sure gcd of all entries is 1.2177g = gram_matrix.gcd()2178if g != 1:2179gram_matrix = gram_matrix / g2180# now get the quadratic form2181Q = QuadraticForm(gram_matrix)2182self.__quadratic_form = Q2183return Q21842185def theta_series(self, B, var='q'):2186"""2187Return normalized theta series of self, as a power series over2188`\\ZZ` in the variable ``var``, which is 'q' by default.21892190The normalized theta series is by definition21912192.. math::21932194\\theta_I(q)=\\sum_{x \\in I} q^{\\frac{N(x)}{N(I)}}21952196INPUT:21972198- ``B`` -- positive integer2199- ``var`` -- string (default: 'q')22002201OUTPUT: power series22022203EXAMPLES::22042205sage: I = BrandtModule(11).right_ideals()[1]; I2206Fractional ideal (2 + 6*j + 4*k, 2*i + 4*j + 2*k, 8*j, 8*k)2207sage: I.norm()2208322209sage: I.theta_series(5)22101 + 12*q^2 + 12*q^3 + 12*q^4 + O(q^5)2211sage: I.theta_series(5,'T')22121 + 12*T^2 + 12*T^3 + 12*T^4 + O(T^5)2213sage: I.theta_series(3)22141 + 12*q^2 + O(q^3)2215"""2216try:2217if self.__theta_series.prec() >= B:2218if var == self.__theta_series.variable():2219return self.__theta_series.add_bigoh(B)2220else:2221ZZ[[var]](self.__theta_series.list()[:B+1])2222except AttributeError: pass2223v = self.theta_series_vector(B)2224theta = ZZ[[var]](v.list()).add_bigoh(B)2225self.__theta_series = theta2226return theta22272228def gram_matrix(self):2229"""2230Return the Gram matrix of this fractional ideal.22312232OUTPUT: 4x4 matrix over `\\QQ`.22332234EXAMPLES::22352236sage: I = BrandtModule(3,5).right_ideals()[1]; I2237Fractional ideal (2 + 6*j + 4*k, 2*i + 4*j + 34*k, 8*j + 32*k, 40*k)2238sage: I.gram_matrix()2239[ 640 1920 2112 1920]2240[ 1920 14080 13440 16320]2241[ 2112 13440 13056 15360]2242[ 1920 16320 15360 19200]2243"""2244try: return self.__gram_matrix2245except AttributeError: pass2246M = []2247A = self.__basis2248B = [z.conjugate() for z in self.__basis]2249two = QQ(2)2250m = [two*(a*b).reduced_trace() for b in B for a in A]2251M44 = MatrixSpace(QQ, 4)2252G = M44(m,coerce=False)2253self.__gram_matrix = G2254return G22552256def norm(self):2257"""2258Return the reduced norm of this fractional ideal.22592260OUTPUT: rational number22612262EXAMPLES::22632264sage: M = BrandtModule(37)2265sage: C = M.right_ideals()2266sage: [I.norm() for I in C]2267[16, 32, 32]22682269sage: (a,b) = M.quaternion_algebra().invariants() # optional - magma2270sage: magma.eval('A<i,j,k> := QuaternionAlgebra<Rationals() | %s, %s>' % (a,b)) # optional - magma2271''2272sage: magma.eval('O := QuaternionOrder(%s)' % str(list(C[0].right_order().basis()))) # optional - magma2273''2274sage: [ magma('rideal<O | %s>' % str(list(I.basis()))).Norm() for I in C] # optional - magma2275[16, 32, 32]22762277sage: A.<i,j,k> = QuaternionAlgebra(-1,-1)2278sage: R = A.ideal([i,j,k,1/2 + 1/2*i + 1/2*j + 1/2*k]) # this is actually an order, so has reduced norm 12279sage: R.norm()228012281sage: [ J.norm() for J in R.cyclic_right_subideals(3) ] # enumerate maximal right R-ideals of reduced norm 3, verify their norms2282[3, 3, 3, 3]2283"""2284G = self.gram_matrix() / QQ(2)2285r = G.det().abs()2286assert r.is_square(), "first is bad!"2287r = r.sqrt()2288# If we know either the left- or the right order, use that one to compute the norm.2289# Otherwise quaternion_order() will raise a RuntimeError and we compute the left order2290try:2291R = self.quaternion_order()2292except RuntimeError:2293R = self.left_order()2294r/= R.discriminant()2295assert r.is_square(), "second is bad!"2296return r.sqrt()22972298def conjugate(self):2299"""2300Return the ideal with generators the conjugates of the generators for self.23012302OUTPUT: a quaternionic fractional ideal23032304EXAMPLES::23052306sage: I = BrandtModule(3,5).right_ideals()[1]; I2307Fractional ideal (2 + 6*j + 4*k, 2*i + 4*j + 34*k, 8*j + 32*k, 40*k)2308sage: I.conjugate()2309Fractional ideal (2 + 2*j + 28*k, 2*i + 4*j + 34*k, 8*j + 32*k, 40*k)2310"""2311return self.quaternion_algebra().ideal([b.conjugate() for b in self.basis()],2312left_order=self.__right_order,2313right_order=self.__left_order)23142315def __mul__(self, right):2316"""2317Return the product of the fractional ideals ``self`` and ``right``.23182319.. note::23202321We do not keep track of left or right order structure.23222323EXAMPLES::23242325sage: I = BrandtModule(3,5).right_ideals()[1]; I2326Fractional ideal (2 + 6*j + 4*k, 2*i + 4*j + 34*k, 8*j + 32*k, 40*k)2327sage: I*I2328Fractional ideal (8 + 24*j + 16*k, 8*i + 16*j + 136*k, 32*j + 128*k, 160*k)2329sage: I*I.conjugate()2330Fractional ideal (16 + 16*j + 224*k, 8*i + 16*j + 136*k, 32*j + 128*k, 320*k)2331sage: I.multiply_by_conjugate(I)2332Fractional ideal (16 + 16*j + 224*k, 8*i + 16*j + 136*k, 32*j + 128*k, 320*k)2333"""2334if not isinstance(right, QuaternionFractionalIdeal_rational):2335return self._scale(right, left=False)2336gens = [a*b for a in self.basis() for b in right.basis()]2337#if self.__right_order == right.__left_order:2338# left_order = self.__left_order2339# right_order = right.__right_order2340basis = tuple(basis_for_quaternion_lattice(gens))2341A = self.quaternion_algebra()2342return A.ideal(basis, check=False)23432344@cached_method2345def free_module(self):2346r"""2347Return the underlying free `\ZZ`-module corresponding to this ideal.23482349EXAMPLES::23502351sage: X = BrandtModule(3,5).right_ideals()2352sage: X[0]2353Fractional ideal (2 + 2*j + 8*k, 2*i + 18*k, 4*j + 16*k, 20*k)2354sage: X[0].free_module()2355Free module of degree 4 and rank 4 over Integer Ring2356Echelon basis matrix:2357[ 2 0 2 8]2358[ 0 2 0 18]2359[ 0 0 4 16]2360[ 0 0 0 20]2361sage: X[0].scale(1/7).free_module()2362Free module of degree 4 and rank 4 over Integer Ring2363Echelon basis matrix:2364[ 2/7 0 2/7 8/7]2365[ 0 2/7 0 18/7]2366[ 0 0 4/7 16/7]2367[ 0 0 0 20/7]23682369The free module method is also useful since it allows for checking if one ideal2370is contained in another, computing quotients I/J, etc.::23712372sage: X = BrandtModule(3,17).right_ideals()2373sage: I = X[0].intersection(X[2]); I2374Fractional ideal (2 + 2*j + 164*k, 2*i + 4*j + 46*k, 16*j + 224*k, 272*k)2375sage: I.free_module().is_submodule(X[3].free_module())2376False2377sage: I.free_module().is_submodule(X[1].free_module())2378True2379sage: X[0].free_module() / I.free_module()2380Finitely generated module V/W over Integer Ring with invariants (4, 4)2381"""2382return self.basis_matrix().row_module(ZZ)23832384def intersection(self, J):2385"""2386Return the intersection of the ideals self and `J`.23872388EXAMPLES::23892390sage: X = BrandtModule(3,5).right_ideals()2391sage: I = X[0].intersection(X[1]); I2392Fractional ideal (2 + 6*j + 4*k, 2*i + 4*j + 34*k, 8*j + 32*k, 40*k)23932394"""2395V = self.free_module().intersection(J.free_module())2396H,d = V.basis_matrix()._clear_denom()2397A = self.quaternion_algebra()2398gens = quaternion_algebra_cython.rational_quaternions_from_integral_matrix_and_denom(A, H, d)2399return A.ideal(gens)24002401def multiply_by_conjugate(self, J):2402"""2403Return product of self and the conjugate Jbar of `J`.24042405INPUT:24062407- ``J`` -- a quaternion ideal.24082409OUTPUT: a quaternionic fractional ideal.24102411EXAMPLES::24122413sage: R = BrandtModule(3,5).right_ideals()2414sage: R[0].multiply_by_conjugate(R[1])2415Fractional ideal (8 + 8*j + 112*k, 8*i + 16*j + 136*k, 32*j + 128*k, 160*k)2416sage: R[0]*R[1].conjugate()2417Fractional ideal (8 + 8*j + 112*k, 8*i + 16*j + 136*k, 32*j + 128*k, 160*k)2418"""2419Jbar = [b.conjugate() for b in J.basis()]2420gens = [a*b for a in self.basis() for b in Jbar]2421basis = tuple(basis_for_quaternion_lattice(gens))2422R = self.quaternion_algebra()2423return R.ideal(basis, check=False)24242425def is_equivalent(I, J, B=10):2426"""2427Return ``True`` if ``I`` and ``J`` are equivalent as right ideals.24282429INPUT:24302431- ``I`` -- a fractional quaternion ideal (self)24322433- ``J`` -- a fractional quaternion ideal with same order as ``I``24342435- ``B`` -- a bound to compute and compare theta series before2436doing the full equivalence test24372438OUTPUT: bool24392440EXAMPLES::24412442sage: R = BrandtModule(3,5).right_ideals(); len(R)244322444sage: R[0].is_equivalent(R[1])2445False2446sage: R[0].is_equivalent(R[0])2447True2448sage: OO = R[0].quaternion_order()2449sage: S = OO.right_ideal([3*a for a in R[0].basis()])2450sage: R[0].is_equivalent(S)2451True2452"""2453if not isinstance(I, QuaternionFractionalIdeal_rational):2454return False24552456if I.right_order() != J.right_order():2457raise ValueError("I and J must be right ideals")24582459# Just test theta series first. If the theta series are2460# different, the ideals are definitely not equivalent.2461if B > 0 and I.theta_series_vector(B) != J.theta_series_vector(B):2462return False24632464# The theta series are the same, so perhaps the ideals are2465# equivalent. We use Prop 1.18 of [Pizer, 1980] to decide.2466# 1. Compute I * Jbar2467# see Prop. 1.17 in Pizer. Note that we use IJbar instead of2468# JbarI since we work with right ideals2469IJbar = I.multiply_by_conjugate(J)24702471# 2. Determine if there is alpha in K such2472# that N(alpha) = N(I)*N(J) as explained by Pizer.2473c = IJbar.theta_series_vector(2)[1]2474return c != 024752476def __contains__(self, x):2477"""2478Returns whether x is in self.24792480EXAMPLES::2481sage: R.<i,j,k> = QuaternionAlgebra(-3, -13)2482sage: I = R.ideal([2+i, 3*i, 5*j, j+k])2483sage: 2+i in I2484True2485sage: 2+i+j+k in I2486True2487sage: 1+i in I2488False2489sage: 101*j + k in I2490True2491"""2492try:2493x = self.quaternion_algebra()(x)2494return self.basis_matrix().transpose().solve_right(vector(x)) in ZZ**42495except (ValueError, TypeError):2496return False24972498@cached_method2499def cyclic_right_subideals(self, p, alpha=None):2500r"""2501Let `I` = ``self``. This function returns the right subideals2502`J` of `I` such that `I/J` is an `\GF{p}`-vector space of2503dimension 2.25042505INPUT:25062507- `p` -- prime number (see below)25082509- `alpha` -- (default: None) element of quaternion algebra,2510which can be used to parameterize the order of the2511ideals `J`. More precisely the `J`'s are the right annihilators2512of `(1,0) \alpha^i` for `i=0,1,2,...,p`25132514OUTPUT:25152516- list of right ideals25172518.. note::25192520Currently, `p` must satisfy a bunch of conditions, or a2521NotImplementedError is raised. In particular, `p` must be2522odd and unramified in the quaternion algebra, must be2523coprime to the index of the right order in the maximal2524order, and also coprime to the normal of self. (The Brandt2525modules code has a more general algorithm in some cases.)25262527EXAMPLES::25282529sage: B = BrandtModule(2,37); I = B.right_ideals()[0]2530sage: I.cyclic_right_subideals(3)2531[Fractional ideal (2 + 2*i + 10*j + 90*k, 4*i + 4*j + 152*k, 12*j + 132*k, 444*k), Fractional ideal (2 + 2*i + 2*j + 150*k, 4*i + 8*j + 196*k, 12*j + 132*k, 444*k), Fractional ideal (2 + 2*i + 6*j + 194*k, 4*i + 8*j + 344*k, 12*j + 132*k, 444*k), Fractional ideal (2 + 2*i + 6*j + 46*k, 4*i + 4*j + 4*k, 12*j + 132*k, 444*k)]25322533sage: B = BrandtModule(5,389); I = B.right_ideals()[0]2534sage: C = I.cyclic_right_subideals(3); C2535[Fractional ideal (2 + 10*j + 546*k, i + 6*j + 133*k, 12*j + 3456*k, 4668*k), Fractional ideal (2 + 2*j + 2910*k, i + 6*j + 3245*k, 12*j + 3456*k, 4668*k), Fractional ideal (2 + i + 2295*k, 3*i + 2*j + 3571*k, 4*j + 2708*k, 4668*k), Fractional ideal (2 + 2*i + 2*j + 4388*k, 3*i + 2*j + 2015*k, 4*j + 4264*k, 4668*k)]2536sage: [(I.free_module()/J.free_module()).invariants() for J in C]2537[(3, 3), (3, 3), (3, 3), (3, 3)]2538sage: I.scale(3).cyclic_right_subideals(3)2539[Fractional ideal (6 + 30*j + 1638*k, 3*i + 18*j + 399*k, 36*j + 10368*k, 14004*k), Fractional ideal (6 + 6*j + 8730*k, 3*i + 18*j + 9735*k, 36*j + 10368*k, 14004*k), Fractional ideal (6 + 3*i + 6885*k, 9*i + 6*j + 10713*k, 12*j + 8124*k, 14004*k), Fractional ideal (6 + 6*i + 6*j + 13164*k, 9*i + 6*j + 6045*k, 12*j + 12792*k, 14004*k)]2540sage: C = I.scale(1/9).cyclic_right_subideals(3); C2541[Fractional ideal (2/9 + 10/9*j + 182/3*k, 1/9*i + 2/3*j + 133/9*k, 4/3*j + 384*k, 1556/3*k), Fractional ideal (2/9 + 2/9*j + 970/3*k, 1/9*i + 2/3*j + 3245/9*k, 4/3*j + 384*k, 1556/3*k), Fractional ideal (2/9 + 1/9*i + 255*k, 1/3*i + 2/9*j + 3571/9*k, 4/9*j + 2708/9*k, 1556/3*k), Fractional ideal (2/9 + 2/9*i + 2/9*j + 4388/9*k, 1/3*i + 2/9*j + 2015/9*k, 4/9*j + 4264/9*k, 1556/3*k)]2542sage: [(I.scale(1/9).free_module()/J.free_module()).invariants() for J in C]2543[(3, 3), (3, 3), (3, 3), (3, 3)]25442545sage: Q.<i,j,k> = QuaternionAlgebra(-2,-5)2546sage: I = Q.ideal([Q(1),i,j,k])2547sage: I.cyclic_right_subideals(3)2548[Fractional ideal (1 + 2*j, i + k, 3*j, 3*k), Fractional ideal (1 + j, i + 2*k, 3*j, 3*k), Fractional ideal (1 + 2*i, 3*i, j + 2*k, 3*k), Fractional ideal (1 + i, 3*i, j + k, 3*k)]25492550The general algorithm is not yet implemented here::25512552sage: I.cyclic_right_subideals(3)[0].cyclic_right_subideals(3)2553Traceback (most recent call last):2554...2555NotImplementedError: general algorithm not implemented (The given basis vectors must be linearly independent.)2556"""2557R = self.right_order()2558Q = self.quaternion_algebra()2559f = Q.modp_splitting_map(p)2560if alpha is not None:2561alpha = f(alpha)2562W = GF(p)**42563try:2564A = W.span_of_basis([W(f(a).list()) for a in self.basis()])2565scale = 12566IB = self.basis_matrix()2567except (ValueError, ZeroDivisionError):2568# try rescaling the ideal.2569B, d = self.basis_matrix()._clear_denom()2570g = gcd(B.list())2571IB = B/g2572scale = g/d2573try:2574A = W.span_of_basis([W(f(Q(a.list())).list()) for a in IB.rows()])2575except (ValueError, ZeroDivisionError) as msg:2576# Here we could replace the ideal by an *equivalent*2577# ideal that works. This is always possible.2578# However, I haven't implemented that algorithm yet.2579raise NotImplementedError("general algorithm not implemented (%s)"%msg)25802581Ai = A.basis_matrix()**(-1)2582AiB = Ai.change_ring(QQ) * IB25832584# Do not care about the denominator since we're really working in I/p*I.2585AiB, _ = AiB._clear_denom()25862587pB = p*IB2588pB, d = pB._clear_denom()25892590ans = []2591Z = matrix(ZZ,2,4)25922593P1 = P1List(p)2594if alpha is None:2595lines = P12596else:2597x = alpha2598lines = []2599for i in range(p+1):2600lines.append(P1.normalize(x[0,0], x[0,1]))2601x *= alpha26022603for u,v in lines:2604# The following does:2605# z = matrix(QQ,2,4,[0,-v,0,u, -v,0,u,0],check=False) * AiB2606Z[0,1]=-v; Z[0,3]=u; Z[1,0]=-v; Z[1,2]=u2607z = Z * AiB2608# Now construct submodule of the ideal I spanned by the2609# linear combinations given by z of the basis for J along2610# with p*I.2611G = (d*z).stack(pB) # have to multiply by d since we divide by it below in the "gens = " line.2612H = G._hnf_pari(0, include_zero_rows=False)2613gens = tuple(quaternion_algebra_cython.rational_quaternions_from_integral_matrix_and_denom(Q, H, d))2614if scale != 1:2615gens = tuple([scale*g for g in gens])2616J = R.right_ideal(gens, check=False)2617ans.append(J)2618return ans26192620#######################################################################2621# Some utility functions that are needed here and are too2622# specialized to go elsewhere.2623#######################################################################26242625def basis_for_quaternion_lattice(gens, reverse = False):2626"""2627Return a basis for the `\\ZZ`-lattice in a quaternion algebra2628spanned by the given gens.26292630INPUT:26312632- ``gens`` -- list of elements of a single quaternion algebra26332634- ``reverse`` -- when computing the HNF do it on the basis2635(k,j,i,1) instead of (1,i,j,k). This ensures2636that if ``gens`` are the generators for an order,2637the first returned basis vector is 1.26382639EXAMPLES::26402641sage: from sage.algebras.quatalg.quaternion_algebra import basis_for_quaternion_lattice2642sage: A.<i,j,k> = QuaternionAlgebra(-1,-7)2643sage: basis_for_quaternion_lattice([i+j, i-j, 2*k, A(1/3)])2644[1/3, i + j, 2*j, 2*k]26452646sage: basis_for_quaternion_lattice([A(1),i,j,k])2647[1, i, j, k]26482649"""2650if len(gens) == 0: return []2651Z, d = quaternion_algebra_cython.integral_matrix_and_denom_from_rational_quaternions(gens, reverse)2652H = Z._hnf_pari(0, include_zero_rows=False)2653A = gens[0].parent()2654return quaternion_algebra_cython.rational_quaternions_from_integral_matrix_and_denom(A, H, d, reverse)265526562657def intersection_of_row_modules_over_ZZ(v):2658"""2659Intersects the `\ZZ`-modules with basis matrices the full rank 4x42660`\QQ`-matrices in the list v. The returned intersection is2661represented by a 4x4 matrix over `\QQ`. This can also be done using2662modules and intersection, but that would take over twice as long2663because of overhead, hence this function.26642665EXAMPLES::26662667sage: a = matrix(QQ,4,[-2, 0, 0, 0, 0, -1, -1, 1, 2, -1/2, 0, 0, 1, 1, -1, 0])2668sage: b = matrix(QQ,4,[0, -1/2, 0, -1/2, 2, 1/2, -1, -1/2, 1, 2, 1, -2, 0, -1/2, -2, 0])2669sage: c = matrix(QQ,4,[0, 1, 0, -1/2, 0, 0, 2, 2, 0, -1/2, 1/2, -1, 1, -1, -1/2, 0])2670sage: v = [a,b,c]2671sage: from sage.algebras.quatalg.quaternion_algebra import intersection_of_row_modules_over_ZZ2672sage: M = intersection_of_row_modules_over_ZZ(v); M2673[ -2 0 1 1]2674[ -4 1 1 -3]2675[ -3 19/2 -1 -4]2676[ 2 -3 -8 4]2677sage: M2 = a.row_module(ZZ).intersection(b.row_module(ZZ)).intersection(c.row_module(ZZ))2678sage: M.row_module(ZZ) == M22679True2680"""2681if len(v) <= 0:2682raise ValueError("v must have positive length")2683if len(v) == 1:2684return v[0]2685elif len(v) == 2:2686# real work - the base case2687a, b = v2688s,_ = a.stack(b)._clear_denom()2689s = s.transpose()2690K = s.right_kernel_matrix(algorithm='pari', basis='computed')2691n = a.nrows()2692return K.matrix_from_columns(range(n)) * a2693else:2694# induct2695w = intersection_of_row_modules_over_ZZ(v[:2])2696return intersection_of_row_modules_over_ZZ([w] + v[2:])269726982699def normalize_basis_at_p(e, p, B = lambda x,y: (x*y.conjugate()).reduced_trace()):2700"""2701Computes a (at ``p``) normalized basis from the given basis ``e``2702of a `\\ZZ`-module.27032704The returned basis is (at ``p``) a `\\ZZ_p` basis for the same2705module, and has the property that with respect to it the quadratic2706form induced by the bilinear form B is represented as a orthogonal2707sum of atomic forms multiplied by p-powers.27082709If `p \\ne 2` this means that the form is diagonal with respect to2710this basis.27112712If `p = 2` there may be additional 2-dimensional subspaces on which2713the form is represented as `2^e (ax^2 + bxy + cx^2)` with2714`0 = v_2(b) = v_2(a) \\le v_2(c)`.27152716INPUT:2717- ``e`` -- list; basis of a `\\ZZ` module.2718WARNING: will be modified!27192720- ``p`` -- prime for at which the basis should be normalized27212722- ``B`` -- (default: lambda x,y: ((x*y).conjugate()).reduced_trace())2723A bilinear form with respect to which to normalize.27242725OUTPUT:27262727- A list containing two-element tuples: The first element of2728each tuple is a basis element, the second the valuation of2729the orthogonal summand to which it belongs. The list is sorted2730by ascending valuation.27312732EXAMPLES::27332734sage: from sage.algebras.quatalg.quaternion_algebra import normalize_basis_at_p2735sage: A.<i,j,k> = QuaternionAlgebra(-1, -1)2736sage: e = [A(1), i, j, k]2737sage: normalize_basis_at_p(e, 2)2738[(1, 0), (i, 0), (j, 0), (k, 0)]27392740sage: A.<i,j,k> = QuaternionAlgebra(210)2741sage: e = [A(1), i, j, k]2742sage: normalize_basis_at_p(e, 2)2743[(1, 0), (i, 1), (j, 1), (k, 2)]27442745sage: A.<i,j,k> = QuaternionAlgebra(286)2746sage: e = [A(1), k, 1/2*j + 1/2*k, 1/2 + 1/2*i + 1/2*k]2747sage: normalize_basis_at_p(e, 5)2748[(1, 0), (1/2*j + 1/2*k, 0), (-5/6*j + 1/6*k, 1), (1/2*i, 1)]27492750sage: A.<i,j,k> = QuaternionAlgebra(-1,-7)2751sage: e = [A(1), k, j, 1/2 + 1/2*i + 1/2*j + 1/2*k]2752sage: normalize_basis_at_p(e, 2)2753[(1, 0), (1/2 + 1/2*i + 1/2*j + 1/2*k, 0), (-34/105*i - 463/735*j + 71/105*k, 1), (-34/105*i - 463/735*j + 71/105*k, 1)]2754"""27552756N = len(e)2757if N == 0:2758return []2759else:2760min_m, min_n, min_v = 0, 0, infinity27612762# Find two basis vector on which the bilinear form has minimal2763# p-valuation. If there is more than one such pair, always2764# prefer diagonal entries over any other and (secondary) take2765# min_m and then min_n as small as possible2766for m in range(N):2767for n in range(m, N):2768v = B(e[m], e[n]).valuation(p)2769if v < min_v or (v == min_v and (min_m != min_n) and (m == n)):2770min_m, min_n, min_v = m, n, v277127722773if (min_m == min_n) or p != 2: # In this case we can diagonalize2774if min_m == min_n: # Diagonal entry has minimal valuation2775f0 = e[min_m]2776else:2777f0 = e[min_m] + e[min_n] # Only off-diagonal entries have min. val., but p!=227782779# Swap with first vector2780e[0], e[min_m] = e[min_m], e[0]27812782# Orthogonalize remaining vectors with respect to f2783c = B(f0, f0)2784for l in range(1, N):2785e[l] = e[l] - B(e[l],f0)/c * f027862787# Recursively normalize remaining vectors2788f = normalize_basis_at_p(e[1:], p)2789f.insert(0, (f0, min_v - valuation(p, 2)))2790return f27912792else: # p = 2 and only off-diagonal entries have min. val., gives 2-dim. block2793# first diagonal entry should have smaller valuation2794if B(e[min_m],e[min_m]).valuation(p) > B(e[min_n],e[min_n]).valuation(p):2795e[min_m], e[min_n] = e[min_n], e[min_m]27962797f0 = p**min_v / B(e[min_m],e[min_n]) * e[min_m]2798f1 = e[min_n]27992800# Ensures that (B(f0,f0)/2).valuation(p) <= B(f0,f1).valuation(p)2801if B(f0,f1).valuation(p) + 1 < B(f0,f0).valuation(p):2802f0 = f0 + f12803f1 = f028042805# Make remaining vectors orthogonal to span of f0, f12806e[min_m] = e[0]2807e[min_n] = e[1]28082809B00 = B(f0,f0)2810B11 = B(f1,f1)2811B01 = B(f0,f1)2812d = B00*B11 - B01**22813tu = [ (B01 * B(f1,e[l]) - B11 * B(f0,e[l]),2814B01 * B(f0,e[l]) - B00 * B(f1,e[l])) for l in range(2,N) ]28152816e[2:n] = [ e[l] + tu[l-2][0]/d * f0 + tu[l-2][1]/d * f1 for l in range(2,N) ]28172818# Recursively normalize remaining vectors2819f = normalize_basis_at_p(e[2:N], p)2820return [(f0, min_v), (f1, min_v)] + f28212822def maxord_solve_aux_eq(a, b, p):2823"""2824Given ``a`` and ``b`` and an even prime ideal ``p`` find2825(y,z,w) with y a unit mod `p^{2e}` such that28262827.. MATH::28281 - ay^2 - bz^2 + abw^2 \equiv 0 mod p^{2e},28292830where `e` is the ramification index of `p`.28312832Currently only `p=2` is implemented by hardcoding solutions.28332834INPUT:28352836- ``a`` -- integer with `v_p(a) = 0`28372838- ``b`` -- integer with `v_p(b) \in \{0,1\}`28392840- ``p`` -- even prime ideal (actually only ``p=ZZ(2)`` is implemented)28412842OUTPUT:28432844- A tuple (y,z,w)28452846EXAMPLES::28472848sage: from sage.algebras.quatalg.quaternion_algebra import maxord_solve_aux_eq2849sage: for a in [1,3]:2850... for b in [1,2,3]:2851... (y,z,w) = maxord_solve_aux_eq(a, b, 2)2852... assert mod(y, 4) == 1 or mod(y, 4) == 32853... assert mod(1 - a*y^2 - b*z^2 + a*b*w^2, 4) == 02854"""2855if p != ZZ(2):2856raise NotImplementedError, "Algorithm only implemented over ZZ at the moment"28572858v_a = a.valuation(p)2859v_b = b.valuation(p)28602861if v_a != 0:2862raise RuntimeError, "a must have v_p(a)=0"2863if v_b != 0 and v_b != 1:2864raise RuntimeError, "b must have v_p(b) in {0,1}"28652866R = ZZ.quo(ZZ(4))2867lut = {2868(R(1), R(1)) : (1,1,1),2869(R(1), R(2)) : (1,0,0),2870(R(1), R(3)) : (1,0,0),2871(R(3), R(1)) : (1,1,1),2872(R(3), R(2)) : (1,0,1),2873(R(3), R(3)) : (1,1,1), }28742875return lut[ (R(a), R(b)) ]287628772878