Quadratic Forms Overview

AUTHORS:

sage.quadratic_forms.quadratic_form.DiagonalQuadraticForm(R, diag)

Return a quadratic form over \(R\) which is a sum of squares.

INPUT:

  • \(R\) – ring

  • diag – list/tuple of elements coercible to R

OUTPUT:

quadratic form

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 0 0 0 ]
[ * 3 0 0 ]
[ * * 5 0 ]
[ * * * 7 ]
class sage.quadratic_forms.quadratic_form.QuadraticForm(R, n=None, entries=None, unsafe_initialization=False, number_of_automorphisms=None, determinant=None)

Bases: sage.structure.sage_object.SageObject

The QuadraticForm class represents a quadratic form in n variables with coefficients in the ring R.

INPUT:

The constructor may be called in any of the following ways.

  1. QuadraticForm(R, n, entries), where

    • \(R\) – ring for which the quadratic form is defined

    • \(n\) – an integer >= 0

    • entries – a list of \(n(n+1)/2\) coefficients of the quadratic form in \(R\) (given lexicographically, or equivalently, by rows of the matrix)

  2. QuadraticForm(R, n), where

    • \(R\) – a ring

    • \(n\) – a symmetric \(n \times n\) matrix with even diagonal (relative to \(R\))

  3. QuadraticForm(R), where

    • \(R\) – a symmetric \(n \times n\) matrix with even diagonal (relative to its base ring)

If the keyword argument unsafe_initialize is True, then the subsequent fields may by used to force the external initialization of various fields of the quadratic form. Currently the only fields which can be set are:

  • number_of_automorphisms

  • determinant

OUTPUT:

quadratic form

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, [1,2,3,4,5,6])
sage: Q
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 1 2 3 ]
[ * 4 5 ]
[ * * 6 ]
sage: Q = QuadraticForm(QQ, 3, [1,2,3,4/3 ,5,6])
sage: Q
Quadratic form in 3 variables over Rational Field with coefficients:
[ 1 2 3 ]
[ * 4/3 5 ]
[ * * 6 ]
sage: Q[0,0]
1
sage: Q[0,0].parent()
Rational Field
sage: Q = QuadraticForm(QQ, 7, range(28))
sage: Q
Quadratic form in 7 variables over Rational Field with coefficients:
[ 0 1 2 3 4 5 6 ]
[ * 7 8 9 10 11 12 ]
[ * * 13 14 15 16 17 ]
[ * * * 18 19 20 21 ]
[ * * * * 22 23 24 ]
[ * * * * * 25 26 ]
[ * * * * * * 27 ]
sage: Q = QuadraticForm(QQ, 2, range(1,4))
sage: A = Matrix(ZZ,2,2,[-1,0,0,1])
sage: Q(A)
Quadratic form in 2 variables over Rational Field with coefficients:
[ 1 -2 ]
[ * 3 ]
sage: m = matrix(2,2,[1,2,3,4])
sage: m + m.transpose()
[2 5]
[5 8]
sage: QuadraticForm(m + m.transpose())
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 1 5 ]
[ * 4 ]
sage: QuadraticForm(ZZ, m + m.transpose())
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 1 5 ]
[ * 4 ]
sage: QuadraticForm(QQ, m + m.transpose())
Quadratic form in 2 variables over Rational Field with coefficients:
[ 1 5 ]
[ * 4 ]
CS_genus_symbol_list(force_recomputation=False)

Returns the list of Conway-Sloane genus symbols in increasing order of primes dividing 2*det.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3,4])
sage: Q.CS_genus_symbol_list()
[Genus symbol at 2:    [2^-2 4^1 8^1]_6, Genus symbol at 3:     1^3 3^-1]
GHY_mass__maximal()

Use the GHY formula to compute the mass of a (maximal?) quadratic lattice. This works for any number field.

Reference: See [GHY, Prop 7.4 and 7.5, p121] and [GY, Thrm 10.20, p25].

INPUT:

none

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.GHY_mass__maximal()
Gram_det()

Gives the determinant of the Gram matrix of Q.

(Note: This is defined over the fraction field of the ring of the quadratic form, but is often not defined over the same ring as the quadratic form.)

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,2,3])
sage: Q.Gram_det()
2
Gram_matrix()

Return a (symmetric) Gram matrix A for the quadratic form Q, meaning that

\[Q(x) = x^t * A * x,\]

defined over the base ring of Q. If this is not possible, then a TypeError is raised.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: A = Q.Gram_matrix(); A
[1 0 0 0]
[0 3 0 0]
[0 0 5 0]
[0 0 0 7]
sage: A.base_ring()
Integer Ring
Gram_matrix_rational()

Return a (symmetric) Gram matrix A for the quadratic form Q, meaning that

\[Q(x) = x^t * A * x,\]

defined over the fraction field of the base ring.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: A = Q.Gram_matrix_rational(); A
[1 0 0 0]
[0 3 0 0]
[0 0 5 0]
[0 0 0 7]
sage: A.base_ring()
Rational Field
Hessian_matrix()

Return the Hessian matrix A for which Q(X) = \((1/2) * X^t * A * X\).

EXAMPLES:

sage: Q = QuadraticForm(QQ, 2, range(1,4))
sage: Q
Quadratic form in 2 variables over Rational Field with coefficients:
[ 1 2 ]
[ * 3 ]
sage: Q.Hessian_matrix()
[2 2]
[2 6]
sage: Q.matrix().base_ring()
Rational Field
Kitaoka_mass_at_2()

Returns the local mass of the quadratic form when \(p=2\), according to Theorem 5.6.3 on pp108–9 of Kitaoka’s Book “The Arithmetic of Quadratic Forms”.

INPUT:

none

OUTPUT:

a rational number > 0

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.Kitaoka_mass_at_2()   ## WARNING:  WE NEED TO CHECK THIS CAREFULLY!
1/2
Pall_mass_density_at_odd_prime(p)

Returns the local representation density of a form (for representing itself) defined over \(ZZ\), at some prime \(p>2\).

REFERENCES:

Pall’s article “The Weight of a Genus of Positive n-ary Quadratic Forms” appearing in Proc. Symp. Pure Math. VIII (1965), pp95–105.

INPUT:

\(p\) – a prime number > 2.

OUTPUT:

a rational number.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, [1,0,0,1,0,1])
sage: Q.Pall_mass_density_at_odd_prime(3)
[(0, Quadratic form in 3 variables over Integer Ring with coefficients:
[ 1 0 0 ]
[ * 1 0 ]
[ * * 1 ])] [(0, 3, 8)] [8/9] 8/9
8/9
Watson_mass_at_2()

Returns the local mass of the quadratic form when \(p=2\), according to Watson’s Theorem 1 of “The 2-adic density of a quadratic form” in Mathematika 23 (1976), pp 94–106.

INPUT:

none

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.Watson_mass_at_2()               ## WARNING:  WE NEED TO CHECK THIS CAREFULLY!
384
add_symmetric(c, i, j, in_place=False)

Performs the substitution \(x_j --> x_j + c*x_i\), which has the effect (on associated matrices) of symmetrically adding \(c * j\)-th row/column to the \(i\)-th row/column.

NOTE: This is meant for compatibility with previous code, which implemented a matrix model for this class. It is used in the local_normal_form() method.

INPUT:

\(c\) – an element of Q.base_ring()

\(i\), \(j\) – integers >= 0

OUTPUT:

a QuadraticForm (by default, otherwise none)

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, range(1,7)); Q
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 1 2 3 ]
[ * 4 5 ]
[ * * 6 ]
sage: Q.add_symmetric(-1, 1, 0)
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 1 0 3 ]
[ * 3 2 ]
[ * * 6 ]
sage: Q.add_symmetric(-3/2, 2, 0)     ## ERROR: -3/2 isn't in the base ring ZZ
Traceback (most recent call last):
...
RuntimeError: Oops!  This coefficient can...t be coerced to an element of the base ring for the quadratic form.
sage: Q = QuadraticForm(QQ, 3, range(1,7)); Q
Quadratic form in 3 variables over Rational Field with coefficients:
[ 1 2 3 ]
[ * 4 5 ]
[ * * 6 ]
sage: Q.add_symmetric(-3/2, 2, 0)
Quadratic form in 3 variables over Rational Field with coefficients:
[ 1 2 0 ]
[ * 4 2 ]
[ * * 15/4 ]
adjoint()

This gives the adjoint (integral) quadratic form associated to the given form, essentially defined by taking the adjoint of the matrix.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,2,5])
sage: Q.adjoint()
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 5 -2 ]
[ * 1 ]
sage: Q = QuadraticForm(ZZ, 3, [1, 0, -1, 2, -1, 5])
sage: Q.adjoint()
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 39 2 8 ]
[ * 19 4 ]
[ * * 8 ]
adjoint_primitive()

Return the primitive adjoint of the quadratic form, which is the smallest discriminant integer-valued quadratic form whose matrix is a scalar multiple of the inverse of the matrix of the given quadratic form.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,2,3])
sage: Q.adjoint_primitive()
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 3 -2 ]
[ *  1 ]
anisotropic_primes()

Return a list with all of the anisotropic primes of the quadratic form.

The infinite place is denoted by \(-1\).

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.anisotropic_primes()
[2, -1]

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.anisotropic_primes()
[2, -1]

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1,1])
sage: Q.anisotropic_primes()
[-1]
antiadjoint()

This gives an (integral) form such that its adjoint is the given form.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, [1, 0, -1, 2, -1, 5])
sage: Q.adjoint().antiadjoint()
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 1 0 -1 ]
[ * 2 -1 ]
[ * * 5 ]
sage: Q.antiadjoint()
Traceback (most recent call last):
...
ValueError: not an adjoint
automorphism_group()

Return the group of automorphisms of the quadratic form.

OUTPUT: a MatrixGroup

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.automorphism_group()
Matrix group over Rational Field with 3 generators (
[-1  0  0]  [0 0 1]  [ 0  0  1]
[ 0 -1  0]  [0 1 0]  [-1  0  0]
[ 0  0 -1], [1 0 0], [ 0  1  0]
)
sage: DiagonalQuadraticForm(ZZ, [1,3,5,7]).automorphism_group()
Matrix group over Rational Field with 4 generators (
[-1  0  0  0]  [ 1  0  0  0]  [ 1  0  0  0]  [ 1  0  0  0]
[ 0 -1  0  0]  [ 0 -1  0  0]  [ 0  1  0  0]  [ 0  1  0  0]
[ 0  0 -1  0]  [ 0  0  1  0]  [ 0  0 -1  0]  [ 0  0  1  0]
[ 0  0  0 -1], [ 0  0  0  1], [ 0  0  0  1], [ 0  0  0 -1]
)

The smallest possible automorphism group has order two, since we can always change all signs:

sage: Q = QuadraticForm(ZZ, 3, [2, 1, 2, 2, 1, 3])
sage: Q.automorphism_group()
Matrix group over Rational Field with 1 generators (
[-1  0  0]
[ 0 -1  0]
[ 0  0 -1]
)
automorphisms()

Return the list of the automorphisms of the quadratic form.

OUTPUT: a list of matrices

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.number_of_automorphisms()
48
sage: 2^3 * factorial(3)
48
sage: len(Q.automorphisms())
48
sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.number_of_automorphisms()
16
sage: aut = Q.automorphisms()
sage: len(aut)
16
sage: all(Q(M) == Q for M in aut)
True

sage: Q = QuadraticForm(ZZ, 3, [2, 1, 2, 2, 1, 3])
sage: sorted(Q.automorphisms())
[
[-1  0  0]  [1 0 0]
[ 0 -1  0]  [0 1 0]
[ 0  0 -1], [0 0 1]
]
base_change_to(R)

Alters the quadratic form to have all coefficients defined over the new base_ring R. Here R must be coercible to from the current base ring.

Note: This is preferable to performing an explicit coercion through the base_ring() method, which does not affect the individual coefficients. This is particularly useful for performing fast modular arithmetic evaluations.

INPUT:

R – a ring

OUTPUT:

quadratic form

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ,[1,1]); Q
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 1 0 ]
[ * 1 ]
sage: Q1 = Q.base_change_to(IntegerModRing(5)); Q1
Quadratic form in 2 variables over Ring of integers modulo 5 with coefficients:
[ 1 0 ]
[ * 1 ]

sage: Q1([35,11])
1
base_ring()

Gives the ring over which the quadratic form is defined.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,2,3])
sage: Q.base_ring()
Integer Ring
basiclemma(M)

Finds a number represented by self and coprime to M.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [2, 1, 3])
sage: Q.basiclemma(6)
71
basiclemmavec(M)

Finds a vector where the value of the quadratic form is coprime to M.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [2, 1, 5])
sage: Q.basiclemmavec(10)
(6, 5)
sage: Q(_)
227
basis_of_short_vectors(show_lengths=False)

Return a basis for \(\ZZ^n\) made of vectors with minimal lengths Q(\(v\)).

OUTPUT:

a tuple of vectors, and optionally a tuple of values for each vector.

This uses pari:qfminim.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.basis_of_short_vectors()
((1, 0, 0, 0), (0, 1, 0, 0), (0, 0, 1, 0), (0, 0, 0, 1))
sage: Q.basis_of_short_vectors(True)
(((1, 0, 0, 0), (0, 1, 0, 0), (0, 0, 1, 0), (0, 0, 0, 1)), (1, 3, 5, 7))

The returned vectors are immutable:

sage: v = Q.basis_of_short_vectors()[0]
sage: v
(1, 0, 0, 0)
sage: v[0] = 0
Traceback (most recent call last):
...
ValueError: vector is immutable; please change a copy instead (use copy())
bilinear_map(v, w)

Return the value of the associated bilinear map on two vectors

Given a quadratic form \(Q\) over some base ring \(R\) with characteristic not equal to 2, this gives the image of two vectors with coefficients in \(R\) under the associated bilinear map \(B\), given by the relation \(2 B(v,w) = Q(v) + Q(w) - Q(v+w)\).

INPUT:

\(v, w\) – two vectors

OUTPUT:

an element of the base ring \(R\).

EXAMPLES:

First, an example over \(\ZZ\):

sage: Q = QuadraticForm(ZZ,3,[1,4,0,1,4,1])
sage: v = vector(ZZ,(1,2,0))
sage: w = vector(ZZ,(0,1,1))
sage: Q.bilinear_map(v,w)
8

This also works over \(\QQ\):

sage: Q = QuadraticForm(QQ,2,[1/2,2,1])
sage: v = vector(QQ,(1,1))
sage: w = vector(QQ,(1/2,2))
sage: Q.bilinear_map(v,w)
19/4

The vectors must have the correct length:

sage: Q = DiagonalQuadraticForm(ZZ,[1,7,7])
sage: v = vector((1,2))
sage: w = vector((1,1,1))
sage: Q.bilinear_map(v,w)
Traceback (most recent call last):
...
TypeError: vectors must have length 3

This does not work if the characteristic is 2:

sage: Q = DiagonalQuadraticForm(GF(2),[1,1,1])
sage: v = vector((1,1,1))
sage: w = vector((1,1,1))
sage: Q.bilinear_map(v,w)
Traceback (most recent call last):
...
TypeError: not defined for rings of characteristic 2
cholesky_decomposition(bit_prec=53)

Give the Cholesky decomposition of this quadratic form \(Q\) as a real matrix of precision bit_prec.

RESTRICTIONS:

Q must be given as a QuadraticForm defined over \(\ZZ\), \(\QQ\), or some real field. If it is over some real field, then an error is raised if the precision given is not less than the defined precision of the real field defining the quadratic form!

REFERENCE:

From Cohen’s “A Course in Computational Algebraic Number Theory” book, p 103.

INPUT:

bit_prec – a natural number (default 53).

OUTPUT:

an upper triangular real matrix of precision bit_prec.

TO DO:

If we only care about working over the real double field (RDF), then we can use the cholesky() method present for square matrices over that.

Note

There is a note in the original code reading

##/////////////////////////////////////////////////////////////////////////////////////////////////
##/// Finds the Cholesky decomposition of a quadratic form -- as an upper-triangular matrix!
##/// (It's assumed to be global, hence twice the form it refers to.)  <-- Python revision asks:  Is this true?!? =|
##/////////////////////////////////////////////////////////////////////////////////////////////////

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.cholesky_decomposition()
[ 1.00000000000000 0.000000000000000 0.000000000000000]
[0.000000000000000  1.00000000000000 0.000000000000000]
[0.000000000000000 0.000000000000000  1.00000000000000]
sage: Q = QuadraticForm(QQ, 3, range(1,7)); Q
Quadratic form in 3 variables over Rational Field with coefficients:
[ 1 2 3 ]
[ * 4 5 ]
[ * * 6 ]
sage: Q.cholesky_decomposition()
[ 1.00000000000000  1.00000000000000  1.50000000000000]
[0.000000000000000  3.00000000000000 0.333333333333333]
[0.000000000000000 0.000000000000000  3.41666666666667]
clifford_conductor()

This is the product of all primes where the Clifford invariant is -1

Note: For ternary forms, this is the discriminant of the quaternion algebra associated to the quadratic space (i.e. the even Clifford algebra)

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, [1, 0, -1, 2, -1, 5])
sage: Q.clifford_invariant(2)
1
sage: Q.clifford_invariant(37)
-1
sage: Q.clifford_conductor()
37
sage: DiagonalQuadraticForm(ZZ, [1, 1, 1]).clifford_conductor()
2
sage: QuadraticForm(ZZ, 3, [2, -2, 0, 2, 0, 5]).clifford_conductor()
30

For hyperbolic spaces, the clifford conductor is 1:

sage: H = QuadraticForm(ZZ, 2, [0, 1, 0])
sage: H.clifford_conductor()
1
sage: (H + H).clifford_conductor()
1
sage: (H + H + H).clifford_conductor()
1
sage: (H + H + H + H).clifford_conductor()
1
clifford_invariant(p)

This is the Clifford invariant, i.e. the class in the Brauer group of the Clifford algebra for even dimension, of the even Clifford Algebra for odd dimension.

See Lam (AMS GSM 67) p. 117 for the definition, and p. 119 for the formula relating it to the Hasse invariant.

EXAMPLES:

For hyperbolic spaces, the clifford invariant is +1:

sage: H = QuadraticForm(ZZ, 2, [0, 1, 0])
sage: H.clifford_invariant(2)
1
sage: (H + H).clifford_invariant(2)
1
sage: (H + H + H).clifford_invariant(2)
1
sage: (H + H + H + H).clifford_invariant(2)
1
coefficients()

Gives the matrix of upper triangular coefficients, by reading across the rows from the main diagonal.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,2,3])
sage: Q.coefficients()
[1, 2, 3]
complementary_subform_to_vector(v)

Finds the \((n-1)\)-dim’l quadratic form orthogonal to the vector \(v\).

Note: This is usually not a direct summand!

Technical Notes: There is a minor difference in the cancellation code here (form the C++ version) since the notation Q \([i,j]\) indexes coefficients of the quadratic polynomial here, not the symmetric matrix. Also, it produces a better splitting now, for the full lattice (as opposed to a sublattice in the C++ code) since we now extend \(v\) to a unimodular matrix.

INPUT:

\(v\) – a list of self.dim() integers

OUTPUT:

a QuadraticForm over \(ZZ\)

EXAMPLES:

sage: Q1 = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q1.complementary_subform_to_vector([1,0,0,0])
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 3 0 0 ]
[ * 5 0 ]
[ * * 7 ]
sage: Q1.complementary_subform_to_vector([1,1,0,0])
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 12 0 0 ]
[ * 5 0 ]
[ * * 7 ]
sage: Q1.complementary_subform_to_vector([1,1,1,1])
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 624 -480 -672 ]
[ * 880 -1120 ]
[ * * 1008 ]
compute_definiteness()

Computes whether the given quadratic form is positive-definite, negative-definite, indefinite, degenerate, or the zero form.

This caches one of the following strings in self.__definiteness_string: “pos_def”, “neg_def”, “indef”, “zero”, “degenerate”. It is called from all routines like:

is_positive_definite(), is_negative_definite(), is_indefinite(), etc.

Note: A degenerate form is considered neither definite nor indefinite. Note: The zero-dim’l form is considered both positive definite and negative definite.

INPUT:

QuadraticForm

OUTPUT:

boolean

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1,1])
sage: Q.compute_definiteness()
sage: Q.is_positive_definite()
True
sage: Q.is_negative_definite()
False
sage: Q.is_indefinite()
False
sage: Q.is_definite()
True
sage: Q = DiagonalQuadraticForm(ZZ, [])
sage: Q.compute_definiteness()
sage: Q.is_positive_definite()
True
sage: Q.is_negative_definite()
True
sage: Q.is_indefinite()
False
sage: Q.is_definite()
True
sage: Q = DiagonalQuadraticForm(ZZ, [1,0,-1])
sage: Q.compute_definiteness()
sage: Q.is_positive_definite()
False
sage: Q.is_negative_definite()
False
sage: Q.is_indefinite()
False
sage: Q.is_definite()
False
compute_definiteness_string_by_determinants()

Compute the (positive) definiteness of a quadratic form by looking at the signs of all of its upper-left subdeterminants. See also self.compute_definiteness() for more documentation.

INPUT:

None

OUTPUT:

string describing the definiteness

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1,1])
sage: Q.compute_definiteness_string_by_determinants()
'pos_def'
sage: Q = DiagonalQuadraticForm(ZZ, [])
sage: Q.compute_definiteness_string_by_determinants()
'zero'
sage: Q = DiagonalQuadraticForm(ZZ, [1,0,-1])
sage: Q.compute_definiteness_string_by_determinants()
'degenerate'
sage: Q = DiagonalQuadraticForm(ZZ, [1,-1])
sage: Q.compute_definiteness_string_by_determinants()
'indefinite'
sage: Q = DiagonalQuadraticForm(ZZ, [-1,-1])
sage: Q.compute_definiteness_string_by_determinants()
'neg_def'
content()

Return the GCD of the coefficients of the quadratic form.

Warning

Only works over Euclidean domains (probably just \(\ZZ\)).

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1, 1])
sage: Q.matrix().gcd()
2
sage: Q.content()
1
sage: DiagonalQuadraticForm(ZZ, [1, 1]).is_primitive()
True
sage: DiagonalQuadraticForm(ZZ, [2, 4]).is_primitive()
False
sage: DiagonalQuadraticForm(ZZ, [2, 4]).primitive()
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 1 0 ]
[ * 2 ]
conway_cross_product_doubled_power(p)

Computes twice the power of p which evaluates the ‘cross product’ term in Conway’s mass formula.

INPUT:

\(p\) – a prime number > 0

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, range(1,8))
sage: Q.conway_cross_product_doubled_power(2)
18
sage: Q.conway_cross_product_doubled_power(3)
10
sage: Q.conway_cross_product_doubled_power(5)
6
sage: Q.conway_cross_product_doubled_power(7)
6
sage: Q.conway_cross_product_doubled_power(11)
0
sage: Q.conway_cross_product_doubled_power(13)
0
conway_diagonal_factor(p)

Computes the diagonal factor of Conway’s \(p\)-mass.

INPUT:

\(p\) – a prime number > 0

OUTPUT:

a rational number > 0

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, range(1,6))
sage: Q.conway_diagonal_factor(3)
81/256
conway_mass()

Compute the mass by using the Conway-Sloane mass formula.

OUTPUT:

a rational number > 0

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.conway_mass()
1/48

sage: Q = DiagonalQuadraticForm(ZZ, [7,1,1])
sage: Q.conway_mass()
3/16

sage: Q = QuadraticForm(ZZ, 3, [7, 2, 2, 2, 0, 2]) + DiagonalQuadraticForm(ZZ, [1])
sage: Q.conway_mass()
3/32

sage: Q = QuadraticForm(Matrix(ZZ,2,[2,1,1,2]))
sage: Q.conway_mass()
1/12
conway_octane_of_this_unimodular_Jordan_block_at_2()

Determines the ‘octane’ of this full unimodular Jordan block at the prime \(p=2\). This is an invariant defined \((mod 8)\), ad.

This assumes that the form is given as a block diagonal form with unimodular blocks of size <= 2 and the 1x1 blocks are all in the upper leftmost position.

INPUT:

none

OUTPUT:

an integer 0 <= x <= 7

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.conway_octane_of_this_unimodular_Jordan_block_at_2()
0
sage: Q = DiagonalQuadraticForm(ZZ, [1,5,13])
sage: Q.conway_octane_of_this_unimodular_Jordan_block_at_2()
3
sage: Q = DiagonalQuadraticForm(ZZ, [3,7,13])
sage: Q.conway_octane_of_this_unimodular_Jordan_block_at_2()
7
conway_p_mass(p)

Computes Conway’s \(p\)-mass.

INPUT:

\(p\) – a prime number > 0

OUTPUT:

a rational number > 0

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, range(1, 6))
sage: Q.conway_p_mass(2)
16/3
sage: Q.conway_p_mass(3)
729/256
conway_species_list_at_2()

Returns an integer called the ‘species’ which determines the type of the orthogonal group over the finite field \(F_p\).

This assumes that the given quadratic form is a unimodular Jordan block at an odd prime \(p\). When the dimension is odd then this number is always positive, otherwise it may be positive or negative.

Note: The species of a zero dim’l form is always 0+, so we interpret the return value of zero as positive here! =)

OUTPUT:

a list of integers

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, range(1,10))
sage: Q.conway_species_list_at_2()
[1, 5, 1, 1, 1, 1]
sage: Q = DiagonalQuadraticForm(ZZ, range(1,8))
sage: Q.conway_species_list_at_2()
[1, 3, 1, 1, 1]
conway_species_list_at_odd_prime(p)

Returns an integer called the ‘species’ which determines the type of the orthogonal group over the finite field \(F_p\).

This assumes that the given quadratic form is a unimodular Jordan block at an odd prime \(p\). When the dimension is odd then this number is always positive, otherwise it may be positive or negative (or zero, but that is considered positive by convention).

Note: The species of a zero dim’l form is always 0+, so we interpret the return value of zero as positive here! =)

INPUT:

a positive prime number

OUTPUT:

a list of integers

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, range(1,10))
sage: Q.conway_species_list_at_odd_prime(3)
[6, 2, 1]
sage: Q = DiagonalQuadraticForm(ZZ, range(1,8))
sage: Q.conway_species_list_at_odd_prime(3)
[5, 2]
sage: Q.conway_species_list_at_odd_prime(5)
[-6, 1]
conway_standard_mass()

Returns the infinite product of the standard mass factors.

INPUT:

none

OUTPUT:

a rational number > 0

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, [2, -2, 0, 3, -5, 4])
sage: Q.conway_standard_mass()
1/6
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.conway_standard_mass()
1/6
conway_standard_p_mass(p)

Computes the standard (generic) Conway-Sloane \(p\)-mass.

INPUT:

\(p\) – a prime number > 0

OUTPUT:

a rational number > 0

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.conway_standard_p_mass(2)
2/3
conway_type_factor()

This is a special factor only present in the mass formula when \(p=2\).

INPUT:

none

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, range(1,8))
sage: Q.conway_type_factor()
4
count_congruence_solutions(p, k, m, zvec, nzvec)

Counts all solutions of Q(\(x\)) \(= m (mod p^k)\) satisfying the additional congruence conditions described in QuadraticForm.count_congruence_solutions_as_vector().

INPUT:

\(p\) – prime number > 0

\(k\) – an integer > 0

\(m\) – an integer (depending only on mod \(p^k\))

zvec, nzvec – a list of integers in range(self.dim()), or None

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.count_congruence_solutions(3, 1, 0, None, None)
15
count_congruence_solutions__bad_type(p, k, m, zvec, nzvec)

Counts the bad-type solutions of Q(\(x\)) \(= m (mod p^k)\) satisfying the additional congruence conditions described in QuadraticForm.count_congruence_solutions_as_vector().

INPUT:

\(p\) – prime number > 0

\(k\) – an integer > 0

\(m\) – an integer (depending only on mod \(p^k\))

zvec, nzvec – a list of integers up to dim(Q)

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.count_congruence_solutions__bad_type(3, 1, 0, None, None)
2
count_congruence_solutions__bad_type_I(p, k, m, zvec, nzvec)

Counts the bad-typeI solutions of Q(\(x\)) = \(m (mod p^k)\) satisfying the additional congruence conditions described in QuadraticForm.count_congruence_solutions_as_vector().

INPUT:

\(p\) – prime number > 0

\(k\) – an integer > 0

\(m\) – an integer (depending only on mod \(p^k\))

zvec, nzvec – a list of integers up to dim(Q)

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.count_congruence_solutions__bad_type_I(3, 1, 0, None, None)
0
count_congruence_solutions__bad_type_II(p, k, m, zvec, nzvec)

Counts the bad-typeII solutions of Q(\(x\)) \(= m (mod p^k)\) satisfying the additional congruence conditions described in QuadraticForm.count_congruence_solutions_as_vector().

INPUT:

\(p\) – prime number > 0

\(k\) – an integer > 0

\(m\) – an integer (depending only on mod \(p^k\))

zvec, nzvec – a list of integers up to dim(Q)

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.count_congruence_solutions__bad_type_II(3, 1, 0, None, None)
2
count_congruence_solutions__good_type(p, k, m, zvec, nzvec)

Counts the good-type solutions of Q(x) = m (mod p^k) satisfying the additional congruence conditions described in QuadraticForm.count_congruence_solutions_as_vector().

INPUT:

\(p\) – prime number > 0

\(k\) – an integer > 0

\(m\) – an integer (depending only on mod \(p^k\))

zvec, nzvec – a list of integers up to dim(Q)

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.count_congruence_solutions__good_type(3, 1, 0, None, None)
12
count_congruence_solutions__zero_type(p, k, m, zvec, nzvec)

Counts the zero-type solutions of Q(\(x\)) = \(m (mod p^k)\) satisfying the additional congruence conditions described in QuadraticForm.count_congruence_solutions_as_vector().

INPUT:

\(p\) – prime number > 0

\(k\) – an integer > 0

\(m\) – an integer (depending only on mod \(p^k\))

zvec, nzvec – a list of integers up to dim(Q)

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.count_congruence_solutions__zero_type(3, 1, 0, None, None)
1
count_congruence_solutions_as_vector(p, k, m, zvec, nzvec)

Gives the number of integer solution vectors \(x\) satisfying the congruence Q(\(x\)) \(= m (mod p^k)\) of each solution type (i.e. All, Good, Zero, Bad, BadI, BadII) which satisfy the additional congruence conditions of having certain coefficients = 0 (mod \(p\)) and certain collections of coefficients not congruent to the zero vector (mod \(p\)).

A solution vector \(x\) satisfies the additional congruence conditions specified by zvec and nzvec (and therefore is counted) iff both of the following conditions hold:

  1. \(x[i] == 0 (mod p)\) for all \(i\) in zvec

  2. \(x[i] != 0 (mod p) for all i\) in nzvec

REFERENCES:

See Hanke’s (????) paper “Local Densities and explicit bounds…”, p??? for the definitions of the solution types and congruence conditions.

INPUT:

  • \(p\) – prime number > 0

  • \(k\) – an integer > 0

  • \(m\) – an integer (depending only on mod \(p^k\))

  • zvec, nzvec – a list of integers in range(self.dim()), or None

OUTPUT:

a list of six integers >= 0 representing the solution types: [All, Good, Zero, Bad, BadI, BadII]

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.count_congruence_solutions_as_vector(3, 1, 1, [], [])
[0, 0, 0, 0, 0, 0]
sage: Q.count_congruence_solutions_as_vector(3, 1, 1, None, [])
[0, 0, 0, 0, 0, 0]
sage: Q.count_congruence_solutions_as_vector(3, 1, 1, [], None)
[6, 6, 0, 0, 0, 0]
sage: Q.count_congruence_solutions_as_vector(3, 1, 1, None, None)
[6, 6, 0, 0, 0, 0]
sage: Q.count_congruence_solutions_as_vector(3, 1, 2, None, None)
[6, 6, 0, 0, 0, 0]
sage: Q.count_congruence_solutions_as_vector(3, 1, 0, None, None)
[15, 12, 1, 2, 0, 2]
count_modp_solutions__by_Gauss_sum(p, m)

Returns the number of solutions of \(Q(x) = m (mod p)\) of a non-degenerate quadratic form over the finite field \(Z/pZ\), where \(p\) is a prime number > 2.

Note: We adopt the useful convention that a zero-dimensional quadratic form has exactly one solution always (i.e. the empty vector).

These are defined in Table 1 on p363 of Hanke’s “Local Densities…” paper.

INPUT:

\(p\) – a prime number > 2

\(m\) – an integer

OUTPUT:

an integer >= 0

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: [Q.count_modp_solutions__by_Gauss_sum(3, m)  for m in range(3)]
[9, 6, 12]
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,2])
sage: [Q.count_modp_solutions__by_Gauss_sum(3, m)  for m in range(3)]
[9, 12, 6]
delta()

This is the omega of the adjoint form, which is the same as the omega of the reciprocal form.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,37])
sage: Q.delta()
148
det()

Gives the determinant of the Gram matrix of 2*Q, or equivalently the determinant of the Hessian matrix of Q.

(Note: This is always defined over the same ring as the quadratic form.)

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,2,3])
sage: Q.det()
8
dim()

Gives the number of variables of the quadratic form.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,2,3])
sage: Q.dim()
2
sage: parent(Q.dim())
Integer Ring
sage: Q = QuadraticForm(Q.matrix())
sage: Q.dim()
2
sage: parent(Q.dim())
Integer Ring
disc()

Return the discriminant of the quadratic form, defined as

  • \((-1)^n {\rm det}(B)\) for even dimension \(2n\)

  • \({\rm det}(B)/2\) for odd dimension

where \(2Q(x) = x^t B x\).

This agrees with the usual discriminant for binary and ternary quadratic forms.

EXAMPLES:

sage: DiagonalQuadraticForm(ZZ, [1]).disc()
1
sage: DiagonalQuadraticForm(ZZ, [1,1]).disc()
-4
sage: DiagonalQuadraticForm(ZZ, [1,1,1]).disc()
4
sage: DiagonalQuadraticForm(ZZ, [1,1,1,1]).disc()
16
discrec()

Return the discriminant of the reciprocal form.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,37])
sage: Q.disc()
148
sage: Q.discrec()
5476
sage: [4 * 37, 4 * 37^2]
[148, 5476]
divide_variable(c, i, in_place=False)

Replace the variables \(x_i\) by \((x_i)/c\) in the quadratic form (replacing the original form if the in_place flag is True).

Here \(c\) must be an element of the base_ring defining the quadratic form, and the division must be defined in the base ring.

INPUT:

\(c\) – an element of Q.base_ring()

\(i\) – an integer >= 0

OUTPUT:

a QuadraticForm (by default, otherwise none)

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,9,5,7])
sage: Q.divide_variable(3,1)
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 0 0 0 ]
[ * 1 0 0 ]
[ * * 5 0 ]
[ * * * 7 ]
elementary_substitution(c, i, j, in_place=False)

Perform the substitution \(x_i --> x_i + c*x_j\) (replacing the original form if the in_place flag is True).

INPUT:

\(c\) – an element of Q.base_ring()

\(i\), \(j\) – integers >= 0

OUTPUT:

a QuadraticForm (by default, otherwise none)

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 4, range(1,11))
sage: Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 2 3 4 ]
[ * 5 6 7 ]
[ * * 8 9 ]
[ * * * 10 ]

sage: Q.elementary_substitution(c=1, i=0, j=3)
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 2 3 6 ]
[ * 5 6 9 ]
[ * * 8 12 ]
[ * * * 15 ]
sage: R = QuadraticForm(ZZ, 4, range(1,11))
sage: R
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 2 3 4 ]
[ * 5 6 7 ]
[ * * 8 9 ]
[ * * * 10 ]
sage: M = Matrix(ZZ, 4, 4, [1,0,0,1,0,1,0,0,0,0,1,0,0,0,0,1])
sage: M
[1 0 0 1]
[0 1 0 0]
[0 0 1 0]
[0 0 0 1]
sage: R(M)
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 2 3 6 ]
[ * 5 6 9 ]
[ * * 8 12 ]
[ * * * 15 ]
extract_variables(var_indices)

Extract the variables (in order) whose indices are listed in var_indices, to give a new quadratic form.

INPUT:

var_indices – a list of integers >= 0

OUTPUT:

a QuadraticForm

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 4, range(10)); Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 0 1 2 3 ]
[ * 4 5 6 ]
[ * * 7 8 ]
[ * * * 9 ]
sage: Q.extract_variables([1,3])
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 4 6 ]
[ * 9 ]
find_entry_with_minimal_scale_at_prime(p)

Finds the entry of the quadratic form with minimal scale at the prime p, preferring diagonal entries in case of a tie. (I.e. If we write the quadratic form as a symmetric matrix M, then this entry M[i,j] has the minimal valuation at the prime p.)

Note: This answer is independent of the kind of matrix (Gram or Hessian) associated to the form.

INPUT:

\(p\) – a prime number > 0

OUTPUT:

a pair of integers >= 0

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [6, 2, 20]); Q
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 6 2 ]
[ * 20 ]
sage: Q.find_entry_with_minimal_scale_at_prime(2)
(0, 1)
sage: Q.find_entry_with_minimal_scale_at_prime(3)
(1, 1)
sage: Q.find_entry_with_minimal_scale_at_prime(5)
(0, 0)
find_p_neighbor_from_vec(p, y)

Return the \(p\)-neighbor of self defined by y.

Let \((L,q)\) be a lattice with \(b(L,L) \subseteq \ZZ\) which is maximal at \(p\). Let \(y \in L\) with \(b(y,y) \in p^2\ZZ\) then the \(p\)-neighbor of \(L\) at \(y\) is given by \(\ZZ y/p + L_y\) where \(L_y = \{x \in L | b(x,y) \in p \ZZ \}\) and \(b(x,y) = q(x+y)-q(x)-q(y)\) is the bilinear form associated to \(q\).

INPUT:

  • p – a prime number

  • y – a vector with \(q(y) \in p \ZZ\).

  • odd – (default=``False``) if \(p=2\) return also odd neighbors

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ,[1,1,1,1])
sage: v = vector([0,2,1,1])
sage: X = Q.find_p_neighbor_from_vec(3,v); X
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 0 0 0 ]
[ * 1 4 4 ]
[ * * 5 12 ]
[ * * * 9 ]

Since the base ring and the domain are not yet separate, for rational, half integral forms we just pretend the base ring is \(ZZ\):

sage: Q = QuadraticForm(QQ,matrix.diagonal([1,1,1,1]))
sage: v = vector([1,1,1,1])
sage: Q.find_p_neighbor_from_vec(2,v)
Quadratic form in 4 variables over Rational Field with coefficients:
[ 1/2 1 1 1 ]
[ * 1 1 2 ]
[ * * 1 2 ]
[ * * * 2 ]
find_primitive_p_divisible_vector__next(p, v=None)

Find the next \(p\)-primitive vector (up to scaling) in \(L/pL\) whose value is \(p\)-divisible, where the last vector returned was \(v\). For an initial call, no \(v\) needs to be passed.

Returns vectors whose last non-zero entry is normalized to 0 or 1 (so no lines are counted repeatedly). The ordering is by increasing the first non-normalized entry. If we have tested all (lines of) vectors, then return None.

OUTPUT:

vector or None

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [10,1,4])
sage: v = Q.find_primitive_p_divisible_vector__next(5); v
(1, 1)
sage: v = Q.find_primitive_p_divisible_vector__next(5, v); v
(1, 0)
sage: v = Q.find_primitive_p_divisible_vector__next(5, v); v
sage: Q = QuadraticForm(QQ,matrix.diagonal([1,1,1,1]))
sage: v = Q.find_primitive_p_divisible_vector__next(2)
sage: Q(v)
2
find_primitive_p_divisible_vector__random(p)

Find a random \(p\)-primitive vector in \(L/pL\) whose value is \(p\)-divisible.

Note

Since there are about \(p^{(n-2)}\) of these lines, we have a \(1/p\) chance of randomly finding an appropriate vector.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [10,1,4])
sage: v = Q.find_primitive_p_divisible_vector__random(5)    # random
sage: v
(3, 3)
sage: 5.divides(Q(v))
True
sage: Q = QuadraticForm(QQ,matrix.diagonal([1,1,1,1]))
sage: v = Q.find_primitive_p_divisible_vector__random(2)
sage: Q(v)
2
gcd()

Return the greatest common divisor of the coefficients of the quadratic form (as a polynomial).

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 4, range(1, 21, 2))
sage: Q.gcd()
1
sage: Q = QuadraticForm(ZZ, 4, range(0, 20, 2))
sage: Q.gcd()
2
static genera(sig_pair, determinant, max_scale=None, even=False)

Return a list of all global genera with the given conditions.

Here a genus is called global if it is non-empty.

INPUT:

  • sig_pair – a pair of non-negative integers giving the signature

  • determinant – an integer; the sign is ignored

  • max_scale – (default: None) an integer; the maximum scale of a jordan block

  • even – boolean (default: False)

OUTPUT:

A list of all (non-empty) global genera with the given conditions.

EXAMPLES:

sage: QuadraticForm.genera((4,0), 125, even=True)
[Genus of
None
Signature:  (4, 0)
Genus symbol at 2:    1^-4
Genus symbol at 5:     1^1 5^3, Genus of
None
Signature:  (4, 0)
Genus symbol at 2:    1^-4
Genus symbol at 5:     1^-2 5^1 25^-1, Genus of
None
Signature:  (4, 0)
Genus symbol at 2:    1^-4
Genus symbol at 5:     1^2 5^1 25^1, Genus of
None
Signature:  (4, 0)
Genus symbol at 2:    1^-4
Genus symbol at 5:     1^3 125^1]
global_genus_symbol()

Returns the genus of a two times a quadratic form over ZZ. These are defined by a collection of local genus symbols (a la Chapter 15 of Conway-Sloane), and a signature.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3,4])
sage: Q.global_genus_symbol()
Genus of
[2 0 0 0]
[0 4 0 0]
[0 0 6 0]
[0 0 0 8]
Signature:  (4, 0)
Genus symbol at 2:    [2^-2 4^1 8^1]_6
Genus symbol at 3:     1^3 3^-1
sage: Q = QuadraticForm(ZZ, 4, range(10))
sage: Q.global_genus_symbol()
Genus of
[ 0  1  2  3]
[ 1  8  5  6]
[ 2  5 14  8]
[ 3  6  8 18]
Signature:  (3, 1)
Genus symbol at 2:    1^-4
Genus symbol at 563:     1^3 563^-1
has_equivalent_Jordan_decomposition_at_prime(other, p)

Determines if the given quadratic form has a Jordan decomposition equivalent to that of self.

INPUT:

a QuadraticForm

OUTPUT:

boolean

EXAMPLES:

sage: Q1 = QuadraticForm(ZZ, 3, [1, 0, -1, 1, 0, 3])
sage: Q2 = QuadraticForm(ZZ, 3, [1, 0, 0, 2, -2, 6])
sage: Q3 = QuadraticForm(ZZ, 3, [1, 0, 0, 1, 0, 11])
sage: [Q1.level(), Q2.level(), Q3.level()]
[44, 44, 44]
sage: Q1.has_equivalent_Jordan_decomposition_at_prime(Q2,2)
False
sage: Q1.has_equivalent_Jordan_decomposition_at_prime(Q2,11)
False
sage: Q1.has_equivalent_Jordan_decomposition_at_prime(Q3,2)
False
sage: Q1.has_equivalent_Jordan_decomposition_at_prime(Q3,11)
True
sage: Q2.has_equivalent_Jordan_decomposition_at_prime(Q3,2)
True
sage: Q2.has_equivalent_Jordan_decomposition_at_prime(Q3,11)
False
has_integral_Gram_matrix()

Return whether the quadratic form has an integral Gram matrix (with respect to its base ring).

A warning is issued if the form is defined over a field, since in that case the return is trivially true.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [7,8,9])
sage: Q.has_integral_Gram_matrix()
True
sage: Q = QuadraticForm(ZZ, 2, [4,5,6])
sage: Q.has_integral_Gram_matrix()
False
hasse_conductor()

This is the product of all primes where the Hasse invariant equals -1

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, [1, 0, -1, 2, -1, 5])
sage: Q.hasse_invariant(2)
-1
sage: Q.hasse_invariant(37)
-1
sage: Q.hasse_conductor()
74
sage: DiagonalQuadraticForm(ZZ, [1, 1, 1]).hasse_conductor()
1
sage: QuadraticForm(ZZ, 3, [2, -2, 0, 2, 0, 5]).hasse_conductor()
10
hasse_invariant(p)

Computes the Hasse invariant at a prime \(p\) or at infinity, as given on p55 of Cassels’s book. If Q is diagonal with coefficients \(a_i\), then the (Cassels) Hasse invariant is given by

\[c_p = \prod_{i < j} (a_i, a_j)_p\]

where \((a,b)_p\) is the Hilbert symbol at \(p\). The underlying quadratic form must be non-degenerate over \(Q_p\) for this to make sense.

Warning

This is different from the O’Meara Hasse invariant, which allows \(i <= j\) in the product. That is given by the method hasse_invariant__OMeara(p).

Note

We should really rename this hasse_invariant__Cassels(), and set hasse_invariant() as a front-end to it.

INPUT:

  • \(p\) – a prime number > 0 or \(-1\) for the infinite place

OUTPUT:

1 or -1

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,2,3])
sage: Q.rational_diagonal_form()
Quadratic form in 2 variables over Rational Field with coefficients:
[ 1 0 ]
[ * 2 ]
sage: [Q.hasse_invariant(p) for p in prime_range(20)]
[1, 1, 1, 1, 1, 1, 1, 1]
sage: [Q.hasse_invariant__OMeara(p) for p in prime_range(20)]
[1, 1, 1, 1, 1, 1, 1, 1]
sage: Q = DiagonalQuadraticForm(ZZ, [1,-1])
sage: [Q.hasse_invariant(p) for p in prime_range(20)]
[1, 1, 1, 1, 1, 1, 1, 1]
sage: [Q.hasse_invariant__OMeara(p) for p in prime_range(20)]
[-1, 1, 1, 1, 1, 1, 1, 1]
sage: Q = DiagonalQuadraticForm(ZZ, [1,-1,5])
sage: [Q.hasse_invariant(p) for p in prime_range(20)]
[1, 1, 1, 1, 1, 1, 1, 1]
sage: [Q.hasse_invariant__OMeara(p) for p in prime_range(20)]
[-1, 1, 1, 1, 1, 1, 1, 1]
sage: K.<a>=NumberField(x^2-23)
sage: Q=DiagonalQuadraticForm(K,[-a,a+2])
sage: [Q.hasse_invariant(p) for p in K.primes_above(19)]
[-1, 1]
hasse_invariant__OMeara(p)

Compute the O’Meara Hasse invariant at a prime \(p\).

This is defined on p167 of O’Meara’s book. If Q is diagonal with coefficients \(a_i\), then the (Cassels) Hasse invariant is given by

\[c_p = \prod_{i <= j} (a_i, a_j)_p\]

where \((a,b)_p\) is the Hilbert symbol at \(p\).

Warning

This is different from the (Cassels) Hasse invariant, which only allows \(i < j\) in the product. That is given by the method hasse_invariant(p).

INPUT:

  • \(p\) – a prime number > 0 or \(-1\) for the infinite place

OUTPUT:

1 or -1

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,2,3])
sage: Q.rational_diagonal_form()
Quadratic form in 2 variables over Rational Field with coefficients:
[ 1 0 ]
[ * 2 ]
sage: [Q.hasse_invariant(p) for p in prime_range(20)]
[1, 1, 1, 1, 1, 1, 1, 1]
sage: [Q.hasse_invariant__OMeara(p) for p in prime_range(20)]
[1, 1, 1, 1, 1, 1, 1, 1]
sage: Q = DiagonalQuadraticForm(ZZ, [1,-1])
sage: [Q.hasse_invariant(p) for p in prime_range(20)]
[1, 1, 1, 1, 1, 1, 1, 1]
sage: [Q.hasse_invariant__OMeara(p) for p in prime_range(20)]
[-1, 1, 1, 1, 1, 1, 1, 1]
sage: Q = DiagonalQuadraticForm(ZZ,[1,-1,-1])
sage: [Q.hasse_invariant(p) for p in prime_range(20)]
[-1, 1, 1, 1, 1, 1, 1, 1]
sage: [Q.hasse_invariant__OMeara(p) for p in prime_range(20)]
[-1, 1, 1, 1, 1, 1, 1, 1]
sage: K.<a>=NumberField(x^2-23)
sage: Q = DiagonalQuadraticForm(K,[-a,a+2])
sage: [Q.hasse_invariant__OMeara(p) for p in K.primes_above(19)]
[1, 1]
is_adjoint()

Determines if the given form is the adjoint of another form

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, [1, 0, -1, 2, -1, 5])
sage: Q.is_adjoint()
False
sage: Q.adjoint().is_adjoint()
True
is_anisotropic(p)

Check if the quadratic form is anisotropic over the p-adic numbers \(\QQ_p\) or \(\RR\).

INPUT:

  • \(p\) – a prime number > 0 or \(-1\) for the infinite place

OUTPUT:

boolean

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1])
sage: Q.is_anisotropic(2)
True
sage: Q.is_anisotropic(3)
True
sage: Q.is_anisotropic(5)
False
sage: Q = DiagonalQuadraticForm(ZZ, [1,-1])
sage: Q.is_anisotropic(2)
False
sage: Q.is_anisotropic(3)
False
sage: Q.is_anisotropic(5)
False
sage: [DiagonalQuadraticForm(ZZ, [1, -least_quadratic_nonresidue(p)]).is_anisotropic(p)  for p in prime_range(3, 30)]
[True, True, True, True, True, True, True, True, True]
sage: [DiagonalQuadraticForm(ZZ, [1, -least_quadratic_nonresidue(p), p, -p*least_quadratic_nonresidue(p)]).is_anisotropic(p)  for p in prime_range(3, 30)]
[True, True, True, True, True, True, True, True, True]
is_definite()

Determines if the given quadratic form is (positive or negative) definite.

Note: A degenerate form is considered neither definite nor indefinite. Note: The zero-dim’l form is considered indefinite.

INPUT:

None

OUTPUT:

boolean – True or False

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [-1,-3,-5])
sage: Q.is_definite()
True
sage: Q = DiagonalQuadraticForm(ZZ, [1,-3,5])
sage: Q.is_definite()
False
is_even(allow_rescaling_flag=True)

Returns true iff after rescaling by some appropriate factor, the form represents no odd integers. For more details, see parity().

Requires that Q is defined over \(ZZ\).

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1, 0, 1])
sage: Q.is_even()
False
sage: Q = QuadraticForm(ZZ, 2, [1, 1, 1])
sage: Q.is_even()
True
is_globally_equivalent_to(other, return_matrix=False)

Determine if the current quadratic form is equivalent to the given form over ZZ.

If return_matrix is True, then we return the transformation matrix \(M\) so that self(M) == other.

INPUT:

  • self, other – positive definite integral quadratic forms

  • return_matrix – (boolean, default False) return the transformation matrix instead of a boolean

OUTPUT:

  • if return_matrix is False: a boolean

  • if return_matrix is True: either False or the transformation matrix

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: M = Matrix(ZZ, 4, 4, [1,2,0,0, 0,1,0,0, 0,0,1,0, 0,0,0,1])
sage: Q1 = Q(M)
sage: Q.is_globally_equivalent_to(Q1)
True
sage: MM = Q.is_globally_equivalent_to(Q1, return_matrix=True)
sage: Q(MM) == Q1
True
sage: Q1 = QuadraticForm(ZZ, 3, [1, 0, -1, 2, -1, 5])
sage: Q2 = QuadraticForm(ZZ, 3, [2, 1, 2, 2, 1, 3])
sage: Q3 = QuadraticForm(ZZ, 3, [8, 6, 5, 3, 4, 2])
sage: Q1.is_globally_equivalent_to(Q2)
False
sage: Q1.is_globally_equivalent_to(Q2, return_matrix=True)
False
sage: Q1.is_globally_equivalent_to(Q3)
True
sage: M = Q1.is_globally_equivalent_to(Q3, True); M
[-1 -1  0]
[ 1  1  1]
[-1  0  0]
sage: Q1(M) == Q3
True
sage: Q = DiagonalQuadraticForm(ZZ, [1, -1])
sage: Q.is_globally_equivalent_to(Q)
Traceback (most recent call last):
...
ValueError: not a definite form in QuadraticForm.is_globally_equivalent_to()

ALGORITHM: this uses the PARI function pari:qfisom, implementing an algorithm by Plesken and Souvignier.

is_hyperbolic(p)

Check if the quadratic form is a sum of hyperbolic planes over the \(p\)-adic numbers \(\QQ_p\) or over the real numbers \(\RR\).

REFERENCES:

This criteria follows from Cassels’s “Rational Quadratic Forms”:

  • local invariants for hyperbolic plane (Lemma 2.4, p58)

  • direct sum formulas (Lemma 2.3, p58)

INPUT:

  • \(p\) – a prime number > 0 or \(-1\) for the infinite place

OUTPUT:

boolean

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1])
sage: Q.is_hyperbolic(-1)
False
sage: Q.is_hyperbolic(2)
False
sage: Q.is_hyperbolic(3)
False
sage: Q.is_hyperbolic(5)     ## Here -1 is a square, so it's true.
True
sage: Q.is_hyperbolic(7)
False
sage: Q.is_hyperbolic(13)    ## Here -1 is a square, so it's true.
True
is_indefinite()

Determines if the given quadratic form is indefinite.

Note: A degenerate form is considered neither definite nor indefinite. Note: The zero-dim’l form is not considered indefinite.

INPUT:

None

OUTPUT:

boolean – True or False

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [-1,-3,-5])
sage: Q.is_indefinite()
False
sage: Q = DiagonalQuadraticForm(ZZ, [1,-3,5])
sage: Q.is_indefinite()
True
is_isotropic(p)

Checks if Q is isotropic over the p-adic numbers \(Q_p\) or \(RR\).

INPUT:

  • \(p\) – a prime number > 0 or \(-1\) for the infinite place

OUTPUT:

boolean

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1])
sage: Q.is_isotropic(2)
False
sage: Q.is_isotropic(3)
False
sage: Q.is_isotropic(5)
True
sage: Q = DiagonalQuadraticForm(ZZ, [1,-1])
sage: Q.is_isotropic(2)
True
sage: Q.is_isotropic(3)
True
sage: Q.is_isotropic(5)
True
sage: [DiagonalQuadraticForm(ZZ, [1, -least_quadratic_nonresidue(p)]).is_isotropic(p)  for p in prime_range(3, 30)]
[False, False, False, False, False, False, False, False, False]
sage: [DiagonalQuadraticForm(ZZ, [1, -least_quadratic_nonresidue(p), p, -p*least_quadratic_nonresidue(p)]).is_isotropic(p)  for p in prime_range(3, 30)]
[False, False, False, False, False, False, False, False, False]
is_locally_equivalent_to(other, check_primes_only=False, force_jordan_equivalence_test=False)

Determine if the current quadratic form (defined over ZZ) is locally equivalent to the given form over the real numbers and the \(p\)-adic integers for every prime p.

This works by comparing the local Jordan decompositions at every prime, and the dimension and signature at the real place.

INPUT:

a QuadraticForm

OUTPUT:

boolean

EXAMPLES:

sage: Q1 = QuadraticForm(ZZ, 3, [1, 0, -1, 2, -1, 5])
sage: Q2 = QuadraticForm(ZZ, 3, [2, 1, 2, 2, 1, 3])
sage: Q1.is_globally_equivalent_to(Q2)
False
sage: Q1.is_locally_equivalent_to(Q2)
True
is_locally_represented_number(m)

Determines if the rational number m is locally represented by the quadratic form.

INPUT:

\(m\) – an integer

OUTPUT:

boolean

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.is_locally_represented_number(2)
True
sage: Q.is_locally_represented_number(7)
False
sage: Q.is_locally_represented_number(-1)
False
sage: Q.is_locally_represented_number(28)
False
sage: Q.is_locally_represented_number(0)
True
is_locally_represented_number_at_place(m, p)

Determines if the rational number m is locally represented by the quadratic form at the (possibly infinite) prime \(p\).

INPUT:

\(m\) – an integer

\(p\) – a prime number > 0 or ‘infinity’

OUTPUT:

boolean

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.is_locally_represented_number_at_place(7, infinity)
True
sage: Q.is_locally_represented_number_at_place(7, 2)
False
sage: Q.is_locally_represented_number_at_place(7, 3)
True
sage: Q.is_locally_represented_number_at_place(7, 5)
True
sage: Q.is_locally_represented_number_at_place(-1, infinity)
False
sage: Q.is_locally_represented_number_at_place(-1, 2)
False
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1,-1])
sage: Q.is_locally_represented_number_at_place(7, infinity)     # long time (8.5 s)
True
sage: Q.is_locally_represented_number_at_place(7, 2)            # long time
True
sage: Q.is_locally_represented_number_at_place(7, 3)            # long time
True
sage: Q.is_locally_represented_number_at_place(7, 5)            # long time
True
is_locally_universal_at_all_places()

Determines if the quadratic form represents \(Z_p\) for all finite/non-archimedean primes, and represents all real numbers.

INPUT:

none

OUTPUT:

boolean

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.is_locally_universal_at_all_places()
False
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.is_locally_universal_at_all_places()
False
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1,-1])
sage: Q.is_locally_universal_at_all_places()        # long time (8.5 s)
True
is_locally_universal_at_all_primes()

Determines if the quadratic form represents \(Z_p\) for all finite/non-archimedean primes.

INPUT:

none

OUTPUT:

boolean

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.is_locally_universal_at_all_primes()
True
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.is_locally_universal_at_all_primes()
True
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.is_locally_universal_at_all_primes()
False
is_locally_universal_at_prime(p)

Determines if the (integer-valued/rational) quadratic form represents all of \(Z_p\).

INPUT:

\(p\) – a positive prime number or “infinity”.

OUTPUT:

boolean

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.is_locally_universal_at_prime(2)
True
sage: Q.is_locally_universal_at_prime(3)
True
sage: Q.is_locally_universal_at_prime(5)
True
sage: Q.is_locally_universal_at_prime(infinity)
False
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.is_locally_universal_at_prime(2)
False
sage: Q.is_locally_universal_at_prime(3)
True
sage: Q.is_locally_universal_at_prime(5)
True
sage: Q.is_locally_universal_at_prime(infinity)
False
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,-1])
sage: Q.is_locally_universal_at_prime(infinity)
True
is_negative_definite()

Determines if the given quadratic form is negative-definite.

Note: A degenerate form is considered neither definite nor indefinite. Note: The zero-dim’l form is considered both positive definite and negative definite.

INPUT:

None

OUTPUT:

boolean – True or False

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [-1,-3,-5])
sage: Q.is_negative_definite()
True
sage: Q = DiagonalQuadraticForm(ZZ, [1,-3,5])
sage: Q.is_negative_definite()
False
is_odd(allow_rescaling_flag=True)

Returns true iff after rescaling by some appropriate factor, the form represents some odd integers. For more details, see parity().

Requires that Q is defined over \(ZZ\).

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1, 0, 1])
sage: Q.is_odd()
True
sage: Q = QuadraticForm(ZZ, 2, [1, 1, 1])
sage: Q.is_odd()
False
is_positive_definite()

Determines if the given quadratic form is positive-definite.

Note: A degenerate form is considered neither definite nor indefinite. Note: The zero-dim’l form is considered both positive definite and negative definite.

INPUT:

None

OUTPUT:

boolean – True or False

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5])
sage: Q.is_positive_definite()
True
sage: Q = DiagonalQuadraticForm(ZZ, [1,-3,5])
sage: Q.is_positive_definite()
False
is_primitive()

Determines if the given integer-valued form is primitive (i.e. not an integer (>1) multiple of another integer-valued quadratic form).

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [2,3,4])
sage: Q.is_primitive()
True
sage: Q = QuadraticForm(ZZ, 2, [2,4,8])
sage: Q.is_primitive()
False
is_rationally_isometric(other, return_matrix=False)

Determine if two regular quadratic forms over a number field are isometric.

INPUT:

  • other – a quadratic form over a number field

  • return_matrix – (boolean, default False) return the transformation matrix instead of a boolean; this is currently only implemented for forms over QQ

OUTPUT:

  • if return_matrix is False: a boolean

  • if return_matrix is True: either False or the transformation matrix

EXAMPLES:

sage: V = DiagonalQuadraticForm(QQ, [1, 1, 2])
sage: W = DiagonalQuadraticForm(QQ, [2, 2, 2])
sage: V.is_rationally_isometric(W)
True
sage: K.<a> = NumberField(x^2-3)
sage: V = QuadraticForm(K, 4, [1, 0, 0, 0, 2*a, 0, 0, a, 0, 2]); V
Quadratic form in 4 variables over Number Field in a with defining polynomial x^2 - 3 with coefficients:
[ 1 0 0 0 ]
[ * 2*a 0 0 ]
[ * * a 0 ]
[ * * * 2 ]
sage: W = QuadraticForm(K, 4, [1, 2*a, 4, 6, 3, 10, 2, 1, 2, 5]); W
Quadratic form in 4 variables over Number Field in a with defining polynomial x^2 - 3 with coefficients:
[ 1 2*a 4 6 ]
[ * 3 10 2 ]
[ * * 1 2 ]
[ * * * 5 ]
sage: V.is_rationally_isometric(W)
False
sage: K.<a> = NumberField(x^4 + 2*x + 6)
sage: V = DiagonalQuadraticForm(K, [a, 2, 3, 2, 1]); V
Quadratic form in 5 variables over Number Field in a with defining polynomial x^4 + 2*x + 6 with coefficients:
[ a 0 0 0 0 ]
[ * 2 0 0 0 ]
[ * * 3 0 0 ]
[ * * * 2 0 ]
[ * * * * 1 ]
sage: W = DiagonalQuadraticForm(K, [a, a, a, 2, 1]); W
Quadratic form in 5 variables over Number Field in a with defining polynomial x^4 + 2*x + 6 with   coefficients:
[ a 0 0 0 0 ]
[ * a 0 0 0 ]
[ * * a 0 0 ]
[ * * * 2 0 ]
[ * * * * 1 ]
sage: V.is_rationally_isometric(W)
False
sage: K.<a> = NumberField(x^2 - 3)
sage: V = DiagonalQuadraticForm(K, [-1, a, -2*a])
sage: W = DiagonalQuadraticForm(K, [-1, -a, 2*a])
sage: V.is_rationally_isometric(W)
True

sage: V = DiagonalQuadraticForm(QQ, [1, 1, 2])
sage: W = DiagonalQuadraticForm(QQ, [2, 2, 2])
sage: T = V.is_rationally_isometric(W, True); T
[   0    0    1]
[-1/2 -1/2    0]
[ 1/2 -1/2    0]
sage: V.Gram_matrix() == T.transpose() * W.Gram_matrix() * T
True

sage: T = W.is_rationally_isometric(V, True); T
[ 0 -1  1]
[ 0 -1 -1]
[ 1  0  0]
sage: W.Gram_matrix() == T.T * V.Gram_matrix() * T
True
sage: L = QuadraticForm(QQ, 3, [2, 2, 0, 2, 2, 5])
sage: M = QuadraticForm(QQ, 3, [2, 2, 0, 3, 2, 3])
sage: L.is_rationally_isometric(M, True)
False
sage: A = DiagonalQuadraticForm(QQ, [1, 5])
sage: B = QuadraticForm(QQ, 2, [1, 12, 81])
sage: T = A.is_rationally_isometric(B, True); T
[  1  -2]
[  0 1/3]
sage: A.Gram_matrix() == T.T * B.Gram_matrix() * T
True
sage: C = DiagonalQuadraticForm(QQ, [1, 5, 9])
sage: D = DiagonalQuadraticForm(QQ, [6, 30, 1])
sage: T = C.is_rationally_isometric(D, True); T
[   0 -5/6  1/2]
[   0  1/6  1/2]
[  -1    0    0]
sage: C.Gram_matrix() == T.T * D.Gram_matrix() * T
True
sage: E = DiagonalQuadraticForm(QQ, [1, 1])
sage: F = QuadraticForm(QQ, 2, [17, 94, 130])
sage: T = F.is_rationally_isometric(E, True); T
[     -4 -189/17]
[     -1  -43/17]
sage: F.Gram_matrix() == T.T * E.Gram_matrix() * T
True
is_zero(v, p=0)

Determines if the vector v is on the conic Q(x) = 0 (mod p).

EXAMPLES:

sage: Q1 = QuadraticForm(ZZ, 3, [1, 0, -1, 2, -1, 5])
sage: Q1.is_zero([0,1,0], 2)
True
sage: Q1.is_zero([1,1,1], 2)
True
sage: Q1.is_zero([1,1,0], 2)
False
is_zero_nonsingular(v, p=0)

Determines if the vector \(v\) is on the conic Q(\(x\)) = 0 (mod \(p\)), and that this point is non-singular point of the conic.

EXAMPLES:

sage: Q1 = QuadraticForm(ZZ, 3, [1, 0, -1, 2, -1, 5])
sage: Q1.is_zero_nonsingular([1,1,1], 2)
True
sage: Q1.is_zero([1, 19, 2], 37)
True
sage: Q1.is_zero_nonsingular([1, 19, 2], 37)
False
is_zero_singular(v, p=0)

Determines if the vector \(v\) is on the conic Q(\(x\)) = 0 (mod \(p\)), and that this point is singular point of the conic.

EXAMPLES:

sage: Q1 = QuadraticForm(ZZ, 3, [1, 0, -1, 2, -1, 5])
sage: Q1.is_zero([1,1,1], 2)
True
sage: Q1.is_zero_singular([1,1,1], 2)
False
sage: Q1.is_zero_singular([1, 19, 2], 37)
True
jordan_blocks_by_scale_and_unimodular(p, safe_flag=True)

Return a list of pairs \((s_i, L_i)\) where \(L_i\) is a maximal \(p^{s_i}\)-unimodular Jordan component which is further decomposed into block diagonals of block size \(\le 2\).

For each \(L_i\) the 2x2 blocks are listed after the 1x1 blocks (which follows from the convention of the local_normal_form() method).

Note

The decomposition of each \(L_i\) into smaller blocks is not unique!

The safe_flag argument allows us to select whether we want a copy of the output, or the original output. By default safe_flag = True, so we return a copy of the cached information. If this is set to False, then the routine is much faster but the return values are vulnerable to being corrupted by the user.

INPUT:

  • \(p\) – a prime number > 0.

OUTPUT:

A list of pairs \((s_i, L_i)\) where:

  • \(s_i\) is an integer,

  • \(L_i\) is a block-diagonal unimodular quadratic form over \(\ZZ_p\).

Note

These forms \(L_i\) are defined over the \(p\)-adic integers, but by a matrix over \(\ZZ\) (or \(\QQ\)?).

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,9,5,7])
sage: Q.jordan_blocks_by_scale_and_unimodular(3)
[(0, Quadratic form in 3 variables over Integer Ring with coefficients:
[ 1 0 0 ]
[ * 5 0 ]
[ * * 7 ]), (2, Quadratic form in 1 variables over Integer Ring with coefficients:
[ 1 ])]
sage: Q2 = QuadraticForm(ZZ, 2, [1,1,1])
sage: Q2.jordan_blocks_by_scale_and_unimodular(2)
[(-1, Quadratic form in 2 variables over Integer Ring with coefficients:
[ 2 2 ]
[ * 2 ])]
sage: Q = Q2 + Q2.scale_by_factor(2)
sage: Q.jordan_blocks_by_scale_and_unimodular(2)
[(-1, Quadratic form in 2 variables over Integer Ring with coefficients:
[ 2 2 ]
[ * 2 ]), (0, Quadratic form in 2 variables over Integer Ring with coefficients:
[ 2 2 ]
[ * 2 ])]
jordan_blocks_in_unimodular_list_by_scale_power(p)

Returns a list of Jordan components, whose component at index i should be scaled by the factor p^i.

This is only defined for integer-valued quadratic forms (i.e. forms with base_ring ZZ), and the indexing only works correctly for p=2 when the form has an integer Gram matrix.

INPUT:

self – a quadratic form over ZZ, which has integer Gram matrix if p == 2 \(p\) – a prime number > 0

OUTPUT:

a list of p-unimodular quadratic forms

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, [2, -2, 0, 3, -5, 4])
sage: Q.jordan_blocks_in_unimodular_list_by_scale_power(2)
Traceback (most recent call last):
...
TypeError: Oops!  The given quadratic form has a Jordan component with a negative scale exponent!
This routine requires an integer-matrix quadratic form for the output indexing to work properly!

sage: Q.scale_by_factor(2).jordan_blocks_in_unimodular_list_by_scale_power(2)
[Quadratic form in 2 variables over Integer Ring with coefficients:
[ 0 2 ]
[ * 0 ], Quadratic form in 0 variables over Integer Ring with coefficients:
, Quadratic form in 1 variables over Integer Ring with coefficients:
[ 345 ]]

sage: Q.jordan_blocks_in_unimodular_list_by_scale_power(3)
[Quadratic form in 2 variables over Integer Ring with coefficients:
[ 2 0 ]
[ * 10 ], Quadratic form in 1 variables over Integer Ring with coefficients:
[ 2 ]]
level()

Determines the level of the quadratic form over a PID, which is a generator for the smallest ideal \(N\) of \(R\) such that N * (the matrix of 2*Q)^(-1) is in R with diagonal in 2*R.

Over \(\ZZ\) this returns a non-negative number.

(Caveat: This always returns the unit ideal when working over a field!)

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, range(1,4))
sage: Q.level()
8

sage: Q1 = QuadraticForm(QQ, 2, range(1,4))
sage: Q1.level()      # random
UserWarning: Warning -- The level of a quadratic form over a field is always 1.  Do you really want to do this?!?
1

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.level()
420
level__Tornaria()

Return the level of the quadratic form, defined as

level(B) for even dimension level(B)/4 for odd dimension

where 2Q(\(x\)) \(= x^t * B * x\).

This agrees with the usual level for even dimension…

EXAMPLES:

sage: DiagonalQuadraticForm(ZZ, [1]).level__Tornaria()
1
sage: DiagonalQuadraticForm(ZZ, [1,1]).level__Tornaria()
4
sage: DiagonalQuadraticForm(ZZ, [1,1,1]).level__Tornaria()
1
sage: DiagonalQuadraticForm(ZZ, [1,1,1,1]).level__Tornaria()
4
level_ideal()

Determines the level of the quadratic form (over R), which is the smallest ideal N of R such that N * (the matrix of 2*Q)^(-1) is in R with diagonal in 2*R. (Caveat: This always returns the principal ideal when working over a field!)

WARNING: THIS ONLY WORKS OVER A PID RING OF INTEGERS FOR NOW!

(Waiting for Sage fractional ideal support.)

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, range(1,4))
sage: Q.level_ideal()
Principal ideal (8) of Integer Ring
sage: Q1 = QuadraticForm(QQ, 2, range(1,4))
sage: Q1.level_ideal()
Principal ideal (1) of Rational Field
sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.level_ideal()
Principal ideal (420) of Integer Ring
list_external_initializations()

Return a list of the fields which were set externally at creation, and not created through the usual QuadraticForm methods. These fields are as good as the external process that made them, and are thus not guaranteed to be correct.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,0,5])
sage: Q.list_external_initializations()
[]
sage: T = Q.theta_series()
sage: Q.list_external_initializations()
[]
sage: Q = QuadraticForm(ZZ, 2, [1,0,5], unsafe_initialization=False, number_of_automorphisms=3, determinant=0)
sage: Q.list_external_initializations()
[]
sage: Q = QuadraticForm(ZZ, 2, [1,0,5], unsafe_initialization=False, number_of_automorphisms=3, determinant=0)
sage: Q.list_external_initializations()
[]
sage: Q = QuadraticForm(ZZ, 2, [1,0,5], unsafe_initialization=True, number_of_automorphisms=3, determinant=0)
sage: Q.list_external_initializations()
['number_of_automorphisms', 'determinant']
lll()

Return an LLL-reduced form of Q (using Pari).

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 4, range(1,11))
sage: Q.is_definite()
True
sage: Q.lll()
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 0 1 0 ]
[ * 4 3 3 ]
[ * * 6 3 ]
[ * * * 6 ]
local_badII_density_congruence(p, m, Zvec=None, NZvec=None)

Finds the Bad-type II local density of Q representing \(m\) at \(p\). (Assuming that \(p\) > 2 and Q is given in local diagonal form.)

INPUT:

Q – quadratic form assumed to be block diagonal and p-integral

\(p\) – a prime number

\(m\) – an integer

Zvec, NZvec – non-repeating lists of integers in range(self.dim()) or None

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.local_badII_density_congruence(2, 1, None, None)
0
sage: Q.local_badII_density_congruence(2, 2, None, None)
0
sage: Q.local_badII_density_congruence(2, 4, None, None)
0
sage: Q.local_badII_density_congruence(3, 1, None, None)
0
sage: Q.local_badII_density_congruence(3, 6, None, None)
0
sage: Q.local_badII_density_congruence(3, 9, None, None)
0
sage: Q.local_badII_density_congruence(3, 27, None, None)
0
sage: Q = DiagonalQuadraticForm(ZZ, [1,3,3,9,9])
sage: Q.local_badII_density_congruence(3, 1, None, None)
0
sage: Q.local_badII_density_congruence(3, 3, None, None)
0
sage: Q.local_badII_density_congruence(3, 6, None, None)
0
sage: Q.local_badII_density_congruence(3, 9, None, None)
4/27
sage: Q.local_badII_density_congruence(3, 18, None, None)
4/9
local_badI_density_congruence(p, m, Zvec=None, NZvec=None)

Finds the Bad-type I local density of Q representing \(m\) at \(p\). (Assuming that p > 2 and Q is given in local diagonal form.)

INPUT:

Q – quadratic form assumed to be block diagonal and \(p\)-integral

\(p\) – a prime number

\(m\) – an integer

Zvec, NZvec – non-repeating lists of integers in range(self.dim()) or None

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.local_badI_density_congruence(2, 1, None, None)
0
sage: Q.local_badI_density_congruence(2, 2, None, None)
1
sage: Q.local_badI_density_congruence(2, 4, None, None)
0
sage: Q.local_badI_density_congruence(3, 1, None, None)
0
sage: Q.local_badI_density_congruence(3, 6, None, None)
0
sage: Q.local_badI_density_congruence(3, 9, None, None)
0
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.local_badI_density_congruence(2, 1, None, None)
0
sage: Q.local_badI_density_congruence(2, 2, None, None)
0
sage: Q.local_badI_density_congruence(2, 4, None, None)
0
sage: Q.local_badI_density_congruence(3, 2, None, None)
0
sage: Q.local_badI_density_congruence(3, 6, None, None)
0
sage: Q.local_badI_density_congruence(3, 9, None, None)
0
sage: Q = DiagonalQuadraticForm(ZZ, [1,3,3,9])
sage: Q.local_badI_density_congruence(3, 1, None, None)
0
sage: Q.local_badI_density_congruence(3, 3, None, None)
4/3
sage: Q.local_badI_density_congruence(3, 6, None, None)
4/3
sage: Q.local_badI_density_congruence(3, 9, None, None)
0
sage: Q.local_badI_density_congruence(3, 18, None, None)
0
local_bad_density_congruence(p, m, Zvec=None, NZvec=None)

Finds the Bad-type local density of Q representing \(m\) at \(p\), allowing certain congruence conditions mod \(p\).

INPUT:

Q – quadratic form assumed to be block diagonal and p-integral

\(p\) – a prime number

\(m\) – an integer

Zvec, NZvec – non-repeating lists of integers in range(self.dim()) or None

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.local_bad_density_congruence(2, 1, None, None)
0
sage: Q.local_bad_density_congruence(2, 2, None, None)
1
sage: Q.local_bad_density_congruence(2, 4, None, None)
0
sage: Q.local_bad_density_congruence(3, 1, None, None)
0
sage: Q.local_bad_density_congruence(3, 6, None, None)
0
sage: Q.local_bad_density_congruence(3, 9, None, None)
0
sage: Q.local_bad_density_congruence(3, 27, None, None)
0
sage: Q = DiagonalQuadraticForm(ZZ, [1,3,3,9,9])
sage: Q.local_bad_density_congruence(3, 1, None, None)
0
sage: Q.local_bad_density_congruence(3, 3, None, None)
4/3
sage: Q.local_bad_density_congruence(3, 6, None, None)
4/3
sage: Q.local_bad_density_congruence(3, 9, None, None)
4/27
sage: Q.local_bad_density_congruence(3, 18, None, None)
4/9
sage: Q.local_bad_density_congruence(3, 27, None, None)
8/27
local_density(p, m)

Gives the local density – should be called by the user. =)

NOTE: This screens for imprimitive forms, and puts the quadratic form in local normal form, which is a requirement of the routines performing the computations!

INPUT:

\(p\) – a prime number > 0 \(m\) – an integer

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])   ## NOTE: This is already in local normal form for *all* primes p!
sage: Q.local_density(p=2, m=1)
1
sage: Q.local_density(p=3, m=1)
8/9
sage: Q.local_density(p=5, m=1)
24/25
sage: Q.local_density(p=7, m=1)
48/49
sage: Q.local_density(p=11, m=1)
120/121
local_density_congruence(p, m, Zvec=None, NZvec=None)

Finds the local density of Q representing \(m\) at \(p\), allowing certain congruence conditions mod \(p\).

INPUT:

Q – quadratic form assumed to be block diagonal and p-integral

\(p\) – a prime number

\(m\) – an integer

Zvec, NZvec – non-repeating lists of integers in range(self.dim()) or None

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.local_density_congruence(p=2, m=1, Zvec=None, NZvec=None)
1
sage: Q.local_density_congruence(p=3, m=1, Zvec=None, NZvec=None)
8/9
sage: Q.local_density_congruence(p=5, m=1, Zvec=None, NZvec=None)
24/25
sage: Q.local_density_congruence(p=7, m=1, Zvec=None, NZvec=None)
48/49
sage: Q.local_density_congruence(p=11, m=1, Zvec=None, NZvec=None)
120/121
sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.local_density_congruence(2, 1, None, None)
1
sage: Q.local_density_congruence(2, 2, None, None)
1
sage: Q.local_density_congruence(2, 4, None, None)
3/2
sage: Q.local_density_congruence(3, 1, None, None)
2/3
sage: Q.local_density_congruence(3, 6, None, None)
4/3
sage: Q.local_density_congruence(3, 9, None, None)
14/9
sage: Q.local_density_congruence(3, 27, None, None)
2
sage: Q = DiagonalQuadraticForm(ZZ, [1,3,3,9,9])
sage: Q.local_density_congruence(3, 1, None, None)
2
sage: Q.local_density_congruence(3, 3, None, None)
4/3
sage: Q.local_density_congruence(3, 6, None, None)
4/3
sage: Q.local_density_congruence(3, 9, None, None)
2/9
sage: Q.local_density_congruence(3, 18, None, None)
4/9
local_genus_symbol(p)

Returns the Conway-Sloane genus symbol of 2 times a quadratic form defined over ZZ at a prime number p. This is defined (in the Genus_Symbol_p_adic_ring() class in the quadratic_forms/genera subfolder) to be a list of tuples (one for each Jordan component p^m*A at p, where A is a unimodular symmetric matrix with coefficients the p-adic integers) of the following form:

  1. If p>2 then return triples of the form [\(m\), \(n\), \(d\)] where

    \(m\) = valuation of the component

    \(n\) = rank of A

    \(d\) = det(A) in {1,u} for normalized quadratic non-residue u.

  2. If p=2 then return quintuples of the form [\(m\),`n`,`s`, \(d\), \(o\)] where

    \(m\) = valuation of the component

    \(n\) = rank of A

    \(d\) = det(A) in {1,3,5,7}

    \(s\) = 0 (or 1) if A is even (or odd)

    \(o\) = oddity of A (= 0 if s = 0) in Z/8Z

    = the trace of the diagonalization of A

NOTE: The Conway-Sloane convention for describing the prime ‘p = -1’ is not supported here, and neither is the convention for including the ‘prime’ Infinity. See note on p370 of Conway-Sloane (3rd ed) for a discussion of this convention.

INPUT:

-\(p\) – a prime number > 0

OUTPUT:

Returns a Conway-Sloane genus symbol at p, which is an instance of the Genus_Symbol_p_adic_ring class.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3,4])
sage: Q.local_genus_symbol(2)
Genus symbol at 2:    [2^-2 4^1 8^1]_6
sage: Q.local_genus_symbol(3)
Genus symbol at 3:     1^3 3^-1
sage: Q.local_genus_symbol(5)
Genus symbol at 5:     1^4
local_good_density_congruence(p, m, Zvec=None, NZvec=None)

Finds the Good-type local density of Q representing \(m\) at \(p\). (Front end routine for parity specific routines for p.)

TO DO: Add Documentation about the additional congruence conditions Zvec and NZvec.

INPUT:

Q – quadratic form assumed to be block diagonal and p-integral

\(p\) – a prime number

\(m\) – an integer

Zvec, NZvec – non-repeating lists of integers in range(self.dim()) or None

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.local_good_density_congruence(2, 1, None, None)
1
sage: Q.local_good_density_congruence(3, 1, None, None)
2/3
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.local_good_density_congruence(2, 1, None, None)
1
sage: Q.local_good_density_congruence(3, 1, None, None)
8/9
local_good_density_congruence_even(m, Zvec, NZvec)

Finds the Good-type local density of Q representing \(m\) at \(p=2\). (Assuming Q is given in local diagonal form.)

The additional congruence condition arguments Zvec and NZvec can be either a list of indices or None. Zvec = [] is equivalent to Zvec = None which both impose no additional conditions, but NZvec = [] returns no solutions always while NZvec = None imposes no additional condition.

WARNING: Here the indices passed in Zvec and NZvec represent indices of the solution vector \(x\) of Q(\(x\)) = \(m (mod p^k)\), and not the Jordan components of Q. They therefore are required (and assumed) to include either all or none of the indices of a given Jordan component of Q. This is only important when \(p=2\) since otherwise all Jordan blocks are 1x1, and so there the indices and Jordan blocks coincide.

TO DO: Add type checking for Zvec, NZvec, and that Q is in local normal form.

INPUT:

Q – quadratic form assumed to be block diagonal and 2-integral

\(p\) – a prime number

\(m\) – an integer

Zvec, NZvec – non-repeating lists of integers in range(self.dim()) or None

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.local_good_density_congruence_even(1, None, None)
1
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.local_good_density_congruence_even(1, None, None)
1
sage: Q.local_good_density_congruence_even(2, None, None)
3/2
sage: Q.local_good_density_congruence_even(3, None, None)
1
sage: Q.local_good_density_congruence_even(4, None, None)
1/2
sage: Q = QuadraticForm(ZZ, 4, range(10))
sage: Q[0,0] = 5
sage: Q[1,1] = 10
sage: Q[2,2] = 15
sage: Q[3,3] = 20
sage: Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 5 1 2 3 ]
[ * 10 5 6 ]
[ * * 15 8 ]
[ * * * 20 ]
sage: Q.theta_series(20)
1 + 2*q^5 + 2*q^10 + 2*q^14 + 2*q^15 + 2*q^16 + 2*q^18 + O(q^20)
sage: Q.local_normal_form(2)
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 0 1 0 0 ]
[ * 0 0 0 ]
[ * * 0 1 ]
[ * * * 0 ]
sage: Q.local_good_density_congruence_even(1, None, None)
3/4
sage: Q.local_good_density_congruence_even(2, None, None)
1
sage: Q.local_good_density_congruence_even(5, None, None)
3/4
local_good_density_congruence_odd(p, m, Zvec, NZvec)

Finds the Good-type local density of Q representing \(m\) at \(p\). (Assuming that \(p\) > 2 and Q is given in local diagonal form.)

The additional congruence condition arguments Zvec and NZvec can be either a list of indices or None. Zvec = [] is equivalent to Zvec = None which both impose no additional conditions, but NZvec = [] returns no solutions always while NZvec = None imposes no additional condition.

TO DO: Add type checking for Zvec, NZvec, and that Q is in local normal form.

INPUT:

Q – quadratic form assumed to be diagonal and p-integral

\(p\) – a prime number

\(m\) – an integer

Zvec, NZvec – non-repeating lists of integers in range(self.dim()) or None

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.local_good_density_congruence_odd(3, 1, None, None)
2/3
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.local_good_density_congruence_odd(3, 1, None, None)
8/9
local_normal_form(p)

Return a locally integrally equivalent quadratic form over the \(p\)-adic integers \(\ZZ_p\) which gives the Jordan decomposition.

The Jordan components are written as sums of blocks of size <= 2 and are arranged by increasing scale, and then by increasing norm. This is equivalent to saying that we put the 1x1 blocks before the 2x2 blocks in each Jordan component.

INPUT:

  • \(p\) – a positive prime number.

OUTPUT:

a quadratic form over \(\ZZ\)

Warning

Currently this only works for quadratic forms defined over \(\ZZ\).

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [10,4,1])
sage: Q.local_normal_form(5)
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 1 0 ]
[ * 6 ]
sage: Q.local_normal_form(3)
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 10 0 ]
[ * 15 ]

sage: Q.local_normal_form(2)
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 1 0 ]
[ * 6 ]
local_primitive_density(p, m)

Gives the local primitive density – should be called by the user. =)

NOTE: This screens for imprimitive forms, and puts the quadratic form in local normal form, which is a requirement of the routines performing the computations!

INPUT:

\(p\) – a prime number > 0 \(m\) – an integer

OUTPUT:

a rational number

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 4, range(10))
sage: Q[0,0] = 5
sage: Q[1,1] = 10
sage: Q[2,2] = 15
sage: Q[3,3] = 20
sage: Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 5 1 2 3 ]
[ * 10 5 6 ]
[ * * 15 8 ]
[ * * * 20 ]
sage: Q.theta_series(20)
1 + 2*q^5 + 2*q^10 + 2*q^14 + 2*q^15 + 2*q^16 + 2*q^18 + O(q^20)
sage: Q.local_normal_form(2)
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 0 1 0 0 ]
[ * 0 0 0 ]
[ * * 0 1 ]
[ * * * 0 ]

sage: Q.local_primitive_density(2, 1)
3/4
sage: Q.local_primitive_density(5, 1)
24/25

sage: Q.local_primitive_density(2, 5)
3/4
sage: Q.local_density(2, 5)
3/4
local_primitive_density_congruence(p, m, Zvec=None, NZvec=None)

Finds the primitive local density of Q representing \(m\) at \(p\), allowing certain congruence conditions mod \(p\).

Note: The following routine is not used internally, but is included for consistency.

INPUT:

Q – quadratic form assumed to be block diagonal and p-integral

\(p\) – a prime number

\(m\) – an integer

Zvec, NZvec – non-repeating lists of integers in range(self.dim()) or None

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.local_primitive_density_congruence(p=2, m=1, Zvec=None, NZvec=None)
1
sage: Q.local_primitive_density_congruence(p=3, m=1, Zvec=None, NZvec=None)
8/9
sage: Q.local_primitive_density_congruence(p=5, m=1, Zvec=None, NZvec=None)
24/25
sage: Q.local_primitive_density_congruence(p=7, m=1, Zvec=None, NZvec=None)
48/49
sage: Q.local_primitive_density_congruence(p=11, m=1, Zvec=None, NZvec=None)
120/121
sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.local_primitive_density_congruence(2, 1, None, None)
1
sage: Q.local_primitive_density_congruence(2, 2, None, None)
1
sage: Q.local_primitive_density_congruence(2, 4, None, None)
1
sage: Q.local_primitive_density_congruence(3, 1, None, None)
2/3
sage: Q.local_primitive_density_congruence(3, 6, None, None)
4/3
sage: Q.local_primitive_density_congruence(3, 9, None, None)
4/3
sage: Q.local_primitive_density_congruence(3, 27, None, None)
4/3
sage: Q = DiagonalQuadraticForm(ZZ, [1,3,3,9,9])
sage: Q.local_primitive_density_congruence(3, 1, None, None)
2
sage: Q.local_primitive_density_congruence(3, 3, None, None)
4/3
sage: Q.local_primitive_density_congruence(3, 6, None, None)
4/3
sage: Q.local_primitive_density_congruence(3, 9, None, None)
4/27
sage: Q.local_primitive_density_congruence(3, 18, None, None)
4/9
sage: Q.local_primitive_density_congruence(3, 27, None, None)
8/27
sage: Q.local_primitive_density_congruence(3, 81, None, None)
8/27
sage: Q.local_primitive_density_congruence(3, 243, None, None)
8/27
local_representation_conditions(recompute_flag=False, silent_flag=False)
WARNING: THIS ONLY WORKS CORRECTLY FOR FORMS IN >=3 VARIABLES,

WHICH ARE LOCALLY UNIVERSAL AT ALMOST ALL PRIMES!

This class finds the local conditions for a number to be integrally represented by an integer-valued quadratic form. These conditions are stored in “self.__local_representability_conditions” and consist of a list of 9 element vectors, with one for each prime with a local obstruction (though only the first 5 are meaningful unless \(p=2\) ). The first element is always the prime \(p\) where the local obstruction occurs, and the next 8 (or 4) entries represent square-classes in the \(p\)-adic integers \(Z_p\), and are labeled by the \(Q_p\) square-classes \(t*(Q_p)^2\) with \(t\) given as follows:

\(p > 2\) ==> [ * 1 u p u p * * * * ]

\(p = 2\) ==> [ * 1 3 5 7 2 6 10 14 ]

The integer appearing in each place tells us how \(p\)-divisible a number needs to be in that square-class in order to be locally represented by Q. A negative number indicates that the entire \(Q_p\) square-class is not represented, while a positive number \(x\) indicates that \(t*p^{(2*x)} (Z_p)^2\) is locally represented but \(t*p^{(2*(x-1))}\) \((Z_p)^2\) is not.

As an example, the vector

[2 3 0 0 0 0 2 0 infinity]

tells us that all positive integers are locally represented at p=2 except those of the forms:

\(2^6 * u * r^2\) with \(u = 1 (mod 8)\)

\(2^5 * u * r^2\) with \(u = 3 (mod 8)\)

\(2 * u * r^2\) with \(u = 7 (mod 8)\)

At the real numbers, the vector which looks like

[infinity, 0, infinity, None, None, None, None, None, None]

means that Q is negative definite (i.e. the 0 tells us all positive reals are represented). The real vector always appears, and is listed before the other ones.

INPUT:

none

OUTPUT:

A list of 9-element vectors describing the representation obstructions at primes dividing the level.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [])
sage: Q.local_representation_conditions()
This 0-dimensional form only represents zero.

sage: Q = DiagonalQuadraticForm(ZZ, [5])
sage: Q.local_representation_conditions()
This 1-dimensional form only represents square multiples of 5.

sage: Q1 = DiagonalQuadraticForm(ZZ, [1,1])
sage: Q1.local_representation_conditions()
This 2-dimensional form represents the p-adic integers of even
valuation for all primes p except [2].
For these and the reals, we have:
 Reals:   [0, +Infinity]
 p = 2:   [0, +Infinity, 0, +Infinity, 0, +Infinity, 0, +Infinity]


sage: Q1 = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q1.local_representation_conditions()
This form represents the p-adic integers Z_p for all primes p except
[2].  For these and the reals, we have:
 Reals:   [0, +Infinity]
 p = 2:   [0, 0, 0, +Infinity, 0, 0, 0, 0]

sage: Q1 = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q1.local_representation_conditions()
This form represents the p-adic integers Z_p for all primes p except
[].  For these and the reals, we have:
 Reals:   [0, +Infinity]

sage: Q1 = DiagonalQuadraticForm(ZZ, [1,3,3,3])
sage: Q1.local_representation_conditions()
This form represents the p-adic integers Z_p for all primes p except
[3].  For these and the reals, we have:
 Reals:   [0, +Infinity]
 p = 3:   [0, 1, 0, 0]

sage: Q2 = DiagonalQuadraticForm(ZZ, [2,3,3,3])
sage: Q2.local_representation_conditions()
This form represents the p-adic integers Z_p for all primes p except
[3].  For these and the reals, we have:
 Reals:   [0, +Infinity]
 p = 3:   [1, 0, 0, 0]

sage: Q3 = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q3.local_representation_conditions()
This form represents the p-adic integers Z_p for all primes p except
[].  For these and the reals, we have:
 Reals:   [0, +Infinity]
local_zero_density_congruence(p, m, Zvec=None, NZvec=None)

Finds the Zero-type local density of Q representing \(m\) at \(p\), allowing certain congruence conditions mod p.

INPUT:

Q – quadratic form assumed to be block diagonal and \(p\)-integral

\(p\) – a prime number

\(m\) – an integer

Zvec, NZvec – non-repeating lists of integers in range(self.dim()) or None

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,2,3])
sage: Q.local_zero_density_congruence(2, 2, None, None)
0
sage: Q.local_zero_density_congruence(2, 4, None, None)
1/2
sage: Q.local_zero_density_congruence(3, 6, None, None)
0
sage: Q.local_zero_density_congruence(3, 9, None, None)
2/9
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.local_zero_density_congruence(2, 2, None, None)
0
sage: Q.local_zero_density_congruence(2, 4, None, None)
1/4
sage: Q.local_zero_density_congruence(3, 6, None, None)
0
sage: Q.local_zero_density_congruence(3, 9, None, None)
8/81
mass__by_Siegel_densities(odd_algorithm='Pall', even_algorithm='Watson')

Gives the mass of transformations (det 1 and -1).

WARNING: THIS IS BROKEN RIGHT NOW… =(

Optional Arguments:

  • When p > 2 – odd_algorithm = “Pall” (only one choice for now)

  • When p = 2 – even_algorithm = “Kitaoka” or “Watson”

REFERENCES:

  • Nipp’s Book “Tables of Quaternary Quadratic Forms”.

  • Papers of Pall (only for p>2) and Watson (for \(p=2\) – tricky!).

  • Siegel, Milnor-Hussemoller, Conway-Sloane Paper IV, Kitoaka (all of which have problems…)

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: m = Q.mass__by_Siegel_densities(); m
1/384
sage: m - (2^Q.dim() * factorial(Q.dim()))^(-1)
0
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: m = Q.mass__by_Siegel_densities(); m
1/48
sage: m - (2^Q.dim() * factorial(Q.dim()))^(-1)
0
mass_at_two_by_counting_mod_power(k)

Computes the local mass at \(p=2\) assuming that it’s stable \((mod 2^k)\).

Note: This is way too slow to be useful, even when k=1!!!

TO DO: Remove this routine, or try to compile it!

INPUT:

k – an integer >= 1

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.mass_at_two_by_counting_mod_power(1)
4
matrix()

Return the Hessian matrix A for which Q(X) = \((1/2) * X^t * A * X\).

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, range(6))
sage: Q.matrix()
[ 0  1  2]
[ 1  6  4]
[ 2  4 10]
minkowski_reduction()

Find a Minkowski-reduced form equivalent to the given one. This means that

\[Q(v_k) <= Q(s_1 * v_1 + ... + s_n * v_n)\]

for all \(s_i\) where GCD`(s_k, … s_n) = 1`.

Note: When Q has dim <= 4 we can take all \(s_i\) in {1, 0, -1}.

References:
Schulze-Pillot’s paper on “An algorithm for computing genera

of ternary and quaternary quadratic forms”, p138.

Donaldson’s 1979 paper “Minkowski Reduction of Integral

Matrices”, p203.

EXAMPLES:

sage: Q = QuadraticForm(ZZ,4,[30, 17, 11, 12, 29, 25, 62, 64, 25, 110])
sage: Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 30 17 11 12 ]
[ * 29 25 62 ]
[ * * 64 25 ]
[ * * * 110 ]
sage: Q.minkowski_reduction()
(
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 30 17 11 -5 ]
[ * 29 25 4 ]
[ * * 64 0 ]
[ * * * 77 ]                                                       ,

[ 1  0  0  0]
[ 0  1  0 -1]
[ 0  0  1  0]
[ 0  0  0  1]
)
sage: Q=QuadraticForm(ZZ,4,[1, -2, 0, 0, 2, 0, 0, 2, 0, 2])
sage: Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 -2 0 0 ]
[ * 2 0 0 ]
[ * * 2 0 ]
[ * * * 2 ]
sage: Q.minkowski_reduction()
(
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 0 0 0 ]
[ * 1 0 0 ]
[ * * 2 0 ]
[ * * * 2 ]                                                        ,

[1 1 0 0]
[0 1 0 0]
[0 0 1 0]
[0 0 0 1]
)
sage: Q=QuadraticForm(ZZ,5,[2,2,0,0,0,2,2,0,0,2,2,0,2,2,2])
sage: Q.Gram_matrix()
[2 1 0 0 0]
[1 2 1 0 0]
[0 1 2 1 0]
[0 0 1 2 1]
[0 0 0 1 2]
sage: Q.minkowski_reduction()
Traceback (most recent call last):
...
NotImplementedError: This algorithm is only for dimensions less than 5
minkowski_reduction_for_4vars__SP()

Find a Minkowski-reduced form equivalent to the given one. This means that

Q(\(v_k\)) <= Q(\(s_1 * v_1 + ... + s_n * v_n\))

for all \(s_i\) where GCD(\(s_k, ... s_n\)) = 1.

Note: When Q has dim <= 4 we can take all \(s_i\) in {1, 0, -1}.

References:
Schulze-Pillot’s paper on “An algorithm for computing genera

of ternary and quaternary quadratic forms”, p138.

Donaldson’s 1979 paper “Minkowski Reduction of Integral

Matrices”, p203.

EXAMPLES:

sage: Q = QuadraticForm(ZZ,4,[30,17,11,12,29,25,62,64,25,110])
sage: Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 30 17 11 12 ]
[ * 29 25 62 ]
[ * * 64 25 ]
[ * * * 110 ]
sage: Q.minkowski_reduction_for_4vars__SP()
(
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 29 -17 25 4 ]
[ * 30 -11 5 ]
[ * * 64 0 ]
[ * * * 77 ]                                                       ,

[ 0  1  0  0]
[ 1  0  0 -1]
[ 0  0  1  0]
[ 0  0  0  1]
)
multiply_variable(c, i, in_place=False)

Replace the variables \(x_i\) by \(c*x_i\) in the quadratic form (replacing the original form if the in_place flag is True).

Here \(c\) must be an element of the base_ring defining the quadratic form.

INPUT:

\(c\) – an element of Q.base_ring()

\(i\) – an integer >= 0

OUTPUT:

a QuadraticForm (by default, otherwise none)

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,9,5,7])
sage: Q.multiply_variable(5,0)
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 25 0 0 0 ]
[ * 9 0 0 ]
[ * * 5 0 ]
[ * * * 7 ]
neighbor_iteration(seeds, p, mass=None, max_classes=1000, algorithm=None, max_neighbors=1000, verbose=False)

Return all classes in the \(p\)-neighbor graph of self.

Starting from the given seeds, this function successively finds p-neighbors until no new quadratic form (class) is obtained.

INPUT:

  • seeds – a list of quadratic forms in the same genus

  • p – a prime number

  • mass – (optional) a rational number; the mass of this genus

  • max_classes – (default: 1000) break the computation when max_classes are found

  • algorithm – (optional) one of ‘orbits’, ‘random’, ‘exaustion’

  • max_random_trys – (default: 1000) the maximum number of neighbors

    computed for a single lattice

OUTPUT:

  • a list of quadratic forms

EXAMPLES:

sage: from sage.quadratic_forms.quadratic_form__neighbors import neighbor_iteration
sage: Q = QuadraticForm(ZZ, 3, [1, 0, 0, 2, 1, 3])
sage: Q.det()
46
sage: mass = Q.conway_mass()
sage: g1 = neighbor_iteration([Q],3, mass=mass, algorithm = 'random') # long time
sage: g2 = neighbor_iteration([Q],3, algorithm = 'exaustion') # long time
sage: g3 = neighbor_iteration([Q],3, algorithm = 'orbits')
sage: mass == sum(1/q.number_of_automorphisms() for q in g1) # long time
True
sage: mass == sum(1/q.number_of_automorphisms() for q in g2) # long time
True
sage: mass == sum(1/q.number_of_automorphisms() for q in g3)
True
number_of_automorphisms()

Return the number of automorphisms (of det 1 and -1) of the quadratic form.

OUTPUT:

an integer >= 2.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, [1, 0, 0, 1, 0, 1], unsafe_initialization=True)
sage: Q.number_of_automorphisms()
48
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.number_of_automorphisms()
384
sage: 2^4 * factorial(4)
384
omega()

This is the content of the adjoint of the primitive associated quadratic form.

Ref: See Dickson’s “Studies in Number Theory”.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,37])
sage: Q.omega()
4
orbits_lines_mod_p(p)

Let \((L, q)\) be a lattice. This returns representatives of the orbits of lines in \(L/pL\) under the orthogonal group of \(q\).

INPUT:

  • p – a prime number

OUTPUT:

  • a list of vectors over GF(p)

EXAMPLES:

sage: from sage.quadratic_forms.quadratic_form__neighbors import orbits_lines_mod_p
sage: Q = QuadraticForm(ZZ, 3, [1, 0, 0, 2, 1, 3])
sage: Q.orbits_lines_mod_p(2)
[(0, 0, 1),
(0, 1, 0),
(0, 1, 1),
(1, 0, 0),
(1, 0, 1),
(1, 1, 0),
(1, 1, 1)]
parity(allow_rescaling_flag=True)

Return the parity (“even” or “odd”) of an integer-valued quadratic form over \(ZZ\), defined up to similitude/rescaling of the form so that its Jordan component of smallest scale is unimodular. After this rescaling, we say a form is even if it only represents even numbers, and odd if it represents some odd number.

If the ‘allow_rescaling_flag’ is set to False, then we require that the quadratic form have a Gram matrix with coefficients in \(ZZ\), and look at the unimodular Jordan block to determine its parity. This returns an error if the form is not integer-matrix, meaning that it has Jordan components at \(p=2\) which do not have an integer scale.

We determine the parity by looking for a 1x1 block in the 0-th Jordan component, after a possible rescaling.

INPUT:

self – a quadratic form with base_ring \(ZZ\), which we may

require to have integer Gram matrix.

OUTPUT:

One of the strings: “even” or “odd”

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 3, [4, -2, 0, 2, 3, 2]); Q
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 4 -2 0 ]
[ * 2 3 ]
[ * * 2 ]
sage: Q.parity()
'even'
sage: Q = QuadraticForm(ZZ, 3, [4, -2, 0, 2, 3, 1]); Q
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 4 -2 0 ]
[ * 2 3 ]
[ * * 1 ]
sage: Q.parity()
'even'
sage: Q = QuadraticForm(ZZ, 3, [4, -2, 0, 2, 2, 2]); Q
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 4 -2 0 ]
[ * 2 2 ]
[ * * 2 ]
sage: Q.parity()
'even'
sage: Q = QuadraticForm(ZZ, 3, [4, -2, 0, 2, 2, 1]); Q
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 4 -2 0 ]
[ * 2 2 ]
[ * * 1 ]
sage: Q.parity()
'odd'
polynomial(names='x')

Return the polynomial in ‘n’ variables of the quadratic form in the ring ‘R[names].’

INPUT:

-‘self’ - a quadratic form over a commutative ring. -‘names’ - the name of the variables. Digits will be appended to the name for each different canonical variable e.g x1, x2, x3 etc.

OUTPUT:

The polynomial form of the quadratic form.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(QQ,[1, 3, 5, 7])
sage: P = Q.polynomial(); P
x0^2 + 3*x1^2 + 5*x2^2 + 7*x3^2
sage: F.<a> = NumberField(x^2 - 5)
sage: Z = F.ring_of_integers()
sage: Q = QuadraticForm(Z,3,[2*a, 3*a, 0 , 1 - a, 0, 2*a + 4])
sage: P = Q.polynomial(names='y'); P
2*a*y0^2 + 3*a*y0*y1 + (-a + 1)*y1^2 + (2*a + 4)*y2^2
sage: Q = QuadraticForm(F,4,[a, 3*a, 0, 1 - a, a - 3, 0, 2*a + 4, 4 + a, 0, 1])
sage: Q.polynomial(names='z')
(a)*z0^2 + (3*a)*z0*z1 + (a - 3)*z1^2 + (a + 4)*z2^2 + (-a + 1)*z0*z3 + (2*a + 4)*z1*z3 + z3^2
sage: B.<i,j,k> = QuaternionAlgebra(F,-1,-1)
sage: Q = QuadraticForm(B, 3, [2*a, 3*a, i, 1 - a, 0, 2*a + 4])
sage: Q.polynomial()
Traceback (most recent call last):
...
ValueError: Can only create polynomial rings over commutative rings.
primitive()

Return a primitive version of an integer-valued quadratic form, defined over \(ZZ\).

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [2,3,4])
sage: Q.primitive()
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 2 3 ]
[ * 4 ]
sage: Q = QuadraticForm(ZZ, 2, [2,4,8])
sage: Q.primitive()
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 1 2 ]
[ * 4 ]
rational_diagonal_form(return_matrix=False)

Return a diagonal form equivalent to the given quadratic from over the fraction field of its defining ring.

INPUT:

  • return_matrix – (boolean, default: False) also return the transformation matrix

OUTPUT: either D (if return_matrix is false) or (D,T) (if return_matrix is true) where

  • D – the diagonalized form of this quadratic form.

  • T – transformation matrix. This is such that T.transpose() * self.matrix() * T gives D.matrix().

Both D and T are defined over the fraction field of the base ring of the given form.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [0,1,-1])
sage: Q
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 0 1 ]
[ * -1 ]
sage: Q.rational_diagonal_form()
Quadratic form in 2 variables over Rational Field with coefficients:
[ 1/4 0 ]
[ * -1 ]

If we start with a diagonal form, we get back the same form defined over the fraction field:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.rational_diagonal_form()
Quadratic form in 4 variables over Rational Field with coefficients:
[ 1 0 0 0 ]
[ * 3 0 0 ]
[ * * 5 0 ]
[ * * * 7 ]

In the following example, we check the consistency of the transformation matrix:

sage: Q = QuadraticForm(ZZ, 4, range(10))
sage: D, T = Q.rational_diagonal_form(return_matrix=True)
sage: D
Quadratic form in 4 variables over Rational Field with coefficients:
[ -1/16 0 0 0 ]
[ * 4 0 0 ]
[ * * 13 0 ]
[ * * * 563/52 ]
sage: T
[     1      0     11 149/26]
[  -1/8      1     -2 -10/13]
[     0      0      1 -29/26]
[     0      0      0      1]
sage: T.transpose() * Q.matrix() * T
[  -1/8      0      0      0]
[     0      8      0      0]
[     0      0     26      0]
[     0      0      0 563/26]
sage: D.matrix()
[  -1/8      0      0      0]
[     0      8      0      0]
[     0      0     26      0]
[     0      0      0 563/26]
sage: Q1 = QuadraticForm(ZZ, 4, [1, 1, 0, 0, 1, 0, 0, 1, 0, 18])
sage: Q1
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 1 0 0 ]
[ * 1 0 0 ]
[ * * 1 0 ]
[ * * * 18 ]
sage: Q1.rational_diagonal_form(return_matrix=True)
(
Quadratic form in 4 variables over Rational Field with coefficients:
[ 1 0 0 0 ]
[ * 3/4 0 0 ]
[ * * 1 0 ]
[ * * * 18 ]                                                         ,

[   1 -1/2    0    0]
[   0    1    0    0]
[   0    0    1    0]
[   0    0    0    1]
)

PARI returns a singular transformation matrix for this case:

sage: Q = QuadraticForm(QQ, 2, [1/2, 1, 1/2])
sage: Q.rational_diagonal_form()
Quadratic form in 2 variables over Rational Field with coefficients:
[ 1/2 0 ]
[ * 0 ]

This example cannot be computed by PARI:

sage: Q = QuadraticForm(RIF, 4, range(10))
sage: Q.__pari__()
Traceback (most recent call last):
...
TypeError
sage: Q.rational_diagonal_form()
Quadratic form in 4 variables over Real Interval Field with 53 bits of precision with coefficients:
[ 5 0.?e-14 0.?e-13 0.?e-13 ]
[ * -0.05000000000000? 0.?e-12 0.?e-12 ]
[ * * 13.00000000000? 0.?e-10 ]
[ * * * 10.8269230769? ]
reciprocal()

This gives the reciprocal quadratic form associated to the given form. This is defined as the multiple of the primitive adjoint with the same content as the given form.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,37])
sage: Q.reciprocal()
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 37 0 0 ]
[ * 37 0 ]
[ * * 1 ]
sage: Q.reciprocal().reciprocal()
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 1 0 0 ]
[ * 1 0 ]
[ * * 37 ]
sage: Q.reciprocal().reciprocal() == Q
True
reduced_binary_form()

Find a form which is reduced in the sense that no further binary form reductions can be done to reduce the original form.

EXAMPLES:

sage: QuadraticForm(ZZ,2,[5,5,2]).reduced_binary_form()
(
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 2 -1 ]
[ * 2 ]                                                            ,

[ 0 -1]
[ 1  1]
)
reduced_binary_form1()

Reduce the form \(ax^2 + bxy+cy^2\) to satisfy the reduced condition \(|b| \le a \le c\), with \(b \ge 0\) if \(a = c\). This reduction occurs within the proper class, so all transformations are taken to have determinant 1.

EXAMPLES:

sage: QuadraticForm(ZZ,2,[5,5,2]).reduced_binary_form1()
(
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 2 -1 ]
[ * 2 ]                                                            ,

[ 0 -1]
[ 1  1]
)
reduced_ternary_form__Dickson()

Find the unique reduced ternary form according to the conditions of Dickson’s “Studies in the Theory of Numbers”, pp164-171.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1, 1, 1])
sage: Q.reduced_ternary_form__Dickson()
Traceback (most recent call last):
...
NotImplementedError: TO DO
representation_number_list(B)

Return the vector of representation numbers < B.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ,[1,1,1,1,1,1,1,1])
sage: Q.representation_number_list(10)
[1, 16, 112, 448, 1136, 2016, 3136, 5504, 9328, 12112]
representation_vector_list(B, maxvectors=100000000)

Find all vectors \(v\) where \(Q(v) < B\).

This only works for positive definite quadratic forms.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1, 1])
sage: Q.representation_vector_list(10)
[[(0, 0)],
 [(0, 1), (0, -1), (1, 0), (-1, 0)],
 [(1, 1), (-1, -1), (1, -1), (-1, 1)],
 [],
 [(0, 2), (0, -2), (2, 0), (-2, 0)],
 [(1, 2), (-1, -2), (1, -2), (-1, 2), (2, 1), (-2, -1), (2, -1), (-2, 1)],
 [],
 [],
 [(2, 2), (-2, -2), (2, -2), (-2, 2)],
 [(0, 3), (0, -3), (3, 0), (-3, 0)]]
sage: list(map(len, _))
[1, 4, 4, 0, 4, 8, 0, 0, 4, 4]
sage: Q.representation_number_list(10)
[1, 4, 4, 0, 4, 8, 0, 0, 4, 4]
scale_by_factor(c, change_value_ring_flag=False)

Scale the values of the quadratic form by the number \(c\), if this is possible while still being defined over its base ring.

If the flag is set to true, then this will alter the value ring to be the field of fractions of the original ring (if necessary).

INPUT:

\(c\) – a scalar in the fraction field of the value ring of the form.

OUTPUT:

A quadratic form of the same dimension

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [3,9,18,27])
sage: Q.scale_by_factor(3)
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 9 0 0 0 ]
[ * 27 0 0 ]
[ * * 54 0 ]
[ * * * 81 ]

sage: Q.scale_by_factor(1/3)
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 0 0 0 ]
[ * 3 0 0 ]
[ * * 6 0 ]
[ * * * 9 ]
set_number_of_automorphisms(num_autos)

Set the number of automorphisms to be the value given. No error checking is performed, to this may lead to erroneous results.

The fact that this result was set externally is recorded in the internal list of external initializations, accessible by the method list_external_initializations().

OUTPUT: None

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1, 1, 1])
sage: Q.list_external_initializations()
[]
sage: Q.set_number_of_automorphisms(-3)
sage: Q.number_of_automorphisms()
-3
sage: Q.list_external_initializations()
['number_of_automorphisms']
shimura_mass__maximal()

Use Shimura’s exact mass formula to compute the mass of a maximal quadratic lattice. This works for any totally real number field, but has a small technical restriction when \(n\) is odd.

INPUT:

none

OUTPUT:

a rational number

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1])
sage: Q.shimura_mass__maximal()
short_primitive_vector_list_up_to_length(len_bound, up_to_sign_flag=False)

Return a list of lists of short primitive vectors \(v\), sorted by length, with Q(\(v\)) < len_bound. The list in output \([i]\) indexes all vectors of length \(i\). If the up_to_sign_flag is set to True, then only one of the vectors of the pair \([v, -v]\) is listed.

Note

This processes the PARI/GP output to always give elements of type \(\ZZ\).

OUTPUT: a list of lists of vectors.

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.short_vector_list_up_to_length(5, True)
[[(0, 0, 0, 0)],
 [(1, 0, 0, 0)],
 [],
 [(0, 1, 0, 0)],
 [(1, 1, 0, 0), (1, -1, 0, 0), (2, 0, 0, 0)]]
sage: Q.short_primitive_vector_list_up_to_length(5, True)
[[], [(1, 0, 0, 0)], [], [(0, 1, 0, 0)], [(1, 1, 0, 0), (1, -1, 0, 0)]]
short_vector_list_up_to_length(len_bound, up_to_sign_flag=False)

Return a list of lists of short vectors \(v\), sorted by length, with Q(\(v\)) < len_bound.

INPUT:

  • len_bound – bound for the length of the vectors.

  • up_to_sign_flag – (default: False) if set to True, then only one of the vectors of the pair \([v, -v]\) is listed.

OUTPUT:

A list of lists of vectors such that entry \([i]\) contains all vectors of length \(i\).

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.short_vector_list_up_to_length(3)
[[(0, 0, 0, 0)], [(1, 0, 0, 0), (-1, 0, 0, 0)], []]
sage: Q.short_vector_list_up_to_length(4)
[[(0, 0, 0, 0)],
 [(1, 0, 0, 0), (-1, 0, 0, 0)],
 [],
 [(0, 1, 0, 0), (0, -1, 0, 0)]]
sage: Q.short_vector_list_up_to_length(5)
[[(0, 0, 0, 0)],
 [(1, 0, 0, 0), (-1, 0, 0, 0)],
 [],
 [(0, 1, 0, 0), (0, -1, 0, 0)],
 [(1, 1, 0, 0),
  (-1, -1, 0, 0),
  (1, -1, 0, 0),
  (-1, 1, 0, 0),
  (2, 0, 0, 0),
  (-2, 0, 0, 0)]]
sage: Q.short_vector_list_up_to_length(5, True)
[[(0, 0, 0, 0)],
 [(1, 0, 0, 0)],
 [],
 [(0, 1, 0, 0)],
 [(1, 1, 0, 0), (1, -1, 0, 0), (2, 0, 0, 0)]]
sage: Q = QuadraticForm(matrix(6, [2, 1, 1, 1, -1, -1, 1, 2, 1, 1, -1, -1, 1, 1, 2, 0, -1, -1, 1, 1, 0, 2, 0, -1, -1, -1, -1, 0, 2, 1, -1, -1, -1, -1, 1, 2]))
sage: vs = Q.short_vector_list_up_to_length(8)
sage: [len(vs[i]) for i in range(len(vs))]
[1, 72, 270, 720, 936, 2160, 2214, 3600]
sage: vs = Q.short_vector_list_up_to_length(30)  # long time (28s on sage.math, 2014)
sage: [len(vs[i]) for i in range(len(vs))]       # long time
[1, 72, 270, 720, 936, 2160, 2214, 3600, 4590, 6552, 5184, 10800, 9360, 12240, 13500, 17712, 14760, 25920, 19710, 26064, 28080, 36000, 25920, 47520, 37638, 43272, 45900, 59040, 46800, 75600]

The cases of len_bound < 2 led to exception or infinite runtime before.

sage: Q.short_vector_list_up_to_length(-1)
[]
sage: Q.short_vector_list_up_to_length(0)
[]
sage: Q.short_vector_list_up_to_length(1)
[[(0, 0, 0, 0, 0, 0)]]

In the case of quadratic forms that are not positive definite an error is raised.

sage: QuadraticForm(matrix(2, [2, 0, 0, -2])).short_vector_list_up_to_length(3)
Traceback (most recent call last):
...
ValueError: Quadratic form must be positive definite in order to enumerate short vectors

Check that PARI does not return vectors which are too long:

sage: Q = QuadraticForm(matrix(2, [72, 12, 12, 120]))
sage: len_bound_pari = 2*22953421 - 2; len_bound_pari
45906840
sage: vs = list(Q.__pari__().qfminim(len_bound_pari)[2])  # long time (18s on sage.math, 2014)
sage: v = vs[0]; v  # long time
[66, -623]~
sage: v.Vec() * Q.__pari__() * v  # long time
45902280
siegel_product(u)

Computes the infinite product of local densities of the quadratic form for the number \(u\).

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Q.theta_series(11)
1 + 8*q + 24*q^2 + 32*q^3 + 24*q^4 + 48*q^5 + 96*q^6 + 64*q^7 + 24*q^8 + 104*q^9 + 144*q^10 + O(q^11)

sage: Q.siegel_product(1)
8
sage: Q.siegel_product(2)      ## This one is wrong -- expect 24, and the higher powers of 2 don't work... =(
24
sage: Q.siegel_product(3)
32
sage: Q.siegel_product(5)
48
sage: Q.siegel_product(6)
96
sage: Q.siegel_product(7)
64
sage: Q.siegel_product(9)
104

sage: Q.local_density(2,1)
1
sage: M = 4; len([v  for v in mrange([M,M,M,M])  if Q(v) % M == 1]) / M^3
1
sage: M = 16; len([v  for v in mrange([M,M,M,M])  if Q(v) % M == 1]) / M^3  # long time (2s on sage.math, 2014)
1

sage: Q.local_density(2,2)
3/2
sage: M = 4; len([v  for v in mrange([M,M,M,M])  if Q(v) % M == 2]) / M^3
3/2
sage: M = 16; len([v  for v in mrange([M,M,M,M])  if Q(v) % M == 2]) / M^3  # long time (2s on sage.math, 2014)
3/2
signature()

Returns the signature of the quadratic form, defined as:

number of positive eigenvalues - number of negative eigenvalues

of the matrix of the quadratic form.

INPUT:

None

OUTPUT:

an integer

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,0,0,-4,3,11,3])
sage: Q.signature()
3
sage: Q = DiagonalQuadraticForm(ZZ, [1,2,-3,-4])
sage: Q.signature()
0
sage: Q = QuadraticForm(ZZ, 4, range(10)); Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 0 1 2 3 ]
[ * 4 5 6 ]
[ * * 7 8 ]
[ * * * 9 ]
sage: Q.signature()
2
signature_vector()

Returns the triple \((p, n, z)\) of integers where

  • \(p\) = number of positive eigenvalues

  • \(n\) = number of negative eigenvalues

  • \(z\) = number of zero eigenvalues

for the symmetric matrix associated to Q.

INPUT:

(none)

OUTPUT:

a triple of integers >= 0

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,0,0,-4])
sage: Q.signature_vector()
(1, 1, 2)
sage: Q = DiagonalQuadraticForm(ZZ, [1,2,-3,-4])
sage: Q.signature_vector()
(2, 2, 0)
sage: Q = QuadraticForm(ZZ, 4, range(10)); Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 0 1 2 3 ]
[ * 4 5 6 ]
[ * * 7 8 ]
[ * * * 9 ]
sage: Q.signature_vector()
(3, 1, 0)
solve(c=0)

Return a vector \(x\) such that self(x) == c.

INPUT:

  • c – (default: 0) a rational number.

OUTPUT:

  • A non-zero vector \(x\) satisfying self(x) == c.

ALGORITHM:

Uses PARI’s qfsolve(). Algorithm described by Jeroen Demeyer; see comments on trac ticket #19112

EXAMPLES:

sage: F = DiagonalQuadraticForm(QQ, [1, -1]); F
Quadratic form in 2 variables over Rational Field with coefficients:
[ 1 0 ]
[ * -1 ]
sage: F.solve()
(1, 1)
sage: F.solve(1)
(1, 0)
sage: F.solve(2)
(3/2, -1/2)
sage: F.solve(3)
(2, -1)
sage: F = DiagonalQuadraticForm(QQ, [1, 1, 1, 1])
sage: F.solve(7)
(1, 2, -1, -1)
sage: F.solve()
Traceback (most recent call last):
...
ArithmeticError: no solution found (local obstruction at -1)
sage: Q = QuadraticForm(QQ, 2, [17, 94, 130])
sage: x = Q.solve(5); x
(17, -6)
sage: Q(x)
5

sage: Q.solve(6)
Traceback (most recent call last):
...
ArithmeticError: no solution found (local obstruction at 3)

sage: G = DiagonalQuadraticForm(QQ, [5, -3, -2])
sage: x = G.solve(10); x
(3/2, -1/2, 1/2)
sage: G(x)
10

sage: F = DiagonalQuadraticForm(QQ, [1, -4])
sage: x = F.solve(); x
(2, 1)
sage: F(x)
0
sage: F = QuadraticForm(QQ, 4, [0, 0, 1, 0, 0, 0, 1, 0, 0, 0]); F
Quadratic form in 4 variables over Rational Field with coefficients:
[ 0 0 1 0 ]
[ * 0 0 1 ]
[ * * 0 0 ]
[ * * * 0 ]
sage: F.solve(23)
(23, 0, 1, 0)

Other fields besides the rationals are currently not supported:

sage: F = DiagonalQuadraticForm(GF(11), [1, 1])
sage: F.solve()
Traceback (most recent call last):
...
TypeError: solving quadratic forms is only implemented over QQ
split_local_cover()

Tries to find subform of the given (positive definite quaternary) quadratic form Q of the form

\[d*x^2 + T(y,z,w)\]

where \(d > 0\) is as small as possible.

This is done by exhaustive search on small vectors, and then comparing the local conditions of its sum with it’s complementary lattice and the original quadratic form Q.

INPUT:

none

OUTPUT:

a QuadraticForm over ZZ

EXAMPLES:

sage: Q1 = DiagonalQuadraticForm(ZZ, [7,5,3])
sage: Q1.split_local_cover()
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 3 0 0 ]
[ * 7 0 ]
[ * * 5 ]
sum_by_coefficients_with(right)

Return the sum (on coefficients) of two quadratic forms of the same size.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,4,10])
sage: Q
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 1 4 ]
[ * 10 ]
sage: Q+Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 4 0 0 ]
[ * 10 0 0 ]
[ * * 1 4 ]
[ * * * 10 ]

sage: Q2 = QuadraticForm(ZZ, 2, [1,4,-10])
sage: Q.sum_by_coefficients_with(Q2)
Quadratic form in 2 variables over Integer Ring with coefficients:
[ 2 8 ]
[ * 0 ]
swap_variables(r, s, in_place=False)

Switch the variables \(x_r\) and \(x_s\) in the quadratic form (replacing the original form if the in_place flag is True).

INPUT:

\(r\), \(s\) – integers >= 0

OUTPUT:

a QuadraticForm (by default, otherwise none)

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 4, range(1,11))
sage: Q
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 2 3 4 ]
[ * 5 6 7 ]
[ * * 8 9 ]
[ * * * 10 ]


sage: Q.swap_variables(0,2)
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 8 6 3 9 ]
[ * 5 2 7 ]
[ * * 1 4 ]
[ * * * 10 ]


sage: Q.swap_variables(0,2).swap_variables(0,2)
Quadratic form in 4 variables over Integer Ring with coefficients:
[ 1 2 3 4 ]
[ * 5 6 7 ]
[ * * 8 9 ]
[ * * * 10 ]
theta_by_cholesky(q_prec)

Uses the real Cholesky decomposition to compute (the \(q\)-expansion of) the theta function of the quadratic form as a power series in \(q\) with terms correct up to the power \(q^{\text{q\_prec}}\). (So its error is \(O(q^ {\text{q\_prec} + 1})\).)

REFERENCE:

From Cohen’s “A Course in Computational Algebraic Number Theory” book, p 102.

EXAMPLES:

## Check the sum of 4 squares form against Jacobi's formula
sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Theta = Q.theta_by_cholesky(10)
sage: Theta
1 + 8*q + 24*q^2 + 32*q^3 + 24*q^4 + 48*q^5 + 96*q^6 + 64*q^7 + 24*q^8 + 104*q^9 + 144*q^10
sage: Expected =  [1] + [8*sum([d for d in divisors(n) if d%4 != 0])  for n in range(1,11)]
sage: Expected
[1, 8, 24, 32, 24, 48, 96, 64, 24, 104, 144]
sage: Theta.list() == Expected
True
## Check the form x^2 + 3y^2 + 5z^2 + 7w^2 represents everything except 2 and 22.
sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Theta = Q.theta_by_cholesky(50)
sage: Theta_list = Theta.list()
sage: [m  for m in range(len(Theta_list))  if Theta_list[m] == 0]
[2, 22]
theta_by_pari(Max, var_str='q', safe_flag=True)

Use PARI/GP to compute the theta function as a power series (or vector) up to the precision \(O(q^Max)\). This also caches the result for future computations.

If var_str = ‘’, then we return a vector \(v\) where \(v[i]\) counts the number of vectors of length \(i\).

The safe_flag allows us to select whether we want a copy of the output, or the original output. It is only meaningful when a vector is returned, otherwise a copy is automatically made in creating the power series. By default safe_flag = True, so we return a copy of the cached information. If this is set to False, then the routine is much faster but the return values are vulnerable to being corrupted by the user.

INPUT:

Max – an integer >=0 var_str – a string

OUTPUT:

a power series or a vector

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1,1,1])
sage: Prec = 100
sage: compute = Q.theta_by_pari(Prec, '')
sage: exact = [1] + [8 * sum([d  for d in divisors(i)  if d % 4 != 0])  for i in range(1, Prec)]
sage: compute == exact
True
theta_series(Max=10, var_str='q', safe_flag=True)

Compute the theta series as a power series in the variable given in var_str (which defaults to ‘\(q\)’), up to the specified precision \(O(q^max)\).

This uses the PARI/GP function qfrep, wrapped by the theta_by_pari() method. This caches the result for future computations.

The safe_flag allows us to select whether we want a copy of the output, or the original output. It is only meaningful when a vector is returned, otherwise a copy is automatically made in creating the power series. By default safe_flag = True, so we return a copy of the cached information. If this is set to False, then the routine is much faster but the return values are vulnerable to being corrupted by the user.

TO DO: Allow the option Max=’mod_form’ to give enough coefficients to ensure we determine the theta series as a modular form. This is related to the Sturm bound, but we’ll need to be careful about this (particularly for half-integral weights!).

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q.theta_series()
1 + 2*q + 2*q^3 + 6*q^4 + 2*q^5 + 4*q^6 + 6*q^7 + 8*q^8 + 14*q^9 + O(q^10)

sage: Q.theta_series(25)
1 + 2*q + 2*q^3 + 6*q^4 + 2*q^5 + 4*q^6 + 6*q^7 + 8*q^8 + 14*q^9 + 4*q^10 + 12*q^11 + 18*q^12 + 12*q^13 + 12*q^14 + 8*q^15 + 34*q^16 + 12*q^17 + 8*q^18 + 32*q^19 + 10*q^20 + 28*q^21 + 16*q^23 + 44*q^24 + O(q^25)
theta_series_degree_2(Q, prec)

Compute the theta series of degree 2 for the quadratic form Q.

INPUT:

  • prec – an integer.

OUTPUT:

dictionary, where:

  • keys are \({\rm GL}_2(\ZZ)\)-reduced binary quadratic forms (given as triples of coefficients)

  • values are coefficients

EXAMPLES:

sage: Q2 = QuadraticForm(ZZ, 4, [1,1,1,1, 1,0,0, 1,0, 1])
sage: S = Q2.theta_series_degree_2(10)
sage: S[(0,0,2)]
24
sage: S[(1,0,1)]
144
sage: S[(1,1,1)]
192

AUTHORS:

  • Gonzalo Tornaria (2010-03-23)

REFERENCE:

  • Raum, Ryan, Skoruppa, Tornaria, ‘On Formal Siegel Modular Forms’ (preprint)

vectors_by_length(bound)

Returns a list of short vectors together with their values.

This is a naive algorithm which uses the Cholesky decomposition, but does not use the LLL-reduction algorithm.

INPUT:

bound – an integer >= 0

OUTPUT:

A list L of length (bound + 1) whose entry L \([i]\) is a list of all vectors of length \(i\).

Reference: This is a slightly modified version of Cohn’s Algorithm 2.7.5 in “A Course in Computational Number Theory”, with the increment step moved around and slightly re-indexed to allow clean looping.

Note: We could speed this up for very skew matrices by using LLL first, and then changing coordinates back, but for our purposes the simpler method is efficient enough. =)

EXAMPLES:

sage: Q = DiagonalQuadraticForm(ZZ, [1,1])
sage: Q.vectors_by_length(5)
[[[0, 0]],
 [[0, -1], [-1, 0]],
 [[-1, -1], [1, -1]],
 [],
 [[0, -2], [-2, 0]],
 [[-1, -2], [1, -2], [-2, -1], [2, -1]]]
sage: Q1 = DiagonalQuadraticForm(ZZ, [1,3,5,7])
sage: Q1.vectors_by_length(5)
[[[0, 0, 0, 0]],
 [[-1, 0, 0, 0]],
 [],
 [[0, -1, 0, 0]],
 [[-1, -1, 0, 0], [1, -1, 0, 0], [-2, 0, 0, 0]],
 [[0, 0, -1, 0]]]
sage: Q = QuadraticForm(ZZ, 4, [1,1,1,1, 1,0,0, 1,0, 1])
sage: list(map(len, Q.vectors_by_length(2)))
[1, 12, 12]
sage: Q = QuadraticForm(ZZ, 4, [1,-1,-1,-1, 1,0,0, 4,-3, 4])
sage: list(map(len, Q.vectors_by_length(3)))
[1, 3, 0, 3]
xi(p)

Return the value of the genus characters Xi_p… which may be missing one character. We allow -1 as a prime.

Reference: Dickson’s “Studies in the Theory of Numbers”

EXAMPLES:

sage: Q1 = QuadraticForm(ZZ, 3, [1, 1, 1, 14, 3, 14])
sage: Q2 = QuadraticForm(ZZ, 3, [2, -1, 0, 2, 0, 50])
sage: [Q1.omega(), Q2.omega()]
[5, 5]
sage: [Q1.hasse_invariant(5), Q2.hasse_invariant(5)]    # equivalent over Q_5
[1, 1]
sage: [Q1.xi(5), Q2.xi(5)]                              # not equivalent over Z_5
[1, -1]
xi_rec(p)

Return Xi(\(p\)) for the reciprocal form.

EXAMPLES:

sage: Q1 = QuadraticForm(ZZ, 3, [1, 1, 1, 14, 3, 14])
sage: Q2 = QuadraticForm(ZZ, 3, [2, -1, 0, 2, 0, 50])
sage: [Q1.clifford_conductor(), Q2.clifford_conductor()]   # equivalent over Q
[3, 3]
sage: Q1.is_locally_equivalent_to(Q2)                # not in the same genus
False
sage: [Q1.delta(), Q2.delta()]
[480, 480]
sage: factor(480)
2^5 * 3 * 5
sage: list(map(Q1.xi_rec, [-1,2,3,5]))
[-1, -1, -1, 1]
sage: list(map(Q2.xi_rec, [-1,2,3,5]))
[-1, -1, -1, -1]
sage.quadratic_forms.quadratic_form.QuadraticForm__constructor(R, n=None, entries=None)

Wrapper for the QuadraticForm class constructor. This is meant for internal use within the QuadraticForm class code only. You should instead directly call QuadraticForm().

EXAMPLES:

sage: from sage.quadratic_forms.quadratic_form import QuadraticForm__constructor
sage: QuadraticForm__constructor(ZZ, 3)   # Makes a generic quadratic form over the integers
Quadratic form in 3 variables over Integer Ring with coefficients:
[ 0 0 0 ]
[ * 0 0 ]
[ * * 0 ]
sage.quadratic_forms.quadratic_form.is_QuadraticForm(Q)

Determine if the object Q is an element of the QuadraticForm class.

EXAMPLES:

sage: Q = QuadraticForm(ZZ, 2, [1,2,3])
sage: from sage.quadratic_forms.quadratic_form import is_QuadraticForm
sage: is_QuadraticForm(Q)  ##random
True
sage: is_QuadraticForm(2)  ##random
False
sage.quadratic_forms.quadratic_form.quadratic_form_from_invariants(F, rk, det, P, sminus)

Return a rational quadratic form with given invariants.

INPUT:

  • F – the base field; currently only QQ is allowed

  • rk – integer; the rank

  • det – rational; the determinant

  • P – a list of primes where Cassel’s Hasse invariant is negative

  • sminus – integer; the number of negative eigenvalues of any Gram matrix

OUTPUT:

  • a quadratic form with the specified invariants

Let \((a_1, \ldots, a_n)\) be the gram marix of a regular quadratic space. Then Cassel’s Hasse invariant is defined as

\[\prod_{i<j} (a_i,a_j),\]

where \((a_i,a_j)\) denotes the Hilbert symbol.

ALGORITHM:

We follow [Kir2016].

EXAMPLES:

sage: P = [3,5]
sage: q = quadratic_form_from_invariants(QQ,2,-15,P,1)
sage: q
Quadratic form in 2 variables over Rational Field with coefficients:
[ 5 0 ]
[ * -3 ]
sage: all(q.hasse_invariant(p)==-1 for p in P)
True