Path: blob/master/src/sage/quadratic_forms/quadratic_form.py
8817 views
"""1Quadratic Forms Overview23AUTHORS:45- Jon Hanke (2007-06-19)6- Anna Haensch (2010-07-01): Formatting and ReSTification7"""89#*****************************************************************************10# Copyright (C) 2007 William Stein and Jonathan Hanke11#12# Distributed under the terms of the GNU General Public License (GPL)13#14# This code is distributed in the hope that it will be useful,15# but WITHOUT ANY WARRANTY; without even the implied warranty of16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU17# General Public License for more details.18#19# The full text of the GPL is available at:20#21# http://www.gnu.org/licenses/22#*****************************************************************************2324from warnings import warn25from copy import deepcopy2627from sage.matrix.constructor import matrix28from sage.matrix.matrix_space import MatrixSpace29from sage.matrix.matrix import is_Matrix30from sage.rings.integer_ring import IntegerRing, ZZ31from sage.rings.ring import Ring32from sage.misc.functional import ideal, denominator, is_even, is_field33from sage.rings.arith import GCD, LCM34from sage.rings.principal_ideal_domain import is_PrincipalIdealDomain35from sage.rings.ring import is_Ring36from sage.matrix.matrix import is_Matrix37from sage.structure.sage_object import SageObject38from sage.structure.element import is_Vector39from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing40from sage.modules.free_module_element import vector4142from sage.quadratic_forms.quadratic_form__evaluate import QFEvaluateVector, QFEvaluateMatrix43444546def QuadraticForm__constructor(R, n=None, entries=None):47"""48Wrapper for the QuadraticForm class constructor. This is meant49for internal use within the QuadraticForm class code only. You50should instead directly call QuadraticForm().5152EXAMPLES::5354sage: from sage.quadratic_forms.quadratic_form import QuadraticForm__constructor55sage: QuadraticForm__constructor(ZZ, 3) ## Makes a generic quadratic form over the integers56Quadratic form in 3 variables over Integer Ring with coefficients:57[ 0 0 0 ]58[ * 0 0 ]59[ * * 0 ]6061"""62return QuadraticForm(R, n, entries)636465def is_QuadraticForm(Q):66"""67Determines if the object Q is an element of the QuadraticForm class.6869EXAMPLES::7071sage: Q = QuadraticForm(ZZ, 2, [1,2,3])72sage: is_QuadraticForm(Q) ##random -- deprecated73True74sage: is_QuadraticForm(2) ##random -- deprecated75False7677"""78return isinstance(Q, QuadraticForm)79808182class QuadraticForm(SageObject):83r"""84The ``QuadraticForm`` class represents a quadratic form in n variables with85coefficients in the ring R.8687INPUT:8889The constructor may be called in any of the following ways.9091#. ``QuadraticForm(R, n, entries)``, where9293- `R` -- ring for which the quadratic form is defined94- `n` -- an integer >= 095- ``entries`` -- a list of `n(n+1)/2` coefficients of the quadratic form96in `R` (given lexographically, or equivalently, by rows of the matrix)9798#. ``QuadraticForm(R, n)``, where99100- `R` -- a ring101- `n` -- a symmetric `n \times n` matrix with even diagonal (relative to102`R`)103104#. ``QuadraticForm(R)``, where105106- `R` -- a symmetric `n \times n` matrix with even diagonal (relative to107its base ring)108109If the keyword argument ``unsafe_initialize`` is True, then the subsequent110fields may by used to force the external initialization of various fields111of the quadratic form. Currently the only fields which can be set are:112113- ``number_of_automorphisms``114- ``determinant``115116117OUTPUT:118119quadratic form120121EXAMPLES::122123sage: Q = QuadraticForm(ZZ, 3, [1,2,3,4,5,6])124sage: Q125Quadratic form in 3 variables over Integer Ring with coefficients:126[ 1 2 3 ]127[ * 4 5 ]128[ * * 6 ]129130::131132sage: Q = QuadraticForm(QQ, 3, [1,2,3,4/3 ,5,6])133sage: Q134Quadratic form in 3 variables over Rational Field with coefficients:135[ 1 2 3 ]136[ * 4/3 5 ]137[ * * 6 ]138sage: Q[0,0]1391140sage: Q[0,0].parent()141Rational Field142143::144145sage: Q = QuadraticForm(QQ, 7, range(28))146sage: Q147Quadratic form in 7 variables over Rational Field with coefficients:148[ 0 1 2 3 4 5 6 ]149[ * 7 8 9 10 11 12 ]150[ * * 13 14 15 16 17 ]151[ * * * 18 19 20 21 ]152[ * * * * 22 23 24 ]153[ * * * * * 25 26 ]154[ * * * * * * 27 ]155156::157158sage: Q = QuadraticForm(QQ, 2, range(1,4))159sage: A = Matrix(ZZ,2,2,[-1,0,0,1])160sage: Q(A)161Quadratic form in 2 variables over Rational Field with coefficients:162[ 1 -2 ]163[ * 3 ]164165::166167sage: m = matrix(2,2,[1,2,3,4])168sage: m + m.transpose()169[2 5]170[5 8]171sage: QuadraticForm(m + m.transpose())172Quadratic form in 2 variables over Integer Ring with coefficients:173[ 1 5 ]174[ * 4 ]175176::177178sage: QuadraticForm(ZZ, m + m.transpose())179Quadratic form in 2 variables over Integer Ring with coefficients:180[ 1 5 ]181[ * 4 ]182183::184185sage: QuadraticForm(QQ, m + m.transpose())186Quadratic form in 2 variables over Rational Field with coefficients:187[ 1 5 ]188[ * 4 ]189"""190191## Import specialized methods:192## ---------------------------193194## Routines to compute the p-adic local normal form195from sage.quadratic_forms.quadratic_form__local_normal_form import \196find_entry_with_minimal_scale_at_prime, \197local_normal_form, \198jordan_blocks_by_scale_and_unimodular, \199jordan_blocks_in_unimodular_list_by_scale_power200201## Routines to perform elementary variable substitutions202from sage.quadratic_forms.quadratic_form__variable_substitutions import \203swap_variables, \204multiply_variable, \205divide_variable, \206scale_by_factor, \207extract_variables, \208elementary_substitution, \209add_symmetric210211## Routines to compute p-adic field invariants212from sage.quadratic_forms.quadratic_form__local_field_invariants import \213rational_diagonal_form, \214signature_vector, \215signature, \216hasse_invariant, \217hasse_invariant__OMeara, \218is_hyperbolic, \219is_anisotropic, \220is_isotropic, \221anisotropic_primes, \222compute_definiteness, \223compute_definiteness_string_by_determinants, \224is_positive_definite, \225is_negative_definite, \226is_indefinite, \227is_definite228229## Routines to compute local densities by the reduction procedure230from sage.quadratic_forms.quadratic_form__local_density_congruence import \231count_modp_solutions__by_Gauss_sum, \232local_good_density_congruence_odd, \233local_good_density_congruence_even, \234local_good_density_congruence, \235local_zero_density_congruence, \236local_badI_density_congruence, \237local_badII_density_congruence, \238local_bad_density_congruence, \239local_density_congruence, \240local_primitive_density_congruence241242## Routines to compute local densities by counting solutions of various types243from sage.quadratic_forms.quadratic_form__count_local_2 import \244count_congruence_solutions_as_vector, \245count_congruence_solutions, \246count_congruence_solutions__good_type, \247count_congruence_solutions__zero_type, \248count_congruence_solutions__bad_type, \249count_congruence_solutions__bad_type_I, \250count_congruence_solutions__bad_type_II251252## Routines to be called by the user to compute local densities253from sage.quadratic_forms.quadratic_form__local_density_interfaces import \254local_density, \255local_primitive_density256257## Routines for computing with ternary forms258from sage.quadratic_forms.quadratic_form__ternary_Tornaria import \259disc, \260content, \261adjoint, \262antiadjoint, \263is_adjoint, \264reciprocal, \265omega, \266delta, \267level__Tornaria, \268discrec, \269hasse_conductor, \270clifford_invariant, \271clifford_conductor, \272basiclemma, \273basiclemmavec, \274xi, \275xi_rec, \276lll, \277representation_number_list, \278representation_vector_list, \279is_zero, \280is_zero_nonsingular, \281is_zero_singular282283## Routines to compute the theta function284from sage.quadratic_forms.quadratic_form__theta import \285theta_series, \286theta_series_degree_2, \287theta_by_pari, \288theta_by_cholesky289290## Routines to compute the product of all local densities291from sage.quadratic_forms.quadratic_form__siegel_product import \292siegel_product293294## Routines to compute p-neighbors295from sage.quadratic_forms.quadratic_form__neighbors import \296find_primitive_p_divisible_vector__random, \297find_primitive_p_divisible_vector__next, \298find_p_neighbor_from_vec299300## Routines to reduce a given quadratic form301from sage.quadratic_forms.quadratic_form__reduction_theory import \302reduced_binary_form1, \303reduced_ternary_form__Dickson, \304reduced_binary_form, \305minkowski_reduction, \306minkowski_reduction_for_4vars__SP307## Wrappers for Conway-Sloane genus routines (in ./genera/)308from sage.quadratic_forms.quadratic_form__genus import \309global_genus_symbol, \310local_genus_symbol, \311CS_genus_symbol_list312313314## Routines to compute local masses for ZZ.315from sage.quadratic_forms.quadratic_form__mass import \316shimura_mass__maximal, \317GHY_mass__maximal318from sage.quadratic_forms.quadratic_form__mass__Siegel_densities import \319mass__by_Siegel_densities, \320Pall_mass_density_at_odd_prime, \321Watson_mass_at_2, \322Kitaoka_mass_at_2, \323mass_at_two_by_counting_mod_power324from sage.quadratic_forms.quadratic_form__mass__Conway_Sloane_masses import \325parity, \326is_even, \327is_odd, \328conway_species_list_at_odd_prime, \329conway_species_list_at_2, \330conway_octane_of_this_unimodular_Jordan_block_at_2, \331conway_diagonal_factor, \332conway_cross_product_doubled_power, \333conway_type_factor, \334conway_p_mass, \335conway_standard_p_mass, \336conway_standard_mass, \337conway_mass338# conway_generic_mass, \339# conway_p_mass_adjustment340341## Routines to check local representability of numbers342from sage.quadratic_forms.quadratic_form__local_representation_conditions import \343local_representation_conditions, \344is_locally_universal_at_prime, \345is_locally_universal_at_all_primes, \346is_locally_universal_at_all_places, \347is_locally_represented_number_at_place, \348is_locally_represented_number349350## Routines to make a split local covering of the given quadratic form.351from sage.quadratic_forms.quadratic_form__split_local_covering import \352cholesky_decomposition, \353vectors_by_length, \354complementary_subform_to_vector, \355split_local_cover356357## Routines to make automorphisms of the given quadratic form.358from sage.quadratic_forms.quadratic_form__automorphisms import \359basis_of_short_vectors, \360short_vector_list_up_to_length, \361short_primitive_vector_list_up_to_length, \362automorphisms, \363number_of_automorphisms, \364number_of_automorphisms__souvigner, \365set_number_of_automorphisms366367## Routines to test the local and global equivalence/isometry of two quadratic forms.368from sage.quadratic_forms.quadratic_form__equivalence_testing import \369is_globally_equivalent__souvigner, \370is_globally_equivalent_to, \371is_locally_equivalent_to, \372has_equivalent_Jordan_decomposition_at_prime373374def __init__(self, R, n=None, entries=None, unsafe_initialization=False, number_of_automorphisms=None, determinant=None):375"""376EXAMPLES::377378sage: s = QuadraticForm(ZZ, 4, range(10))379sage: s == loads(dumps(s))380True381"""382## Deal with: QuadraticForm(ring, matrix)383matrix_init_flag = False384if isinstance(R, Ring):385if is_Matrix(n):386## Test if n is symmetric and has even diagonal387if not self._is_even_symmetric_matrix_(n, R):388raise TypeError, "Oops! The matrix is not a symmetric with even diagonal defined over R."389390## Rename the matrix and ring391M = n392M_ring = R393matrix_init_flag = True394395396## Deal with: QuadraticForm(matrix)397if is_Matrix(R) and (n == None):398399## Test if R is symmetric and has even diagonal400if not self._is_even_symmetric_matrix_(R):401raise TypeError, "Oops! The matrix is not a symmetric with even diagonal."402403## Rename the matrix and ring404M = R405M_ring = R.base_ring()406matrix_init_flag = True407408## Perform the quadratic form initialization409if matrix_init_flag == True:410self.__n = M.nrows()411self.__base_ring = M_ring412self.__coeffs = []413for i in range(M.nrows()):414for j in range(i, M.nrows()):415if (i == j):416self.__coeffs += [ M_ring(M[i,j] / 2) ]417else:418self.__coeffs += [ M_ring(M[i,j]) ]419420return421422## -----------------------------------------------------------423424## Verify the size of the matrix is an integer >= 0425try:426n = int(n)427except Exception:428raise TypeError, "Oops! The size " + str(n) + " must be an integer."429if (n < 0):430raise TypeError, "Oops! The size " + str(n) + " must be a non-negative integer."431432## TODO: Verify that R is a ring...433434## Store the relevant variables435N = int(n*(n+1))/2436self.__n = int(n)437self.__base_ring = R438self.__coeffs = [self.__base_ring(0) for i in range(N)]439440## Check if entries is a list for the current size, and if so, write the upper-triangular matrix441if isinstance(entries, list) and (len(entries) == N):442for i in range(N):443self.__coeffs[i] = self.__base_ring(entries[i])444elif (entries != None):445raise TypeError, "Oops! The entries " + str(entries) + "must be a list of size n(n+1)/2."446447## -----------------------------------------------------------448449## Process possible forced initialization of various fields450self._external_initialization_list = []451if unsafe_initialization:452453## Set the number of automorphisms454if number_of_automorphisms != None:455self.set_number_of_automorphisms(number_of_automorphisms)456#self.__number_of_automorphisms = number_of_automorphisms457#self.__external_initialization_list.append('number_of_automorphisms')458459## Set the determinant460if determinant != None:461self.__det = determinant462self._external_initialization_list.append('determinant')463464465def list_external_initializations(self):466"""467Returns a list of the fields which were set externally at468creation, and not created through the usual QuadraticForm469methods. These fields are as good as the external process470that made them, and are thus not guaranteed to be correct.471472EXAMPLES::473474sage: Q = QuadraticForm(ZZ, 2, [1,0,5])475sage: Q.list_external_initializations()476[]477sage: T = Q.theta_series()478sage: Q.list_external_initializations()479[]480sage: Q = QuadraticForm(ZZ, 2, [1,0,5], unsafe_initialization=False, number_of_automorphisms=3, determinant=0)481sage: Q.list_external_initializations()482[]483484::485486sage: Q = QuadraticForm(ZZ, 2, [1,0,5], unsafe_initialization=False, number_of_automorphisms=3, determinant=0)487sage: Q.list_external_initializations()488[]489sage: Q = QuadraticForm(ZZ, 2, [1,0,5], unsafe_initialization=True, number_of_automorphisms=3, determinant=0)490sage: Q.list_external_initializations()491['number_of_automorphisms', 'determinant']492"""493return deepcopy(self._external_initialization_list)494495496def _pari_(self):497"""498Return a PARI-formatted Hessian matrix for Q.499500EXAMPLES::501502sage: Q = QuadraticForm(ZZ, 2, [1,0,5])503sage: Q._pari_()504[2, 0; 0, 10]505506"""507return self.matrix()._pari_()508509def _pari_init_(self):510"""511Return a PARI-formatted Hessian matrix for Q, as string.512513EXAMPLES::514515sage: Q = QuadraticForm(ZZ, 2, [1,0,5])516sage: Q._pari_init_()517'Mat([2,0;0,10])'518519"""520return self.matrix()._pari_init_()521522523def _repr_(self):524"""525Give a text representation for the quadratic form given as an upper-triangular matrix of coefficients.526527EXAMPLES::528529sage: QuadraticForm(ZZ, 2, [1,3,5])530Quadratic form in 2 variables over Integer Ring with coefficients:531[ 1 3 ]532[ * 5 ]533534"""535n = self.dim()536out_str = "Quadratic form in " + str(n) + " variables over " + str(self.base_ring()) + " with coefficients: \n"537for i in range(n):538if i > 0:539out_str += '\n'540out_str += "[ "541for j in range(n):542if (i > j):543out_str += "* "544else:545out_str += str(self[i,j]) + " "546out_str += "]"547return out_str548549550def _latex_(self):551"""552Give a LaTeX representation for the quadratic form given as an upper-triangular matrix of coefficients.553554EXAMPLES::555556sage: Q = QuadraticForm(ZZ, 2, [2,3,5])557sage: Q._latex_()558'Quadratic form in 2 variables over Integer Ring with coefficients: \\newline\\left[ \\begin{array}{cc}2 & 3 & * & 5 & \\end{array} \\right]'559560"""561n = self.dim()562out_str = ""563out_str += "Quadratic form in " + str(n) + " variables over " + str(self.base_ring())564out_str += " with coefficients: \\newline"565out_str += "\\left[ \\begin{array}{" + n * "c" + "}"566for i in range(n):567for j in range(n):568if (i > j):569out_str += " * & "570else:571out_str += str(self[i,j]) + " & "572# if i < (n-1):573# out_str += "\\"574out_str += "\\end{array} \\right]"575return out_str576577578579def __getitem__(self, ij):580"""581Return the coefficient `a_{ij}` of `x_i * x_j`.582583EXAMPLES::584585sage: Q = QuadraticForm(ZZ, 3, [1,2,3,4,5,6])586sage: matrix(ZZ, 3, 3, [Q[i,j] for i in range(3) for j in range(3)])587[1 2 3]588[2 4 5]589[3 5 6]590591"""592## Unpack the list of indices593i, j = ij594i = int(i)595j = int(j)596597## Ensure we're using upper-triangular coordinates598if i > j:599tmp = i600i = j601j = tmp602603return self.__coeffs[i*self.__n - i*(i-1)/2 + j - i]604605606def __setitem__(self, ij, coeff):607"""608Set the coefficient `a_{ij}` in front of `x_i * x_j`.609610EXAMPLES::611612sage: Q = QuadraticForm(ZZ, 3, [1,2,3,4,5,6])613sage: Q614Quadratic form in 3 variables over Integer Ring with coefficients:615[ 1 2 3 ]616[ * 4 5 ]617[ * * 6 ]618sage: Q[2,1] = 17619sage: Q620Quadratic form in 3 variables over Integer Ring with coefficients:621[ 1 2 3 ]622[ * 4 17 ]623[ * * 6 ]624625"""626## Unpack the list of indices627i, j = ij628i = int(i)629j = int(j)630631## TO DO: Verify that 0 <= i, j <= (n-1)632633## Ensure we're using upper-triangular coordinates634if i > j:635tmp = i636i = j637j = tmp638639## Set the entry640try:641self.__coeffs[i*self.__n - i*(i-1)/2 + j -i] = self.__base_ring(coeff)642except Exception:643raise RuntimeError, "Oops! This coefficient can't be coerced to an element of the base ring for the quadratic form."644645646######################################647# TO DO: def __cmp__(self, other):648######################################649650def __eq__(self, right):651"""652Determines if two quadratic forms are equal.653654EXAMPLES::655656sage: Q = QuadraticForm(ZZ, 2, [1,4,10])657sage: Q == Q658True659660sage: Q1 = QuadraticForm(QQ, 2, [1,4,10])661sage: Q == Q1662False663664sage: Q2 = QuadraticForm(ZZ, 2, [1,4,-10])665sage: Q == Q1666False667sage: Q == Q2668False669sage: Q1 == Q2670False671672"""673if not isinstance(right, QuadraticForm):674return False675return (self.__base_ring == right.__base_ring) and (self.__coeffs == right.__coeffs)676677678def __add__(self, right):679"""680Returns the direct sum of two quadratic forms.681682EXAMPLES::683sage: Q = QuadraticForm(ZZ, 2, [1,4,10])684sage: Q685Quadratic form in 2 variables over Integer Ring with coefficients:686[ 1 4 ]687[ * 10 ]688sage: Q2 = QuadraticForm(ZZ, 2, [1,4,-10])689sage: Q + Q2690Quadratic form in 4 variables over Integer Ring with coefficients:691[ 1 4 0 0 ]692[ * 10 0 0 ]693[ * * 1 4 ]694[ * * * -10 ]695696"""697if not isinstance(right, QuadraticForm):698raise TypeError, "Oops! Can't add these objects since they're not both quadratic forms. =("699elif (self.base_ring() != right.base_ring()):700raise TypeError, "Oops! Can't add these since the quadratic forms don't have the same base rings... =("701else:702Q = QuadraticForm(self.base_ring(), self.dim() + right.dim())703n = self.dim()704m = right.dim()705706for i in range(n):707for j in range(i,n):708Q[i,j] = self[i,j]709710for i in range(m):711for j in range(i,m):712Q[n+i,n+j] = right[i,j]713714return Q715716717def sum_by_coefficients_with(self, right):718"""719Returns the sum (on coefficients) of two quadratic forms of the same size.720721EXAMPLES::722723sage: Q = QuadraticForm(ZZ, 2, [1,4,10])724sage: Q725Quadratic form in 2 variables over Integer Ring with coefficients:726[ 1 4 ]727[ * 10 ]728sage: Q+Q729Quadratic form in 4 variables over Integer Ring with coefficients:730[ 1 4 0 0 ]731[ * 10 0 0 ]732[ * * 1 4 ]733[ * * * 10 ]734735sage: Q2 = QuadraticForm(ZZ, 2, [1,4,-10])736sage: Q.sum_by_coefficients_with(Q2)737Quadratic form in 2 variables over Integer Ring with coefficients:738[ 2 8 ]739[ * 0 ]740741"""742if not isinstance(right, QuadraticForm):743raise TypeError, "Oops! Can't add these objects since they're not both quadratic forms. =("744elif (self.__n != right.__n):745raise TypeError, "Oops! Can't add these since the quadratic forms don't have the same sizes... =("746elif (self.__base_ring != right.__base_ring):747raise TypeError, "Oops! Can't add these since the quadratic forms don't have the same base rings... =("748else:749return QuadraticForm(self.__base_ring, self.__n, [self.__coeffs[i] + right.__coeffs[i] for i in range(len(self.__coeffs))])750751752## ======================== CHANGE THIS TO A TENSOR PRODUCT?!? Even in Characteristic 2?!? =======================753# def __mul__(self, right):754# """755# Multiply (on the right) the quadratic form Q by an element of the ring that Q is defined over.756#757# EXAMPLES::758# sage: Q = QuadraticForm(ZZ, 2, [1,4,10])759# sage: Q*2760# Quadratic form in 2 variables over Integer Ring with coefficients:761# [ 2 8 ]762# [ * 20 ]763#764# sage: Q+Q == Q*2765# True766#767# """768# try:769# c = self.base_ring()(right)770# except Exception:771# raise TypeError, "Oh no! The multiplier cannot be coerced into the base ring of the quadratic form. =("772#773# return QuadraticForm(self.base_ring(), self.dim(), [c * self.__coeffs[i] for i in range(len(self.__coeffs))])774# =========================================================================================================================775776777778def __call__(self, v):779"""780Evaluate this quadratic form Q on a vector or matrix of elements781coercible to the base ring of the quadratic form. If a vector782is given then the output will be the ring element Q(`v`), but if a783matrix is given then the output will be the quadratic form Q'784which in matrix notation is given by:785786.. math::787Q' = v^t * Q * v.788789790EXAMPLES::791792## Evaluate a quadratic form at a vector:793## --------------------------------------794sage: Q = QuadraticForm(QQ, 3, range(6))795sage: Q796Quadratic form in 3 variables over Rational Field with coefficients:797[ 0 1 2 ]798[ * 3 4 ]799[ * * 5 ]800sage: Q([1,2,3])80189802sage: Q([1,0,0])8030804sage: Q([1,1,1])80515806807::808809## Evaluate a quadratic form using a column matrix:810## ------------------------------------------------811sage: Q = QuadraticForm(QQ, 2, range(1,4))812sage: A = Matrix(ZZ,2,2,[-1,0,0,1])813sage: Q(A)814Quadratic form in 2 variables over Rational Field with coefficients:815[ 1 -2 ]816[ * 3 ]817sage: Q([1,0])8181819sage: type(Q([1,0]))820<type 'sage.rings.rational.Rational'>821sage: Q = QuadraticForm(QQ, 2, range(1,4))822sage: Q(matrix(2, [1,0]))823Quadratic form in 1 variables over Rational Field with coefficients:824[ 1 ]825826::827828## Simple 2x2 change of variables:829## -------------------------------830sage: Q = QuadraticForm(ZZ, 2, [1,0,1])831sage: Q832Quadratic form in 2 variables over Integer Ring with coefficients:833[ 1 0 ]834[ * 1 ]835sage: M = Matrix(ZZ, 2, 2, [1,1,0,1])836sage: M837[1 1]838[0 1]839sage: Q(M)840Quadratic form in 2 variables over Integer Ring with coefficients:841[ 1 2 ]842[ * 2 ]843844::845846## Some more tests:847## ----------------848sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])849sage: Q([1,2,3])85014851sage: v = vector([1,2,3])852sage: Q(v)85314854sage: t = tuple([1,2,3])855sage: Q(v)85614857sage: M = Matrix(ZZ, 3, [1,2,3])858sage: Q(M)859Quadratic form in 1 variables over Integer Ring with coefficients:860[ 14 ]861862"""863## If we are passed a matrix A, return the quadratic form Q(A(x))864## (In matrix notation: A^t * Q * A)865n = self.dim()866867if is_Matrix(v):868## Check that v has the correct number of rows869if v.nrows() != n:870raise TypeError, "Oops! The matrix must have " + str(n) + " rows. =("871872## Create the new quadratic form873m = v.ncols()874Q2 = QuadraticForm(self.base_ring(), m)875return QFEvaluateMatrix(self, v, Q2)876877elif (is_Vector(v) or isinstance(v, (list, tuple))):878## Check the vector/tuple/list has the correct length879if not (len(v) == n):880raise TypeError, "Oops! Your vector needs to have length " + str(n) + " ."881882## TO DO: Check that the elements can be coerced into the base ring of Q -- on first elt.883if len(v) > 0:884try:885x = self.base_ring()(v[0])886except Exception:887raise TypeError, "Oops! Your vector is not coercible to the base ring of the quadratic form... =("888889## Attempt to evaluate Q[v]890return QFEvaluateVector(self, v)891892else:893raise(TypeError, "Oops! Presently we can only evaluate a quadratic form on a list, tuple, vector or matrix.")894895896897898## =====================================================================================================899900def _is_even_symmetric_matrix_(self, A, R=None):901"""902Tests if a matrix is symmetric, defined over R, and has even diagonal in R.903904INPUT:905A -- matrix906907R -- ring908909EXAMPLES::910911sage: Q = QuadraticForm(ZZ, 2, [2,3,5])912sage: A = Q.matrix()913sage: A914[ 4 3]915[ 3 10]916sage: Q._is_even_symmetric_matrix_(A)917True918sage: A[0,0] = 1919sage: Q._is_even_symmetric_matrix_(A)920False921922"""923if not is_Matrix(A):924raise TypeError, "A is not a matrix."925926ring_coerce_test = True927if R == None: ## This allows us to omit the ring from the variables, and take it from the matrix928R = A.base_ring()929ring_coerce_test = False930931if not isinstance(R, Ring):932raise TypeError, "R is not a ring."933934if not A.is_square():935return False936937## Test that the matrix is symmetric938n = A.nrows()939for i in range(n):940for j in range(i+1, n):941if A[i,j] != A[j,i]:942return False943944## Test that all entries coerce to R945if not ((A.base_ring() == R) or (ring_coerce_test == True)):946try:947for i in range(n):948for j in range(i, n):949x = R(A[i,j])950except Exception:951return False952953## Test that the diagonal is even (if 1/2 isn't in R)954if not R(2).is_unit():955for i in range(n):956if not is_even(R(A[i,i])):957return False958959return True960961962## =====================================================================================================963964def matrix(self):965"""966Returns the Hessian matrix A for which Q(X) = `(1/2) * X^t * A * X`.967968EXAMPLES::969970sage: Q = QuadraticForm(ZZ, 3, range(6))971sage: Q.matrix()972[ 0 1 2]973[ 1 6 4]974[ 2 4 10]975976"""977return self.Hessian_matrix()978979980def Hessian_matrix(self):981"""982Returns the Hessian matrix A for which Q(X) = `(1/2) * X^t * A * X`.983984EXAMPLES::985986sage: Q = QuadraticForm(QQ, 2, range(1,4))987sage: Q988Quadratic form in 2 variables over Rational Field with coefficients:989[ 1 2 ]990[ * 3 ]991sage: Q.Hessian_matrix()992[2 2]993[2 6]994sage: Q.matrix().base_ring()995Rational Field996997"""998mat_entries = []999for i in range(self.dim()):1000for j in range(self.dim()):1001if (i == j):1002mat_entries += [ 2 * self[i,j] ]1003else:1004mat_entries += [ self[i,j] ]10051006return matrix(self.base_ring(), self.dim(), self.dim(), mat_entries)100710081009def Gram_matrix_rational(self):1010"""1011Returns a (symmetric) Gram matrix A for the quadratic form Q,1012meaning that10131014.. math::10151016Q(x) = x^t * A * x,10171018defined over the fraction field of the base ring.10191020EXAMPLES::10211022sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])1023sage: A = Q.Gram_matrix_rational(); A1024[1 0 0 0]1025[0 3 0 0]1026[0 0 5 0]1027[0 0 0 7]1028sage: A.base_ring()1029Rational Field10301031"""1032return (ZZ(1) / ZZ(2)) * self.matrix()103310341035def Gram_matrix(self):1036"""1037Returns a (symmetric) Gram matrix A for the quadratic form Q,1038meaning that10391040.. math::10411042Q(x) = x^t * A * x,10431044defined over the base ring of Q. If this is not possible,1045then a TypeError is raised.10461047EXAMPLES::10481049sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])1050sage: A = Q.Gram_matrix(); A1051[1 0 0 0]1052[0 3 0 0]1053[0 0 5 0]1054[0 0 0 7]1055sage: A.base_ring()1056Integer Ring10571058"""1059A = (ZZ(1) / ZZ(2)) * self.matrix()1060n = self.dim()10611062## Test to see if it has an integral Gram matrix1063Int_flag = True1064for i in range(n):1065for j in range(i,n):1066Int_flag = Int_flag and A[i,j] in self.base_ring()10671068## Return the Gram matrix, or an error1069if Int_flag:1070return MatrixSpace(self.base_ring(), n, n)(A)1071else:1072raise TypeError, "Oops! This form does not have an integral Gram matrix. =("107310741075def has_integral_Gram_matrix(self):1076"""1077Returns whether the quadratic form has an integral Gram matrix (with respect to its base ring).10781079A warning is issued if the form is defined over a field, since in that case the return is trivially true.10801081EXAMPLES::10821083sage: Q = QuadraticForm(ZZ, 2, [7,8,9])1084sage: Q.has_integral_Gram_matrix()1085True10861087::10881089sage: Q = QuadraticForm(ZZ, 2, [4,5,6])1090sage: Q.has_integral_Gram_matrix()1091False10921093"""1094## Warning over fields1095if is_field(self.base_ring()):1096warn("Warning -- A quadratic form over a field always has integral Gram matrix. Do you really want to do this?!?")10971098## Determine integrality of the Gram matrix1099flag = True1100try:1101self.Gram_matrix()1102except Exception:1103flag = False11041105return flag110611071108def gcd(self):1109"""1110Returns the greatest common divisor of the coefficients of the1111quadratic form (as a polynomial).11121113EXAMPLES::11141115sage: Q = QuadraticForm(ZZ, 4, range(1, 21, 2))1116sage: Q.gcd()1117111181119::11201121sage: Q = QuadraticForm(ZZ, 4, range(0, 20, 2))1122sage: Q.gcd()112321124"""1125if self.base_ring() != ZZ:1126raise TypeError, "Oops! The given quadratic form must be defined over ZZ."11271128return GCD(self.coefficients())112911301131def polynomial(self,names='x'):1132r"""1133Returns the polynomial in 'n' variables of the quadratic form in the ring 'R[names].'11341135INPUT:11361137-'self' - a quadratic form over a commatitive ring.1138-'names' - the name of the variables. Digits will be appended to the name for each different canonical1139variable e.g x1, x2, x3 etc.11401141OUTPUT:11421143The polynomial form of the quadratic form.11441145EXAMPLES::11461147sage: Q = DiagonalQuadraticForm(QQ,[1, 3, 5, 7])1148sage: P = Q.polynomial(); P11492*x0^2 + 6*x1^2 + 10*x2^2 + 14*x3^211501151::11521153sage: F.<a> = NumberField(x^2 - 5)1154sage: Z = F.ring_of_integers()1155sage: Q = QuadraticForm(Z,3,[2*a, 3*a, 0 , 1 - a, 0, 2*a + 4])1156sage: P = Q.polynomial(names='y'); P11574*a*y0^2 + 6*a*y0*y1 + (-2*a + 2)*y1^2 + (4*a + 8)*y2^21158sage: Q = QuadraticForm(F,4,[a, 3*a, 0, 1 - a, a - 3, 0, 2*a + 4, 4 + a, 0, 1])1159sage: Q.polynomial(names='z')1160(2*a)*z0^2 + (6*a)*z0*z1 + (2*a - 6)*z1^2 + (2*a + 8)*z2^2 + (-2*a + 2)*z0*z3 + (4*a + 8)*z1*z3 + 2*z3^21161sage: B.<i,j,k> = QuaternionAlgebra(F,-1,-1)1162sage: Q = QuadraticForm(B, 3, [2*a, 3*a, i, 1 - a, 0, 2*a + 4])1163sage: Q.polynomial()1164Traceback (most recent call last):1165...1166ValueError: Can only create polynomial rings over commutative rings.11671168"""1169M = self.matrix()1170n = self.dim()1171B = self.base_ring()1172try:1173R = PolynomialRing(self.base_ring(),names,n)1174except Exception:1175raise ValueError, 'Can only create polynomial rings over commutative rings.'1176V = vector(R.gens())1177P = (V*M).dot_product(V)1178return P11791180118111821183def is_primitive(self):1184"""1185Determines if the given integer-valued form is primitive1186(i.e. not an integer (>1) multiple of another integer-valued1187quadratic form).11881189EXAMPLES::11901191sage: Q = QuadraticForm(ZZ, 2, [2,3,4])1192sage: Q.is_primitive()1193True1194sage: Q = QuadraticForm(ZZ, 2, [2,4,8])1195sage: Q.is_primitive()1196False11971198"""1199return (self.gcd() == 1)120012011202def primitive(self):1203"""1204Returns a primitive version of an integer-valued quadratic form, defined over `ZZ`.12051206EXAMPLES::12071208sage: Q = QuadraticForm(ZZ, 2, [2,3,4])1209sage: Q.primitive()1210Quadratic form in 2 variables over Integer Ring with coefficients:1211[ 2 3 ]1212[ * 4 ]1213sage: Q = QuadraticForm(ZZ, 2, [2,4,8])1214sage: Q.primitive()1215Quadratic form in 2 variables over Integer Ring with coefficients:1216[ 1 2 ]1217[ * 4 ]12181219"""1220if self.base_ring() != ZZ:1221raise TypeError, "Oops! The given quadratic form must be defined over ZZ."12221223g = self.gcd()1224return QuadraticForm(self.base_ring(), self.dim(), [ZZ(x/g) for x in self.coefficients()])1225122612271228def adjoint_primitive(self):1229"""1230Returns the primitive adjoint of the quadratic form, which is1231the smallest discriminant integer-valued quadratic form whose1232matrix is a scalar multiple of the inverse of the matrix of1233the given quadratic form.12341235EXAMPLES::12361237sage: Q = QuadraticForm(ZZ, 2, [1,2,3])1238sage: Q.adjoint_primitive()1239Quadratic form in 2 variables over Integer Ring with coefficients:1240[ 3 -2 ]1241[ * 1 ]12421243"""1244return QuadraticForm(self.Hessian_matrix().adjoint()).primitive()1245124612471248def dim(self):1249"""1250Gives the number of variables of the quadratic form.12511252EXAMPLES::12531254sage: Q = QuadraticForm(ZZ, 2, [1,2,3])1255sage: Q.dim()1256212571258"""1259return self.__n126012611262def base_ring(self):1263"""1264Gives the ring over which the quadratic form is defined.12651266EXAMPLES::12671268sage: Q = QuadraticForm(ZZ, 2, [1,2,3])1269sage: Q.base_ring()1270Integer Ring12711272"""1273return self.__base_ring127412751276def coefficients(self):1277"""1278Gives the matrix of upper triangular coefficients,1279by reading across the rows from the main diagonal.12801281EXAMPLES::12821283sage: Q = QuadraticForm(ZZ, 2, [1,2,3])1284sage: Q.coefficients()1285[1, 2, 3]12861287"""1288return self.__coeffs128912901291def det(self):1292"""1293Gives the determinant of the Gram matrix of 2*Q, or1294equivalently the determinant of the Hessian matrix of Q.12951296(Note: This is always defined over the same ring as the1297quadratic form.)12981299EXAMPLES::13001301sage: Q = QuadraticForm(ZZ, 2, [1,2,3])1302sage: Q.det()1303813041305"""1306try:1307return self.__det1308except AttributeError:1309## Compute the determinant1310if self.dim() == 0:1311new_det = self.base_ring()(1)1312else:1313new_det = self.matrix().det()13141315## Cache and return the determinant1316self.__det = new_det1317return new_det131813191320def Gram_det(self):1321"""1322Gives the determinant of the Gram matrix of Q.13231324(Note: This is defined over the fraction field of the ring of1325the quadratic form, but is often not defined over the same1326ring as the quadratic form.)13271328EXAMPLES::13291330sage: Q = QuadraticForm(ZZ, 2, [1,2,3])1331sage: Q.Gram_det()1332213331334"""1335return self.det() / ZZ(2**self.dim())133613371338def base_change_to(self, R):1339"""1340Alters the quadratic form to have all coefficients1341defined over the new base_ring R. Here R must be1342coercible to from the current base ring.13431344Note: This is preferable to performing an explicit1345coercion through the base_ring() method, which does1346not affect the individual coefficients. This is1347particularly useful for performing fast modular1348arithmetic evaluations.13491350INPUT:1351R -- a ring13521353OUTPUT:1354quadratic form13551356EXAMPLES::13571358sage: Q = DiagonalQuadraticForm(ZZ,[1,1]); Q1359Quadratic form in 2 variables over Integer Ring with coefficients:1360[ 1 0 ]1361[ * 1 ]13621363::13641365sage: Q1 = Q.base_change_to(IntegerModRing(5)); Q11366Quadratic form in 2 variables over Ring of integers modulo 5 with coefficients:1367[ 1 0 ]1368[ * 1 ]13691370sage: Q1([35,11])1371113721373"""1374## Check that a canonical coercion is possible1375if not is_Ring(R):1376raise TypeError, "Oops! R is not a ring. =("1377if not R.has_coerce_map_from(self.base_ring()):1378raise TypeError, "Oops! There is no canonical coercion from " + str(self.base_ring()) + " to R."13791380## Return the coerced form1381return QuadraticForm(R, self.dim(), [R(x) for x in self.coefficients()])138213831384def level(self):1385r"""1386Determines the level of the quadratic form over a PID, which is a1387generator for the smallest ideal `N` of `R` such that N * (the matrix of13882*Q)^(-1) is in R with diagonal in 2*R.13891390Over `\ZZ` this returns a non-negative number.13911392(Caveat: This always returns the unit ideal when working over a field!)13931394EXAMPLES::13951396sage: Q = QuadraticForm(ZZ, 2, range(1,4))1397sage: Q.level()1398813991400sage: Q1 = QuadraticForm(QQ, 2, range(1,4))1401sage: Q1.level() # random1402UserWarning: Warning -- The level of a quadratic form over a field is always 1. Do you really want to do this?!?1403114041405sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])1406sage: Q.level()140742014081409"""1410## Try to return the cached level1411try:1412return self.__level1413except AttributeError:14141415## Check that the base ring is a PID1416if not is_PrincipalIdealDomain(self.base_ring()):1417raise TypeError, "Oops! The level (as a number) is only defined over a Principal Ideal Domain. Try using level_ideal()."141814191420## Warn the user if the form is defined over a field!1421if self.base_ring().is_field():1422warn("Warning -- The level of a quadratic form over a field is always 1. Do you really want to do this?!?")1423#raise RuntimeError, "Warning -- The level of a quadratic form over a field is always 1. Do you really want to do this?!?"142414251426## Check invertibility and find the inverse1427try:1428mat_inv = self.matrix()**(-1)1429except ZeroDivisionError:1430raise TypeError, "Oops! The quadratic form is degenerate (i.e. det = 0). =("14311432## Compute the level1433inv_denoms = []1434for i in range(self.dim()):1435for j in range(i, self.dim()):1436if (i == j):1437inv_denoms += [denominator(mat_inv[i,j] / 2)]1438else:1439inv_denoms += [denominator(mat_inv[i,j])]1440lvl = LCM(inv_denoms)1441lvl = ideal(self.base_ring()(lvl)).gen()1442##############################################################1443## To do this properly, the level should be the inverse of the1444## fractional ideal (over R) generated by the entries whose1445## denominators we take above. =)1446##############################################################14471448## Normalize the result over ZZ1449if self.base_ring() == IntegerRing():1450lvl = abs(lvl)14511452## Cache and return the level1453self.__level = lvl1454return lvl1455145614571458def level_ideal(self):1459"""1460Determines the level of the quadratic form (over R), which is the1461smallest ideal N of R such that N * (the matrix of 2*Q)^(-1) is1462in R with diagonal in 2*R.1463(Caveat: This always returns the principal ideal when working over a field!)14641465WARNING: THIS ONLY WORKS OVER A PID RING OF INTEGERS FOR NOW!1466(Waiting for Sage fractional ideal support.)14671468EXAMPLES::14691470sage: Q = QuadraticForm(ZZ, 2, range(1,4))1471sage: Q.level_ideal()1472Principal ideal (8) of Integer Ring14731474::14751476sage: Q1 = QuadraticForm(QQ, 2, range(1,4))1477sage: Q1.level_ideal()1478Principal ideal (1) of Rational Field14791480::14811482sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])1483sage: Q.level_ideal()1484Principal ideal (420) of Integer Ring14851486"""1487##############################################################1488## To do this properly, the level should be the inverse of the1489## fractional ideal (over R) generated by the entries whose1490## denominators we take above. =)1491##############################################################14921493return ideal(self.base_ring()(self.level()))14941495def bilinear_map(self,v,w):1496r"""1497Returns the value of the associated bilinear map on two vectors14981499Given a quadratic form `Q` over some base ring `R` with1500characteristic not equal to 2, this gives the image of two1501vectors with coefficients in `R` under the associated bilinear1502map `B`, given by the relation `2 B(v,w) = Q(v) + Q(w) - Q(v+w)`.15031504INPUT:15051506`v, w` -- two vectors15071508OUTPUT:15091510an element of the base ring `R`.15111512EXAMPLES:15131514First, an example over `\ZZ`::15151516sage: Q = QuadraticForm(ZZ,3,[1,4,0,1,4,1])1517sage: v = vector(ZZ,(1,2,0))1518sage: w = vector(ZZ,(0,1,1))1519sage: Q.bilinear_map(v,w)1520815211522This also works over `\QQ`::15231524sage: Q = QuadraticForm(QQ,2,[1/2,2,1])1525sage: v = vector(QQ,(1,1))1526sage: w = vector(QQ,(1/2,2))1527sage: Q.bilinear_map(v,w)152819/415291530The vectors must have the correct length::15311532sage: Q = DiagonalQuadraticForm(ZZ,[1,7,7])1533sage: v = vector((1,2))1534sage: w = vector((1,1,1))1535sage: Q.bilinear_map(v,w)1536Traceback (most recent call last):1537...1538TypeError: vectors must have length 315391540This does not work if the characteristic is 2::15411542sage: Q = DiagonalQuadraticForm(GF(2),[1,1,1])1543sage: v = vector((1,1,1))1544sage: w = vector((1,1,1))1545sage: Q.bilinear_map(v,w)1546Traceback (most recent call last):1547...1548TypeError: not defined for rings of characteristic 21549"""1550if len(v) != self.dim() or len(w) != self.dim():1551raise TypeError("vectors must have length " + str(self.dim()))1552if self.base_ring().characteristic() == 2:1553raise TypeError("not defined for rings of characteristic 2")1554return (self(v+w) - self(v) - self(w))/2155515561557## =====================================================================================================15581559def DiagonalQuadraticForm(R, diag):1560"""1561Returns a quadratic form over `R` which is a sum of squares.15621563INPUT:15641565- `R` -- ring1566- ``diag`` -- list/tuple of elements coercible to R15671568OUTPUT:15691570quadratic form15711572EXAMPLES::15731574sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])1575sage: Q1576Quadratic form in 4 variables over Integer Ring with coefficients:1577[ 1 0 0 0 ]1578[ * 3 0 0 ]1579[ * * 5 0 ]1580[ * * * 7 ]1581"""1582Q = QuadraticForm(R, len(diag))1583for i in range(len(diag)):1584Q[i,i] = diag[i]1585return Q158615871588