Matrix groups
Introduction
A matrix group is a group that consists of invertible square matrices over a common ring, the base ring of the group (see base_ring(G::MatrixGroup{RE}) where RE <: RingElem
).
Matrix groups in Oscar have the type MatrixGroup{RE<:RingElem, T<:MatElem{RE}}
, their elements have the type MatrixGroupElem{RE<:RingElem, T<:MatElem{RE}}
.
Basic Creation
In order to create a matrix group, one first creates matrices that generate the group, and then calls matrix_group
with these generators.
julia> mats = [[0 -1; 1 -1], [0 1; 1 0]]
2-element Vector{Matrix{Int64}}:
[0 -1; 1 -1]
[0 1; 1 0]
julia> matelms = map(m -> matrix(ZZ, m), mats)
2-element Vector{ZZMatrix}:
[0 -1; 1 -1]
[0 1; 1 0]
julia> g = matrix_group(matelms)
Matrix group of degree 2
over integer ring
julia> describe(g)
"S3"
julia> t = matrix_group(ZZ, 2, typeof(matelms[1])[])
Matrix group of degree 2
over integer ring
julia> t == trivial_subgroup(g)[1]
true
julia> F = GF(3); matelms = map(m -> matrix(F, m), mats)
2-element Vector{FqMatrix}:
[0 2; 1 2]
[0 1; 1 0]
julia> g = matrix_group(matelms)
Matrix group of degree 2
over prime field of characteristic 3
julia> describe(g)
"S3"
Alternatively, one can start with a classical matrix group.
julia> g = GL(3, GF(2))
GL(3,2)
julia> F = GF(4)
Finite field of degree 2 and characteristic 2
julia> G = GL(3, F)
GL(3,4)
julia> is_subset(g, G)
false
julia> h = change_base_ring(F, g)
Matrix group of degree 3
over finite field of degree 2 and characteristic 2
julia> flag, mp = is_subgroup(h, G)
(true, Hom: h -> G)
julia> mp(gen(h, 1)) in G
true
Functions for matrix groups
matrix_group
— Methodmatrix_group(R::Ring, m::Int, V::T...) where T<:Union{MatElem,MatrixGroupElem}
matrix_group(R::Ring, m::Int, V::AbstractVector{T}) where T<:Union{MatElem,MatrixGroupElem}
matrix_group(V::T...) where T<:Union{MatElem,MatrixGroupElem}
matrix_group(V::AbstractVector{T}) where T<:Union{MatElem,MatrixGroupElem}
Return the matrix group generated by matrices in V
. If the degree m
and coefficient ring R
are not given, then V
must be non-empty
base_ring
— Methodbase_ring(G::MatrixGroup)
Return the base ring of the matrix group G
.
Examples
julia> F = GF(4); g = general_linear_group(2, F);
julia> base_ring(g) == F
true
degree
— Methoddegree(G::MatrixGroup)
Return the degree of G
, i.e., the number of rows of its matrices.
Examples
julia> degree(GL(4, 2))
4
centralizer
— Methodcentralizer(G::MatrixGroup{T}, x::MatrixGroupElem{T})
Return (C
,f
), where C
is the centralizer of x
in C
and f
is the embedding of C
into G
. If G
= GL(n,F)
or SL(n,F)
, then f
= nothing
. In this case, use is_subgroup(C, G)[2]
to get the embedding homomorphism of C
into G
.
Examples
julia> g = Sp(4, 2); x = gen(g, 1);
julia> C, emb = centralizer(g, x)
(Matrix group of degree 4 over GF(2), Hom: C -> g)
julia> order(C)
8
map_entries
— Methodmap_entries(f, G::MatrixGroup)
Return the matrix group obtained by applying f
element-wise to each generator of G
.
f
can be a ring or a field, a suitable map, or a Julia function.
Examples
julia> mat = matrix(ZZ, 2, 2, [1, 1, 0, 1]);
julia> G = matrix_group(mat);
julia> G2 = map_entries(x -> -x, G)
Matrix group of degree 2
over integer ring
julia> is_finite(G2)
false
julia> order(map_entries(GF(3), G))
3
Elements of matrix groups
The following functions delegate to the underlying matrix of the given matrix group element.
matrix
— Methodmatrix(x::MatrixGroupElem)
Return the underlying matrix of x
.
Examples
julia> F = GF(4); g = general_linear_group(2, F);
julia> x = gen(g, 1)
[o 0]
[0 1]
julia> m = matrix(x)
[o 0]
[0 1]
julia> x == m
false
julia> x == g(m)
true
base_ring
— Methodbase_ring(x::MatrixGroupElem)
Return the base ring of the matrix group to which x
belongs. This is also the base ring of the underlying matrix of x
.
Examples
julia> F = GF(4); g = general_linear_group(2, F);
julia> x = gen(g, 1)
[o 0]
[0 1]
julia> base_ring(x) == F
true
julia> base_ring(x) == base_ring(matrix(x))
true
nrows
— Methodnumber_of_rows(x::MatrixGroupElem)
Return the number of rows of the underlying matrix of x
.
det
— Methoddet(x::MatrixGroupElem)
Return the determinant of the underlying matrix of x
.
tr
— Methodtr(x::MatrixGroupElem)
Return the trace of the underlying matrix of x
.
Examples
julia> F = GF(4); g = general_linear_group(2, F);
julia> x = gen(g, 1)
[o 0]
[0 1]
julia> t = tr(x)
o + 1
julia> t in F
true
multiplicative_jordan_decomposition
— Methodmultiplicative_jordan_decomposition(M::MatrixGroupElem)
Return S
and U
in the group G = parent(M)
such that S
is semisimple, U
is unipotent and M = SU = US
.
this is NOT, in general, equal to multiplicative_jordan_decomposition(matrix(M))
.
is_semisimple
— Methodis_semisimple(x::MatrixGroupElem{T}) where T <: FinFieldElem
Return whether x
is semisimple, i.e. has order coprime with the characteristic of its base ring.
is_unipotent
— Methodis_unipotent(x::MatrixGroupElem{T}) where T <: FinFieldElem
Return whether x
is unipotent, i.e. its order is a power of the characteristic of its base ring.
Sesquilinear forms
Sesquilinear forms are alternating and symmetric bilinear forms, hermitian forms, and quadratic forms.
Sesquilinear forms in Oscar have the type SesquilinearForm{T<:RingElem}
.
A sesquilinear form can be created from its Gram matrix or as an invariant form of a matrix group.
Create sesquilinear forms from Gram matrices
alternating_form
— Methodalternating_form(B::MatElem{T})
Return the alternating form with Gram matrix B
. An exception is thrown if B
is not square or does not have zeros on the diagonal or does not satisfy B = -transpose(B)
.
Examples
julia> f = alternating_form(matrix(GF(3), 2, 2, [0, 1, -1, 0]))
alternating form with Gram matrix
[0 1]
[2 0]
julia> describe(isometry_group(f))
"SL(2,3)"
symmetric_form
— Methodsymmetric_form(B::MatElem{T})
Return the symmetric form with Gram matrix B
. An exception is thrown if B
is not square or not symmetric.
Examples
julia> f = symmetric_form(matrix(GF(3), 2, 2, [0, 1, 1, 0]))
symmetric form with Gram matrix
[0 1]
[1 0]
julia> describe(isometry_group(f))
"C2 x C2"
hermitian_form
— Methodhermitian_form(B::MatElem{T})
Return the hermitian form with Gram matrix B
. An exception is thrown if B
is not square or does not satisfy B = conjugate_transpose(B)
, see conjugate_transpose
.
Examples
julia> F = GF(4); z = gen(F);
julia> f = hermitian_form(matrix(F, 2, 2, [0, z, z^2, 0]))
hermitian form with Gram matrix
[ 0 o]
[o + 1 0]
julia> describe(isometry_group(f))
"C3 x S3"
quadratic_form
— Methodquadratic_form(B::MatElem{T})
Return the quadratic form with Gram matrix B
.
Examples
julia> f = quadratic_form(matrix(GF(3), 2, 2, [2, 2, 0, 1]))
quadratic form with Gram matrix
[2 2]
[0 1]
julia> describe(isometry_group(f))
"D8"
quadratic_form
— Methodquadratic_form(f::MPolyRingElem{T}; check=true)
Return the quadratic form described by the polynomial f
. Here, f
must be a homogeneous polynomial of degree 2. If check
is set to false
, it is not checked whether f
is homogeneous of degree 2. To define quadratic forms of dimension 1, f
can also have type PolyRingElem{T}
.
Examples
julia> _, (x1, x2) = polynomial_ring(GF(3), [:x1, :x2]);
julia> f = quadratic_form(2*x1^2 + 2*x1*x2 + x2^2)
quadratic form with Gram matrix
[2 2]
[0 1]
julia> describe(isometry_group(f))
"D8"
Functions for sesquilinear forms
is_alternating
— Methodis_alternating(f::SesquilinearForm)
Return whether the form f
is an alternating form, see is_alternating(B::MatElem)
.
is_hermitian
— Methodis_hermitian(f::SesquilinearForm)
Return whether the form f
is a hermitian form, see is_hermitian(B::MatElem{T}) where T <: FinFieldElem
.
is_quadratic
— Methodis_quadratic(f::SesquilinearForm)
Return whether the form f
is a quadratic form.
is_symmetric
— Methodis_symmetric(f::SesquilinearForm)
Return whether the form f
is a symmetric form.
corresponding_bilinear_form
— Methodcorresponding_bilinear_form(Q::SesquilinearForm)
Given a quadratic form Q
, return the symmetric form with Gram matrix B
defined by B(u,v) = Q(u+v)-Q(u)-Q(v)
.
Examples
Q = quadratic_form(invariant_quadratic_form(GO(3,3)))
corresponding_bilinear_form(Q)
corresponding_quadratic_form
— Methodcorresponding_quadratic_form(f::SesquilinearForm)
Given a symmetric form f
, return the quadratic form Q
defined by Q(v) = f(v,v)/2
. It is defined only in odd characteristic.
Examples
f = symmetric_form(invariant_bilinear_form(GO(3, 3)))
corresponding_quadratic_form(f)
gram_matrix
— Methodgram_matrix(f::SesquilinearForm)
Return the Gram matrix that defines f
.
defining_polynomial
— Methoddefining_polynomial(f::SesquilinearForm)
Return the polynomial that defines f
.
radical
— Methodradical(f::SesquilinearForm{T})
Return (U, emb)
where U
is the radical of f
and emb
is the embedding of U
into the full vector space on which f
is defined.
For a quadratic form f
, the radical is defined as the set of vectors v
such that f(v) = 0
and v
lies in the radical of the corresponding bilinear form.
Otherwise the radical of f
is defined as the subspace of all v
such that f(u, v) = 0
for all u
.
Examples
julia> g = GO(7, 2);
julia> f = symmetric_form(invariant_symmetric_forms(g)[1]);
julia> U, emb = radical(f);
julia> dim(U)
1
witt_index
— Methodwitt_index(f::SesquilinearForm{T})
Return the Witt index of the form induced by f
on V/Rad(f)
. The Witt index is the dimension of a maximal totally isotropic (singular for quadratic forms) subspace.
Examples
julia> g = GO(1, 6, 2);
julia> witt_index(quadratic_form(invariant_quadratic_forms(g)[1]))
3
julia> g = GO(-1, 6, 2);
julia> witt_index(quadratic_form(invariant_quadratic_forms(g)[1]))
2
is_degenerate
— Methodis_degenerate(f::SesquilinearForm{T})
Return whether f
is degenerate, i.e. f
has nonzero radical. A quadratic form is degenerate if the corresponding bilinear form is.
is_singular
— Methodis_singular(Q::SesquilinearForm{T})
For a quadratic form Q
, return whether Q
is singular, i.e. Q
has nonzero radical.
is_congruent
— Methodis_congruent(f::SesquilinearForm{T}, g::SesquilinearForm{T}) where T <: RingElem
If f
and g
are quadratic forms, return (true
, C
) if there exists a matrix C
such that f^C = ag
for some scalar a
. If such C
does not exist, then return (false
, nothing
).
Otherwise return (true
, C
) if there exists a matrix C
such that f^C = g
, or equivalently, CBC* = A
, where A
and B
are the Gram matrices of f
and g
respectively, and C*
is the transpose-conjugate matrix of C
. If such C
does not exist, then return (false
, nothing
).
isometry_group
— Methodisometry_group(f::SesquilinearForm{T})
Return the group of isometries for the sesquilinear form f
.
isometry_group(L::AbstractLat; depth::Int = -1, bacher_depth::Int = 0) -> MatrixGroup
Return the group of isometries of the lattice L
.
The transformations are represented with respect to the ambient space 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
.
Invariant forms
The following functions compute (Gram matrices of) sesquilinear forms that are invariant under the given matrix group.
invariant_bilinear_forms
— Methodinvariant_bilinear_forms(G::MatrixGroup)
Return a generating set for the vector spaces of bilinear forms preserved by the group G
.
At the moment, elements of the generating set are returned of type mat_elem_type(G)
.
invariant_sesquilinear_forms
— Methodinvariant_sesquilinear_forms(G::MatrixGroup)
Return a generating set for the vector spaces of sesquilinear non-bilinear forms preserved by the group G
. An exception is thrown if base_ring(G)
is not a finite field with even degree over its prime subfield.
At the moment, elements of the generating set are returned of type mat_elem_type(G)
.
invariant_quadratic_forms
— Methodinvariant_quadratic_forms(G::MatrixGroup)
Return a generating set for the vector spaces of quadratic forms preserved by the group G
.
At the moment, elements of the generating set are returned of type mat_elem_type(G)
.
invariant_symmetric_forms
— Methodinvariant_symmetric_forms(G::MatrixGroup)
Return a generating set for the vector spaces of symmetric forms preserved by the group G
.
At the moment, elements of the generating set are returned of type mat_elem_type(G)
.
Work properly only in odd characteristic. In even characteristic, only alternating forms are found.
invariant_alternating_forms
— Methodinvariant_alternating_forms(G::MatrixGroup)
Return a generating set for the vector spaces of alternating forms preserved by the group G
.
At the moment, elements of the generating set are returned of type mat_elem_type(G)
.
invariant_hermitian_forms
— Methodinvariant_hermitian_forms(G::MatrixGroup)
Return a generating set for the vector spaces of hermitian forms preserved by the group G
. An exception is thrown if base_ring(G)
is not a finite field with even degree over its prime subfield.
At the moment, elements of the generating set are returned of type mat_elem_type(G)
.
invariant_bilinear_form
— Methodinvariant_bilinear_form(G::MatrixGroup)
Return the Gram matrix of an invariant bilinear form for G
. An exception is thrown if the module induced by the action of G
is not absolutely irreducible.
At the moment, the output is returned of type mat_elem_type(G)
.
Examples
julia> invariant_bilinear_form(Sp(4, 2))
[0 0 0 1]
[0 0 1 0]
[0 1 0 0]
[1 0 0 0]
invariant_sesquilinear_form
— Methodinvariant_sesquilinear_form(G::MatrixGroup)
Return the Gram matrix of an invariant sesquilinear (non bilinear) form for G
.
An exception is thrown if the module induced by the action of G
is not absolutely irreducible or if G
is defined over a finite field of odd degree over the prime field.
At the moment, the output is returned of type mat_elem_type(G)
.
Examples
julia> invariant_sesquilinear_form(GU(4, 2))
[0 0 0 1]
[0 0 1 0]
[0 1 0 0]
[1 0 0 0]
invariant_quadratic_form
— Methodinvariant_quadratic_form(G::MatrixGroup)
Return the Gram matrix of an invariant quadratic form for G
. An exception is thrown if the module induced by the action of G
is not absolutely irreducible.
At the moment, the output is returned of type mat_elem_type(G)
.
Examples
julia> invariant_quadratic_form(GO(1, 4, 2))
[0 1 0 0]
[0 0 0 0]
[0 0 0 1]
[0 0 0 0]
preserved_quadratic_forms
— Methodpreserved_quadratic_forms(G::MatrixGroup)
Uses random methods to find all of the quadratic forms preserved by G
up to a scalar (i.e. such that G
is a group of similarities for the forms). Since the procedure relies on a pseudo-random generator, the user may need to execute the operation more than once to find all invariant quadratic forms.
preserved_sesquilinear_forms
— Methodpreserved_sesquilinear_forms(G::MatrixGroup)
Uses random methods to find all of the sesquilinear forms preserved by G
up to a scalar (i.e. such that G
is a group of similarities for the forms). Since the procedure relies on a pseudo-random generator, the user may need to execute the operation more than once to find all invariant sesquilinear forms.
orthogonal_sign
— Methodorthogonal_sign(G::MatrixGroup)
For absolutely irreducible G
of degree n
and such that base_ring(G)
is a finite field, return
nothing
ifG
does not preserve a nonzero quadratic form,0
ifn
is odd andG
preserves a nonzero quadratic form,1
ifn
is even andG
preserves a nonzero quadratic form of+
type,-1
ifn
is even andG
preserves a nonzero quadratic form of-
type.
Examples
julia> orthogonal_sign(GO(1, 4, 2))
1
julia> orthogonal_sign(GO(-1, 4, 2))
-1
Utilities for matrices
(The inputs/outputs of the functions described in this section are matrices not matrix group elements.)
pol_elementary_divisors
— Methodpol_elementary_divisors(x::MatElem)
pol_elementary_divisors(x::MatrixGroupElem)
Return a list of pairs (f_i,m_i)
, for irreducible polynomials f_i
and positive integers m_i
, where the f_i^m_i
are the elementary divisors of x
.
generalized_jordan_block
— Methodgeneralized_jordan_block(f::T, n::Int) where T<:PolyRingElem
Return the Jordan block of dimension n
corresponding to the polynomial f
.
generalized_jordan_form
— Methodgeneralized_jordan_form(A::MatElem{T}; with_pol::Bool=false) where T
Return (J
,Z
), where Z^-1*J*Z = A
and J
is a diagonal join of Jordan blocks (corresponding to irreducible polynomials).
matrix
— Methodmatrix(A::Vector{AbstractAlgebra.Generic.FreeModuleElem})
Return the matrix whose rows are the vectors in A
. All vectors in A
must have the same length and the same base ring.
upper_triangular_matrix
— Methodupper_triangular_matrix(L::AbstractVector{T}) where {T <: NCRingElement}
Return the $n$ by $n$ matrix whose entries on and above the main diagonal are the elements of L
, and which has zeroes elsewhere. The value of $n$ is determined by the condition that L
has length $n(n+1)/2$.
An exception is thrown if there is no integer $n$ with this property.
Examples
julia> upper_triangular_matrix([1, 2, 3])
[1 2]
[0 3]
lower_triangular_matrix
— Methodlower_triangular_matrix(L::AbstractVector{T}) where {T <: NCRingElement}
Return the $n$ by $n$ matrix whose entries on and below the main diagonal are the elements of L
, and which has zeroes elsewhere. The value of $n$ is determined by the condition that L
has length $n(n+1)/2$.
An exception is thrown if there is no integer $n$ with this property.
Examples
julia> lower_triangular_matrix([1, 2, 3])
[1 0]
[2 3]
conjugate_transpose
— Methodconjugate_transpose(x::MatElem{T}) where T <: FinFieldElem
If the base ring of x
is GF(q^2)
, return the matrix transpose( map ( y -> y^q, x) )
. An exception is thrown if the base ring does not have even degree.
complement
— Methodcomplement(V::AbstractAlgebra.Generic.FreeModule{T}, W::AbstractAlgebra.Generic.Submodule{T}) where T <: FieldElem
Return a complement for W
in V
, i.e. a subspace U
of V
such that V
is the direct sum of U
and W
.
permutation_matrix
— Methodpermutation_matrix(F::Ring, Q::AbstractVector{T}) where T <: Int
permutation_matrix(F::Ring, p::PermGroupElem)
Return the permutation matrix over the ring R
corresponding to the sequence Q
or to the permutation p
. If Q
is a sequence, then Q
must contain exactly once every integer from 1 to some n
.
Examples
julia> s = perm([3,1,2])
(1,3,2)
julia> permutation_matrix(QQ,s)
[0 0 1]
[1 0 0]
[0 1 0]
is_alternating
— Methodis_alternating(B::MatElem)
Return whether the form corresponding to the matrix B
is alternating, i.e. B = -transpose(B)
and B
has zeros on the diagonal. Return false
if B
is not a square matrix.
is_hermitian
— Methodis_hermitian(B::MatElem{T}) where T <: FinFieldElem
Return whether the matrix B
is hermitian, i.e. B = conjugate_transpose(B)
. Return false
if B
is not a square matrix, or the field has not even degree.
Classical groups
general_linear_group
— Methodgeneral_linear_group(n::Int, q::Int)
general_linear_group(n::Int, R::Ring)
GL = general_linear_group
Return the general linear group of dimension n
over the ring R
respectively the field GF(q)
.
Currently R
must be either a finite field or a residue ring or ZZ
.
Examples
julia> F = GF(7)
Prime field of characteristic 7
julia> H = general_linear_group(2, F)
GL(2,7)
julia> gens(H)
2-element Vector{MatrixGroupElem{FqFieldElem, FqMatrix}}:
[3 0; 0 1]
[6 1; 6 0]
julia> order(general_linear_group(2, residue_ring(ZZ, 6)[1]))
288
special_linear_group
— Methodspecial_linear_group(n::Int, q::Int)
special_linear_group(n::Int, R::Ring)
SL = special_linear_group
Return the special linear group of dimension n
over the ring R
respectively the field GF(q)
.
Currently R
must be either a finite field or a residue ring or ZZ
.
Examples
julia> F = GF(7)
Prime field of characteristic 7
julia> H = special_linear_group(2, F)
SL(2,7)
julia> gens(H)
2-element Vector{MatrixGroupElem{FqFieldElem, FqMatrix}}:
[3 0; 0 5]
[6 1; 6 0]
julia> order(special_linear_group(2, residue_ring(ZZ, 6)[1]))
144
symplectic_group
— Methodsymplectic_group(n::Int, q::Int)
symplectic_group(n::Int, R::Ring)
Sp = symplectic_group
Return the symplectic group of dimension n
over the ring R
respectively the field GF(q)
. The dimension n
must be even.
Currently R
must be either a finite field or a residue ring of prime power order.
Examples
julia> F = GF(7)
Prime field of characteristic 7
julia> H = symplectic_group(2, F)
Sp(2,7)
julia> gens(H)
2-element Vector{MatrixGroupElem{FqFieldElem, FqMatrix}}:
[3 0; 0 5]
[6 1; 6 0]
julia> order(symplectic_group(2, residue_ring(ZZ, 4)[1]))
48
orthogonal_group
— Methodorthogonal_group(e::Int, n::Int, R::Ring)
orthogonal_group(e::Int, n::Int, q::Int)
GO = orthogonal_group
Return the orthogonal group of dimension n
over the ring R
respectively the field GF(q)
, and of type e
, where e
in {+1
,-1
} for n
even and e
=0
for n
odd. If n
is odd, e
can be omitted.
Currently R
must be either a finite field or a residue ring of odd prime power order.
Examples
julia> F = GF(7)
Prime field of characteristic 7
julia> H = orthogonal_group(1, 2, F)
GO+(2,7)
julia> gens(H)
2-element Vector{MatrixGroupElem{FqFieldElem, FqMatrix}}:
[3 0; 0 5]
[0 1; 1 0]
julia> order(orthogonal_group(-1, 2, residue_ring(ZZ, 9)[1]))
24
special_orthogonal_group
— Methodspecial_orthogonal_group(e::Int, n::Int, R::Ring)
special_orthogonal_group(e::Int, n::Int, q::Int)
SO = special_orthogonal_group
Return the special orthogonal group of dimension n
over the ring R
respectively the field GF(q)
, and of type e
, where e
in {+1
,-1
} for n
even and e
=0
for n
odd. If n
is odd, e
can be omitted.
Currently R
must be either a finite field or a residue ring of odd prime power order.
Examples
julia> F = GF(7)
Prime field of characteristic 7
julia> H = special_orthogonal_group(1, 2, F)
SO+(2,7)
julia> gens(H)
3-element Vector{MatrixGroupElem{FqFieldElem, FqMatrix}}:
[3 0; 0 5]
[5 0; 0 3]
[1 0; 0 1]
julia> order(special_orthogonal_group(-1, 2, residue_ring(ZZ, 9)[1]))
12
omega_group
— Methodomega_group(e::Int, n::Int, R::Ring)
omega_group(e::Int, n::Int, q::Int)
Return the Omega group of dimension n
over the field GF(q)
of type e
, where e
in {+1
,-1
} for n
even and e
=0
for n
odd. If n
is odd, e
can be omitted.
Currently R
must be either a finite field or a residue ring of odd prime power order.
Examples
julia> F = GF(7)
Prime field of characteristic 7
julia> H = omega_group(1, 2, F)
Omega+(2,7)
julia> gens(H)
1-element Vector{MatrixGroupElem{FqFieldElem, FqMatrix}}:
[2 0; 0 4]
julia> order(omega_group(0, 3, residue_ring(ZZ, 9)[1]))
324
unitary_group
— Methodunitary_group(n::Int, q::Int)
GU = unitary_group
Return the unitary group of dimension n
over the field GF(q^2)
.
Examples
julia> H = unitary_group(2,3)
GU(2,3)
julia> gens(H)
2-element Vector{MatrixGroupElem{FqFieldElem, FqMatrix}}:
[o 0; 0 2*o]
[2 2*o+2; 2*o+2 0]
special_unitary_group
— Methodspecial_unitary_group(n::Int, q::Int)
SU = special_unitary_group
Return the special unitary group of dimension n
over the field with q^2
elements.
Examples
julia> H = special_unitary_group(2,3)
SU(2,3)
julia> gens(H)
2-element Vector{MatrixGroupElem{FqFieldElem, FqMatrix}}:
[1 2*o+2; 0 1]
[0 2*o+2; 2*o+2 0]
Technicalities
MatrixGroup
— TypeMatrixGroup{RE<:RingElem, T<:MatElem{RE}} <: GAPGroup
Type of groups G
of n x n
matrices over the ring R
, where n = degree(G)
and R = base_ring(G)
.
MatrixGroupElem
— TypeMatrixGroupElem{RE<:RingElem, T<:MatElem{RE}} <: AbstractMatrixGroupElem
Elements of a group of type MatrixGroup{RE<:RingElem, T<:MatElem{RE}}
ring_elem_type
— Methodring_elem_type(G::MatrixGroup{S,T}) where {S,T}
ring_elem_type(::Type{MatrixGroup{S,T}}) where {S,T}
Return the type S
of the entries of the elements of G
. One can enter the type of G
instead of G
.
Examples
julia> g = GL(2, 3);
julia> ring_elem_type(typeof(g)) == elem_type(typeof(base_ring(g)))
true
mat_elem_type
— Methodmat_elem_type(G::MatrixGroup{S,T}) where {S,T}
mat_elem_type(::Type{MatrixGroup{S,T}}) where {S,T}
Return the type T
of matrix(x)
, for elements x
of G
. One can enter the type of G
instead of G
.
Examples
julia> g = GL(2, 3);
julia> mat_elem_type(typeof(g)) == typeof(matrix(one(g)))
true
SesquilinearForm
— TypeSesquilinearForm{T<:RingElem}
Type of alternating and symmetric bilinear forms, hermitian forms, and quadratic forms, defined by a Gram matrix (or, in the case of a quadratic form, alternatively by a polynomial).