Introduction
Nemo allow the creation of dense matricses over any computable ring . There are two different kinds of implementation: a generic one for the case where no specific implementation exists, and efficient implementations of matrices over numerous specific rings, usually provided by C/C++ libraries.
The following table shows each of the matrix types available in Nemo, the base ring , and the Julia/Nemo types for that kind of matrix (the type information is mainly of concern to developers).
Base ring | Library | Element type | Parent type |
---|---|---|---|
Generic ring | Nemo | GenMat{T} | GenMatSpace{T} |
Flint | fmpz_mat | FmpzMatSpace | |
(small ) | Flint | nmod_mat | NmodMatSpace |
Flint | fmpq_mat | FmpqMatSpace | |
Arb | arb_mat | ArbMatSpace | |
Arb | acb_mat | AcbMatSpace |
The dimensions and base ring of a generic matrix are stored in its parent object.
All matrix element types belong to the abstract type MatElem
and all of the matrix space types belong to the abstract type MatSpace
. This enables one to write generic functions that can accept any Nemo matrix type.
Matrix space constructors
In Nemo we have the concept of a matrix space. This is the collection of matrices with specified dimensions and base ring.
In order to construct matrices in Nemo, one usually first constructs the matrix space itself. This is accomplished with the following constructor.
#Nemo.MatrixSpace
— Method.
MatrixSpace(R::Ring, r::Int, c::Int, cached::Bool = true)
Return parent object corresponding to the space of matrices over the ring . If
cached == true
(the default), the returned parent object is cached so that it can returned by future calls to the constructor with the same dimensions and base ring.
We also allow matrices over a given base ring to be constructed directly. In such cases, Nemo automatically constructs the matrix space internally. See the matrix element constructors below for examples. However, note that there may be a small peformance disadvantage to doing it that way, since the matrix space needs to be looked up internally every time a matrix is constructed.
Here are some examples of creating matrix spaces and making use of the resulting parent objects to coerce various elements into the matrix space.
R, t = PolynomialRing(QQ, "t")
S = MatrixSpace(R, 3, 3)
A = S()
B = S(12)
C = S(R(11))
Matrix element constructors
Once a matrix space is constructed, there are various ways to construct matrices in that space.
In addition to coercing elements into the matrix space as above, we provide the following functions for constructing certain useful matrices.
#Base.zero
— Method.
zero(a::MatSpace)
Construct the zero matrix in the given matrix space.
#Base.one
— Method.
one(a::MatSpace)
Construct the matrix in the given matrix space with ones down the diagonal and zeroes elsewhere.
In addition, there are various shorthand notations for constructing matrices over a given base ring without first constructing the matrix space parent object.
R[a b c...;...]
Create the matrix over the base ring consisting of the given rows (separated by semicolons). Each entry is coerced into automatically. Note that parentheses may be placed around individual entries if the lists would otherwise be ambiguous, e.g. R[1 2; 2 (-3)]
.
Beware that this syntax does not support the creation of column vectors. See the notation below for creating those.
R[a b c...]
Create the row vector with entries in consisting of the given entries (separated by spaces). Each entry is coerced into automatically. Note that parentheses may be placed around individual entries if the list would otherwise be ambiguous, e.g. R[1 2 (-3)]
.
R[a b c...]'
Create the column vector with entries in consisting of the given entries (separated by spaces). Each entry is coerced into automatically. Observe the dash that is used to transpose the row vector notation (for free) to turn it into a column vector. Note that parentheses may be placed around individual entries if the list would otherwise be ambiguous, e.g. R[1 2 (-3)]'
.
Here are some examples of constructing matrices.
R, t = PolynomialRing(QQ, "t")
S = MatrixSpace(R, 3, 3)
A = zero(S)
B = one(S)
C = R[t + 1 1; t^2 0]
D = R[t + 1 2 t]
F = R[1 2 t]'
Basic functionality
All matric modules in Nemo must provide the functionality listed in this section. (Note that only some of these functions are useful to a user.)
Developers who are writing their own matrix module, whether as an interface to a C library, or as some kind of generic module, must provide all of these functions for custom matrix types in Nemo.
We write U
for the type of the matrices in the matrix space and T
for the type of elements of the base ring.
All of these functions are provided for all existing matrix types in Nemo.
parent_type{U <: MatElem}(::Type{U})
Given the type of matrix elements, should return the type of the corresponding parent object.
elem_type(R::MatSpace)
Given a parent object for the matrix space, return the type of elements of the matrix space.
Base.hash(a::MatElem, h::UInt)
Return a UInt
hexadecimal hash of the matrix . This should be xor'd with a fixed random hexadecimal specific to the matrix type. The hash of each entry should be xor'd with the supplied parameter h
as part of computing the hash.
deepcopy(a::MatElem)
Construct a copy of the given matrix and return it. This function must recursively construct copies of all of the internal data in the given matrix. Nemo matricess are mutable and so returning shallow copies is not sufficient.
To access entries of a Nemo matrix, we overload the square bracket notation.
M[r::Int, c::Int]
One can both assign to and access a given entry at row and column of a matrix with this notation. Note that Julia and Nemo matrices are -indexed, i.e. the first row has index , not , etc. This is in accordance with many papers on matrices and with systems such as Pari/GP.
Given a parent object S
for a matrix space, the following coercion functions are provided to coerce various elements into the matrix space. Developers provide these by overloading the call
operator for the matrix parent objects.
S()
Coerce zero into the space .
S(n::Integer)
S(n::fmpz)
Return the diagonal matrix with the given integer along the diagonal and zeroes elsewhere.
S(n::T)
Coerces an element of the base ring, of type T
into .
S(A::Array{T, 2})
Take a Julia two dimensional array of elements in the base ring, of type T
and construct the matrix with those entries.
S(f::MatElem)
Take a matrix that is already in the space and simply return it. A copy of the original is not made.
S(c::RingElem)
Try to coerce the given ring element into the matrix space (as a diagonal matrix). This only succeeds if can be coerced into the base ring.
In addition to the above, developers of custom matrices must ensure the parent object of a matrix type constains a field base_ring
specifying the base ring, and fields rows
and cols
to specify the dimensions. They must also ensure that each matrix element contains a field parent
specifying the parent object of the matrix, or that there is at least a function parent(a::MatElem)
which returns the parent of the given matrix.
Typically a developer will also overload the MatrixSpace
generic function to create matrices of the custom type they are implementing.
Basic manipulation
Numerous functions are provided to manipulate matricess and to set and retrieve entries and other basic data associated with the matrices. Also see the section on basic functionality above.
base_ring(::MatSpace)
#Nemo.base_ring
— Method.
base_ring(r::MatElem)
Return the base ring of the matrix space that the supplied matrix belongs to.
#Base.parent
— Method.
parent(a::MatElem)
Return the parent object of the given matrix.
#Nemo.rows
— Method.
rows(a::MatElem)
Return the number of rows of the given matrix.
#Nemo.cols
— Method.
cols(a::MatElem)
Return the number of columns of the given matrix.
#Nemo.iszero
— Method.
iszero(a::MatElem)
Return
true
if the supplied matrix is the zero matrix, otherwise returnfalse
.
#Nemo.isone
— Method.
isone(a::MatElem)
Return
true
if the supplied matrix is diagonal with ones along the diagonal, otherwise returnfalse
.
Here are some examples of basic manipulation of matrices.
R, t = PolynomialRing(QQ, "t")
S = MatrixSpace(R, 3, 3)
A = S([t + 1 t R(1); t^2 t t; R(-2) t + 2 t^2 + t + 1])
B = S([R(2) R(3) R(1); t t + 1 t + 2; R(-1) t^2 t^3])
C = zero(S)
D = one(S)
f = iszero(C)
g = isone(D)
r = rows(B)
c = cols(B)
U = base_ring(C)
V = base_ring(S)
W = parent(D)
Arithmetic operators
All the usual arithmetic operators are overloaded for Nemo matrices. Note that Julia uses the single slash for floating point division. Therefore to perform exact division by a constant we use divexact
.
Function | Operation |
---|---|
-(a::MatElem) | unary minus |
+{T <: RingElem}(a::MatElem{T}, b::MatElem{T}) | addition |
-{T <: RingElem}(a::MatElem{T}, b::MatElem{T}) | subtraction |
*{T <: RingElem}(a::MatElem{T}, b::MatElem{T}) | multiplication |
divexact{T <: RingElem}(a::MatElem{T}, b::MatElem{T}) | exact division |
An exception is raised if the matrix dimensions are not compatible for the given operation. The divexact
function computes a*inv(b)
where inv(b)
is the inverse of the matrix . This assumes that can be inverted.
The following ad hoc operators are also provided.
Function | Operation |
---|---|
+(a::Integer, b::MatElem) | addition |
+(a::MatElem, b::Integer) | addition |
+(a::fmpz, b::MatElem) | addition |
+(a::MatElem, b::fmpz) | addition |
+{T <: RingElem}(a::T, b::MatElem{T}) | addition |
+{T <: RingElem}(a::MatElem{T}, b::T) | addition |
-(a::Integer, b::MatElem) | subtraction |
-(a::MatElem, b::Integer) | subtraction |
-(a::fmpz, b::MatElem) | subtraction |
-(a::MatElem, b::fmpz) | subtraction |
-{T <: RingElem}(a::T, b::MatElem{T}) | subtraction |
-{T <: RingElem}(a::MatElem{T}, b::T) | subtraction |
*(a::Integer, b::MatElem) | multiplication |
*(a::MatElem, b::Integer) | multiplication |
*(a::fmpz, b::MatElem) | multiplication |
*(a::MatElem, b::fmpz) | multiplication |
*{T <: RingElem}(a::T, b::MatElem{T}) | multiplication |
*{T <: RingElem}(a::MatElem{T}, b::T) | multiplication |
divexact(a::MatElem, b::Integer) | exact division |
divexact(a::MatElem, b::fmpz) | exact division |
divexact{T <: RingElem}(a::MatElem{T}, b::T) | exact division |
^(a::MatElem, n::Int) | powering |
The following function is also provided.
#Nemo.powers
— Method.
powers{T <: RingElem}(a::MatElem{T}, d::Int)
Return an array of matrices wher for
If the appropriate promote_rule
and coercion exists, these operators can also be used with elements of other rings. Nemo will try to coerce the operands to the dominating type and then apply the operator.
Here are some examples of arithmetic operations on matrices.
R, t = PolynomialRing(QQ, "t")
S = MatrixSpace(R, 3, 3)
A = S([t + 1 t R(1); t^2 t t; R(-2) t + 2 t^2 + t + 1])
B = S([R(2) R(3) R(1); t t + 1 t + 2; R(-1) t^2 t^3])
C = -A
D = A + B
F = A - B
G = A*B
H = 3*A
K = B + 2
M = fmpz(3) - B
N = t - A
P = A^3
Q = powers(A, 3)
R = divexact(A*3, 3)
Comparison operators
The following comparison operators are implemented for matrices in Nemo.
Function
isequal{T <: RingElem}(a::MatElem{T}, b::MatElem{T})
=={T <: RingElem}(a::MatElem{T}, b::MatElem{T})
The isequal
operation returns true
if and only if all the entries of the matrix are precisely equal as compared by isequal
. This is a stronger form of equality, used for comparing inexact coefficients, such as elements of a power series ring, the -adics, or the reals or complex numbers. Two elements are precisely equal only if they have the same precision or bounds in addition to being arithmetically equal.
#Nemo.overlaps
— Method.
overlaps(x::arb_mat, y::arb_mat)
Returns
true
if all entries of overlap with the corresponding entry of , otherwise returnfalse
.
#Nemo.overlaps
— Method.
overlaps(x::acb_mat, y::acb_mat)
Returns
true
if all entries of overlap with the corresponding entry of , otherwise returnfalse
.
#Base.contains
— Method.
contains(x::arb_mat, y::arb_mat)
Returns
true
if all entries of contain the corresponding entry of , otherwise returnfalse
.
#Base.contains
— Method.
contains(x::acb_mat, y::acb_mat)
Returns
true
if all entries of contain the corresponding entry of , otherwise returnfalse
.
In addition we have the following ad hoc comparison operators.
Function
=={T <: RingElem}(a::MatElem{T}, b::T)
=={T <: RingElem}(a::T, b::MatElem{T})
==(a::MatElem, b::Integer)
==(a::Integer, b::MatElem)
==(a::MatElem, b::fmpz)
==(a::fmpz, b::MatElem)
Here are some examples of comparisons.
R, t = PolynomialRing(QQ, "t")
S = MatrixSpace(R, 3, 3)
A = S([t + 1 t R(1); t^2 t t; R(-2) t + 2 t^2 + t + 1])
B = S([R(2) R(3) R(1); t t + 1 t + 2; R(-1) t^2 t^3])
A != B
A == deepcopy(A)
A != 12
fmpz(11) != A
B != t
S(11) == 11
C = RR[1 2; 3 4]
D = RR["1 +/- 0.1" "2 +/- 0.1"; "3 +/- 0.1" "4 +/- 0.1"]
overlaps(C, D)
contains(D, C)
Scaling
#Base.:<<
— Method.
<<(x::fmpz_mat, y::Int)
Return .
>>(x::fmpz_mat, y::Int)
Return where rounding is towards zero.
Here are some examples of scaling matrices.
S = MatrixSpace(ZZ, 3, 3)
A = S([fmpz(2) 3 5; 1 4 7; 9 6 3])
B = A<<5
C = B>>2
Transpose
#Base.transpose
— Method.
transpose(x::MatElem)
Return the transpose of the given matrix.
Here is an example of transposing a matrix.
R, t = PolynomialRing(QQ, "t")
S = MatrixSpace(R, 3, 3)
A = S([t + 1 t R(1); t^2 t t; R(-2) t + 2 t^2 + t + 1])
B = transpose(A)
Gram matrix
#Nemo.gram
— Method.
gram(x::MatElem)
Return the Gram matrix of , i.e. if is an matrix return the matrix whose entries are the dot products of the -th and -th rows, respectively.
Here is an example of computing the Gram matrix.
R, t = PolynomialRing(QQ, "t")
S = MatrixSpace(R, 3, 3)
A = S([t + 1 t R(1); t^2 t t; R(-2) t + 2 t^2 + t + 1])
B = gram(A)
Trace
#Base.LinAlg.trace
— Method.
trace(x::MatElem)
Return the trace of the matrix , i.e. the sum of the diagonal elements. We require the matrix to be square.
Here is an example of computing the trace.
R, t = PolynomialRing(QQ, "t")
S = MatrixSpace(R, 3, 3)
A = S([t + 1 t R(1); t^2 t t; R(-2) t + 2 t^2 + t + 1])
b = trace(A)
Content
#Nemo.content
— Method.
content(x::MatElem)
Return the content of the matrix , i.e. the greatest common divisor of all its entries, assuming it exists.
Here is an example of computing the content of a matrix.
R, t = PolynomialRing(QQ, "t")
S = MatrixSpace(R, 3, 3)
A = S([t + 1 t R(1); t^2 t t; R(-2) t + 2 t^2 + t + 1])
b = content(A)
Concatenation
#Base.hcat
— Method.
hcat(a::MatElem, b::MatElem)
Return the horizontal concatenation of and . Assumes that the number of rows is the same in and .
#Base.vcat
— Method.
vcat(a::MatElem, b::MatElem)
Return the vertical concatenation of and . Assumes that the number of columns is the same in and .
Here are some examples of concatenation.
R, t = PolynomialRing(QQ, "t")
S = MatrixSpace(R, 3, 3)
A = S([t + 1 t R(1); t^2 t t; R(-2) t + 2 t^2 + t + 1])
B = S([R(2) R(3) R(1); t t + 1 t + 2; R(-1) t^2 t^3])
hcat(A, B)
vcat(A, B)
Permutation
#Base.:*
— Method.
*(P::perm, x::MatElem)
Apply the pemutation to the rows of the matrix and return the result.
Here is an example of applying a permutation to a matrix.
R, t = PolynomialRing(QQ, "t")
S = MatrixSpace(R, 3, 3)
G = PermGroup(3)
A = S([t + 1 t R(1); t^2 t t; R(-2) t + 2 t^2 + t + 1])
P = G([1, 3, 2])
B = P*A
LU factorisation
#Base.LinAlg.lufact
— Method.
lufact{T <: FieldElem}(A::MatElem{T}, P = PermGroup(rows(A)))
Return a tuple consisting of the rank of , a permutation of belonging to , a lower triangular matrix and an upper triangular matrix such that , where stands for the matrix whose rows are the given permutation of the rows of .
#Nemo.fflu
— Method.
fflu{T <: RingElem}(A::MatElem{T}, P = PermGroup(rows(A)))
Return a tuple consisting of the rank of , a denominator , a permutation of belonging to , a lower triangular matrix and an upper triangular matrix such that , where stands for the matrix whose rows are the given permutation of the rows of and such that is the diagonal matrix diag where the are the inverses of the diagonal entries of . The denominator is set to where is an appropriate submatrix of ( if is square) and the sign is decided by the parity of the permutation.
Here are some examples of LU factorisation.
R, x = PolynomialRing(QQ, "x")
K, a = NumberField(x^3 + 3x + 1, "a")
S = MatrixSpace(K, 3, 3)
A = S([K(0) 2a + 3 a^2 + 1; a^2 - 2 a - 1 2a; a^2 - 2 a - 1 2a])
r, P, L, U = lufact(A)
r, d, P, L, U = fflu(A)
Reduced row-echelon form
#Nemo.rref
— Method.
rref{T <: RingElem}(M::MatElem{T})
Returns a tuple consisting of the rank of and a denominator in the base ring of and a matrix such that is the reduced row echelon form of . Note that the denominator is not usually minimal.
#Nemo.rref
— Method.
rref{T <: FieldElem}(M::MatElem{T})
Returns a tuple consisting of the rank of and a reduced row echelon form of .
#Nemo.isrref
— Method.
isrref{T <: RingElem}(M::MatElem{T})
Return
true
if is in reduced row echelon form, otherwise returnfalse
.
#Nemo.isrref
— Method.
isrref{T <: FieldElem}(M::MatElem{T})
Return
true
if is in reduced row echelon form, otherwise returnfalse
.
Here are some examples of computing reduced row echelon form.
R, x = PolynomialRing(QQ, "x")
K, a = NumberField(x^3 + 3x + 1, "a")
S = MatrixSpace(K, 3, 3)
M = S([K(0) 2a + 3 a^2 + 1; a^2 - 2 a - 1 2a; a^2 + 3a + 1 2a K(1)])
r, d, A = rref(M)
isrref(A)
R, x = PolynomialRing(ZZ, "x")
S = MatrixSpace(R, 3, 3)
U = MatrixSpace(R, 3, 2)
M = S([R(0) 2x + 3 x^2 + 1; x^2 - 2 x - 1 2x; x^2 + 3x + 1 2x R(1)])
r, A = rref(M)
isrref(A)
Hermite normal form
#Nemo.hnf
— Method.
hnf{T <: RingElem}(A::GenMat{T}) -> GenMat{T}
Return the upper right row Hermite normal form of .
#Nemo.hnf_with_trafo
— Method.
hnf{T <: RingElem}(A::GenMat{T}) -> GenMaT{T}, GenMat{T}
Return the upper right row Hermite normal form of together with invertible matrix such that .
Determinant
#Base.LinAlg.det
— Method.
det{T <: RingElem}(M::MatElem{T})
Return the determinant of the matrix . We assume is square.
#Base.LinAlg.det
— Method.
det{T <: FieldElem}(M::MatElem{T})
Return the determinant of the matrix . We assume is square.
#Nemo.det_divisor
— Method.
det_divisor(x::fmpz_mat)
Return some positive divisor of the determinant of , if the determinant is nonzero, otherwise return zero.
#Nemo.det_given_divisor
— Method.
det_given_divisor(x::fmpz_mat, d::Integer, proved=true)
Return the determinant of given a positive divisor of its determinant. If
proved == true
(the default), the output is guaranteed to be correct, otherwise a heuristic algorithm is used.
#Nemo.det_given_divisor
— Method.
det_given_divisor(x::fmpz_mat, d::fmpz, proved=true)
Return the determinant of given a positive divisor of its determinant. If
proved == true
(the default), the output is guaranteed to be correct, otherwise a heuristic algorithm is used.
Here are some examples of computing the determinant.
R, x = PolynomialRing(QQ, "x")
K, a = NumberField(x^3 + 3x + 1, "a")
S = MatrixSpace(K, 3, 3)
A = S([K(0) 2a + 3 a^2 + 1; a^2 - 2 a - 1 2a; a^2 + 3a + 1 2a K(1)])
d = det(A)
S = MatrixSpace(ZZ, 3, 3)
A = S([fmpz(2) 3 5; 1 4 7; 9 6 3])
c = det_divisor(A)
d = det_given_divisor(A, c)
Rank
#Base.LinAlg.rank
— Method.
rank{T <: RingElem}(M::MatElem{T})
Return the rank of the matrix .
#Base.LinAlg.rank
— Method.
rank{T <: FieldElem}(M::MatElem{T})
Return the rank of the matrix .
Here are some examples of computing the rank of a matrix.
R, x = PolynomialRing(QQ, "x")
K, a = NumberField(x^3 + 3x + 1, "a")
S = MatrixSpace(K, 3, 3)
A = S([K(0) 2a + 3 a^2 + 1; a^2 - 2 a - 1 2a; a^2 + 3a + 1 2a K(1)])
d = rank(A)
Linear solving
solve{T <: RingElem}(::MatElem{T}, ::MatElem{T})
#Nemo.cansolve
— Method.
cansolve(a::fmpz_mat, b::fmpz_mat) -> Bool, fmpz_mat
Return true and a matrix such that , or false and some matrix in case does not exist.
#Nemo.solve_rational
— Method.
solve_rational{T <: RingElem}(M::MatElem{T}, b::MatElem{T})
Given a non-singular matrix over a ring and an matrix over the same ring, return a tuple consisting of an matrix and a denominator such that . The denominator will be the determinant of up to sign. If is singular an exception is raised.
#Nemo.solve_triu
— Method.
solve_triu{T <: FieldElem}(U::MatElem{T}, b::MatElem{T}, unit=false)
Given a non-singular matrix over a field which is upper triangular, and an matrix over the same field, return an matrix such that . If is singular an exception is raised. If unit is true then is assumed to have ones on its diagonal, and the diagonal will not be read.
#Nemo.solve_dixon
— Method.
solve_dixon(a::fmpz_mat, b::fmpz_mat)
Return a tuple consisting of a column vector such that . The element must be a column vector with the same number > of rows as and must be a square matrix. If these conditions are not met or does not exist, an exception is raised.
#Nemo.solve_dixon
— Method.
solve_dixon(a::fmpq_mat, b::fmpq_mat)
Solve by clearing denominators and using Dixon's algorithm. This is usually faster for large systems.
Here are some examples of linear solving.
R, x = PolynomialRing(QQ, "x")
K, a = NumberField(x^3 + 3x + 1, "a")
S = MatrixSpace(K, 3, 3)
U = MatrixSpace(K, 3, 1)
A = S([K(0) 2a + 3 a^2 + 1; a^2 - 2 a - 1 2a; a^2 + 3a + 1 2a K(1)])
b = U([2a a + 1 (-a - 1)]')
x = solve(A, b)
A = S([a + 1 2a + 3 a^2 + 1; K(0) a^2 - 1 2a; K(0) K(0) a])
b = U([2a a + 1 (-a - 1)]')
x = solve_triu(A, b, false)
R, x = PolynomialRing(ZZ, "x")
S = MatrixSpace(R, 3, 3)
U = MatrixSpace(R, 3, 2)
A = S([R(0) 2x + 3 x^2 + 1; x^2 - 2 x - 1 2x; x^2 + 3x + 1 2x R(1)])
b = U([2x x + 1 (-x - 1); x + 1 (-x) x^2]')
x, d = solve_rational(A, b)
S = MatrixSpace(ZZ, 3, 3)
T = MatrixSpace(ZZ, 3, 1)
A = S([fmpz(2) 3 5; 1 4 7; 9 2 2])
B = T([fmpz(4), 5, 7])
X, d = solve_rational(A, B)
X, m = solve_dixon(A, B)
Inverse
#Base.inv
— Method.
inv{T <: RingElem}(M::MatElem{T})
Given a non-singular matrix over a ring the tuple consisting of an matrix and a denominator such that , where is the identity matrix. The denominator will be the determinant of up to sign. If is singular an exception is raised.
#Base.inv
— Method.
inv{T <: FieldElem}(M::MatElem{T})
Given a non-singular matrix over a field, return an matrix such that where is the identity matrix. If is singular an exception is raised.
#Base.inv
— Method.
inv(M::arb_mat)
Given a matrix of type
arb_mat
, return an matrix such that contains the identity matrix. If cannot be inverted numerically an exception is raised.
#Base.inv
— Method.
inv(M::acb_mat)
Given a matrix of type
acb_mat
, return an matrix such that contains the identity matrix. If cannot be inverted numerically an exception is raised.
#Nemo.pseudo_inv
— Method.
pseudo_inv(x::fmpz_mat)
Return a tuple consisting of a matrix and denominator such that is the inverse of .
Here are some examples of taking the inverse of a matrix.
R, x = PolynomialRing(QQ, "x")
K, a = NumberField(x^3 + 3x + 1, "a")
S = MatrixSpace(K, 3, 3)
A = S([K(0) 2a + 3 a^2 + 1; a^2 - 2 a - 1 2a; a^2 + 3a + 1 2a K(1)])
X = inv(A)
R, x = PolynomialRing(ZZ, "x")
S = MatrixSpace(R, 3, 3)
A = S([R(0) 2x + 3 x^2 + 1; x^2 - 2 x - 1 2x; x^2 + 3x + 1 2x R(1)])
X, d = inv(A)
S = MatrixSpace(ZZ, 3, 3)
A = S([1 0 1; 2 3 1; 5 6 7])
B, d = pseudo_inv(A)
A = RR[1 0 1; 2 3 1; 5 6 7]
X = inv(A)
Nullspace
#Base.LinAlg.nullspace
— Method.
nullspace{T <: RingElem}(M::MatElem{T})
Returns a tuple consisting of the nullity of and a basis (consisting of column vectors) for the right nullspace of , i.e. such that is the zero matrix. If is an matrix will be an matrix. Note that the nullspace is taken to be the vector space kernel over the fraction field of the base ring if the latter is not a field. In Nemo we use the name ``kernel'' for a function to compute an integral kernel.
#Base.LinAlg.nullspace
— Method.
nullspace{T <: FieldElem}(M::MatElem{T})
Returns a tuple consisting of the nullity of and a basis (consisting of column vectors) for the right nullspace of , i.e. such that is the zero matrix. If is an matrix will be an matrix. Note that the nullspace is taken to be the vector space kernel over the fraction field of the base ring if the latter is not a field. In Nemo we use the name ``kernel'' for a function to compute an integral kernel.
#Nemo.nullspace_right_rational
— Method.
nullspace_right_rational(x::fmpz_mat)
Return the right rational nullspace of , i.e. a set of vectors over giving a -basis for the nullspace of considered as a matrix over .
Here are some examples of computing the nullspace of a matrix.
R, x = PolynomialRing(ZZ, "x")
S = MatrixSpace(R, 4, 4)
M = S([-6*x^2+6*x+12 -12*x^2-21*x-15 -15*x^2+21*x+33 -21*x^2-9*x-9;
-8*x^2+8*x+16 -16*x^2+38*x-20 90*x^2-82*x-44 60*x^2+54*x-34;
-4*x^2+4*x+8 -8*x^2+13*x-10 35*x^2-31*x-14 22*x^2+21*x-15;
-10*x^2+10*x+20 -20*x^2+70*x-25 150*x^2-140*x-85 105*x^2+90*x-50])
n, N = nullspace(M)
Hessenberg form
#Nemo.hessenberg
— Method.
hessenberg{T <: RingElem}(A::MatElem{T})
Returns the Hessenberg form of , i.e. an upper Hessenberg matrix which is similar to . The upper Hessenberg form has nonzero entries above and on the diagonal and in the diagonal line immediately below the diagonal.
#Nemo.ishessenberg
— Method.
ishessenberg{T <: RingElem}(A::MatElem{T})
Returns
true
if is in Hessenberg form, otherwise returnsfalse
.
Here are some examples of computing the Hessenberg form.
R = ResidueRing(ZZ, 7)
S = MatrixSpace(R, 4, 4)
M = S([R(1) R(2) R(4) R(3); R(2) R(5) R(1) R(0);
R(6) R(1) R(3) R(2); R(1) R(1) R(3) R(5)])
A = hessenberg(M)
ishessenberg(A) == true
Characteristic polynomial
#Nemo.charpoly
— Method.
charpoly{T <: RingElem}(V::Ring, Y::MatElem{T})
Returns the characteristic polynomial of the matrix . The polynomial ring of the resulting polynomial must be supplied and the matrix is assumed to be square.
Here are some examples of computing the characteristic polynomial.
R = ResidueRing(ZZ, 7)
S = MatrixSpace(R, 4, 4)
T, x = PolynomialRing(R, "x")
M = S([R(1) R(2) R(4) R(3); R(2) R(5) R(1) R(0);
R(6) R(1) R(3) R(2); R(1) R(1) R(3) R(5)])
A = charpoly(T, M)
Minimal polynomial
#Nemo.minpoly
— Method.
minpoly{T <: RingElem}(S::Ring, M::MatElem{T}, charpoly_only = false)
Returns the minimal polynomial of the matrix . The polynomial ring of the resulting polynomial must be supplied and the matrix must be square.
#Nemo.minpoly
— Method.
minpoly{T <: FieldElem}(S::Ring, M::MatElem{T}, charpoly_only = false)
Returns the minimal polynomial of the matrix . The polynomial ring of the resulting polynomial must be supplied and the matrix must be square.
minpoly{T <: RingElem}(S::Ring, M::MatElem{T}, charpoly_only = false)
Returns the minimal polynomial of the matrix . The polynomial ring of the resulting polynomial must be supplied and the matrix must be square.
Here are some examples of computing the minimal polynomial of a matrix.
R, x = FiniteField(13, 1, "x")
T, y = PolynomialRing(R, "y")
M = R[7 6 1;
7 7 5;
8 12 5]
A = minpoly(T, M)
Transforms
#Nemo.similarity!
— Method.
similarity!{T <: RingElem}(A::MatElem{T}, r::Int, d::T)
Applies a similarity transform to the matrix in-place. Let be the identity matrix that has had all zero entries of row replaced with , then the transform applied is equivalent to . We require to be a square matrix. A similarity transform preserves the minimal and characteristic polynomials of a matrix.
Here is an example of applying a similarity transform to a matrix.
R = ResidueRing(ZZ, 7)
S = MatrixSpace(R, 4, 4)
M = S([R(1) R(2) R(4) R(3); R(2) R(5) R(1) R(0);
R(6) R(1) R(3) R(2); R(1) R(1) R(3) R(5)])
similarity!(M, 1, R(3))
Modular reduction
#Nemo.reduce_mod
— Method.
reduce_mod(x::fmpz_mat, y::Integer)
Reduce the entries of modulo and return the result.
#Nemo.reduce_mod
— Method.
reduce_mod(x::fmpz_mat, y::fmpz)
Reduce the entries of modulo and return the result.
Here are some examples of modular reduction.
S = MatrixSpace(ZZ, 3, 3)
A = S([fmpz(2) 3 5; 1 4 7; 9 2 2])
reduce_mod(A, ZZ(5))
reduce_mod(A, 2)
Lifting
#Nemo.lift
— Method.
lift(a::nmod_mat)
Return a lift of the matrix to a matrix over , i.e. where the entries of the returned matrix are those of lifted to .
Here are some examples of lifting.
R = ResidueRing(ZZ, 7)
S = MatrixSpace(R, 3, 3)
a = S([4 5 6; 7 3 2; 1 4 5])
b = lift(a)
Special matrices
#Nemo.hadamard
— Method.
hadamard(R::FmpzMatSpace)
Return the Hadamard matrix for the given matrix space. The number of rows and columns must be equal.
#Nemo.ishadamard
— Method.
ishadamard(x::fmpz_mat)
Return
true
if the given matrix is Hadamard, otherwise returnfalse
.
#Nemo.hilbert
— Method.
hilbert(R::FmpqMatSpace)
Return the Hilbert matrix in the given matrix space. This is the matrix with entries .
Here are some examples of computing special matrices.
R = MatrixSpace(ZZ, 3, 3)
S = MatrixSpace(QQ, 3, 3)
A = hadamard(R)
ishadamard(A)
B = hilbert(R)
Hermite Normal Form
#Nemo.hnf
— Method.
hnf(x::fmpz_mat)
Return the Hermite Normal Form of .
#Nemo.hnf_with_transform
— Method.
hnf_with_transform(x::fmpz_mat)
Compute a tuple where is the Hermite normal form of and is a transformation matrix so that .
#Nemo.hnf_modular
— Method.
hnf_modular(x::fmpz_mat, d::fmpz)
Compute the Hermite normal form of given that is a multiple of the determinant of the nonzero rows of .
#Nemo.hnf_modular_eldiv
— Method.
hnf_modular_eldiv(x::fmpz_mat, d::fmpz)
Compute the Hermite normal form of given that is a multiple of the largest elementary divisor of . The matrix must have full rank.
#Nemo.ishnf
— Method.
ishnf(x::fmpz_mat)
Return
true
if the given matrix is in Hermite Normal Form, otherwise returnfalse
.
Here are some examples of computing the Hermite Normal Form.
S = MatrixSpace(ZZ, 3, 3)
A = S([fmpz(2) 3 5; 1 4 7; 19 3 7])
B = hnf(A)
H, T = hnf_with_transform(A)
M = hnf_modular(A, fmpz(27))
N = hnf_modular_eldiv(A, fmpz(27))
ishnf(M)
(Weak) Popov form
Nemo provides algorithms for computing the (weak) Popov of a matrix with entries in a univariate polynomial ring over a field.
weak_popov{T <: PolyElem}(::GenMat{T})
weak_popov_with_trafo{T <: PolyElem}(::GenMat{T})
popov{T <: PolyElem}(::GenMat{T})
popov_with_trafo{T <: PolyElem}(::GenMat{T})
det_popov{T <: PolyElem}(::GenMat{T})
Lattice basis reduction
Nemo provides LLL lattice basis reduction. Optionally one can specify the setup using a context object created by the following function.
lll_ctx(delta::Float64, eta::Float64, rep=:zbasis, gram=:approx)
Return a LLL context object specifying LLL parameters and and specifying the representation as either :zbasis
or :gram
and the Gram type as either :approx
or :exact
.
#Nemo.lll
— Method.
lll(x::fmpz_mat, ctx=lll_ctx(0.99, 0.51))
Return the LLL reduction of the matrix . By default the matrix is a -basis and the Gram matrix is maintained throughout in approximate form. The LLL is performed with reduction parameters and . All of these defaults can be overridden by specifying an optional context object.
#Nemo.lll_with_transform
— Method.
Compute a tuple where is the LLL reduction of and is a transformation matrix so that . All the default parameters can be overridden by supplying an optional context object.
#Nemo.lll_gram
— Method.
lll_gram(x::fmpz_mat, ctx=lll_ctx(0.99, 0.51, :gram))
Given the Gram matrix of a matrix, compute the Gram matrix of its LLL reduction.
#Nemo.lll_gram_with_transform
— Method.
lll_gram_with_transform(x::fmpz_mat, ctx=lll_ctx(0.99, 0.51, :gram))
Given the Gram matrix of a matrix , compute a tuple where is the gram matrix of the LLL reduction of the matrix and is a transformation matrix so that .
#Nemo.lll_with_removal
— Method.
lll_with_removal(x::fmpz_mat, b::fmpz, ctx=lll_ctx(0.99, 0.51))
Compute the LLL reduction of and throw away rows whose norm exceeds the given bound . Return a tuple where the first rows of are the rows remaining after removal.
#Nemo.lll_with_removal_transform
— Method.
lll_with_removal_transform(x::fmpz_mat, b::fmpz, ctx=lll_ctx(0.99, 0.51))
Compute a tuple where the first rows of are those remaining from the LLL reduction after removal of vectors with norm exceeding the bound and is a transformation matrix so that .
Here are some examples of lattice basis reduction.
S = MatrixSpace(ZZ, 3, 3)
A = S([fmpz(2) 3 5; 1 4 7; 19 3 7])
L = lll(A, lll_ctx(0.95, 0.55, :zbasis, :approx)
L, T = lll_with_transform(A)
G == lll_gram(gram(A))
G, T = lll_gram_with_transform(gram(A))
r, L = lll_with_removal(A, fmpz(100))
r, L, T = lll_with_removal_transform(A, fmpz(100))
Smith Normal Form
#Nemo.snf
— Method.
snf(x::fmpz_mat)
Compute the Smith normal form of .
#Nemo.snf_diagonal
— Method.
snf_diagonal(x::fmpz_mat)
Given a diagonal matrix compute the Smith normal form of .
#Nemo.issnf
— Method.
issnf(x::fmpz_mat)
Return
true
if is in Smith normal form, otherwise returnfalse
.
Here are some examples of computing the Smith Normal Form.
S = MatrixSpace(ZZ, 3, 3)
A = S([fmpz(2) 3 5; 1 4 7; 19 3 7])
B = snf(A)
issnf(B) == true
B = S([fmpz(2) 0 0; 0 4 0; 0 0 7])
C = snf_diagonal(B)
Strong Echelon Form
#Nemo.strong_echelon_form
— Method.
strong_echelon_form(a::nmod_mat)
Return the strong echeleon form of . The matrix must have at least as many rows as columns.
Here is an example of computing the strong echelon form.
R = ResidueRing(ZZ, 12)
S = MatrixSpace(R, 3, 3)
A = S([4 1 0; 0 0 5; 0 0 0 ])
B = strong_echelon_form(A)
Howell Form
#Nemo.howell_form
— Method.
howell_form(a::nmod_mat)
Return the Howell normal form of . The matrix must have at least as many rows as columns.
Here is an example of computing the Howell form.
R = ResidueRing(ZZ, 12)
S = MatrixSpace(R, 3, 3)
A = S([4 1 0; 0 0 5; 0 0 0 ])
B = howell_form(A)
Gram-Schmidt Orthogonalisation
#Nemo.gso
— Method.
gso(x::fmpq_mat)
Return the Gram-Schmidt Orthogonalisation of the matrix .
Here are some examples of computing the Gram-Schmidt Orthogonalisation.
S = MatrixSpace(QQ, 3, 3)
A = S([4 7 3; 2 9 1; 0 5 3])
B = gso(A)
Exponential
#Base.exp
— Method.
exp(x::arb_mat)
Returns the exponential of the matrix .
#Base.exp
— Method.
exp(x::acb_mat)
Returns the exponential of the matrix .
Here are some examples of computing the exponential function of matrix.
A = RR[2 0 0; 0 3 0; 0 0 1]
B = exp(A)
Norm
#Nemo.bound_inf_norm
— Method.
bound_inf_norm(x::arb_mat)
Returns a nonnegative element of type
arb
, such that is an upper bound for the infinity norm for every matrix in
#Nemo.bound_inf_norm
— Method.
bound_inf_norm(x::acb_mat)
Returns a nonnegative element of type
acb
, such that is an upper bound for the infinity norm for every matrix in
Here are some examples of computing bounds on the infinity norm of a matrix.
A = RR[1 2 3; 4 5 6; 7 8 9]
d = bound_inf_norm(A)
Shifting
#Base.Math.ldexp
— Method.
ldexp(x::acb_mat, y::Int)
Return . Note that can be positive, zero or negative.
#Base.Math.ldexp
— Method.
ldexp(x::acb_mat, y::Int)
Return . Note that can be positive, zero or negative.
Here are some examples of shifting.
A = RR[1 2 3; 4 5 6; 7 8 9]
B = ldexp(A, 4)
overlaps(16*A, B)
Predicates
#Base.isreal
— Method.
isreal(M::acb_mat)
Returns whether every entry of has vanishing imaginary part.
Here are some examples for predicates.
A = CC[1 2 3; 4 5 6; 7 8 9]
isreal(A)
isreal(onei(CC)*A)
Conversion to Julia matrices
Julia matrices use a different data structure than Nemo matrices. Conversion to Julia matrices is usually only required for interfacing with other packages. It isn't necessary to convert Nemo matrices to Julia matrices in order to manipulate them.
This conversion can be performed with standard Julia syntax, such as the following, where A
is a fmpz_mat
:
Matrix{Int}(A)
Matrix{BigInt}(A)
In case the matrix cannot be converted without loss, an InexactError
is thrown: in this case, cast to a matrix of BigInt
s rather than Int
s.