Spaces

Creation of spaces

quadratic_spaceMethod
quadratic_space(K::NumField, n::Int; cached::Bool = true) -> QuadSpace

Create the quadratic space over K with dimension n and Gram matrix equals to the identity matrix.

source
hermitian_spaceMethod
hermitian_space(E::NumField, n::Int; cached::Bool = true) -> HermSpace

Create the hermitian space over E with dimension n and Gram matrix equals to the identity matrix. The number field E must be a quadratic extension, that is, $degree(E) == 2$ must hold.

source
quadratic_spaceMethod
quadratic_space(K::NumField, G::MatElem; cached::Bool = true) -> QuadSpace

Create the quadratic space over K with Gram matrix G. The matrix G must be square and symmetric.

source
hermitian_spaceMethod
hermitian_space(E::NumField, gram::MatElem; cached::Bool = true) -> HermSpace

Create the hermitian space over E with Gram matrix equals to gram. The matrix gram must be square and hermitian with respect to the non-trivial automorphism of E. The number field E must be a quadratic extension, that is, $degree(E) == 2$ must hold.

source

Examples

Here are easy examples to see how these constructors work. We will keep the two following spaces for the rest of this section:


julia> K, a = cyclotomic_real_subfield(7);
julia> Kt, t = K["t"];
julia> E, b = number_field(t^2-a*t+1, "b");
julia> Q = quadratic_space(K, K[0 1; 1 0])Quadratic space of dimension 2 over maximal real subfield of cyclotomic field of order 7 with gram matrix [0 1] [1 0]
julia> H = hermitian_space(E, 3)Hermitian space of dimension 3 over relative number field with defining polynomial t^2 - (z_7 + 1/z_7)*t + 1 over number field with defining polynomial $^3 + $^2 - 2*$ - 1 over rational field with gram matrix [1 0 0] [0 1 0] [0 0 1]

Attributes

Let $(V, \Phi)$ be a space over $E/K$. We define its dimension to be its dimension as a vector space over its base ring $E$ and its rank to be the rank of its Gram matrix. If these two invariants agree, the space is said to be regular.

While dealing with lattices, one always works with regular ambient spaces.

The determinant $\text{det}(V, \Phi)$ of $(V, \Phi)$ is defined to be the class of the determinant of its Gram matrix in $K^{\times}/N(E^{\times})$ (which is similar to $K^{\times}/(K^{\times})^2$ in the quadratic case). The discriminant $\text{disc}(V, \Phi)$ of $(V, \Phi)$ is defined to be $(-1)^{(m(m-1)/2)}\text{det}(V, \Phi)$, where $m$ is the rank of $(V, \Phi)$.

rankMethod
rank(V::AbstractSpace) -> Int

Return the rank of the space V.

source
dimMethod
dim(V::AbstractSpace) -> Int

Return the dimension of the space V.

source
gram_matrixMethod
gram_matrix(V::AbstractSpace) -> MatElem

Return the Gram matrix of the space V.

source
involutionMethod
involution(V::AbstractSpace) -> NumFieldHom

Return the involution of the space V.

source
base_ringMethod
base_ring(V::AbstractSpace) -> NumField

Return the algebra over which the space V is defined.

source
fixed_fieldMethod
fixed_field(V::AbstractSpace) -> NumField

Return the fixed field of the space V.

source
detMethod
det(V::AbstractSpace) -> FieldElem

Return the determinant of the space V as an element of its fixed field.

source
discriminantMethod
discriminant(V::AbstractSpace) -> FieldElem

Return the discriminant of the space V as an element of its fixed field.

source

Examples

So for instance, one could get the following information about the hermitian space $H$:


julia> K, a = cyclotomic_real_subfield(7);
julia> Kt, t = K["t"];
julia> E, b = number_field(t^2-a*t+1, "b");
julia> H = hermitian_space(E, 3);
julia> rank(H), dim(H)(3, 3)
julia> gram_matrix(H)[1 0 0] [0 1 0] [0 0 1]
julia> involution(H)Map from relative number field of degree 2 over maximal real subfield of cyclotomic field of order 7 to relative number field of degree 2 over maximal real subfield of cyclotomic field of order 7
julia> base_ring(H)Relative number field with defining polynomial t^2 - (z_7 + 1/z_7)*t + 1 over number field with defining polynomial $^3 + $^2 - 2*$ - 1 over rational field
julia> fixed_field(H)Number field with defining polynomial $^3 + $^2 - 2*$ - 1 over rational field
julia> det(H), discriminant(H)(1, -1)

Predicates

Let $(V, \Phi)$ be a hermitian space over $E/K$ (resp. quadratic space $K$). We say that $(V, \Phi)$ is definite if $E/K$ is CM (resp. $K$ is totally real) and if there exists an orthogonal basis of $V$ for which the diagonal elements of the associated Gram matrix of $(V, \Phi)$ are either all totally positive or all totally negative. In the former case, $V$ is said to be positive definite, while in the latter case it is negative definite. In all the other cases, we say that $V$ is indefinite.

is_regularMethod
is_regular(V::AbstractSpace) -> Bool

Return whether the space V is regular, that is, if the Gram matrix has full rank.

source
is_quadraticMethod
is_quadratic(V::AbstractSpace) -> Bool

Return whether the space V is quadratic.

source
is_hermitianMethod
is_hermitian(V::AbstractSpace) -> Bool

Return whether the space V is hermitian.

source
is_positive_definiteMethod
is_positive_definite(V::AbstractSpace) -> Bool

Return whether the space V is positive definite.

source
is_negative_definiteMethod
is_negative_definite(V::AbstractSpace) -> Bool

Return whether the space V is negative definite.

source
is_definiteMethod
is_definite(V::AbstractSpace) -> Bool

Return whether the space V is definite.

source

Note that the is_hermitian function tests whether the space is non-quadratic.

Examples


julia> K, a = cyclotomic_real_subfield(7);
julia> Kt, t = K["t"];
julia> E, b = number_field(t^2-a*t+1, "b");
julia> Q = quadratic_space(K, K[0 1; 1 0]);
julia> H = hermitian_space(E, 3);
julia> is_regular(Q), is_regular(H)(true, true)
julia> is_quadratic(Q), is_hermitian(H)(true, true)
julia> is_definite(Q), is_positive_definite(H)(false, true)

Inner products and diagonalization

gram_matrixMethod
gram_matrix(V::AbstractSpace, M::MatElem) -> MatElem

Return the Gram matrix of the rows of M with respect to the Gram matrix of the space V.

source
gram_matrixMethod
gram_matrix(V::AbstractSpace, S::Vector{Vector}) -> MatElem

Return the Gram matrix of the sequence S with respect to the Gram matrix of the space V.

source
inner_productMethod
inner_product(V::AbstractSpace, v::Vector, w::Vector) -> FieldElem

Return the inner product of v and w with respect to the bilinear form of the space V.

source
orthogonal_basisMethod
orthogonal_basis(V::AbstractSpace) -> MatElem

Return a matrix M, such that the rows of M form an orthogonal basis of the space V.

source
diagonalMethod
diagonal(V::AbstractSpace) -> Vector{FieldElem}

Return a vector of elements $a_1,\dotsc,a_n$ such that the space V is isometric to the diagonal space $\langle a_1,\dotsc,a_n \rangle$.

The elements are contained in the fixed field of V.

source
diagonal_with_transformMethod
diagonal_with_transform(V::AbstractSpace) -> Vector{FieldElem},
                                                         MatElem{FieldElem}

Return a vector of elements $a_1,\dotsc,a_n$ such that the space V is isometric to the diagonal space $\langle a_1,\dotsc,a_n \rangle$. The second output is a matrix U whose rows span an orthogonal basis of V for which the Gram matrix is given by the diagonal matrix of the $a_i$'s.

The elements are contained in the fixed field of V.

source
restrict_scalarsMethod
restrict_scalars(V::AbstractSpace, K::QQField,
                                   alpha::FieldElem = one(base_ring(V)))
                                                -> QuadSpace, AbstractSpaceRes

Given a space $(V, \Phi)$ and a subfield K of the base algebra E of V, return the quadratic space W obtained by restricting the scalars of $(V, \alpha\Phi)$ to K, together with the map f for extending the scalars back. The form on the restriction is given by $Tr \circ \Phi$ where $Tr: E \to K$ is the trace form. The rescaling factor $\alpha$ is set to 1 by default.

Note that for now one can only restrict scalars to $\mathbb Q$.

source

Examples


julia> K, a = cyclotomic_real_subfield(7);
julia> Kt, t = K["t"];
julia> E, b = number_field(t^2-a*t+1, "b");
julia> Q = quadratic_space(K, K[0 1; 1 0]);
julia> H = hermitian_space(E, 3);
julia> gram_matrix(Q, K[1 1; 2 0])[2 2] [2 0]
julia> gram_matrix(H, E[1 0 0; 0 1 0; 0 0 1])[1 0 0] [0 1 0] [0 0 1]
julia> inner_product(Q, K[1 1], K[0 2])[2]
julia> orthogonal_basis(H)[1 0 0] [0 1 0] [0 0 1]
julia> diagonal(Q), diagonal(H)(AbsSimpleNumFieldElem[1, -1], AbsSimpleNumFieldElem[1, 1, 1])

Equivalence

Let $(V, \Phi)$ and $(V', \Phi')$ be spaces over the same extension $E/K$. A homomorphism of spaces from $V$ to $V'$ is a $E$-linear mapping $f \colon V \to V'$ such that for all $x,y \in V$, one has

\[ \Phi'(f(x), f(y)) = \Phi(x,y).\]

An automorphism of spaces is called an isometry and a monomorphism is called an embedding.

hasse_invariantMethod
hasse_invariant(V::QuadSpace, p::Union{InfPlc, AbsNumFieldOrderIdeal{AbsSimpleNumField, AbsSimpleNumFieldElem}}) -> Int

Returns the Hasse invariant of the quadratic space V at p. This is equal to the product of local Hilbert symbols $(a_i, a_j)_p$, $i < j$, where $V$ is isometric to $\langle a_1, \dotsc, a_n\rangle$. If V is degenerate return the hasse invariant of V/radical(V).

source
witt_invariantMethod
witt_invariant(V::QuadSpace, p::Union{InfPlc, AbsNumFieldOrderIdeal{AbsSimpleNumField, AbsSimpleNumFieldElem}}) -> Int

Returns the Witt invariant of the quadratic space V at p.

See [Definition 3.2.1, Kir16].

source
is_isometricMethod
is_isometric(L::AbstractSpace, M::AbstractSpace) -> Bool

Return whether the spaces L and M are isometric.

source
is_isometricMethod
is_isometric(L::AbstractSpace, M::AbstractSpace, p::Union{InfPlc, AbsNumFieldOrderIdeal{AbsSimpleNumField, AbsSimpleNumFieldElem}}) -> Bool

Return whether the spaces L and M are isometric over the completion at p.

source
invariantsMethod
invariants(M::QuadSpace)
      -> FieldElem, Dict{AbsNumFieldOrderIdeal{AbsSimpleNumField, AbsSimpleNumFieldElem}, Int}, Vector{Tuple{InfPlc, Int}}

Returns a tuple (n, k, d, H, I) of invariants of M, which determine the isometry class completely. Here n is the dimension. The dimension of the kernel is k. The element d is the determinant of a Gram matrix of the non-degenerate part, H contains the non-trivial Hasse invariants and I contains for each real place the negative index of inertia.

Note that d is determined only modulo squares.

source

Examples

For instance, for the case of $Q$ and the totally ramified prime $\mathfrak p$ of $O_K$ above $7$, one can get:


julia> K, a = cyclotomic_real_subfield(7);
julia> Q = quadratic_space(K, K[0 1; 1 0]);
julia> OK = maximal_order(K);
julia> p = prime_decomposition(OK, 7)[1][1];
julia> hasse_invariant(Q, p), witt_invariant(Q, p)(1, 1)
julia> Q2 = quadratic_space(K, K[-1 0; 0 1]);
julia> is_isometric(Q, Q2, p)true
julia> is_isometric(Q, Q2)true
julia> invariants(Q2)(2, 0, -1, Dict{AbsSimpleNumFieldOrderIdeal, Int64}(), Tuple{InfPlc{AbsSimpleNumField, AbsSimpleNumFieldEmbedding}, Int64}[(Infinite place corresponding to (Complex embedding corresponding to -1.80 of maximal real subfield of cyclotomic field of order 7), 1), (Infinite place corresponding to (Complex embedding corresponding to -0.45 of maximal real subfield of cyclotomic field of order 7), 1), (Infinite place corresponding to (Complex embedding corresponding to 1.25 of maximal real subfield of cyclotomic field of order 7), 1)])

Embeddings

Let $(V, \Phi)$ and $(V', \Phi')$ be two spaces over the same extension $E/K$, and let $\sigma \colon V \to V'$ be an $E$-linear morphism. $\sigma$ is called a representation of $V$ into $V'$ if for all $x \in V$

\[ \Phi'(\sigma(x), \sigma(x)) = \Phi(x,x).\]

In such a case, $V$ is said to be represented by $V'$ and $\sigma$ can be seen as an embedding of $V$ into $V'$. This representation property can be also tested locally with respect to the completions at some finite places. Note that in both quadratic and hermitian cases, completions are taken at finite places of the fixed field $K$.

is_locally_represented_byMethod
is_locally_represented_by(U::T, V::T, p::AbsNumFieldOrderIdeal{AbsSimpleNumField, AbsSimpleNumFieldElem}) where T <: AbstractSpace -> Bool

Given two spaces U and V over the same algebra E, and a prime ideal p in the maximal order $\mathcal O_K$ of their fixed field K, return whether U is represented by V locally at p, i.e. whether $U_p$ embeds in $V_p$.

source
is_represented_byMethod
is_represented_by(U::T, V::T) where T <: AbstractSpace -> Bool

Given two spaces U and V over the same algebra E, return whether U is represented by V, i.e. whether U embeds in V.

source

Examples

Still using the spaces $Q$ and $H$, we can decide whether some other spaces embed respectively locally or globally into $Q$ or $H$:


julia> K, a = cyclotomic_real_subfield(7);
julia> Kt, t = K["t"];
julia> E, b = number_field(t^2-a*t+1, "b");
julia> Q = quadratic_space(K, K[0 1; 1 0]);
julia> H = hermitian_space(E, 3);
julia> OK = maximal_order(K);
julia> p = prime_decomposition(OK, 7)[1][1];
julia> Q2 = quadratic_space(K, K[-1 0; 0 1]);
julia> H2 = hermitian_space(E, E[-1 0 0; 0 1 0; 0 0 -1]);
julia> is_locally_represented_by(Q2, Q, p)true
julia> is_represented_by(Q2, Q)true
julia> is_locally_represented_by(H2, H, p)true
julia> is_represented_by(H2, H)false

Categorical constructions

One can construct direct sums of spaces of the same kind. Since those are also direct products, they are called biproducts in this context. Depending on the user usage, one of the following three methods can be called to obtain the direct sum of a finite collection of spaces. Note that the corresponding copies of the original spaces in the direct sum are pairwise orthogonal.

direct_sumMethod
direct_sum(x::Vararg{T}) where T <: AbstractSpace -> T, Vector{AbstractSpaceMor}
direct_sum(x::Vector{T}) where T <: AbstractSpace -> T, Vector{AbstractSpaceMor}

Given a collection of quadratic or hermitian spaces $V_1, \ldots, V_n$, return their direct sum $V := V_1 \oplus \ldots \oplus V_n$, together with the injections $V_i \to V$.

For objects of type AbstractSpace, finite direct sums and finite direct products agree and they are therefore called biproducts. If one wants to obtain V as a direct product with the projections $V \to V_i$, one should call direct_product(x). If one wants to obtain V as a biproduct with the injections $V_i \to V$ and the projections $V \to V_i$, one should call biproduct(x).

source
direct_sum(g1::QuadSpaceCls, g2::QuadSpaceCls) -> QuadSpaceCls

Return the isometry class of the direct sum of two representatives.

source
direct_sum(M::ModuleFP{T}...; task::Symbol = :sum) where T

Given modules $M_1\dots M_n$, say, return the direct sum $\bigoplus_{i=1}^n M_i$.

Additionally, return

  • a vector containing the canonical injections $M_i\to\bigoplus_{i=1}^n M_i$ if task = :sum (default),
  • a vector containing the canonical projections $\bigoplus_{i=1}^n M_i\to M_i$ if task = :prod,
  • two vectors containing the canonical injections and projections, respectively, if task = :both,
  • none of the above maps if task = :none.
source
direct_productMethod
direct_product(algebras::StructureConstantAlgebra...; task::Symbol = :sum)
  -> StructureConstantAlgebra, Vector{AbsAlgAssMor}, Vector{AbsAlgAssMor}
direct_product(algebras::Vector{StructureConstantAlgebra}; task::Symbol = :sum)
  -> StructureConstantAlgebra, Vector{AbsAlgAssMor}, Vector{AbsAlgAssMor}

Returns the algebra $A = A_1 \times \cdots \times A_k$. task can be ":sum", ":prod", ":both" or ":none" and determines which canonical maps are computed as well: ":sum" for the injections, ":prod" for the projections.

source
direct_product(x::Vararg{T}) where T <: AbstractSpace -> T, Vector{AbstractSpaceMor}
direct_product(x::Vector{T}) where T <: AbstractSpace -> T, Vector{AbstractSpaceMor}

Given a collection of quadratic or hermitian spaces $V_1, \ldots, V_n$, return their direct product $V := V_1 \times \ldots \times V_n$, together with the projections $V \to V_i$.

For objects of type AbstractSpace, finite direct sums and finite direct products agree and they are therefore called biproducts. If one wants to obtain V as a direct sum with the injections $V_i \to V$, one should call direct_sum(x). If one wants to obtain V as a biproduct with the injections $V_i \to V$ and the projections $V \to V_i$, one should call biproduct(x).

source
direct_product(F::FreeMod{T}...; task::Symbol = :prod) where T

Given free modules $F_1\dots F_n$, say, return the direct product $\prod_{i=1}^n F_i$.

Additionally, return

  • a vector containing the canonical projections $\prod_{i=1}^n F_i\to F_i$ if task = :prod (default),
  • a vector containing the canonical injections $F_i\to\prod_{i=1}^n F_i$ if task = :sum,
  • two vectors containing the canonical projections and injections, respectively, if task = :both,
  • none of the above maps if task = :none.
source
direct_product(M::ModuleFP{T}...; task::Symbol = :prod) where T

Given modules $M_1\dots M_n$, say, return the direct product $\prod_{i=1}^n M_i$.

Additionally, return

  • a vector containing the canonical projections $\prod_{i=1}^n M_i\to M_i$ if task = :prod (default),
  • a vector containing the canonical injections $M_i\to\prod_{i=1}^n M_i$ if task = :sum,
  • two vectors containing the canonical projections and injections, respectively, if task = :both,
  • none of the above maps if task = :none.
source
biproductMethod
biproduct(x::Vararg{T}) where T <: AbstractSpace -> T, Vector{AbstractSpaceMor}, Vector{AbstractSpaceMor}
biproduct(x::Vector{T}) where T <: AbstractSpace -> T, Vector{AbstractSpaceMor}, Vector{AbstractSpaceMor}

Given a collection of quadratic or hermitian spaces $V_1, \ldots, V_n$, return their biproduct $V := V_1 \oplus \ldots \oplus V_n$, together with the injections $V_i \to V$ and the projections $V \to V_i$.

For objects of type AbstractSpace, finite direct sums and finite direct products agree and they are therefore called biproducts. If one wants to obtain V as a direct sum with the injections $V_i \to V$, one should call direct_sum(x). If one wants to obtain V as a direct product with the projections $V \to V_i$, one should call direct_product(x).

source

Example


julia> E, b = cyclotomix_field_as_cm_extensions(7);ERROR: UndefVarError: `cyclotomix_field_as_cm_extensions` not defined
julia> H = hermitian_space(E, 3);
julia> H2 = hermitian_space(E, E[-1 0 0; 0 1 0; 0 0 -1]);
julia> H3, inj, proj = biproduct(H, H2)(Hermitian space of dimension 6, AbstractSpaceMor[Map: hermitian space -> hermitian space, Map: hermitian space -> hermitian space], AbstractSpaceMor[Map: hermitian space -> hermitian space, Map: hermitian space -> hermitian space])
julia> is_one(matrix(compose(inj[1], proj[1])))true
julia> is_zero(matrix(compose(inj[1], proj[2])))true

Orthogonality operations

orthogonal_complementMethod
orthogonal_complement(V::AbstractSpace, M::T) where T <: MatElem -> T

Given a space V and a subspace W with basis matrix M, return a basis matrix of the orthogonal complement of W inside V.

source
orthogonal_projectionMethod
orthogonal_projection(V::AbstractSpace, M::T) where T <: MatElem -> AbstractSpaceMor

Given a space V and a non-degenerate subspace W with basis matrix M, return the endomorphism of V corresponding to the projection onto the complement of W in V.

source

Example


julia> K, a = cyclotomic_real_subfield(7);
julia> Kt, t = K["t"];
julia> Q = quadratic_space(K, K[0 1; 1 0]);
julia> orthogonal_complement(Q, matrix(K, 1, 2, [1 0]))[1 0]

Isotropic spaces

Let $(V, \Phi)$ be a space over $E/K$ and let $\mathfrak p$ be a place in $K$. $V$ is said to be isotropic locally at $\mathfrak p$ if there exists an element $x \in V_{\mathfrak p}$ such that $\Phi_{\mathfrak p}(x,x) = 0$, where $\Phi_{\mathfrak p}$ is the continuous extension of $\Phi$ to $V_{\mathfrak p} \times V_{\mathfrak p}$.

is_isotropicMethod
is_isotropic(V::AbstractSpace, p::Union{AbsNumFieldOrderIdeal{AbsSimpleNumField, AbsSimpleNumFieldElem}, InfPlc}) -> Bool

Given a space V and a place p in the fixed field K of V, return whether the completion of V at p is isotropic.

source

Example


julia> K, a = cyclotomic_real_subfield(7);
julia> Kt, t = K["t"];
julia> E, b = number_field(t^2-a*t+1, "b");
julia> H = hermitian_space(E, 3);
julia> OK = maximal_order(K);
julia> p = prime_decomposition(OK, 7)[1][1];
julia> is_isotropic(H, p)true

Hyperbolic spaces

Let $(V, \Phi)$ be a space over $E/K$ and let $\mathfrak p$ be a prime ideal of $\mathcal O_K$. $V$ is said to be hyperbolic locally at $\mathfrak p$ if the completion $V_{\mathfrak p}$ of $V$ can be decomposed as an orthogonal sum of hyperbolic planes. The hyperbolic plane is the space $(H, \Psi)$ of rank 2 over $E/K$ such that there exists a basis $e_1, e_2$ of $H$ such that $\Psi(e_1, e_1) = \Psi(e_2, e_2) = 0$ and $\Psi(e_1, e_2) = 1$.

is_locally_hyperbolicMethod
is_locally_hyperbolic(V::Hermspace, p::AbsNumFieldOrderIdeal{AbsSimpleNumField, AbsSimpleNumFieldElem}) -> Bool

Return whether the completion of the hermitian space V over $E/K$ at the prime ideal p of $\mathcal O_K$ is hyperbolic.

source

Example


julia> K, a = cyclotomic_real_subfield(7);
julia> Kt, t = K["t"];
julia> E, b = number_field(t^2-a*t+1, "b");
julia> H = hermitian_space(E, 3);
julia> OK = maximal_order(K);
julia> p = prime_decomposition(OK, 7)[1][1];
julia> is_locally_hyperbolic(H, p)false