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}$FlintZZMatrixZZMatrixSpace
$\mathbb{Z}/n\mathbb{Z}$ (small $n$)FlintzzModMatrixzzModMatrixSpace
$\mathbb{Z}/n\mathbb{Z}$ (large $n$)FlintZZModMatrixZZModMatrixSpace
$\mathbb{Q}$FlintQQMatrixQQMatrixSpace
$\mathbb{Z}/p\mathbb{Z}$ (small $p$)FlintfpMatrixfpMatrixSpace
$\mathbb{F}_{p^n}$ (small $p$)FlintfqPolyRepMatrixfqPolyRepMatrixSpace
$\mathbb{F}_{p^n}$ (large $p$)FlintFqPolyRepMatrixFqPolyRepMatrixSpace
$\mathbb{R}$ (arbitrary precision)ArbRealMatrixRealMatrixSpace
$\mathbb{C}$ (arbitrary precision)ArbComplexMatrixComplexMatrixSpace
$\mathbb{R}$ (fixed precision)ArbArbMatrixArbMatrixSpace
$\mathbb{C}$ (fixed precision)ArbAcbMatrixAcbMatrixSpace

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.

Note that the preferred way to create matrices is not to use the type constructors but to use the matrix function, see also the Matrix element constructors section of the AbstractAlgebra manual.

Matrix functionality

All matrix spaces in Nemo provide the matrix functionality of AbstractAlgebra:

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

Some of this functionality is provided in Nemo by C libraries, such as Flint, for various specific rings.

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

Comparison operators

overlapsMethod
overlaps(x::RealMatrix, y::RealMatrix)

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

source
overlapsMethod
overlaps(x::ComplexMatrix, y::ComplexMatrix)

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

source
containsMethod
contains(x::RealMatrix, y::RealMatrix)

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

source
containsMethod
contains(x::ComplexMatrix, y::ComplexMatrix)

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

source

In addition we have the following ad hoc comparison operators.

Examples

julia> C = RR[1 2; 3 4]
[1.0000000000000000000   2.0000000000000000000]
[3.0000000000000000000   4.0000000000000000000]

julia> D = RR["1 +/- 0.1" "2 +/- 0.1"; "3 +/- 0.1" "4 +/- 0.1"]
[[1e+0 +/- 0.101]   [2e+0 +/- 0.101]]
[[3e+0 +/- 0.101]   [4e+0 +/- 0.101]]

julia> overlaps(C, D)
true

julia> contains(D, C)
true

Scaling

<<Method
<<(x::ZZMatrix, y::Int)

Return $2^yx$.

source
>>Method
>>(x::ZZMatrix, y::Int)

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

source

Examples

julia> A = ZZ[2 3 5; 1 4 7; 9 6 3]
[2   3   5]
[1   4   7]
[9   6   3]

julia> B = A<<5
[ 64    96   160]
[ 32   128   224]
[288   192    96]

julia> C = B>>2
[16   24   40]
[ 8   32   56]
[72   48   24]

Determinant

det_divisorMethod
det_divisor(x::ZZMatrix)

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

source
det_given_divisorMethod
det_given_divisor(x::ZZMatrix, 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.

source
det_given_divisorMethod
det_given_divisor(x::ZZMatrix, d::ZZRingElem, 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.

source

Examples

julia> A = ZZ[2 3 5; 1 4 7; 9 6 3]
[2   3   5]
[1   4   7]
[9   6   3]

julia> c = det_divisor(A)
3

julia> d = det_given_divisor(A, c)
-30

Pseudo inverse

pseudo_invMethod
pseudo_inv(x::ZZMatrix)

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

Examples

julia> A = ZZ[1 0 1; 2 3 1; 5 6 7]
[1   0   1]
[2   3   1]
[5   6   7]

julia> B, d = pseudo_inv(A)
([15 6 -3; -9 2 1; -3 -6 3], 12)
source

Nullspace

nullspace_right_rationalMethod
nullspace_right_rational(x::ZZMatrix)

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}$.

source

Modular reduction

reduce_modMethod
reduce_mod(x::ZZMatrix, y::Integer)

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

source
reduce_modMethod
reduce_mod(x::ZZMatrix, y::ZZRingElem)

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

source

Examples

julia> A = ZZ[2 3 5; 1 4 7; 9 2 2]
[2   3   5]
[1   4   7]
[9   2   2]

julia> reduce_mod(A, ZZ(5))
[2   3   0]
[1   4   2]
[4   2   2]

julia> reduce_mod(A, 2)
[0   1   1]
[1   0   1]
[1   0   0]

Lifting

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}$.

source
lift(M::smodule{T}, SM::smodule{T}) where T

Represents the generators of SM in terms of the generators of M. If SM is in M, rest is the null module, otherwise rest = reduce(SM, std(M)). Returns (result, rest) with for global orderings: Matrix(SM) - Matrix(rest) = Matrix(M)*Matrix(result) for non-global orderings: Matrix(SM)*U - Matrix(rest) = Matrix(M)*Matrix(result) where U is some diagonal matrix of units. To compute this U, see lift(M::smodule, SM::smodule, goodShape::Bool, isSB::Bool, divide::Bool).

source
lift(M::smodule{T}, SM::smodule{T}, goodShape::Bool, isSB::Bool, divide::Bool) where T

Represents the generators of SM in terms of the generators of M. Returns (result, rest, U) with Matrix(SM)*U - Matrix(rest) = Matrix(M)*Matrix(result) If SM is in M, then rest is the null module. Otherwise, rest = SM if !divide, and rest = normalform(SM, std(M)) if divide. U is a diagonal matrix of units, differing from the identity matrix only for local ring orderings.

There are three boolean options. goodShape: maximal non-zero index in generators of SM <= that of M, which should be come from a rank check rank(SM)==rank(M). isSB: generators of M form a Groebner basis. divide: allow SM not to be a submodule of M.

source
lift(M::sideal{T}, SM::sideal{T}) where T

Represents the generators of SM in terms of the generators of M. If SM is in M, rest is the null module, otherwise rest = reduce(SM, std(M)). Returns (result, rest) with for global orderings: Matrix(SM) - Matrix(rest) = Matrix(M)*Matrix(result) for non-global orderings: Matrix(SM)*U - Matrix(rest) = Matrix(M)*Matrix(result) where U is some diagonal matrix of units. To compute this U, see lift(M::sideal, SM::sideal, goodShape::Bool, isSB::Bool, divide::Bool).

source
lift(M::sideal{T}, SM::sideal{T}, goodShape::Bool, isSB::Bool, divide::Bool) where T

Represents the generators of SM in terms of the generators of M. Returns (result, rest, U) with Matrix(SM)*U - Matrix(rest) = Matrix(M)*Matrix(result) If SM is in M, then rest is the null module. Otherwise, rest = SM if !divide, and rest = normalform(SM, std(M)) if divide. U is a diagonal matrix of units, differing from the identity matrix only for local ring orderings.

There are three boolean options. goodShape: maximal non-zero index in generators of SM <= that of M, which should be come from a rank check rank(SM)==rank(M). isSB: generators of M form a Groebner basis divide: allow SM not to be a submodule of M, which is useful for division with remainder.

source
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}$.

source
lift(M::smodule{T}, SM::smodule{T}) where T

Represents the generators of SM in terms of the generators of M. If SM is in M, rest is the null module, otherwise rest = reduce(SM, std(M)). Returns (result, rest) with for global orderings: Matrix(SM) - Matrix(rest) = Matrix(M)*Matrix(result) for non-global orderings: Matrix(SM)*U - Matrix(rest) = Matrix(M)*Matrix(result) where U is some diagonal matrix of units. To compute this U, see lift(M::smodule, SM::smodule, goodShape::Bool, isSB::Bool, divide::Bool).

source
lift(M::smodule{T}, SM::smodule{T}, goodShape::Bool, isSB::Bool, divide::Bool) where T

Represents the generators of SM in terms of the generators of M. Returns (result, rest, U) with Matrix(SM)*U - Matrix(rest) = Matrix(M)*Matrix(result) If SM is in M, then rest is the null module. Otherwise, rest = SM if !divide, and rest = normalform(SM, std(M)) if divide. U is a diagonal matrix of units, differing from the identity matrix only for local ring orderings.

There are three boolean options. goodShape: maximal non-zero index in generators of SM <= that of M, which should be come from a rank check rank(SM)==rank(M). isSB: generators of M form a Groebner basis. divide: allow SM not to be a submodule of M.

source
lift(M::sideal{T}, SM::sideal{T}) where T

Represents the generators of SM in terms of the generators of M. If SM is in M, rest is the null module, otherwise rest = reduce(SM, std(M)). Returns (result, rest) with for global orderings: Matrix(SM) - Matrix(rest) = Matrix(M)*Matrix(result) for non-global orderings: Matrix(SM)*U - Matrix(rest) = Matrix(M)*Matrix(result) where U is some diagonal matrix of units. To compute this U, see lift(M::sideal, SM::sideal, goodShape::Bool, isSB::Bool, divide::Bool).

source
lift(M::sideal{T}, SM::sideal{T}, goodShape::Bool, isSB::Bool, divide::Bool) where T

Represents the generators of SM in terms of the generators of M. Returns (result, rest, U) with Matrix(SM)*U - Matrix(rest) = Matrix(M)*Matrix(result) If SM is in M, then rest is the null module. Otherwise, rest = SM if !divide, and rest = normalform(SM, std(M)) if divide. U is a diagonal matrix of units, differing from the identity matrix only for local ring orderings.

There are three boolean options. goodShape: maximal non-zero index in generators of SM <= that of M, which should be come from a rank check rank(SM)==rank(M). isSB: generators of M form a Groebner basis divide: allow SM not to be a submodule of M, which is useful for division with remainder.

source

Examples

julia> R, = residue_ring(ZZ, 7)
(Integers modulo 7, Map: ZZ -> ZZ/(7))

julia> a = R[4 5 6; 7 3 2; 1 4 5]
[4   5   6]
[0   3   2]
[1   4   5]

julia> b = lift(a)
[-3   -2   -1]
[ 0    3    2]
[ 1   -3   -2]

Special matrices

hadamardMethod
hadamard(R::ZZMatrixSpace)

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

source
is_hadamardMethod
is_hadamard(x::ZZMatrix)

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

source
hilbertMethod
hilbert(R::QQMatrixSpace)

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

source

Examples

julia> hadamard(matrix_space(ZZ, 3, 3))
ERROR: Unable to create Hadamard matrix
[...]

julia> A = hadamard(matrix_space(ZZ, 4, 4))
[1    1    1    1]
[1   -1    1   -1]
[1    1   -1   -1]
[1   -1   -1    1]

julia> is_hadamard(A)
true

julia> B = hilbert(matrix_space(QQ, 3, 3))
[   1   1//2   1//3]
[1//2   1//3   1//4]
[1//3   1//4   1//5]

Hermite Normal Form

hnfMethod
hnf(x::ZZMatrix)

Return the Hermite Normal Form of $x$.

source
hnf_with_transformMethod
hnf_with_transform(x::ZZMatrix)

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

source
hnf_modularMethod
hnf_modular(x::ZZMatrix, d::ZZRingElem)

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

source
hnf_modular_eldivMethod
hnf_modular_eldiv(x::ZZMatrix, d::ZZRingElem)

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.

source
is_hnfMethod
is_hnf(x::ZZMatrix)

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

source

Examples

julia> A = ZZ[2 3 5; 1 4 7; 19 3 7]
[ 2   3   5]
[ 1   4   7]
[19   3   7]

julia> B = hnf(A)
[1   0   16]
[0   1   18]
[0   0   27]

julia> H, T = hnf_with_transform(A)
([1 0 16; 0 1 18; 0 0 27], [-43 30 3; -44 31 3; -73 51 5])

julia> M = hnf_modular(A, ZZ(27))
[1   0   16]
[0   1   18]
[0   0   27]

julia> N = hnf_modular_eldiv(A, ZZ(27))
[1   0   16]
[0   1   18]
[0   0   27]

julia> is_hnf(M)
true

Lattice basis reduction

Nemo provides LLL lattice basis reduction. Optionally one can specify the setup using a context object created by the following function.

LLLContext(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.

lllMethod
lll(x::ZZMatrix, ctx::LLLContext = LLLContext(0.99, 0.51))

Return a matrix $L$ whose rows form an LLL-reduced basis of the $\mathbb{Z}$-lattice generated by the rows of $x$. $L$ may contain additional zero rows.

By default, the LLL is performed with reduction parameters $\delta = 0.99$ and $\eta = 0.51$. These defaults can be overridden by specifying an optional context object.

See lll_gram for a function taking the Gram matrix as input.

source
lll_with_transformMethod
lll_with_transform(x::ZZMatrix, ctx::LLLContext = LLLContext(0.99, 0.51))

Return a tuple $(L, T)$ where the rows of $L$ form an LLL-reduced basis of the $\mathbb{Z}$-lattice generated by the rows of $x$ and $T$ is a transformation matrix so that $L = Tx$. $L$ may contain additional zero rows. See lll for the used default parameters which can be overridden by supplying an optional context object.

See lll_gram_with_transform for a function taking the Gram matrix as input.

source
lll_gramMethod
lll_gram(x::ZZMatrix, ctx::LLLContext = LLLContext(0.99, 0.51, :gram))

Return the Gram matrix $L$ of an LLL-reduced basis of the lattice given by the Gram matrix $x$. The matrix $x$ must be symmetric and non-singular.

By default, the LLL is performed with reduction parameters $\delta = 0.99$ and $\eta = 0.51$. These defaults can be overridden by specifying an optional context object.

source
lll_gram_with_transformMethod
lll_gram_with_transform(x::ZZMatrix, ctx::LLLContext = LLLContext(0.99, 0.51, :gram))

Return a tuple $(L, T)$ where $L$ is the Gram matrix of an LLL-reduced basis of the lattice given by the Gram matrix $x$ and $T$ is a transformation matrix with $L = T^\top x T$. The matrix $x$ must be symmetric and non-singular.

See lll_gram for the used default parameters which can be overridden by supplying an optional context object.

source
lll_with_removalMethod
lll_with_removal(x::ZZMatrix, b::ZZRingElem, ctx::LLLContext = LLLContext(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.

source
lll_with_removal_transformMethod
lll_with_removal_transform(x::ZZMatrix, b::ZZRingElem, ctx::LLLContext = LLLContext(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$.

source
lll!Method
lll!(x::ZZMatrix, ctx::LLLContext = LLLContext(0.99, 0.51))

Compute an LLL-reduced basis of the $\mathbb{Z}$-lattice generated by the rows of $x$ inplace.

By default, the LLL is performed with reduction parameters $\delta = 0.99$ and $\eta = 0.51$. These defaults can be overridden by specifying an optional context object.

See lll_gram! for a function taking the Gram matrix as input.

source
lll_gram!Method
lll_gram!(x::ZZMatrix, ctx::LLLContext = LLLContext(0.99, 0.51, :gram))

Compute the Gram matrix of an LLL-reduced basis of the lattice given by the Gram matrix $x$ inplace. The matrix $x$ must be symmetric and non-singular.

By default, the LLL is performed with reduction parameters $\delta = 0.99$ and $\eta = 0.51$. These defaults can be overridden by specifying an optional context object.

source

Examples

julia> A = ZZ[2 3 5; 1 4 7; 19 3 7]
[ 2   3   5]
[ 1   4   7]
[19   3   7]

julia> L = lll(A, LLLContext(0.95, 0.55, :zbasis, :approx))
[-1    1   2]
[-1   -2   2]
[ 4    1   1]

julia> L, T = lll_with_transform(A)
([-1 1 2; -1 -2 2; 4 1 1], [-1 1 0; -15 10 1; 3 -2 0])

julia> G = lll_gram(gram(A))
[ 6    3   -1]
[ 3    9   -4]
[-1   -4   18]

julia> G, T = lll_gram_with_transform(gram(A))
([6 3 -1; 3 9 -4; -1 -4 18], [-1 1 0; -15 10 1; 3 -2 0])

julia> r, L = lll_with_removal(A, ZZ(100))
(3, [-1 1 2; -1 -2 2; 4 1 1])

julia> r, L, T = lll_with_removal_transform(A, ZZ(100))
(3, [-1 1 2; -1 -2 2; 4 1 1], [-1 1 0; -15 10 1; 3 -2 0])

Smith Normal Form

snfMethod
snf(x::ZZMatrix)

Compute the Smith normal form of $x$.

source
snf_diagonalMethod
snf_diagonal(x::ZZMatrix)

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

source
is_snfMethod
is_snf(x::ZZMatrix)

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

source

Examples

julia> A = ZZ[2 3 5; 1 4 7; 19 3 7]
[ 2   3   5]
[ 1   4   7]
[19   3   7]

julia> B = snf(A)
[1   0    0]
[0   1    0]
[0   0   27]

julia> is_snf(B) == true
true

julia> B = ZZ[2 0 0; 0 4 0; 0 0 7]
[2   0   0]
[0   4   0]
[0   0   7]

julia> C = snf_diagonal(B)
[1   0    0]
[0   2    0]
[0   0   28]

Strong Echelon Form

strong_echelon_formMethod
strong_echelon_form(a::zzModMatrix)

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

source
strong_echelon_formMethod
strong_echelon_form(a::fpMatrix)

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

source

Examples

julia> R, = residue_ring(ZZ, 12);

julia> A = R[4 1 0; 0 0 5; 0 0 0 ]
[4   1   0]
[0   0   5]
[0   0   0]

julia> B = strong_echelon_form(A)
[4   1   0]
[0   3   0]
[0   0   1]

Howell Form

howell_formMethod
howell_form(a::zzModMatrix)

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

source
howell_formMethod
howell_form(a::fpMatrix)

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

source

Examples

julia> R, = residue_ring(ZZ, 12);

julia> A = R[4 1 0; 0 0 5; 0 0 0 ]
[4   1   0]
[0   0   5]
[0   0   0]

julia> B = howell_form(A)
[4   1   0]
[0   3   0]
[0   0   1]

Gram-Schmidt Orthogonalisation

gram_schmidt_orthogonalisationMethod
gram_schmidt_orthogonalisation(x::QQMatrix)

Takes the columns of $x$ as the generators of a subset of $\mathbb{Q}^m$ and returns a matrix whose columns are an orthogonal generating set for the same subspace.

Examples

julia> S = matrix_space(QQ, 3, 3);

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

julia> B = gram_schmidt_orthogonalisation(A)
[4   -11//5     95//123]
[2    22//5   -190//123]
[0        5    209//123]
source

Exponential

Examples

julia> A = RR[2 0 0; 0 3 0; 0 0 1]
[2.0000000000000000000                       0                       0]
[                    0   3.0000000000000000000                       0]
[                    0                       0   1.0000000000000000000]

julia> B = exp(A)
[[7.389056098930650227 +/- 4.72e-19]                                     0                                     0]
[                                  0   [20.08553692318766774 +/- 1.94e-18]                                     0]
[                                  0                                     0   [2.718281828459045235 +/- 4.30e-19]]

Norm

bound_inf_normMethod
bound_inf_norm(x::RealMatrix)

Returns a non-negative element $z$ of type ArbFieldElem, such that $z$ is an upper bound for the infinity norm for every matrix in $x$

source
bound_inf_normMethod
bound_inf_norm(x::ComplexMatrix)

Returns a non-negative element $z$ of type AcbFieldElem, such that $z$ is an upper bound for the infinity norm for every matrix in $x$

source

Examples

julia> A = RR[1 2 3; 4 5 6; 7 8 9]
[1.0000000000000000000   2.0000000000000000000   3.0000000000000000000]
[4.0000000000000000000   5.0000000000000000000   6.0000000000000000000]
[7.0000000000000000000   8.0000000000000000000   9.0000000000000000000]

julia> d = bound_inf_norm(A)
[24.000000059604644775 +/- 3.91e-19]

Shifting

Examples

julia> A = RR[1 2 3; 4 5 6; 7 8 9]
[1.0000000000000000000   2.0000000000000000000   3.0000000000000000000]
[4.0000000000000000000   5.0000000000000000000   6.0000000000000000000]
[7.0000000000000000000   8.0000000000000000000   9.0000000000000000000]

julia> B = ldexp(A, 4)
[16.000000000000000000   32.000000000000000000   48.000000000000000000]
[64.000000000000000000   80.000000000000000000   96.000000000000000000]
[112.00000000000000000   128.00000000000000000   144.00000000000000000]

julia> overlaps(16*A, B)
true

Predicates

Examples

julia> A = CC[1 2 3; 4 5 6; 7 8 9]
[1.0000000000000000000   2.0000000000000000000   3.0000000000000000000]
[4.0000000000000000000   5.0000000000000000000   6.0000000000000000000]
[7.0000000000000000000   8.0000000000000000000   9.0000000000000000000]

julia> isreal(A)
true

julia> isreal(onei(CC)*A)
false

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 ZZMatrix:

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.

Eigenvalues and Eigenvectors (experimental)

eigenvaluesMethod
eigenvalues(A::ComplexMatrix)

Return the eigenvalues of A.

This function is experimental.

source
eigenvalues_with_multiplicitiesMethod
eigenvalues_with_multiplicities(A::ComplexMatrix)

Return the eigenvalues of A with their algebraic multiplicities as a vector of tuples (ComplexFieldElem, Int). Each tuple (z, k) corresponds to a cluster of k eigenvalues of $A$.

This function is experimental.

source
eigenvalues_simpleMethod
eigenvalues_simple(A::ComplexMatrix, algorithm::Symbol = :default)

Returns the eigenvalues of A as a vector of AcbFieldElem. It is assumed that A has only simple eigenvalues.

The algorithm used can be changed by setting the algorithm keyword to :vdhoeven_mourrain or :rump.

This function is experimental.

source
julia> A = CC[1 2 3; 0 4 5; 0 0 6]
[1.0000000000000000000   2.0000000000000000000   3.0000000000000000000]
[                    0   4.0000000000000000000   5.0000000000000000000]
[                    0                       0   6.0000000000000000000]

julia> eigenvalues_simple(A)
3-element Vector{ComplexFieldElem}:
 1.0000000000000000000
 4.0000000000000000000
 6.0000000000000000000

julia> A = CC[2 2 3; 0 2 5; 0 0 2]
[2.0000000000000000000   2.0000000000000000000   3.0000000000000000000]
[                    0   2.0000000000000000000   5.0000000000000000000]
[                    0                       0   2.0000000000000000000]

julia> eigenvalues(A)
1-element Vector{ComplexFieldElem}:
 2.0000000000000000000

julia> eigenvalues_with_multiplicities(A)
1-element Vector{Tuple{ComplexFieldElem, Int64}}:
 (2.0000000000000000000, 3)