Matrices

Matrices

Nemo allow the creation of dense matrices over any computable ring $R$. There are two different kinds of implementation: a generic one for the case where no specific implementation exists (provided by AbstractAlgebra.jl), 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 $R$, and the Julia/Nemo types for that kind of matrix (the type information is mainly of concern to developers).

Base ringLibraryElement typeParent type
Generic ring $R$AbstractAlgebra.jlGeneric.Mat{T}Generic.MatSpace{T}
$\mathbb{Z}$Flintfmpz_matFmpzMatSpace
$\mathbb{Z}/n\mathbb{Z}$ (small $n$)Flintnmod_matNmodMatSpace
$\mathbb{Q}$Flintfmpq_matFmpqMatSpace
$\mathbb{Z}/p\mathbb{Z}$ (small $p$)Flintgfp_matGFPMatSpace
$\mathbb{F}_{p^n}$ (small $p$)Flintfq_nmod_matFqNmodMatSpace
$\mathbb{F}_{p^n}$ (large $p$)Flintfq_mat`FqMatSpace
$\mathbb{R}$Arbarb_matArbMatSpace
$\mathbb{C}$Arbacb_matAcbMatSpace

The dimensions and base ring $R$ 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 functionality

All matrix spaces in Nemo follow the AbstractAlgebra.jl matrix interface:

https://nemocas.github.io/AbstractAlgebra.jl/matrix_spaces.html

In addition, AbstractAlgebra.jl provides a great deal of generic functionality for matrices. Some of this functionality is also provided by C libraries, such as Flint, for various specific rings.

https://nemocas.github.io/AbstractAlgebra.jl/matrix.html

In the following, we list the functionality which is provided in addition to the generic matrix functionality, for specific rings in Nemo.

Comparison operators

Nemo.overlapsMethod.
overlaps(x::arb_mat, y::arb_mat)

Returns true if all entries of $x$ overlap with the corresponding entry of $y$, otherwise return false.

Nemo.overlapsMethod.
overlaps(x::acb_mat, y::acb_mat)

Returns true if all entries of $x$ overlap with the corresponding entry of $y$, otherwise return false.

Nemo.containsMethod.
contains(x::arb_mat, y::arb_mat)

Returns true if all entries of $x$ contain the corresponding entry of $y$, otherwise return false.

Nemo.containsMethod.
contains(x::acb_mat, y::acb_mat)

Returns true if all entries of $x$ contain the corresponding entry of $y$, otherwise return false.

In addition we have the following ad hoc comparison operators.

Examples

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 $2^yx$.

Base.:>>Method.
>>(x::fmpz_mat, y::Int)

Return $x/2^y$ where rounding is towards zero.

Examples

S = MatrixSpace(ZZ, 3, 3)

A = S([fmpz(2) 3 5; 1 4 7; 9 6 3])
 
B = A<<5
C = B>>2

##i# Determinant

Nemo.det_divisorMethod.
det_divisor(x::fmpz_mat)

Return some positive divisor of the determinant of $x$, if the determinant is nonzero, otherwise return zero.

det_given_divisor(x::fmpz_mat, d::Integer, proved=true)

Return the determinant of $x$ 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.

det_given_divisor(x::fmpz_mat, d::fmpz, proved=true)

Return the determinant of $x$ 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.

Examples

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)

Linear solving

Nemo.cansolveMethod.
cansolve(a::fmpz_mat, b::fmpz_mat) -> Bool, fmpz_mat

Return true and a matrix $x$ such that $ax = b$, or false and some matrix in case $x$ does not exist.

Nemo.solve_dixonMethod.
solve_dixon(a::fmpz_mat, b::fmpz_mat)

Return a tuple $(x, m)$ consisting of a column vector $x$ such that $ax = b \pmod{m}$. The element $b$ must be a column vector with the same number > of rows as $a$ and $a$ must be a square matrix. If these conditions are not met or $(x, d)$ does not exist, an exception is raised.

Nemo.solve_dixonMethod.
solve_dixon(a::fmpq_mat, b::fmpq_mat)

Solve $ax = b$ by clearing denominators and using Dixon's algorithm. This is usually faster for large systems.

Examples

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, m = solve_dixon(A, B)

Pseudo inverse

pseudo_inv(x::fmpz_mat)

Return a tuple $(z, d)$ consisting of a matrix $z$ and denominator $d$ such that $z/d$ is the inverse of $x$.

Examples

S = MatrixSpace(ZZ, 3, 3)

A = S([1 0 1; 2 3 1; 5 6 7])
  
B, d = pseudo_inv(A)

Nullspace

nullspace_right_rational(x::fmpz_mat)

Return a tuple $(r, U)$ consisting of a matrix $U$ such that the first $r$ columns form the right rational nullspace of $x$, i.e. a set of vectors over $\mathbb{Z}$ giving a $\mathbb{Q}$-basis for the nullspace of $x$ considered as a matrix over

\[\mathbb{Q}\]

.

Modular reduction

Nemo.reduce_modMethod.
reduce_mod(x::fmpz_mat, y::Integer)

Reduce the entries of $x$ modulo $y$ and return the result.

Nemo.reduce_modMethod.
reduce_mod(x::fmpz_mat, y::fmpz)

Reduce the entries of $x$ modulo $y$ and return the result.

Examples

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.liftMethod.
lift(a::T) where {T <: Zmodn_mat}

Return a lift of the matrix $a$ to a matrix over $\mathbb{Z}$, i.e. where the entries of the returned matrix are those of $a$ lifted to $\mathbb{Z}$.

Nemo.liftMethod.
lift(a::gfp_mat)

Return a lift of the matrix $a$ to a matrix over $\mathbb{Z}$, i.e. where the entries of the returned matrix are those of $a$ lifted to $\mathbb{Z}$.

Examples

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.hadamardMethod.
hadamard(R::FmpzMatSpace)

Return the Hadamard matrix for the given matrix space. The number of rows and columns must be equal.

Nemo.ishadamardMethod.
ishadamard(x::fmpz_mat)

Return true if the given matrix is Hadamard, otherwise return false.

Nemo.hilbertMethod.
hilbert(R::FmpqMatSpace)

Return the Hilbert matrix in the given matrix space. This is the matrix with entries $H_{i,j} = 1/(i + j - 1)$.

Examples

R = MatrixSpace(ZZ, 3, 3)
S = MatrixSpace(QQ, 3, 3)

A = hadamard(R)
ishadamard(A)
B = hilbert(R)

Hermite Normal Form

hnf(x::fmpz_mat)

Return the Hermite Normal Form of $x$.

hnf_with_transform(x::fmpz_mat)

Compute a tuple $(H, T)$ where $H$ is the Hermite normal form of $x$ and $T$ is a transformation matrix so that $H = Tx$.

Nemo.hnf_modularMethod.
hnf_modular(x::fmpz_mat, d::fmpz)

Compute the Hermite normal form of $x$ given that $d$ is a multiple of the determinant of the nonzero rows of $x$.

hnf_modular_eldiv(x::fmpz_mat, d::fmpz)

Compute the Hermite normal form of $x$ given that $d$ is a multiple of the largest elementary divisor of $x$. The matrix $x$ must have full rank.

Nemo.ishnfMethod.
ishnf(x::fmpz_mat)

Return true if the given matrix is in Hermite Normal Form, otherwise return false.

Examples

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)

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 $\delta$ and $\eta$ and specifying the representation as either :zbasis or :gram and the Gram type as either :approx or :exact.

Nemo.lllMethod.
lll(x::fmpz_mat, ctx::lll_ctx = lll_ctx(0.99, 0.51))

Return the LLL reduction of the matrix $x$. By default the matrix $x$ is a $\mathbb{Z}$-basis and the Gram matrix is maintained throughout in approximate form. The LLL is performed with reduction parameters $\delta = 0.99$ and $\eta = 0.51$. All of these defaults can be overridden by specifying an optional context object.

lll_with_transform(x::fmpz_mat, ctx::lll_ctx = lll_ctx(0.99, 0.51))

Compute a tuple $(L, T)$ where $L$ is the LLL reduction of $a$ and $T$ is a transformation matrix so that $L = Ta$. All the default parameters can be overridden by supplying an optional context object.

Nemo.lll_gramMethod.
lll_gram(x::fmpz_mat, ctx::lll_ctx = lll_ctx(0.99, 0.51, :gram))

Given the Gram matrix $x$ of a matrix, compute the Gram matrix of its LLL reduction.

lll_gram_with_transform(x::fmpz_mat, ctx::lll_ctx = lll_ctx(0.99, 0.51, :gram))

Given the Gram matrix $x$ of a matrix $M$, compute a tuple $(L, T)$ where $L$ is the gram matrix of the LLL reduction of the matrix and $T$ is a transformation matrix so that $L = TM$.

lll_with_removal(x::fmpz_mat, b::fmpz, ctx::lll_ctx = lll_ctx(0.99, 0.51))

Compute the LLL reduction of $x$ and throw away rows whose norm exceeds the given bound $b$. Return a tuple $(r, L)$ where the first $r$ rows of $L$ are the rows remaining after removal.

lll_with_removal_transform(x::fmpz_mat, b::fmpz, ctx::lll_ctx = lll_ctx(0.99, 0.51))

Compute a tuple $(r, L, T)$ where the first $r$ rows of $L$ are those remaining from the LLL reduction after removal of vectors with norm exceeding the bound $b$ and $T$ is a transformation matrix so that $L = Tx$.

Nemo.lll!Method.
lll!(x::fmpz_mat, ctx::lll_ctx = lll_ctx(0.99, 0.51))

Perform the LLL reduction of the matrix $x$ inplace. By default the matrix $x$ is a > $\mathbb{Z}$-basis and the Gram matrix is maintained throughout in approximate form. The LLL is performed with reduction parameters $\delta = 0.99$ and $\eta = 0.51$. All of these defaults can be overridden by specifying an optional context object.

Nemo.lll_gram!Method.
lll_gram!(x::fmpz_mat, ctx::lll_ctx = lll_ctx(0.99, 0.51, :gram))

Given the Gram matrix $x$ of a matrix, compute the Gram matrix of its LLL reduction inplace.

Examples

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

snf(x::fmpz_mat)

Compute the Smith normal form of $x$.

Nemo.snf_diagonalMethod.
snf_diagonal(x::fmpz_mat)

Given a diagonal matrix $x$ compute the Smith normal form of $x$.

Nemo.issnfMethod.
issnf(x::fmpz_mat)

Return true if $x$ is in Smith normal form, otherwise return false.

Examples

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

strong_echelon_form(a::nmod_mat)

Return the strong echeleon form of $a$. The matrix $a$ must have at least as many rows as columns.

strong_echelon_form(a::gfp_mat)

Return the strong echeleon form of $a$. The matrix $a$ must have at least as many rows as columns.

Examples

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_formMethod.
howell_form(a::nmod_mat)

Return the Howell normal form of $a$. The matrix $a$ must have at least as many rows as columns.

Nemo.howell_formMethod.
howell_form(a::gfp_mat)

Return the Howell normal form of $a$. The matrix $a$ must have at least as many rows as columns.

Examples

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.gsoMethod.
gso(x::fmpq_mat)

Return the Gram-Schmidt Orthogonalisation of the matrix $x$.

Examples

S = MatrixSpace(QQ, 3, 3)

A = S([4 7 3; 2 9 1; 0 5 3])

B = gso(A)

Exponential

Base.expMethod.
exp(x::arb_mat)

Returns the exponential of the matrix $x$.

Base.expMethod.
exp(x::acb_mat)

Returns the exponential of the matrix $x$.

Examples

A = RR[2 0 0; 0 3 0; 0 0 1]

B = exp(A)

Norm

bound_inf_norm(x::arb_mat)

Returns a nonnegative element $z$ of type arb, such that $z$ is an upper bound for the infinity norm for every matrix in $x$

bound_inf_norm(x::acb_mat)

Returns a nonnegative element $z$ of type acb, such that $z$ is an upper bound for the infinity norm for every matrix in $x$

Examples

A = RR[1 2 3; 4 5 6; 7 8 9]

d = bound_inf_norm(A)

Shifting

Base.Math.ldexpMethod.
ldexp(x::arb_mat, y::Int)

Return $2^yx$. Note that $y$ can be positive, zero or negative.

Base.Math.ldexpMethod.
ldexp(x::acb_mat, y::Int)

Return $2^yx$. Note that $y$ can be positive, zero or negative.

Examples

A = RR[1 2 3; 4 5 6; 7 8 9]

B = ldexp(A, 4)

overlaps(16*A, B)

Predicates

Base.isrealMethod.
isreal(x::acb_mat)

Returns whether every entry of $x$ has vanishing imaginary part.

Examples

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 an 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 BigInts rather than Ints.