Dense matrices over the integer ring¶
AUTHORS:
William Stein
Robert Bradshaw
Marc Masdeu (August 2014). Implemented using FLINT, see Issue #16803.
Jeroen Demeyer (October 2014): lots of fixes, see Issue #17090 and Issue #17094.
Vincent Delecroix (February 2015): make it faster, see Issue #17822.
Vincent Delecroix (May 2017): removed duplication of entries and cleaner linbox interface
EXAMPLES:
sage: a = matrix(ZZ, 3,3, range(9)); a
[0 1 2]
[3 4 5]
[6 7 8]
sage: a.det()
0
sage: a[0,0] = 10; a.det()
-30
sage: a.charpoly()
x^3 - 22*x^2 + 102*x + 30
sage: b = -3*a
sage: a == b
False
sage: b < a
True
- class sage.matrix.matrix_integer_dense.Matrix_integer_dense[source]¶
Bases:
Matrix_denseMatrix over the integers, implemented using FLINT.
On a 32-bit machine, they can have at most
rows or columns. On a 64-bit machine, matrices can have at most rows or columns.EXAMPLES:
sage: a = MatrixSpace(ZZ,3)(2); a [2 0 0] [0 2 0] [0 0 2] sage: a = matrix(ZZ,1,3, [1,2,-3]); a [ 1 2 -3] sage: a = MatrixSpace(ZZ,2,4)(2); a Traceback (most recent call last): ... TypeError: nonzero scalar matrix must be square
- BKZ(delta=None, algorithm='fpLLL', fp=None, block_size=10, prune=0, use_givens=False, precision=0, proof=None, **kwds)[source]¶
Return the result of running Block Korkin-Zolotarev reduction on
selfinterpreted as a lattice.INPUT:
delta– (default:0.99) LLL parameteralgorithm– (default:'fpLLL')'fpLLL'or"NTL"fp– floating point number implementationNone– NTL’s exact reduction or fpLLL’s wrapper (default)'fp'– double precision: NTL’s FP or fpLLL’s double'ld'– long doubles (fpLLL only)'qd'– NTL’s QP'qd1'– quad doubles: Usesquad_floatprecision to compute Gram-Schmidt, but uses double precision in the search phase of the block reduction algorithm. This seems adequate for most purposes, and is faster than'qd', which uses quad_float precision uniformly throughout (NTL only).'xd'– extended exponent: NTL’s XD or fpLLL’s dpe'rr'– arbitrary precision: NTL’RR or fpLLL’s MPFR
block_size– (default:10) specifies the size of the blocks in the reduction. High values yield shorter vectors, but the running time increases double exponentially withblock_size.block_sizeshould be between 2 and the number of rows ofself.proof– (default: same asproof.linear_algebra()) Insist on full BKZ reduction. If disabled and fplll is called, reduction is much faster but the result is not fully BKZ reduced.
NTL SPECIFIC INPUT:
prune– (default:0) the optional parameterprunecan be set to any positive number to invoke the Volume Heuristic from [SH1995]. This can significantly reduce the running time, and hence allow much bigger block size, but the quality of the reduction is of course not as good in general. Higher values ofprunemean better quality, and slower running time. Whenpruneis0, pruning is disabled. Recommended usage: forblock_size==30, set10 <= prune <=15.use_givens– use Givens orthogonalization. Only applies to approximate reduction using NTL. This is a bit slower, but generally much more stable, and is really the preferred orthogonalization strategy. For a nice description of this, see Chapter 5 of [GL1996].
fpLLL SPECIFIC INPUT:
precision– (default:0for automatic choice) bit precision to use iffp='rr'is set**kwds– keywords to be passed tofpylll; seefpylll.BKZ.Paramfor details
Also, if the verbose level is at least
, some output is printed during the computation.EXAMPLES:
sage: A = Matrix(ZZ,3,3,range(1,10)) sage: A.BKZ() [ 0 0 0] [ 2 1 0] [-1 1 3] sage: A = Matrix(ZZ,3,3,range(1,10)) sage: A.BKZ(use_givens=True) [ 0 0 0] [ 2 1 0] [-1 1 3] sage: A = Matrix(ZZ,3,3,range(1,10)) sage: A.BKZ(fp='fp') [ 0 0 0] [ 2 1 0] [-1 1 3]
ALGORITHM:
Calls either NTL or fpLLL.
- LLL(delta=None, eta=None, algorithm='fpLLL:wrapper', fp=None, prec=0, early_red=False, use_givens=False, use_siegel=False, transformation=False, **kwds)[source]¶
Return LLL-reduced or approximated LLL reduced matrix
of the lattice generated by the rows ofself.A set of vectors
is -LLL-reduced if the two following conditions hold:For any
, we have .For any
, we have ,
where
and is the -th vector of the Gram-Schmidt orthogonalisation of .The default reduction parameters are
and . The parameters and must satisfy and . Polynomial time complexity is only guaranteed for . Not every algorithm admits the case .If the matrix has a nonzero kernel, the LLL-reduced matrix will contain zero rows, so that the output has the same dimensions as the input. The transformation matrix is always invertible over the integers.
Also the rank of
selfis cached if it is computed during the reduction. Note that in general this only happens whenself.rank() == self.ncols()and the exact algorithm is used.INPUT:
delta– (default:0.99) parameter as described above, ignored by parieta– (default:0.501) parameter as described above, ignored by NTL and parialgorithm– string; one of the algorithms listed below (default:'fpLLL:wrapper')fp– floating point number implementation, ignored by pari:None– NTL’s exact reduction or fpLLL’s wrapper'fp'– double precision: NTL’s FP or fpLLL’s double'ld'– long doubles (fpLLL only)'qd'– NTL’s QP'xd'– extended exponent: NTL’s XD or fpLLL’s dpe'rr'– arbitrary precision: NTL’s RR or fpLLL’s MPFR
prec– (default: auto choose) precision, ignored by NTL and pariearly_red– boolean (default:False); perform early reduction, ignored by NTL and pariuse_givens– boolean (default:False); use Givens orthogonalization. Only applies to approximate reduction using NTL. This is slower but generally more stable.use_siegel– boolean (default:False); use Siegel’s condition instead of Lovász’s condition, ignored by NTL and paritransformation– boolean (default:False); also return transformation matrix**kwds– keywords to be passed tofpylll; seefpylll.LLL.reduction()for details
Also, if the verbose level is at least
, some output is printed during the computation.AVAILABLE ALGORITHMS:
'NTL:LLL'– NTL’s LLL + choice offp'fpLLL:heuristic'– fpLLL’s heuristic + choice offp'fpLLL:fast'– fpLLL’s fast + choice offp'fpLLL:proved'– fpLLL’s proved + choice offp'fpLLL:wrapper'– fpLLL’s automatic choice (default)'pari'– pari’s qflll'flatter'– external executableflatter, requires manual install from https://github.com/keeganryan/flatter. When the input matrix does not have full row rank, versions before https://github.com/keeganryan/flatter/pull/23 would error out, versions after that might work but remove zero rows, unlike'fplll'algorithm. Note that sufficiently new version ofparialso supports FLATTER algorithm, see pari:qflll.
OUTPUT: a matrix over the integers
EXAMPLES:
sage: A = Matrix(ZZ,3,3,range(1,10)) sage: A.LLL() [ 0 0 0] [ 2 1 0] [-1 1 3]
We compute the extended GCD of a list of integers using LLL, this example is from the Magma handbook:
sage: Q = [ 67015143, 248934363018, 109210, 25590011055, 74631449, ....: 10230248, 709487, 68965012139, 972065, 864972271 ] sage: n = len(Q) sage: S = 100 sage: X = Matrix(ZZ, n, n + 1) sage: for i in range(n): ....: X[i, i + 1] = 1 sage: for i in range(n): ....: X[i, 0] = S * Q[i] sage: L = X.LLL() sage: M = L.row(n-1).list()[1:] sage: M [-3, -1, 13, -1, -4, 2, 3, 4, 5, -1] sage: add(Q[i]*M[i] for i in range(n)) -1
The case
is not always supported:sage: L = X.LLL(delta=2) Traceback (most recent call last): ... TypeError: delta must be <= 1 sage: L = X.LLL(delta=1) # not tested, will eat lots of ram Traceback (most recent call last): ... RuntimeError: infinite loop in LLL sage: L = X.LLL(delta=1, algorithm='NTL:LLL') sage: L[-1] (-100, -3, -1, 13, -1, -4, 2, 3, 4, 5, -1)
We return the transformation matrix:
sage: A = random_matrix(ZZ, 10, 20) sage: R, U = A.LLL(transformation=True) sage: U * A == R True sage: R, U = A.LLL(algorithm='NTL:LLL', transformation=True) sage: U * A == R True sage: R, U = A.LLL(algorithm='pari', transformation=True) sage: U * A == R True
Example with a nonzero kernel:
sage: M = matrix(4,3,[1,2,3,2,4,6,7,0,1,-1,-2,-3]) sage: M.LLL()[0:2] [0 0 0] [0 0 0] sage: M.LLL(algorithm="NTL:LLL")[0:2] [0 0 0] [0 0 0] sage: M.LLL(algorithm='pari')[0:2] [0 0 0] [0 0 0]
Note
See
sage.libs.ntl.ntl_mat_ZZ.ntl_mat_ZZ.LLLandfpylll.fplll.lllfor details on the algorithms used.Although LLL is a deterministic algorithm, the output for different implementations and CPUs (32-bit vs. 64-bit) may vary, while still being correct.
Check
flatter:sage: # needs flatter sage: M = matrix(ZZ, 2, 2, [-1,1,1,1]) sage: L = M.LLL(algorithm="flatter") sage: abs(M.det()) == abs(L.det()) True sage: L = M.LLL(algorithm="flatter", delta=0.99) sage: abs(M.det()) == abs(L.det()) True
In sufficiently new versions of flatter, the following works:
sage: # needs flatter sage: matrix.identity(3).stack(matrix.identity(3)).LLL(algorithm="flatter") [1 0 0] [0 1 0] [0 0 1] sage: matrix.identity(4)[:,:1].LLL(algorithm="flatter") [1] sage: matrix.zero(1, 2).LLL(algorithm="flatter") [] sage: matrix.zero(2, 1).LLL(algorithm="flatter") []
- antitranspose()[source]¶
Return the antitranspose of
self, without changingself.EXAMPLES:
sage: A = matrix(2,3,range(6)) sage: type(A) <class 'sage.matrix.matrix_integer_dense.Matrix_integer_dense'> sage: A.antitranspose() [5 2] [4 1] [3 0] sage: A [0 1 2] [3 4 5] sage: A.subdivide(1,2); A [0 1|2] [---+-] [3 4|5] sage: A.antitranspose() [5|2] [-+-] [4|1] [3|0]
- augment(right, subdivide=False)[source]¶
Return a new matrix formed by appending the matrix (or vector)
righton the right side ofself.INPUT:
right– a matrix, vector or free module element, whose dimensions are compatible withselfsubdivide– (default:False) request the resulting matrix to have a new subdivision, separatingselffromright
OUTPUT:
A new matrix formed by appending
rightonto the right side ofself. Ifrightis a vector (or free module element) then in this context it is appropriate to consider it as a column vector. (The code first converts a vector to a 1-column matrix.)EXAMPLES:
sage: A = matrix(ZZ, 4, 5, range(20)) sage: B = matrix(ZZ, 4, 3, range(12)) sage: A.augment(B) [ 0 1 2 3 4 0 1 2] [ 5 6 7 8 9 3 4 5] [10 11 12 13 14 6 7 8] [15 16 17 18 19 9 10 11]
A vector may be augmented to a matrix.
sage: A = matrix(ZZ, 3, 5, range(15)) sage: v = vector(ZZ, 3, range(3)) sage: A.augment(v) [ 0 1 2 3 4 0] [ 5 6 7 8 9 1] [10 11 12 13 14 2]
The
subdivideoption will add a natural subdivision betweenselfandright. For more details about how subdivisions are managed when augmenting, seesage.matrix.matrix1.Matrix.augment().sage: A = matrix(ZZ, 3, 5, range(15)) sage: B = matrix(ZZ, 3, 3, range(9)) sage: A.augment(B, subdivide=True) [ 0 1 2 3 4| 0 1 2] [ 5 6 7 8 9| 3 4 5] [10 11 12 13 14| 6 7 8]
Errors are raised if the sizes are incompatible.
sage: A = matrix(ZZ, [[1, 2],[3, 4]]) sage: B = matrix(ZZ, [[10, 20], [30, 40], [50, 60]]) sage: A.augment(B) Traceback (most recent call last): ... TypeError: number of rows must be the same, not 2 != 3
- charpoly(var='x', algorithm=None)[source]¶
Note
The characteristic polynomial is defined as
.INPUT:
var– a variable namealgorithm– (default:'linbox') either'generic','flint'or'linbox'
EXAMPLES:
sage: A = matrix(ZZ,6, range(36)) sage: f = A.charpoly(); f x^6 - 105*x^5 - 630*x^4 sage: f(A) == 0 True sage: g = A.charpoly(algorithm='flint') sage: f == g True sage: n=20; A = Mat(ZZ,n)(range(n^2)) sage: A.charpoly() x^20 - 3990*x^19 - 266000*x^18 sage: A.minpoly() x^3 - 3990*x^2 - 266000*x
On non square matrices, this method raises an ArithmeticError:
sage: matrix(ZZ, 2, 1).charpoly() Traceback (most recent call last): ... ArithmeticError: only valid for square matrix
- column(i, from_list=False)[source]¶
Return the
-th column of this matrix as a dense vector.INPUT:
i– integerfrom_list– ignored
EXAMPLES:
sage: m = matrix(ZZ, 3, 2, [1, -2, 3, 4, -1, 0]) sage: m.column(1) (-2, 4, 0) sage: m.column(1, from_list=True) (-2, 4, 0) sage: m.column(-1) (-2, 4, 0) sage: m.column(-2) (1, 3, -1) sage: m.column(2) Traceback (most recent call last): ... IndexError: column index out of range sage: m.column(-3) Traceback (most recent call last): ... IndexError: column index out of range
- decomposition(**kwds)[source]¶
Return the decomposition of the free module on which this matrix A acts from the right (i.e., the action is x goes to x A), along with whether this matrix acts irreducibly on each factor. The factors are guaranteed to be sorted in the same way as the corresponding factors of the characteristic polynomial, and are saturated as ZZ modules.
INPUT:
self– a matrix over the integers**kwds– these are passed onto to the decomposition over QQ command
EXAMPLES:
sage: t = ModularSymbols(11,sign=1).hecke_matrix(2) sage: w = t.change_ring(ZZ) sage: w [ 3 -2] [ 0 -2] sage: w.charpoly().factor() (x - 3) * (x + 2) sage: w.decomposition() [(Free module of degree 2 and rank 1 over Integer Ring Echelon basis matrix: [ 5 -2], True), (Free module of degree 2 and rank 1 over Integer Ring Echelon basis matrix: [0 1], True)]
- determinant(algorithm='default', proof=None, stabilize=2)[source]¶
Return the determinant of this matrix.
INPUT:
algorithm'default'– useflint'flint'– let flint do the determinant'padic'– uses a -adic / multimodular algorithm that relies on code in IML and linbox'linbox'– calls linbox det (you must set proof=False to use this!)'ntl'– calls NTL’s det function'pari'– uses PARI
proof– boolean orNone; ifNoneuse proof.linear_algebra(); only relevant for the padic algorithmNote
It would be VERY VERY hard for det to fail even with proof=False.
stabilize– if proof isFalse, require det to be the same for this many CRT primes in a row. Ignored if proof isTrue.
ALGORITHM: The
-adic algorithm works by first finding a random vector v, then solving and taking the denominator . This gives a divisor of the determinant. Then we compute using a multimodular algorithm and the Hadamard bound, skipping primes that divide .EXAMPLES:
sage: A = matrix(ZZ,8,8,[3..66]) sage: A.determinant() 0
sage: A = random_matrix(ZZ,20,20) sage: D1 = A.determinant() sage: A._clear_cache() sage: D2 = A.determinant(algorithm='ntl') sage: D1 == D2 True
We have a special-case algorithm for 4 x 4 determinants:
sage: A = matrix(ZZ,4,[1,2,3,4,4,3,2,1,0,5,0,1,9,1,2,3]) sage: A.determinant() 270
Next we try the Linbox det. Note that we must have proof=False.
sage: A = matrix(ZZ,5,[1,2,3,4,5,4,6,3,2,1,7,9,7,5,2,1,4,6,7,8,3,2,4,6,7]) sage: A.determinant(algorithm='linbox') Traceback (most recent call last): ... RuntimeError: you must pass the proof=False option to the determinant command to use LinBox's det algorithm sage: A.determinant(algorithm='linbox', proof=False) -21 sage: A._clear_cache() sage: A.determinant() -21
Try the other algorithms on the same example:
sage: A._clear_cache(); A.determinant(algorithm='padic') -21 sage: A._clear_cache(); A.determinant(algorithm='pari') -21 sage: A._clear_cache(); A.determinant(algorithm='ntl') -21 sage: A._clear_cache(); A.determinant(algorithm='padic') -21
A bigger example:
sage: A = random_matrix(ZZ,30) sage: d = A.determinant() sage: A._clear_cache() sage: A.determinant(algorithm='linbox',proof=False) == d True
- echelon_form(algorithm='default', proof=None, include_zero_rows=True, transformation=False, D=None)[source]¶
Return the echelon form of this matrix over the integers, also known as the hermite normal form (HNF).
INPUT:
algorithm– string; the algorithm to use. Valid options are:'default'– let Sage pick an algorithm (default). Up to 75 rows or columns with no transformation matrix, use pari with flag 0; otherwise, use flint.'flint'– use flint'ntl'– use NTL (only works for square matrices of full rank!)'padic'– an asymptotically fast -adic modular algorithm, If your matrix has large coefficients and is small, you may also want to try this.'pari'– use PARI with flag 1'pari0'– use PARI with flag 0'pari1'– use PARI with flag 1'pari4'– use PARI with flag 4 (use heuristic LLL)
proof– (default:True) if proof=False certain determinants are computed using a randomized hybrid -adic multimodular strategy until it stabilizes twice (instead of up to the Hadamard bound). It is incredibly unlikely that one would ever get an incorrect result with proof=False.include_zero_rows– boolean (default:True); ifFalse, don’t include zero rowstransformation– if given, also compute transformation matrix; only valid for flint and padic algorithmD– (default:None) if given and the algorithm is'ntl', then D must be a multiple of the determinant and this function will use that fact
OUTPUT:
The Hermite normal form (=echelon form over
) ofselfas an immutable matrix.EXAMPLES:
sage: A = MatrixSpace(ZZ,2)([1,2,3,4]) sage: A.echelon_form() [1 0] [0 2] sage: A = MatrixSpace(ZZ,5)(range(25)) sage: A.echelon_form() [ 5 0 -5 -10 -15] [ 0 1 2 3 4] [ 0 0 0 0 0] [ 0 0 0 0 0] [ 0 0 0 0 0]
Getting a transformation matrix in the nonsquare case:
sage: A = matrix(ZZ,5,3,[1..15]) sage: H, U = A.hermite_form(transformation=True, include_zero_rows=False) sage: H [1 2 3] [0 3 6] sage: U [ 0 0 0 4 -3] [ 0 0 0 13 -10] sage: U*A == H True
Note
If ‘ntl’ is chosen for a non square matrix this function raises a
ValueError.Special cases: 0 or 1 rows:
sage: a = matrix(ZZ, 1,2,[0,-1]) sage: a.hermite_form() [0 1] sage: a.pivots() (1,) sage: a = matrix(ZZ, 1,2,[0,0]) sage: a.hermite_form() [0 0] sage: a.pivots() () sage: a = matrix(ZZ,1,3); a [0 0 0] sage: a.echelon_form(include_zero_rows=False) [] sage: a.echelon_form(include_zero_rows=True) [0 0 0]
Illustrate using various algorithms.:
sage: matrix(ZZ,3,[1..9]).hermite_form(algorithm='pari') [1 2 3] [0 3 6] [0 0 0] sage: matrix(ZZ,3,[1..9]).hermite_form(algorithm='pari0') [1 2 3] [0 3 6] [0 0 0] sage: matrix(ZZ,3,[1..9]).hermite_form(algorithm='pari4') [1 2 3] [0 3 6] [0 0 0] sage: matrix(ZZ,3,[1..9]).hermite_form(algorithm='padic') [1 2 3] [0 3 6] [0 0 0] sage: matrix(ZZ,3,[1..9]).hermite_form(algorithm='default') [1 2 3] [0 3 6] [0 0 0]
The ‘ntl’ algorithm doesn’t work on matrices that do not have full rank.:
sage: matrix(ZZ,3,[1..9]).hermite_form(algorithm='ntl') Traceback (most recent call last): ... ValueError: ntl only computes HNF for square matrices of full rank. sage: matrix(ZZ,3,[0] +[2..9]).hermite_form(algorithm='ntl') [1 0 0] [0 1 0] [0 0 3]
- elementary_divisors(algorithm='pari')[source]¶
Return the elementary divisors of self, in order.
Warning
This is MUCH faster than the
smith_form()function.The elementary divisors are the invariants of the finite abelian group that is the cokernel of left multiplication of this matrix. They are ordered in reverse by divisibility.
INPUT:
self– matrixalgorithm– (default:'pari')'pari'– works robustly, but is slower.'linbox'– use linbox (currently off, broken)
OUTPUT: list of integers
Note
These are the invariants of the cokernel of left multiplication:
sage: M = Matrix([[3,0,1],[0,1,0]]) sage: M [3 0 1] [0 1 0] sage: M.elementary_divisors() [1, 1] sage: M.transpose().elementary_divisors() [1, 1, 0]
EXAMPLES:
sage: matrix(3, range(9)).elementary_divisors() [1, 3, 0] sage: matrix(3, range(9)).elementary_divisors(algorithm='pari') [1, 3, 0] sage: C = MatrixSpace(ZZ,4)([3,4,5,6,7,3,8,10,14,5,6,7,2,2,10,9]) sage: C.elementary_divisors() [1, 1, 1, 687]
sage: M = matrix(ZZ, 3, [1,5,7, 3,6,9, 0,1,2]) sage: M.elementary_divisors() [1, 1, 6]
This returns a copy, which is safe to change:
sage: edivs = M.elementary_divisors() sage: edivs.pop() 6 sage: M.elementary_divisors() [1, 1, 6]
See also
- frobenius(*args, **kwds)[source]¶
Deprecated: Use
frobenius_form()instead. See Issue #36396 for details.
- frobenius_form(flag=0, var='x')[source]¶
Return the Frobenius form (rational canonical form) of this matrix.
INPUT:
flag– 0 (default), 1 or 2 as follows:0– (default) return the Frobenius form of thismatrix
1– return only the elementary divisorpolynomials, as polynomials in var
2– return a two-components vector [F,B] where Fis the Frobenius form and B is the basis change so that
var– string (default:'x')
ALGORITHM: uses PARI’s pari:matfrobenius
EXAMPLES:
sage: A = MatrixSpace(ZZ, 3)(range(9)) sage: A.frobenius_form(0) [ 0 0 0] [ 1 0 18] [ 0 1 12] sage: A.frobenius_form(1) [x^3 - 12*x^2 - 18*x] sage: A.frobenius_form(1, var='y') [y^3 - 12*y^2 - 18*y] sage: F, B = A.frobenius_form(2) sage: A == B^(-1)*F*B True sage: a=matrix([]) sage: a.frobenius_form(2) ([], []) sage: a.frobenius_form(0) [] sage: a.frobenius_form(1) [] sage: B = random_matrix(ZZ,2,3) sage: B.frobenius_form() Traceback (most recent call last): ... ArithmeticError: frobenius matrix of non-square matrix not defined.
AUTHORS:
Martin Albrecht (2006-04-02)
TODO: - move this to work for more general matrices than just over Z. This will require fixing how PARI polynomials are coerced to Sage polynomials.
- gcd()[source]¶
Return the gcd of all entries of self; very fast.
EXAMPLES:
sage: a = matrix(ZZ,2, [6,15,-6,150]) sage: a.gcd() 3
- height()[source]¶
Return the height of this matrix, i.e., the max absolute value of the entries of the matrix.
OUTPUT: nonnegative integer
EXAMPLES:
sage: a = Mat(ZZ,3)(range(9)) sage: a.height() 8 sage: a = Mat(ZZ,2,3)([-17,3,-389,15,-1,0]); a [ -17 3 -389] [ 15 -1 0] sage: a.height() 389
- hermite_form(algorithm='default', proof=None, include_zero_rows=True, transformation=False, D=None)[source]¶
alias of
echelon_form().
- index_in_saturation(proof=None)[source]¶
Return the index of
selfin its saturation.INPUT:
proof– (default: use proof.linear_algebra()); ifFalse, the determinant calculations are done withproof=False
OUTPUT: positive integer; the index of the row span of this matrix in its saturation
ALGORITHM:
Use Hermite normal form twice to find an invertible matrix whose inverse transforms a matrix with the same row span as
selfto its saturation, then compute the determinant of that matrix.EXAMPLES:
sage: A = matrix(ZZ, 2,3, [1..6]); A [1 2 3] [4 5 6] sage: A.index_in_saturation() 3 sage: A.saturation() [1 2 3] [1 1 1]
- insert_row(index, row)[source]¶
Create a new matrix from
selfwith.INPUT:
index– integerrow– a vector
EXAMPLES:
sage: X = matrix(ZZ,3,range(9)); X [0 1 2] [3 4 5] [6 7 8] sage: X.insert_row(1, [1,5,-10]) [ 0 1 2] [ 1 5 -10] [ 3 4 5] [ 6 7 8] sage: X.insert_row(0, [1,5,-10]) [ 1 5 -10] [ 0 1 2] [ 3 4 5] [ 6 7 8] sage: X.insert_row(3, [1,5,-10]) [ 0 1 2] [ 3 4 5] [ 6 7 8] [ 1 5 -10]
- integer_valued_polynomials_generators()[source]¶
Determine the generators of the ring of integer valued polynomials on this matrix.
OUTPUT:
A pair
(mu_B, P)wherePis a list of polynomials in such thatwhere
is this matrix.EXAMPLES:
sage: B = matrix(ZZ, [[1, 0, 1], [1, -2, -1], [10, 0, 0]]) sage: B.integer_valued_polynomials_generators() (x^3 + x^2 - 12*x - 20, [1, 1/4*x^2 + 3/4*x + 1/2])
- inverse_of_unit()[source]¶
If
selfis a matrix with determinant or return the inverse ofselfas a matrix over .EXAMPLES:
sage: m = matrix(ZZ, 2, [2,1,1,1]).inverse_of_unit() sage: m [ 1 -1] [-1 2] sage: parent(m) Full MatrixSpace of 2 by 2 dense matrices over Integer Ring sage: matrix(2, [2,1,0,1]).inverse_of_unit() Traceback (most recent call last): ... ArithmeticError: non-invertible matrix
- is_LLL_reduced(delta=None, eta=None, algorithm='fpLLL')[source]¶
Return
Trueif this lattice is -LLL reduced. Seeself.LLLfor a definition of LLL reduction.INPUT:
delta– (default: ) parameter as described aboveeta– (default: ) parameter as described abovealgorithm– either'fpLLL'(default) or'sage'
EXAMPLES:
sage: A = random_matrix(ZZ, 10, 10) sage: L = A.LLL() sage: A.is_LLL_reduced() False sage: L.is_LLL_reduced() True
The
'sage'algorithm currently does not work for matrices with linearly dependent rows:sage: A = matrix(ZZ, [[1, 2, 3], [2, 4, 6]]) sage: A.is_LLL_reduced(algorithm='sage') Traceback (most recent call last): ... ValueError: linearly dependent input for module version of Gram-Schmidt
- is_one()[source]¶
Test whether
selfis the identity matrix.EXAMPLES:
sage: matrix(2, [1,0,0,1]).is_one() True sage: matrix(2, [1,1,0,1]).is_one() False sage: matrix(2, 3, [1,0,0,0,1,0]).is_one() False
- is_primitive()[source]¶
Test whether the matrix is primitive.
An integral matrix
is primitive if all its entries are nonnegative and for some positive integer the matrix has all its entries positive.EXAMPLES:
sage: m = matrix(3, [1,1,0,0,0,1,1,0,0]) sage: m.is_primitive() True sage: m**4 [3 2 1] [1 1 1] [2 1 1] sage: m = matrix(4, [[1,1,0,0],[0,0,1,0],[0,0,0,1],[1,0,0,0]]) sage: m.is_primitive() True sage: m**6 [4 3 2 1] [1 1 1 1] [2 1 1 1] [3 2 1 1] sage: m = matrix(4, [[0,1,0,1],[1,0,1,0],[0,1,0,1],[1,0,1,0]]) sage: m.is_primitive() False
Testing extremal matrices:
sage: def matrix1(d): ....: m = matrix(d) ....: m[0,0] = 1 ....: for i in range(d-1): ....: m[i,i+1] = m[i+1,i] = 1 ....: return m sage: all(matrix1(d).is_primitive() for d in range(2,20)) True sage: def matrix2(d): ....: m = matrix(d) ....: for i in range(d-1): ....: m[i,i+1] = 1 ....: m[d-1,0] = m[d-1,1] = 1 ....: return m sage: all(matrix2(d).is_primitive() for d in range(2,20)) True
Non-primitive families:
sage: def matrix3(d): ....: m = matrix(d) ....: m[0,0] = 1 ....: for i in range(d-1): ....: m[i,i+1] = 1 ....: return m sage: any(matrix3(d).is_primitive() for d in range(2,20)) False
- minpoly(var='x', algorithm=None)[source]¶
INPUT:
var– a variable namealgorithm– either'linbox'(default) or'generic'
EXAMPLES:
sage: A = matrix(ZZ, 6, range(36)) sage: A.minpoly() x^3 - 105*x^2 - 630*x sage: A = Mat(ZZ, 6)([k^2 for k in range(36)]) sage: A.minpoly(algorithm='linbox') x^4 - 2695*x^3 - 257964*x^2 + 1693440*x sage: A.minpoly(algorithm='generic') x^4 - 2695*x^3 - 257964*x^2 + 1693440*x
On non square matrices, this method raises an ArithmeticError:
sage: matrix(ZZ, 2, 1).minpoly() Traceback (most recent call last): ... ArithmeticError: only valid for square matrix
- null_ideal(b=0)[source]¶
Return the
-ideal of this matrix.Let
be a matrix. The null ideal modulo , or -ideal, isINPUT:
b– an element of (default: 0)
OUTPUT: an ideal in
EXAMPLES:
sage: B = matrix(ZZ, [[1, 0, 1], [1, -2, -1], [10, 0, 0]]) sage: B.null_ideal() Principal ideal (x^3 + x^2 - 12*x - 20) of Univariate Polynomial Ring in x over Integer Ring sage: B.null_ideal(8) Ideal (8, x^3 + x^2 - 12*x - 20, 2*x^2 + 6*x + 4) of Univariate Polynomial Ring in x over Integer Ring sage: B.null_ideal(6) Ideal (6, 2*x^3 + 2*x^2 - 24*x - 40, 3*x^2 + 3*x) of Univariate Polynomial Ring in x over Integer Ring
See also
- p_minimal_polynomials(p, s_max=None)[source]¶
Compute
-minimal polynomials of this matrix.For
, a -minimal polynomial of a matrix is a monic polynomial of minimal degree such that all entries of are divisible by .Compute a finite subset
of the positive integers and -minimal polynomials for .For
, a -minimal polynomial is given by where . For , the minimal polynomial of is also a -minimal polynomial.INPUT:
p– a prime ins_max– positive integer (default:None); if set, only -minimal polynomials fors <= s_maxare computed (see below for details)
OUTPUT:
A dictionary. Keys are the finite set
, the values are the associated -minimal polynomials , .Setting
s_maxonly affects the output ifs_maxis at most where denotes the full set. In that case, only those withs <= s_maxare returned wheres_maxis always included even if it is not included in the full set .EXAMPLES:
sage: B = matrix(ZZ, [[1, 0, 1], [1, -2, -1], [10, 0, 0]]) sage: B.p_minimal_polynomials(2) {2: x^2 + 3*x + 2}
See also
- pivots()[source]¶
Return the pivot column positions of this matrix.
OUTPUT: a tuple of Python integers: the position of the first nonzero entry in each row of the echelon form.
EXAMPLES:
sage: n = 3; A = matrix(ZZ,n,range(n^2)); A [0 1 2] [3 4 5] [6 7 8] sage: A.pivots() (0, 1) sage: A.echelon_form() [ 3 0 -3] [ 0 1 2] [ 0 0 0]
- prod_of_row_sums(cols)[source]¶
Return the product of the sums of the entries in the submatrix of
selfwith given columns.INPUT:
cols– list (or set) of integers representing columns ofself
OUTPUT: integer
EXAMPLES:
sage: a = matrix(ZZ,2,3,[1..6]); a [1 2 3] [4 5 6] sage: a.prod_of_row_sums([0,2]) 40 sage: (1+3)*(4+6) 40 sage: a.prod_of_row_sums(set([0,2])) 40
- randomize(density=1, x=None, y=None, distribution=None, nonzero=False)[source]¶
Randomize
densityproportion of the entries of this matrix, leaving the rest unchanged.The parameters are the same as the ones for the integer ring’s
random_elementfunction.If
xandyare given, randomized entries of this matrix have to be betweenxandyand have density 1.INPUT:
self– a mutable matrix over ZZdensity– a float between 0 and 1x,y– if notNone, these are passed to theZZ.random_elementfunction as the upper and lower endpoints in the uniform distribution
distribution– would also be passed intoZZ.random_elementif givennonzero– boolean (default:False); whether the new entries are guaranteed to be zero
OUTPUT: None, the matrix is modified in-place
EXAMPLES:
sage: A = matrix(ZZ, 2,3, [1..6]) sage: ranks = [True, True, True] sage: while any(ranks): ....: A.randomize() ....: ranks[A.rank()] = False sage: mini = 0 sage: maxi = 0 sage: while mini != -30 and maxi != 30: ....: A.randomize(x=-30, y=30) ....: mini = min(min(A.list()), mini) ....: maxi = min(min(A.list()), maxi)
- rank(algorithm='modp')[source]¶
Return the rank of this matrix.
INPUT:
algorithm– either'modp'(default) or'flint'or'linbox'
OUTPUT: nonnegative integer – the rank
Note
The rank is cached.
ALGORITHM:
If set to
'modp', first check if the matrix has maximum possible rank by working modulo one random prime. If not, call LinBox’s rank function.EXAMPLES:
sage: a = matrix(ZZ,2,3,[1..6]); a [1 2 3] [4 5 6] sage: a.rank() 2 sage: a = matrix(ZZ,3,3,[1..9]); a [1 2 3] [4 5 6] [7 8 9] sage: a.rank() 2
Here is a bigger example - the rank is of course still 2:
sage: a = matrix(ZZ,100,[1..100^2]); a.rank() 2
- rational_reconstruction(N)[source]¶
Use rational reconstruction to lift
selfto a matrix over the rational numbers (if possible), where we viewselfas a matrix modulo N.INPUT:
N– integer
OUTPUT: matrix over
or raise aValueErrorEXAMPLES: We create a random 4x4 matrix over ZZ.
sage: A = matrix(ZZ, 4, [4, -4, 7, 1, -1, 1, -1, -12, -1, -1, 1, -1, -3, 1, 5, -1])
There isn’t a unique rational reconstruction of it:
sage: A.rational_reconstruction(11) Traceback (most recent call last): ... ValueError: rational reconstruction does not exist
We throw in a denominator and reduce the matrix modulo 389 - it does rationally reconstruct.
sage: B = (A/3 % 389).change_ring(ZZ) sage: B.rational_reconstruction(389) == A/3 True
- row(i, from_list=False)[source]¶
Return the
-th row of this matrix as a dense vector.INPUT:
i– integerfrom_list– ignored
EXAMPLES:
sage: m = matrix(ZZ, 2, [1, -2, 3, 4]) sage: m.row(0) (1, -2) sage: m.row(1) (3, 4) sage: m.row(1, from_list=True) (3, 4) sage: m.row(-2) (1, -2) sage: m.row(2) Traceback (most recent call last): ... IndexError: row index out of range sage: m.row(-3) Traceback (most recent call last): ... IndexError: row index out of range
- saturation(p=0, proof=None, max_dets=5)[source]¶
Return a saturation matrix of self, which is a matrix whose rows span the saturation of the row span of
self. This is not unique.The saturation of a
module embedded in is a module that contains with finite index such that is torsion free. This function takes the row span of self, and finds another matrix of full rank with row span the saturation of .INPUT:
p– (default: 0) if nonzero given, saturate only at the prime , i.e., return a matrix whose row span is a -module that containsselfand such that the index of in its saturation is coprime to . If is None, return full saturation ofself.proof– (default: use proof.linear_algebra()); ifFalse, the determinant calculations are done withproof=Falsemax_dets– (default: 5) technical parameter - max number of determinant to compute when bounding prime divisor ofselfin its saturation.
OUTPUT: matrix over ZZ
Note
The result is not cached.
ALGORITHM: 1. Replace input by a matrix of full rank got from a subset of the rows. 2. Divide out any common factors from rows. 3. Check max_dets random dets of submatrices to see if their GCD (with p) is 1 - if so matrix is saturated and we’re done. 4. Finally, use that if A is a matrix of full rank, then
is a saturation of A.EXAMPLES:
sage: A = matrix(ZZ, 3, 5, [-51, -1509, -71, -109, -593, -19, -341, 4, 86, 98, 0, -246, -11, 65, 217]) sage: A.echelon_form() [ 1 5 2262 20364 56576] [ 0 6 35653 320873 891313] [ 0 0 42993 386937 1074825] sage: S = A.saturation(); S [ -51 -1509 -71 -109 -593] [ -19 -341 4 86 98] [ 35 994 43 51 347]
Notice that the saturation spans a different module than A.
sage: S.echelon_form() [ 1 2 0 8 32] [ 0 3 0 -2 -6] [ 0 0 1 9 25] sage: V = A.row_space(); W = S.row_space() sage: V.is_submodule(W) True sage: V.index_in(W) 85986 sage: V.index_in_saturation() 85986
We illustrate each option:
sage: S = A.saturation(p=2) sage: S = A.saturation(proof=False) sage: S = A.saturation(max_dets=2)
- smith_form(transformation=True, integral=None)[source]¶
Return the smith normal form of this matrix, that is the diagonal matrix
with diagonal entries the ordered elementary divisors of this matrix.INPUT:
transformation– boolean (default:True); whether to return the transformation matrices and such thatintegral– a subring of the base ring orTrue(default:None); ignored for matrices with integer entries
Note
The
elementary_divisors()function, which returns the diagonal entries of , is VASTLY faster than this function.The elementary divisors are the invariants of the finite abelian group that is the cokernel of this matrix. They are ordered in reverse by divisibility.
EXAMPLES:
sage: A = MatrixSpace(IntegerRing(), 3)(range(9)) sage: D, U, V = A.smith_form() sage: D [1 0 0] [0 3 0] [0 0 0] sage: U [ 0 2 -1] [ 0 -1 1] [ 1 -2 1] sage: V [ 0 0 1] [-1 2 -2] [ 1 -1 1] sage: U*A*V [1 0 0] [0 3 0] [0 0 0]
It also makes sense for nonsquare matrices:
sage: A = Matrix(ZZ,3,2,range(6)) sage: D, U, V = A.smith_form() sage: D [1 0] [0 2] [0 0] sage: U [ 0 2 -1] [ 0 -1 1] [ 1 -2 1] sage: V [-1 1] [ 1 0] sage: U * A * V [1 0] [0 2] [0 0]
Empty matrices are handled sensibly (see Issue #3068):
sage: m = MatrixSpace(ZZ, 2,0)(0); d,u,v = m.smith_form(); u*m*v == d True sage: m = MatrixSpace(ZZ, 0,2)(0); d,u,v = m.smith_form(); u*m*v == d True sage: m = MatrixSpace(ZZ, 0,0)(0); d,u,v = m.smith_form(); u*m*v == d True
See also
- symplectic_form()[source]¶
Find a symplectic basis for
selfifselfis an anti-symmetric, alternating matrix.Return a pair (F, C) such that the rows of C form a symplectic basis for
selfandF = C * self * C.transpose().Raise a
ValueErrorifselfis not anti-symmetric, orselfis not alternating.Anti-symmetric means that
. Alternating means that the diagonal of is identically zero.A symplectic basis is a basis of the form
such that for all vectors for all for all for all for all not equal .
The ordering for the factors
and for the placement of zeroes was chosen to agree with the output ofsmith_form().See the example for a pictorial description of such a basis.
EXAMPLES:
sage: E = matrix(ZZ, 5, 5, [0, 14, 0, -8, -2, -14, 0, -3, -11, 4, 0, 3, 0, 0, 0, 8, 11, 0, 0, 8, 2, -4, 0, -8, 0]); E [ 0 14 0 -8 -2] [-14 0 -3 -11 4] [ 0 3 0 0 0] [ 8 11 0 0 8] [ 2 -4 0 -8 0] sage: F, C = E.symplectic_form() sage: F [ 0 0 1 0 0] [ 0 0 0 2 0] [-1 0 0 0 0] [ 0 -2 0 0 0] [ 0 0 0 0 0] sage: F == C * E * C.transpose() True sage: E.smith_form()[0] [1 0 0 0 0] [0 1 0 0 0] [0 0 2 0 0] [0 0 0 2 0] [0 0 0 0 0]
- transpose()[source]¶
Return the transpose of self, without changing
self.EXAMPLES:
We create a matrix, compute its transpose, and note that the original matrix is not changed.
sage: A = matrix(ZZ, 2, 3, range(6)) sage: type(A) <class 'sage.matrix.matrix_integer_dense.Matrix_integer_dense'> sage: B = A.transpose() sage: print(B) [0 3] [1 4] [2 5] sage: print(A) [0 1 2] [3 4 5]
.Tis a convenient shortcut for the transpose:sage: A.T [0 3] [1 4] [2 5]
sage: A.subdivide(None, 1); A [0|1 2] [3|4 5] sage: A.transpose() [0 3] [---] [1 4] [2 5]