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
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}}
base_ring
— Methodbase_ring(G::MatrixGroup)
Return the base ring of the matrix group G
.
degree
— Methoddegree(G::MatrixGroup)
Return the degree of the matrix group G
, i.e. the number of rows of its matrices.
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, to get the embedding homomorphism of C
into G
, use
is_subgroup(C, G)[2]
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
matrix
— Methodmatrix(x::MatrixGroupElem)
Return the underlying matrix of x
.
base_ring
— Methodbase_ring(x::MatrixGroupElem)
Return the base ring of the underlying matrix of x
.
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
.
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, the same output returned when M
has type MatElem
.
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
SesquilinearForm
— TypeSesquilinearForm{T<:RingElem}
Type of groups G
of n x n
matrices over the ring R
, where n = degree(G)
and R = base_ring(G)
. At the moment, only rings of type fqPolyRepField
are supported.
is_alternating
— Methodis_alternating(f::SesquilinearForm)
Return whether the form f
is an alternating form.
is_hermitian
— Methodis_hermitian(f::SesquilinearForm)
Return whether the form f
is a hermitian form.
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.
alternating_form
— Methodalternating_form(B::MatElem{T})
Return the alternating form with Gram matrix B
.
symmetric_form
— Methodsymmetric_form(B::MatElem{T})
Return the symmetric form with Gram matrix B
.
hermitian_form
— Methodhermitian_form(B::MatElem{T})
Return the hermitian form with Gram matrix B
.
quadratic_form
— Methodquadratic_form(B::MatElem{T})
Return the quadratic form with Gram matrix B
.
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 as false
, it does not check whether the polynomial is homogeneous of degree 2. To define quadratic forms of dimension 1, f
can also have type PolyRingElem{T}
.
corresponding_bilinear_form
— Methodcorresponding_bilinear_form(Q::SesquilinearForm)
Given a quadratic form Q
, return the bilinear form B
defined by B(u,v) = Q(u+v)-Q(u)-Q(v)
.
corresponding_quadratic_form
— Methodcorresponding_quadratic_form(Q::SesquilinearForm)
Given a symmetric form f
, returns the quadratic form Q
defined by Q(v) = f(v,v)/2
. It is defined only in odd characteristic.
gram_matrix
— Methodgram_matrix(B::SesquilinearForm)
Return the Gram matrix of a sesquilinear or quadratic form B
.
defining_polynomial
— Methoddefining_polynomial(f::SesquilinearForm)
Return the polynomial that defines the quadratic form f
.
radical
— Methodradical(f::SesquilinearForm{T})
Return the radical of the sesquilinear form f
, i.e. the subspace of all v
such that f(u,v)=0
for all u
. The radical of a quadratic form Q
is the set of vectors v
such that Q(v)=0
and v
lies in the radical of the corresponding bilinear form.
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.
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 sesquilinear forms, 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
).
If f
and g
are quadratic forms, return (true
, C
) if there exists a matrix C
such that f^A = ag
for some scalar a
. If such C
does not exist, then return (false
, nothing
).
Invariant forms
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 an invariant bilinear form for the group 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)
.
invariant_sesquilinear_form
— Methodinvariant_sesquilinear_form(G::MatrixGroup)
Return an invariant sesquilinear (non bilinear) form for the group G
. An exception is thrown if the module induced by the action of G
is not absolutely irreducible or if the group 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)
.
invariant_quadratic_form
— Methodinvariant_quadratic_form(G::MatrixGroup)
Return an invariant quadratic form for the group 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)
.
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.
isometry_group
— Methodisometry_group(f::SesquilinearForm{T})
Return the group of isometries for the sesquilinear form f
.
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.
Utilities for matrices
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 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]