Integer Lattices
An integer lattice
A ZZLat
. It is given in terms of its ambient quadratic space
To access ambient_space(L::ZZLat)
and basis_matrix(L::ZZLat)
.
Creation of integer lattices
From a gram matrix
integer_lattice Method
integer_lattice([B::MatElem]; gram) -> ZZLat
Return the Z-lattice with basis matrix gram
.
If the keyword gram
is not specified, the Gram matrix is the identity matrix. If
Examples
julia> L = integer_lattice(matrix(QQ, 2, 2, [1//2, 0, 0, 2]));
julia> gram_matrix(L) == matrix(QQ, 2, 2, [1//4, 0, 0, 4])
true
julia> L = integer_lattice(gram = matrix(ZZ, [2 -1; -1 2]));
julia> gram_matrix(L) == matrix(ZZ, [2 -1; -1 2])
true
In a quadratic space
lattice Method
lattice(V::AbstractSpace, basis::MatElem ; check::Bool = true) -> AbstractLat
Given an ambient space V
and a matrix basis
, return the lattice spanned by the rows of basis
inside V
. If V
is hermitian (resp. quadratic) then the output is a hermitian (resp. quadratic) lattice.
By default, basis
is checked to be of full rank. This test can be disabled by setting check
to false.
Special lattices
root_lattice Method
root_lattice(R::Symbol, n::Int) -> ZZLat
Return the root lattice of type R
given by :A
, :D
or :E
with parameter n
.
The type :I
with parameter n = 1
is also allowed and denotes the odd unimodular lattice of rank 1.
hyperbolic_plane_lattice Method
hyperbolic_plane_lattice(n::RationalUnion = 1) -> ZZLat
Return the hyperbolic plane with intersection form of scale n
, that is, the unique (up to isometry) even unimodular hyperbolic n
.
Examples
julia> L = hyperbolic_plane_lattice(6);
julia> gram_matrix(L)
[0 6]
[6 0]
julia> L = hyperbolic_plane_lattice(ZZ(-13));
julia> gram_matrix(L)
[ 0 -13]
[-13 0]
integer_lattice Method
integer_lattice(S::Symbol, n::RationalUnion = 1) -> ZZlat
Given S = :H
or S = :U
, return a
leech_lattice Function
leech_lattice() -> ZZLat
Return the Leech lattice.
sourceleech_lattice(niemeier_lattice::ZZLat) -> ZZLat, QQMatrix, Int
Return a triple L, v, h
where L
is the Leech lattice.
L is an h
-neighbor of the Niemeier lattice N
with respect to v
. This means that L / L ∩ N ≅ ℤ / h ℤ
. Here h
is the Coxeter number of the Niemeier lattice.
This implements the 23 holy constructions of the Leech lattice in [5].
Examples
julia> R = integer_lattice(gram=2 * identity_matrix(ZZ, 24));
julia> N = maximal_even_lattice(R) # Some Niemeier lattice
Integer lattice of rank 24 and degree 24
with gram matrix
[2 1 1 1 0 0 0 0 0 0 0 0 0 0 0 0 1 0 1 1 0 0 0 0]
[1 2 1 1 0 0 0 0 0 0 0 0 0 0 0 0 1 1 0 1 0 0 0 0]
[1 1 2 1 0 0 0 0 0 0 0 0 0 0 0 0 1 1 1 0 0 0 0 0]
[1 1 1 2 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0]
[0 0 0 0 2 1 1 1 0 0 0 0 1 0 1 1 0 0 0 0 0 0 0 0]
[0 0 0 0 1 2 1 1 0 0 0 0 1 1 0 1 0 0 0 0 0 0 0 0]
[0 0 0 0 1 1 2 1 0 0 0 0 1 1 1 0 0 0 0 0 0 0 0 0]
[0 0 0 0 1 1 1 2 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0]
[0 0 0 0 0 0 0 0 2 1 1 1 0 0 0 0 0 0 0 0 1 1 1 0]
[0 0 0 0 0 0 0 0 1 2 1 1 0 0 0 0 0 0 0 0 1 0 1 1]
[0 0 0 0 0 0 0 0 1 1 2 1 0 0 0 0 0 0 0 0 1 1 0 1]
[0 0 0 0 0 0 0 0 1 1 1 2 0 0 0 0 0 0 0 0 0 0 0 0]
[0 0 0 0 1 1 1 0 0 0 0 0 2 1 1 1 0 0 0 0 0 0 0 0]
[0 0 0 0 0 1 1 0 0 0 0 0 1 2 0 0 0 0 0 0 0 0 0 0]
[0 0 0 0 1 0 1 0 0 0 0 0 1 0 2 0 0 0 0 0 0 0 0 0]
[0 0 0 0 1 1 0 0 0 0 0 0 1 0 0 2 0 0 0 0 0 0 0 0]
[1 1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 2 1 1 1 0 0 0 0]
[0 1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 1 2 0 0 0 0 0 0]
[1 0 1 0 0 0 0 0 0 0 0 0 0 0 0 0 1 0 2 0 0 0 0 0]
[1 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 1 0 0 2 0 0 0 0]
[0 0 0 0 0 0 0 0 1 1 1 0 0 0 0 0 0 0 0 0 2 1 1 1]
[0 0 0 0 0 0 0 0 1 0 1 0 0 0 0 0 0 0 0 0 1 2 0 0]
[0 0 0 0 0 0 0 0 1 1 0 0 0 0 0 0 0 0 0 0 1 0 2 0]
[0 0 0 0 0 0 0 0 0 1 1 0 0 0 0 0 0 0 0 0 1 0 0 2]
julia> minimum(N)
2
julia> det(N)
1
julia> L, v, h = leech_lattice(N);
julia> minimum(L)
4
julia> det(L)
1
julia> h == index(L, intersect(L, N))
true
We illustrate how the Leech lattice is constructed from N
, h
and v
.
julia> Zmodh, _ = residue_ring(ZZ, h);
julia> V = ambient_space(N);
julia> vG = map_entries(x->Zmodh(ZZ(x)), inner_product(V, v, basis_matrix(N)));
julia> LN = transpose(lift(Hecke.kernel(vG; side = :right)))*basis_matrix(N); # vectors whose inner product with `v` is divisible by `h`.
julia> M = lattice(V, LN) + h*N;
julia> M == intersect(L, N)
true
julia> M + lattice(V, 1//h * v) == L
true
k3_lattice Function
k3_lattice()
Return the integer lattice corresponding to the Beauville-Bogomolov-Fujiki form associated to a K3 surface.
Examples
julia> L = k3_lattice();
julia> is_unimodular(L)
true
julia> signature_tuple(L)
(3, 0, 19)
mukai_lattice Method
mukai_lattice(S::Symbol = :K3; extended::Bool = false)
Return the (extended) Mukai lattice.
If S == :K3
, it returns the (extended) Mukai lattice associated to hyperkaehler manifolds which are deformation equivalent to a moduli space of stable sheaves on a K3 surface.
If S == :Ab
, it returns the (extended) Mukai lattice associated to hyperkaehler manifolds which are deformation equivalent to a moduli space of stable sheaves on an abelian surface.
Examples
julia> L = mukai_lattice();
julia> genus(L)
Genus symbol for integer lattices
Signatures: (4, 0, 20)
Local symbol:
Local genus symbol at 2: 1^24
julia> L = mukai_lattice(; extended = true);
julia> genus(L)
Genus symbol for integer lattices
Signatures: (5, 0, 21)
Local symbol:
Local genus symbol at 2: 1^26
julia> L = mukai_lattice(:Ab);
julia> genus(L)
Genus symbol for integer lattices
Signatures: (4, 0, 4)
Local symbol:
Local genus symbol at 2: 1^8
julia> L = mukai_lattice(:Ab; extended = true);
julia> genus(L)
Genus symbol for integer lattices
Signatures: (5, 0, 5)
Local symbol:
Local genus symbol at 2: 1^10
hyperkaehler_lattice Method
hyperkaehler_lattice(S::Symbol; n::Int = 2)
Return the integer lattice corresponding to the Beauville-Bogomolov-Fujiki form on a hyperkaehler manifold whose deformation type is determined by S
and n
.
If
S == :K3
orS == :Kum
, thenn
must be an integer bigger than 2;If
S == :OG6
orS == :OG10
, the value ofn
has no effect.
Examples
julia> L = hyperkaehler_lattice(:Kum; n = 3)
Integer lattice of rank 7 and degree 7
with gram matrix
[0 1 0 0 0 0 0]
[1 0 0 0 0 0 0]
[0 0 0 1 0 0 0]
[0 0 1 0 0 0 0]
[0 0 0 0 0 1 0]
[0 0 0 0 1 0 0]
[0 0 0 0 0 0 -8]
julia> L = hyperkaehler_lattice(:OG6)
Integer lattice of rank 8 and degree 8
with gram matrix
[0 1 0 0 0 0 0 0]
[1 0 0 0 0 0 0 0]
[0 0 0 1 0 0 0 0]
[0 0 1 0 0 0 0 0]
[0 0 0 0 0 1 0 0]
[0 0 0 0 1 0 0 0]
[0 0 0 0 0 0 -2 0]
[0 0 0 0 0 0 0 -2]
julia> L = hyperkaehler_lattice(:OG10);
julia> genus(L)
Genus symbol for integer lattices
Signatures: (3, 0, 21)
Local symbols:
Local genus symbol at 2: 1^-24
Local genus symbol at 3: 1^-23 3^1
julia> L = hyperkaehler_lattice(:K3; n = 3);
julia> genus(L)
Genus symbol for integer lattices
Signatures: (3, 0, 20)
Local symbol:
Local genus symbol at 2: 1^22 4^1_7
From a genus
Integer lattices can be created as representatives of a genus. See (representative(L::ZZGenus)
)
Rescaling the Quadratic Form
rescale Method
rescale(L::ZZLat, r::RationalUnion) -> ZZLat
Return the lattice L
in the quadratic space with form r \Phi
.
Examples
This can be useful to apply methods intended for positive definite lattices.
julia> L = integer_lattice(gram=ZZ[-1 0; 0 -1])
Integer lattice of rank 2 and degree 2
with gram matrix
[-1 0]
[ 0 -1]
julia> shortest_vectors(rescale(L, -1))
2-element Vector{Vector{ZZRingElem}}:
[0, 1]
[1, 0]
Attributes
ambient_space Method
ambient_space(L::AbstractLat) -> AbstractSpace
Return the ambient space of the lattice L
. If the ambient space is not known, an error is raised.
basis_matrix Method
basis_matrix(L::ZZLat) -> QQMatrix
Return the basis matrix
The lattice is given by the row span of
gram_matrix Method
gram_matrix(L::ZZLat) -> QQMatrix
Return the gram matrix of
Examples
julia> L = integer_lattice(matrix(ZZ, [2 0; -1 2]));
julia> gram_matrix(L)
[ 4 -2]
[-2 5]
rational_span Method
rational_span(L::ZZLat) -> QuadSpace
Return the rational span of gram_matrix(L)
.
Examples
julia> L = integer_lattice(matrix(ZZ, [2 0; -1 2]));
julia> rational_span(L)
Quadratic space of dimension 2
over rational field
with gram matrix
[ 4 -2]
[-2 5]
Invariants
rank Method
rank(L::AbstractLat) -> Int
Return the rank of the underlying module of the lattice L
.
scale Method
scale(L::ZZLat) -> QQFieldElem
Return the scale of L
.
The scale of L
is defined as the positive generator of the
norm Method
norm(L::ZZLat) -> QQFieldElem
Return the norm of L
.
The norm of L
is defined as the positive generator of the
iseven Method
iseven(L::ZZLat) -> Bool
Return whether L
is even.
An integer lattice L
in the rational quadratic space
is_integral Method
is_integral(L::AbstractLat) -> Bool
Return whether the lattice L
is integral.
is_primary_with_prime Method
is_primary_with_prime(L::ZZLat) -> Bool, ZZRingElem
Given a L
, return whether L
is primary, that is whether L
is integral and its discriminant group (see discriminant_group
) is a p
-group for some prime number p
. In case it is, p
is also returned as second output.
Note that for unimodular lattices, this function returns (true, 1)
. If the lattice is not primary, the second return value is -1
by default.
is_primary Method
is_primary(L::ZZLat, p::Union{Integer, ZZRingElem}) -> Bool
Given an integral L
and a prime number p
, return whether L
is p
-primary, that is whether its discriminant group (see discriminant_group
) is a p
-group.
is_elementary_with_prime Method
is_elementary_with_prime(L::ZZLat) -> Bool, ZZRingElem
Given a L
, return whether L
is elementary, that is whether L
is integral and its discriminant group (see discriminant_group
) is an elemenentary p
-group for some prime number p
. In case it is, p
is also returned as second output.
Note that for unimodular lattices, this function returns (true, 1)
. If the lattice is not elementary, the second return value is -1
by default.
is_elementary Method
is_elementary(L::ZZLat, p::Union{Integer, ZZRingElem}) -> Bool
Given an integral L
and a prime number p
, return whether L
is p
-elementary, that is whether its discriminant group (see discriminant_group
) is an elementary p
-group.
The Genus
For an integral lattice The genus of an integer lattice collects its local invariants. genus(::ZZLat)
genus_representatives Method
genus_representatives(L::ZZLat) -> Vector{ZZLat}
Return representatives for the isometry classes in the genus of L
.
Real invariants
signature_tuple Method
signature_tuple(L::ZZLat) -> Tuple{Int,Int,Int}
Return the number of (positive, zero, negative) inertia of L
.
is_positive_definite Method
is_positive_definite(L::AbstractLat) -> Bool
Return whether the rational span of the lattice L
is positive definite.
is_negative_definite Method
is_negative_definite(L::AbstractLat) -> Bool
Return whether the rational span of the lattice L
is negative definite.
is_definite Method
is_definite(L::AbstractLat) -> Bool
Return whether the rational span of the lattice L
is definite.
Isometries
automorphism_group_generators Method
automorphism_group_generators(E::EllipticCurve) -> Vector{EllCrvIso}
Return generators of the automorphism group of
automorphism_group_generators(L::AbstractLat; ambient_representation::Bool = true,
depth::Int = -1, bacher_depth::Int = 0)
-> Vector{MatElem}
Given a definite lattice L
, return generators for the automorphism group of L
. If ambient_representation == true
(the default), the transformations are represented with respect to the ambient space of L
. Otherwise, the transformations are represented with respect to the (pseudo-)basis of L
.
Setting the parameters depth
and bacher_depth
to a positive value may improve performance. If set to -1
(default), the used value of depth
is chosen heuristically depending on the rank of L
. By default, bacher_depth
is set to 0
.
automorphism_group_order Method
automorphism_group_order(L::AbstractLat; depth::Int = -1, bacher_depth::Int = 0) -> Int
Given a definite lattice L
, return the order of the automorphism group of L
.
Setting the parameters depth
and bacher_depth
to a positive value may improve performance. If set to -1
(default), the used value of depth
is chosen heuristically depending on the rank of L
. By default, bacher_depth
is set to 0
.
is_isometric Method
is_isometric(L::AbstractLat, M::AbstractLat; depth::Int = -1, bacher_depth::Int = 0) -> Bool
Return whether the lattices L
and M
are isometric.
Setting the parameters depth
and bacher_depth
to a positive value may improve performance. If set to -1
(default), the used value of depth
is chosen heuristically depending on the rank of L
. By default, bacher_depth
is set to 0
.
is_locally_isometric Method
is_locally_isometric(L::ZZLat, M::ZZLat, p::Int) -> Bool
Return whether L
and M
are isometric over the p
-adic integers.
i.e. whether
Root lattices
root_lattice_recognition Method
root_lattice_recognition(L::ZZLat)
Return the ADE type of the root sublattice of L
.
The root sublattice is the lattice spanned by the vectors of squared length (:I, 1)
.
Input:
L
– a definite and integral
Output:
Two lists, the first one containing the ADE types and the second one the irreducible root sublattices.
For more recognizable gram matrices use root_lattice_recognition_fundamental
.
Examples
julia> L = integer_lattice(; gram=ZZ[4 0 0 0 3 0 3 0;
0 16 8 12 2 12 6 10;
0 8 8 6 2 8 4 5;
0 12 6 10 2 9 5 8;
3 2 2 2 4 2 4 2;
0 12 8 9 2 12 6 9;
3 6 4 5 4 6 6 5;
0 10 5 8 2 9 5 8])
Integer lattice of rank 8 and degree 8
with gram matrix
[4 0 0 0 3 0 3 0]
[0 16 8 12 2 12 6 10]
[0 8 8 6 2 8 4 5]
[0 12 6 10 2 9 5 8]
[3 2 2 2 4 2 4 2]
[0 12 8 9 2 12 6 9]
[3 6 4 5 4 6 6 5]
[0 10 5 8 2 9 5 8]
julia> R = root_lattice_recognition(L)
([(:A, 1), (:D, 6)], ZZLat[Integer lattice of rank 1 and degree 8, Integer lattice of rank 6 and degree 8])
julia> L = integer_lattice(; gram = QQ[1 0 0 0;
0 9 3 3;
0 3 2 1;
0 3 1 11])
Integer lattice of rank 4 and degree 4
with gram matrix
[1 0 0 0]
[0 9 3 3]
[0 3 2 1]
[0 3 1 11]
julia> root_lattice_recognition(L)
([(:A, 1), (:I, 1)], ZZLat[Integer lattice of rank 1 and degree 4, Integer lattice of rank 1 and degree 4])
root_lattice_recognition_fundamental Method
root_lattice_recognition_fundamental(L::ZZLat)
Return the ADE type of the root sublattice of L
as well as the corresponding irreducible root sublattices with basis given by a fundamental root system.
The type (:I, 1)
corresponds to the odd unimodular root lattice of rank 1.
Input:
L
– a definite and integral
Output:
the root sublattice, with basis given by a fundamental root system
the ADE types
a Vector consisting of the irreducible root sublattices.
Examples
julia> L = integer_lattice(gram=ZZ[4 0 0 0 3 0 3 0;
0 16 8 12 2 12 6 10;
0 8 8 6 2 8 4 5;
0 12 6 10 2 9 5 8;
3 2 2 2 4 2 4 2;
0 12 8 9 2 12 6 9;
3 6 4 5 4 6 6 5;
0 10 5 8 2 9 5 8])
Integer lattice of rank 8 and degree 8
with gram matrix
[4 0 0 0 3 0 3 0]
[0 16 8 12 2 12 6 10]
[0 8 8 6 2 8 4 5]
[0 12 6 10 2 9 5 8]
[3 2 2 2 4 2 4 2]
[0 12 8 9 2 12 6 9]
[3 6 4 5 4 6 6 5]
[0 10 5 8 2 9 5 8]
julia> R = root_lattice_recognition_fundamental(L);
julia> gram_matrix(R[1])
[2 0 0 0 0 0 0]
[0 2 0 -1 0 0 0]
[0 0 2 -1 0 0 0]
[0 -1 -1 2 -1 0 0]
[0 0 0 -1 2 -1 0]
[0 0 0 0 -1 2 -1]
[0 0 0 0 0 -1 2]
ADE_type Method
ADE_type(G::MatrixElem) -> Tuple{Symbol,Int64}
Return the type of the irreducible root lattice with gram matrix G
.
See also root_lattice_recognition
.
Examples
julia> Hecke.ADE_type(gram_matrix(root_lattice(:A,3)))
(:A, 3)
coxeter_number Method
coxeter_number(ADE::Symbol, n) -> Int
Return the Coxeter number of the corresponding ADE root lattice.
If n
is the rank of
Examples
julia> coxeter_number(:D, 4)
6
highest_root Method
highest_root(ADE::Symbol, n) -> ZZMatrix
Return coordinates of the highest root of root_lattice(ADE, n)
.
Examples
julia> highest_root(:E, 6)
[1 2 3 2 1 2]
Module operations
Most module operations assume that the lattices live in the same ambient space. For instance only lattices in the same ambient space compare.
== Method
Return true
if both lattices have the same ambient quadratic space and the same underlying module.
is_sublattice Method
is_sublattice(L::AbstractLat, M::AbstractLat) -> Bool
Return whether M
is a sublattice of the lattice L
.
is_sublattice_with_relations Method
is_sublattice_with_relations(M::ZZLat, N::ZZLat) -> Bool, QQMatrix
Returns whether
+ Method
+(L::AbstractLat, M::AbstractLat) -> AbstractLat
Return the sum of the lattices L
and M
.
The lattices L
and M
must have the same ambient space.
* Method
*(a::RationalUnion, L::ZZLat) -> ZZLat
Return the lattice
intersect Method
intersect(L::AbstractLat, M::AbstractLat) -> AbstractLat
Return the intersection of the lattices L
and M
.
The lattices L
and M
must have the same ambient space.
in Method
Base.in(v::Vector, L::ZZLat) -> Bool
Return whether the vector v
lies in the lattice L
.
in Method
Base.in(v::QQMatrix, L::ZZLat) -> Bool
Return whether the row span of v
lies in the lattice L
.
primitive_closure Method
primitive_closure(M::ZZLat, N::ZZLat) -> ZZLat
Given two M
and N
with N
in M
.
Examples
julia> M = root_lattice(:D, 6);
julia> N = lattice_in_same_ambient_space(M, 3*basis_matrix(M)[1:1,:]);
julia> basis_matrix(N)
[3 0 0 0 0 0]
julia> N2 = primitive_closure(M, N)
Integer lattice of rank 1 and degree 6
with gram matrix
[2]
julia> basis_matrix(N2)
[1 0 0 0 0 0]
julia> M2 = primitive_closure(dual(M), M);
julia> is_integral(M2)
false
is_primitive Method
is_primitive(M::ZZLat, N::ZZLat) -> Bool
Given two N
is a primitive sublattice of M
.
Examples
julia> U = hyperbolic_plane_lattice(3);
julia> bU = basis_matrix(U);
julia> e1, e2 = bU[1:1,:], bU[2:2,:]
([1 0], [0 1])
julia> N = lattice_in_same_ambient_space(U, e1 + e2)
Integer lattice of rank 1 and degree 2
with gram matrix
[6]
julia> is_primitive(U, N)
true
julia> M = root_lattice(:A, 3);
julia> f = matrix(QQ, 3, 3, [0 1 1; -1 -1 -1; 1 1 0]);
julia> N = kernel_lattice(M, f+1)
Integer lattice of rank 1 and degree 3
with gram matrix
[4]
julia> is_primitive(M, N)
true
is_primitive Method
is_primitive(L::ZZLat, v::Union{Vector, QQMatrix}) -> Bool
Return whether the vector v
is primitive in L
.
A vector v
in a L
is called primitive if for all w
in L
such that d
, then
divisibility Method
divisibility(L::ZZLat, v::Union{Vector, QQMatrix}) -> QQFieldElem
Return the divisibility of v
with respect to L
.
For a vector v
in the ambient quadratic space L
, we call the divisibility of v
with the respect to L
the non-negative generator of the fractional
Embeddings
Categorical constructions
direct_sum Method
direct_sum(x::Vararg{ZZLat}) -> ZZLat, Vector{AbstractSpaceMor}
direct_sum(x::Vector{ZZLat}) -> ZZLat, Vector{AbstractSpaceMor}
Given a collection of
For objects of type ZZLat
, finite direct sums and finite direct products agree and they are therefore called biproducts. If one wants to obtain L
as a direct product with the projections direct_product(x)
. If one wants to obtain L
as a biproduct with the injections biproduct(x)
.
direct_product Method
direct_product(x::Vararg{ZZLat}) -> ZZLat, Vector{AbstractSpaceMor}
direct_product(x::Vector{ZZLat}) -> ZZLat, Vector{AbstractSpaceMor}
Given a collection of
For objects of type ZZLat
, finite direct sums and finite direct products agree and they are therefore called biproducts. If one wants to obtain L
as a direct sum with the injections direct_sum(x)
. If one wants to obtain L
as a biproduct with the injections biproduct(x)
.
biproduct Method
biproduct(x::Vararg{ZZLat}) -> ZZLat, Vector{AbstractSpaceMor}, Vector{AbstractSpaceMor}
biproduct(x::Vector{ZZLat}) -> ZZLat, Vector{AbstractSpaceMor}, Vector{AbstractSpaceMor}
Given a collection of
For objects of type ZZLat
, finite direct sums and finite direct products agree and they are therefore called biproducts. If one wants to obtain L
as a direct sum with the injections direct_sum(x)
. If one wants to obtain L
as a direct product with the projections direct_product(x)
.
Orthogonal sublattices
orthogonal_submodule Method
orthogonal_submodule(L::ZZLat, S::ZZLat) -> ZZLat
Return the largest submodule of
irreducible_components Method
irreducible_components(L::ZZLat) -> Vector{ZZLat}
Return the irreducible components
This yields a maximal orthogonal splitting of L
as
The decomposition is unique up to ordering of the factors.
Examples
julia> R = integer_lattice(; gram=2 * identity_matrix(ZZ, 24));
julia> N = maximal_even_lattice(R); # Niemeier lattice
julia> irr_comp = irreducible_components(N)
3-element Vector{ZZLat}:
Integer lattice of rank 8 and degree 24
Integer lattice of rank 8 and degree 24
Integer lattice of rank 8 and degree 24
julia> genus.(irr_comp)
3-element Vector{ZZGenus}:
Genus symbol: II_(8, 0)
Genus symbol: II_(8, 0)
Genus symbol: II_(8, 0)
Dual lattice
Discriminant group
See discriminant_group(L::ZZLat)
.
Overlattices
glue_map Method
glue_map(L::ZZLat, S::ZZLat, R::ZZLat; check=true)
-> Tuple{TorQuadModuleMap, TorQuadModuleMap, TorQuadModuleMap}
Given three integral L
, S
and R
, with S
and R
primitive sublattices of L
and such that the sum of the ranks of S
and R
is equal to the rank of L
, return the glue map S
and R
.
Example
julia> M = root_lattice(:E,8);
julia> f = matrix(QQ, 8, 8, [-1 -1 0 0 0 0 0 0;
1 0 0 0 0 0 0 0;
0 1 1 0 0 0 0 0;
0 0 0 1 0 0 0 0;
0 0 0 0 1 0 0 0;
0 0 0 0 0 1 1 0;
-2 -4 -6 -5 -4 -3 -2 -3;
0 0 0 0 0 0 0 1]);
julia> S = kernel_lattice(M ,f-1)
Integer lattice of rank 4 and degree 8
with gram matrix
[12 -3 0 -3]
[-3 2 -1 0]
[ 0 -1 2 0]
[-3 0 0 2]
julia> R = kernel_lattice(M , f^2+f+1)
Integer lattice of rank 4 and degree 8
with gram matrix
[ 2 -1 0 0]
[-1 2 -6 0]
[ 0 -6 30 -3]
[ 0 0 -3 2]
julia> glue, iS, iR = glue_map(M, S, R)
(Map: finite quadratic module -> finite quadratic module, Map: finite quadratic module -> finite quadratic module, Map: finite quadratic module -> finite quadratic module)
julia> is_bijective(glue)
true
overlattice Method
overlattice(glue_map::TorQuadModuleMap) -> ZZLat
Given the glue map of a primitive extension of L
.
Example
julia> M = root_lattice(:E,8);
julia> f = matrix(QQ, 8, 8, [ 1 0 0 0 0 0 0 0;
0 1 0 0 0 0 0 0;
1 2 4 4 3 2 1 2;
-2 -4 -6 -5 -4 -3 -2 -3;
2 4 6 4 3 2 1 3;
-1 -2 -3 -2 -1 0 0 -2;
0 0 0 0 0 -1 0 0;
-1 -2 -3 -3 -2 -1 0 -1]);
julia> S = kernel_lattice(M ,f-1)
Integer lattice of rank 4 and degree 8
with gram matrix
[ 2 -1 0 0]
[-1 2 -1 0]
[ 0 -1 12 -15]
[ 0 0 -15 20]
julia> R = kernel_lattice(M , f^4+f^3+f^2+f+1)
Integer lattice of rank 4 and degree 8
with gram matrix
[10 -4 0 1]
[-4 2 -1 0]
[ 0 -1 4 -3]
[ 1 0 -3 4]
julia> glue, iS, iR = glue_map(M, S, R);
julia> overlattice(glue) == M
true
local_modification Method
local_modification(M::ZZLat, L::ZZLat, p)
Return a local modification of M
that matches L
at p
.
INPUT:
– a \mathbb{Z}_p
-maximal lattice– the a lattice isomorphic to M
over\QQ_p
– a prime number
OUTPUT:
an integral lattice M'
in the ambient space of M
such that M
and M'
are locally equal at all completions except at p
where M'
is locally isometric to the lattice L
.
maximal_integral_lattice Method
maximal_integral_lattice(L::AbstractLat) -> AbstractLat
Given a lattice L
with integral norm, return a maximal integral overlattice M
of L
.
Sublattices defined by endomorphisms
kernel_lattice Method
kernel_lattice(L::ZZLat, f::MatElem;
ambient_representation::Bool = true) -> ZZLat
Given a
If ambient_representation
is true
(the default), the endomorphism is represented with respect to the ambient space of
invariant_lattice Method
invariant_lattice(L::ZZLat, G::Vector{MatElem};
ambient_representation::Bool = true) -> ZZLat
invariant_lattice(L::ZZLat, G::MatElem;
ambient_representation::Bool = true) -> ZZLat
Given a
If ambient_representation
is true
(the default), the endomorphism is represented with respect to the ambient space of
coinvariant_lattice Method
coinvariant_lattice(L::ZZLat, G::Vector{MatElem};
ambient_representation::Bool = true) -> ZZLat
coinvariant_lattice(L::ZZLat, G::MatElem;
ambient_representation::Bool = true) -> ZZLat
Given a invariant_lattice
).
If ambient_representation
is true
(the default), the endomorphism is represented with respect to the ambient space of
Computing embeddings
embed Method
embed(S::ZZLat, G::Genus, primitive::Bool=true) -> Bool, embedding
Return a (primitive) embedding of the integral lattice S
into some lattice in the genus of G
.
julia> G = integer_genera((8,0), 1, even=true)[1];
julia> L, S, i = embed(root_lattice(:A,5), G);
embed_in_unimodular Method
embed_in_unimodular(S::ZZLat, pos::Int, neg::Int, primitive=true, even=true) -> Bool, L, S', iS, iR
Return a (primitive) embedding of the integral lattice S
into some (even) unimodular lattice of signature (pos, neg)
.
For now this works only for even lattices.
julia> NS = direct_sum(integer_lattice(:U), rescale(root_lattice(:A, 16), -1))[1];
julia> LK3, iNS, i = embed_in_unimodular(NS, 3, 19);
julia> genus(LK3)
Genus symbol for integer lattices
Signatures: (3, 0, 19)
Local symbol:
Local genus symbol at 2: 1^22
julia> iNS
Integer lattice of rank 18 and degree 22
with gram matrix
[0 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0]
[1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0]
[0 0 -2 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0]
[0 0 1 -2 1 0 0 0 0 0 0 0 0 0 0 0 0 0]
[0 0 0 1 -2 1 0 0 0 0 0 0 0 0 0 0 0 0]
[0 0 0 0 1 -2 1 0 0 0 0 0 0 0 0 0 0 0]
[0 0 0 0 0 1 -2 1 0 0 0 0 0 0 0 0 0 0]
[0 0 0 0 0 0 1 -2 1 0 0 0 0 0 0 0 0 0]
[0 0 0 0 0 0 0 1 -2 1 0 0 0 0 0 0 0 0]
[0 0 0 0 0 0 0 0 1 -2 1 0 0 0 0 0 0 0]
[0 0 0 0 0 0 0 0 0 1 -2 1 0 0 0 0 0 0]
[0 0 0 0 0 0 0 0 0 0 1 -2 1 0 0 0 0 0]
[0 0 0 0 0 0 0 0 0 0 0 1 -2 1 0 0 0 0]
[0 0 0 0 0 0 0 0 0 0 0 0 1 -2 1 0 0 0]
[0 0 0 0 0 0 0 0 0 0 0 0 0 1 -2 1 0 0]
[0 0 0 0 0 0 0 0 0 0 0 0 0 0 1 -2 1 0]
[0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 1 -2 1]
[0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 1 -2]
julia> is_primitive(LK3, iNS)
true
LLL, Short and Close Vectors
LLL and indefinite LLL
lll Method
lll(L::ZZLat, same_ambient::Bool = true) -> ZZLat
Given an integral L
with basis matrix B
, compute a basis C
of L
such that the gram matrix L
with respect to C
is LLL-reduced.
By default, it creates the lattice in the same ambient space as L
. This can be disabled by setting same_ambient = false
. Works with both definite and indefinite lattices.
Short Vectors
short_vectors Function
short_vectors(L::ZZLat, [lb = 0], ub, [elem_type = ZZRingElem]; check::Bool = true)
-> Vector{Tuple{Vector{elem_type}, QQFieldElem}}
Return all tuples (v, n)
such that lb <= n <= ub
, where G
is the Gram matrix of L
and v
is non-zero.
Note that the vectors are computed up to sign (so only one of v
and -v
appears).
It is assumed and checked that L
is definite.
See also short_vectors_iterator
for an iterator version.
shortest_vectors Function
shortest_vectors(L::ZZLat, [elem_type = ZZRingElem]; check::Bool = true)
-> QQFieldElem, Vector{elem_type}, QQFieldElem}
Return the list of shortest non-zero vectors in absolute value. Note that the vectors are computed up to sign (so only one of v
and -v
appears).
It is assumed and checked that L
is definite.
See also minimum
.
short_vectors_iterator Function
short_vectors_iterator(L::ZZLat, [lb = 0], ub,
[elem_type = ZZRingElem]; check::Bool = true)
-> Tuple{Vector{elem_type}, QQFieldElem} (iterator)
Return an iterator for all tuples (v, n)
such that lb <= n <= ub
, where G
is the Gram matrix of L
and v
is non-zero.
Note that the vectors are computed up to sign (so only one of v
and -v
appears).
It is assumed and checked that L
is definite.
See also short_vectors
.
minimum Method
minimum(L::ZZLat) -> QQFieldElem
Return the minimum absolute squared length among the non-zero vectors in L
.
kissing_number Method
kissing_number(L::ZZLat) -> Int
Return the Kissing number of the sphere packing defined by L
.
This is the number of non-overlapping spheres touching any other given sphere.
sourceVectors with given square and divisibility
vectors_of_square_and_divisibility Function
vectors_of_square_and_divisibility(
L::ZZLat,
S::ZZLat,
n::RationalUnion,
d::RationalUnion = scale(L);
coordinates_representation::Symbol=:S,
check::Bool=true,
) -> Vector{Tuple{QQMatrix}, RationalUnion, RationalUnion}}
Given a nondegenerate
For a vector
The entry n
must be nonzero and d
must be a positive rational number, set to scale(L)
by default.
Note
Alternatively, instead of single values n
and d
one can input:
a list of pairs of rational numbers
(n, d)
wheren
is nonzero andd
is positive;a dictionary whose keys are positive rational numbers
d
and the associated list of numbers consist of nonzero rational numbersn
.
The output consists of a list of triples (v, n', d')
where v
is a vector of
Note
In the case where one wants to choose vectors_of_square_and_divisibility(L, n, d)
.
One can choose in which coordinates system each vector v
in output is represented by changing the symbol coordinates_representation
. There are three possibilities:
coordinates_representation = :L
: the vectorv
is given in terms of its coordinates in the standard basis of the rational span of the lattice; coordinates_representation = :S
(default): the vectorv
is given in terms of its coordinates in the fixed basis of the lattice; coordinates_representation = :ambient
: the vectorv
is given in terms of its coordinates in the standard basis of the ambient space of.
If the keyword argument check
is set to true, the function checks whether
Examples
julia> E6 = root_lattice(:E, 6)
Integer lattice of rank 6 and degree 6
with gram matrix
[ 2 -1 0 0 0 0]
[-1 2 -1 0 0 0]
[ 0 -1 2 -1 0 -1]
[ 0 0 -1 2 -1 0]
[ 0 0 0 -1 2 0]
[ 0 0 -1 0 0 2]
julia> A2 = lattice_in_same_ambient_space(E6, basis_matrix(E6)[1:2, :])
Integer lattice of rank 2 and degree 6
with gram matrix
[ 2 -1]
[-1 2]
julia> C = orthogonal_submodule(E6, A2)
Integer lattice of rank 4 and degree 6
with gram matrix
[12 -3 0 -3]
[-3 2 -1 0]
[ 0 -1 2 0]
[-3 0 0 2]
julia> vectors_of_square_and_divisibility(E6, C, 12, 3; coordinates_representation=:L)
9-element Vector{Tuple{QQMatrix, QQFieldElem, QQFieldElem}}:
([1 2 3 1 -1 3], 12, 3)
([2 4 6 2 1 3], 12, 3)
([1 2 3 1 -1 0], 12, 3)
([2 4 6 5 1 3], 12, 3)
([1 2 3 4 2 0], 12, 3)
([1 2 3 1 2 3], 12, 3)
([2 4 6 5 4 3], 12, 3)
([1 2 3 4 2 3], 12, 3)
([1 2 3 1 2 0], 12, 3)
julia> L = integer_lattice(; gram=matrix(QQ, 2, 2, [2 1; 1 4]))
Integer lattice of rank 2 and degree 2
with gram matrix
[2 1]
[1 4]
julia> vectors_of_square_and_divisibility(L, 8, 2)
1-element Vector{Tuple{QQMatrix, QQFieldElem, QQFieldElem}}:
([2 0], 8, 2)
julia> length(short_vectors(L, 8, 8))
3
Short Vectors affine
enumerate_quadratic_triples Function
enumerate_quadratic_triples(
Q::MatrixElem{T},
b::MatrixElem{T},
c::RationalUnion;
algorithm::Symbol=:short_vectors,
equal::Bool=false
) where T <: Union{ZZRingElem, QQFieldElem}
-> Vector{Tuple{Vector{ZZRingElem}, QQFieldElem}}
Given the Gram matrix
For the moment, only the algorithm :short_vectors
is available. The function uses the data required for the closest vector problem for the quadratic triple [Q, b, c]
; in particular, close_vectors
).
If equal
is set to true
, the function returns only the pairs
short_vectors_affine Method
short_vectors_affine(
S::ZZLat,
v::QQMatrix,
alpha::RationalUnion,
d::RationalUnion
) -> Vector{QQMatrix}
Given a hyperbolic or negative definite
where
The vectors in output are given in terms of their coordinates in the ambient space of
The implementation is based on Algorithm 2.2 in [7]
sourceshort_vectors_affine Method
short_vectors_affine
gram::MatrixElem{T},
v::MatrixElem{T},
alpha::RationalUnion,
d::RationalUnion
) where T <: Union{ZZRingElem, QQFieldElem} -> Vector{MatrixElem{T}}
Given the Gram matrix gram
of a hyperbolic or negative definite
where
The vectors in output are given in terms of their coordinates in the rational span of
The implementation is based on Algorithm 2.2 in [7]
sourceClose Vectors
close_vectors Method
close_vectors(L:ZZLat, v:Vector, [lb,], ub; check::Bool = false)
-> Vector{Tuple{Vector{Int}}, QQFieldElem}
Return all tuples (x, d)
where x
is an element of L
such that d = b(v - x, v - x) <= ub
. If lb
is provided, then also lb <= d
.
If filter
is not nothing
, then only those x
with filter(x)
evaluating to true
are returned.
By default, it will be checked whether L
is positive definite. This can be disabled setting check = false
.
Both input and output are with respect to the basis matrix of L
.
Examples
julia> L = integer_lattice(matrix(QQ, 2, 2, [1, 0, 0, 2]));
julia> close_vectors(L, [1, 1], 1)
3-element Vector{Tuple{Vector{ZZRingElem}, QQFieldElem}}:
([2, 1], 1)
([0, 1], 1)
([1, 1], 0)
julia> close_vectors(L, [1, 1], 1, 1)
2-element Vector{Tuple{Vector{ZZRingElem}, QQFieldElem}}:
([2, 1], 1)
([0, 1], 1)