Dense matrices over GF(2) using the M4RI library#
AUTHOR: Martin Albrecht <malb@informatik.uni-bremen.de>
EXAMPLES:
sage: a = matrix(GF(2),3,range(9),sparse=False); a
[0 1 0]
[1 0 1]
[0 1 0]
sage: a.rank()
2
sage: type(a)
<class 'sage.matrix.matrix_mod2_dense.Matrix_mod2_dense'>
sage: a[0,0] = 1
sage: a.rank()
3
sage: parent(a)
Full MatrixSpace of 3 by 3 dense matrices over Finite Field of size 2
sage: a^2
[0 1 1]
[1 0 0]
[1 0 1]
sage: a+a
[0 0 0]
[0 0 0]
[0 0 0]
sage: b = a.new_matrix(2,3,range(6)); b
[0 1 0]
[1 0 1]
sage: a*b
Traceback (most recent call last):
...
TypeError: unsupported operand parent(s) for *: 'Full MatrixSpace of 3 by 3 dense matrices over Finite Field of size 2' and 'Full MatrixSpace of 2 by 3 dense matrices over Finite Field of size 2'
sage: b*a
[1 0 1]
[1 0 0]
sage: TestSuite(a).run()
sage: TestSuite(b).run()
sage: a.echelonize(); a
[1 0 0]
[0 1 0]
[0 0 1]
sage: b.echelonize(); b
[1 0 1]
[0 1 0]
Todo
make LinBox frontend and use it
charpoly ?
minpoly ?
make Matrix_modn_frontend and use it (?)
- class sage.matrix.matrix_mod2_dense.Matrix_mod2_dense#
Bases:
sage.matrix.matrix_dense.Matrix_dense
Dense matrix over GF(2).
- augment(right, subdivide=False)#
Augments
self
withright
.EXAMPLES:
sage: MS = MatrixSpace(GF(2),3,3) sage: A = MS([0, 1, 0, 1, 1, 0, 1, 1, 1]); A [0 1 0] [1 1 0] [1 1 1] sage: B = A.augment(MS(1)); B [0 1 0 1 0 0] [1 1 0 0 1 0] [1 1 1 0 0 1] sage: B.echelonize(); B [1 0 0 1 1 0] [0 1 0 1 0 0] [0 0 1 0 1 1] sage: C = B.matrix_from_columns([3,4,5]); C [1 1 0] [1 0 0] [0 1 1] sage: C == ~A True sage: C*A == MS(1) True
A vector may be augmented to a matrix.
sage: A = matrix(GF(2), 3, 4, range(12)) sage: v = vector(GF(2), 3, range(3)) sage: A.augment(v) [0 1 0 1 0] [0 1 0 1 1] [0 1 0 1 0]
The
subdivide
option will add a natural subdivision betweenself
andright
. For more details about how subdivisions are managed when augmenting, seesage.matrix.matrix1.Matrix.augment()
.sage: A = matrix(GF(2), 3, 5, range(15)) sage: B = matrix(GF(2), 3, 3, range(9)) sage: A.augment(B, subdivide=True) [0 1 0 1 0|0 1 0] [1 0 1 0 1|1 0 1] [0 1 0 1 0|0 1 0]
- density(approx=False)#
Return the density of this matrix.
By density we understand the ratio of the number of nonzero positions and the self.nrows() * self.ncols(), i.e. the number of possible nonzero positions.
INPUT:
approx – return floating point approximation (default: False)
EXAMPLES:
sage: A = random_matrix(GF(2), 1000, 1000) sage: d = A.density() sage: float(d) == A.density(approx=True) True sage: len(A.nonzero_positions())/1000^2 == d True sage: total = 1.0 sage: density_sum = A.density() sage: while abs(density_sum/total - 0.5) > 0.001: ....: A = random_matrix(GF(2), 1000, 1000) ....: total += 1 ....: density_sum += A.density()
- determinant()#
Return the determinant of this matrix over GF(2).
EXAMPLES:
sage: matrix(GF(2),2,[1,1,0,1]).determinant() 1 sage: matrix(GF(2),2,[1,1,1,1]).determinant() 0
- echelonize(algorithm='heuristic', cutoff=0, reduced=True, **kwds)#
Puts self in (reduced) row echelon form.
INPUT:
self – a mutable matrix
algorithm
‘heuristic’ – uses M4RI and PLUQ (default)
‘m4ri’ – uses M4RI
‘pluq’ – uses PLUQ factorization
‘classical’ – uses classical Gaussian elimination
k – the parameter ‘k’ of the M4RI algorithm. It MUST be between 1 and 16 (inclusive). If it is not specified it will be calculated as 3/4 * log_2( min(nrows, ncols) ) as suggested in the M4RI paper.
reduced – return reduced row echelon form (default:True)
EXAMPLES:
sage: A = random_matrix(GF(2), 10, 10) sage: B = A.__copy__(); B.echelonize() # fastest sage: C = A.__copy__(); C.echelonize(k=2) # force k sage: E = A.__copy__(); E.echelonize(algorithm='classical') # force Gaussian elimination sage: B == C == E True
ALGORITHM:
Uses M4RI library
REFERENCES:
- randomize(density=1, nonzero=False)#
Randomize
density
proportion of the entries of this matrix, leaving the rest unchanged.INPUT:
density
- float; proportion (roughly) to be considered for changesnonzero
- Bool (default:False
); whether the new entries are forced to be non-zero
OUTPUT:
None, the matrix is modified in-space
EXAMPLES:
sage: A = matrix(GF(2), 5, 5, 0) sage: A.randomize(0.5) sage: A.density() < 0.5 True sage: expected = 0.5 sage: A = matrix(GF(2), 5, 5, 0) sage: A.randomize() sage: density_sum = float(A.density()) sage: total = 1 sage: while abs(density_sum/total - expected) > 0.001: ....: A = matrix(GF(2), 5, 5, 0) ....: A.randomize() ....: density_sum += float(A.density()) ....: total += 1
- rank(algorithm='ple')#
Return the rank of this matrix.
On average ‘ple’ should be faster than ‘m4ri’ and hence it is the default choice. However, for small - i.e. quite few thousand rows & columns - and sparse matrices ‘m4ri’ might be a better choice.
INPUT:
algorithm
- either “ple” or “m4ri”
EXAMPLES:
sage: while random_matrix(GF(2), 1000, 1000).rank() != 999: ....: pass sage: A = matrix(GF(2),10, 0) sage: A.rank() 0
- row(i, from_list=False)#
Return the
i
’th row of this matrix as a vector.This row is a dense vector if and only if the matrix is a dense matrix.
INPUT:
i
- integerfrom_list
- bool (default:False
); ifTrue
, returns thei
’th element ofself.rows()
(seerows()
), which may be faster, but requires building a list of all rows the first time it is called after an entry of the matrix is changed.
EXAMPLES:
sage: l = [GF(2).random_element() for _ in range(100)] sage: A = matrix(GF(2), 10, 10 , l) sage: list(A.row(0)) == l[:10] True sage: list(A.row(-1)) == l[-10:] True sage: list(A.row(2, from_list=True)) == l[20:30] True sage: A = Matrix(GF(2),1,0) sage: A.row(0) ()
- str(rep_mapping=None, zero=None, plus_one=None, minus_one=None, unicode=False, shape=None, character_art=False)#
Return a nice string representation of the matrix.
INPUT:
rep_mapping
– a dictionary or callable used to override the usual representation of elements. For a dictionary, keys should be elements of the base ring and values the desired string representation.zero
– string (default:None
); if notNone
use the value ofzero
as the representation of the zero element.plus_one
– string (default:None
); if notNone
use the value ofplus_one
as the representation of the one element.minus_one
– Ignored. Only for compatibility with generic matrices.unicode
– boolean (default:False
). Whether to use Unicode symbols instead of ASCII symbols for brackets and subdivision lines.shape
– one of"square"
or"round"
(default:None
). Switches between round and square brackets. The default depends on the setting of theunicode
keyword argument. For Unicode symbols, the default is round brackets in accordance with the TeX rendering, while the ASCII rendering defaults to square brackets.character_art
– boolean (default:False
); ifTrue
, the result will be of typeAsciiArt
orUnicodeArt
which support line breaking of wide matrices that exceed the window width
EXAMPLES:
sage: B = matrix(GF(2), 3, 3, [0, 1, 0, 0, 1, 1, 0, 0, 0]) sage: B # indirect doctest [0 1 0] [0 1 1] [0 0 0] sage: block_matrix([[B, 1], [0, B]]) [0 1 0|1 0 0] [0 1 1|0 1 0] [0 0 0|0 0 1] [-----+-----] [0 0 0|0 1 0] [0 0 0|0 1 1] [0 0 0|0 0 0] sage: B.str(zero='.') '[. 1 .]\n[. 1 1]\n[. . .]' sage: M = matrix.identity(GF(2), 3) sage: M.subdivide(None, 2) sage: print(M.str(unicode=True, shape='square')) ⎡1 0│0⎤ ⎢0 1│0⎥ ⎣0 0│1⎦ sage: print(unicode_art(M)) # indirect doctest ⎛1 0│0⎞ ⎜0 1│0⎟ ⎝0 0│1⎠
- submatrix(row=0, col=0, nrows=-1, ncols=-1)#
Return submatrix from the index row, col (inclusive) with dimension nrows x ncols.
INPUT:
row – index of start row
col – index of start column
nrows – number of rows of submatrix
ncols – number of columns of submatrix
EXAMPLES:
sage: A = random_matrix(GF(2),200,200) sage: A[0:2,0:2] == A.submatrix(0,0,2,2) True sage: A[0:100,0:100] == A.submatrix(0,0,100,100) True sage: A == A.submatrix(0,0,200,200) True sage: A[1:3,1:3] == A.submatrix(1,1,2,2) True sage: A[1:100,1:100] == A.submatrix(1,1,99,99) True sage: A[1:200,1:200] == A.submatrix(1,1,199,199) True
TESTS for handling of default arguments (trac ticket #18761):
sage: A.submatrix(17,15) == A.submatrix(17,15,183,185) True sage: A.submatrix(row=100,col=37,nrows=1,ncols=3) == A.submatrix(100,37,1,3) True
- transpose()#
Returns transpose of self and leaves self untouched.
EXAMPLES:
sage: A = Matrix(GF(2),3,5,[1, 0, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 1, 0]) sage: A [1 0 1 0 0] [0 1 1 0 0] [1 1 0 1 0] sage: B = A.transpose(); B [1 0 1] [0 1 1] [1 1 0] [0 0 1] [0 0 0] sage: B.transpose() == A True
.T
is a convenient shortcut for the transpose:sage: A.T [1 0 1] [0 1 1] [1 1 0] [0 0 1] [0 0 0]
- sage.matrix.matrix_mod2_dense.from_png(filename)#
Returns a dense matrix over GF(2) from a 1-bit PNG image read from
filename
. No attempt is made to verify that the filename string actually points to a PNG image.INPUT:
filename – a string
EXAMPLES:
sage: from sage.matrix.matrix_mod2_dense import from_png, to_png sage: A = random_matrix(GF(2),10,10) sage: fn = tmp_filename() sage: to_png(A, fn) sage: B = from_png(fn) sage: A == B True
- sage.matrix.matrix_mod2_dense.parity(a)#
Return the parity of the number of bits in a.
EXAMPLES:
sage: from sage.matrix.matrix_mod2_dense import parity sage: parity(1) 1 sage: parity(3) 0 sage: parity(0x10000101011) 1
- sage.matrix.matrix_mod2_dense.ple(A, algorithm='standard', param=0)#
Return PLE factorization of A.
INPUT:
A – matrix
algorithm
‘standard’ asymptotically fast (default)
‘russian’ M4RI inspired
‘naive’ naive cubic
param – either k for ‘mmpf’ is chosen or matrix multiplication cutoff for ‘standard’ (default: 0)
EXAMPLES:
sage: from sage.matrix.matrix_mod2_dense import ple sage: A = matrix(GF(2), 4, 4, [0, 1, 0, 1, 0, 1, 1, 1, 0, 0, 0, 1, 0, 1, 1, 0]) sage: A [0 1 0 1] [0 1 1 1] [0 0 0 1] [0 1 1 0] sage: LU, P, Q = ple(A) sage: LU [1 0 0 1] [1 1 0 0] [0 0 1 0] [1 1 1 0] sage: P [0, 1, 2, 3] sage: Q [1, 2, 3, 3] sage: A = random_matrix(GF(2),1000,1000) sage: ple(A) == ple(A,'russian') == ple(A,'naive') True
- sage.matrix.matrix_mod2_dense.pluq(A, algorithm='standard', param=0)#
Return PLUQ factorization of A.
INPUT:
A – matrix
algorithm
‘standard’ asymptotically fast (default)
‘mmpf’ M4RI inspired
‘naive’ naive cubic
param – either k for ‘mmpf’ is chosen or matrix multiplication cutoff for ‘standard’ (default: 0)
EXAMPLES:
sage: from sage.matrix.matrix_mod2_dense import pluq sage: A = matrix(GF(2), 4, 4, [0, 1, 0, 1, 0, 1, 1, 1, 0, 0, 0, 1, 0, 1, 1, 0]) sage: A [0 1 0 1] [0 1 1 1] [0 0 0 1] [0 1 1 0] sage: LU, P, Q = pluq(A) sage: LU [1 0 1 0] [1 1 0 0] [0 0 1 0] [1 1 1 0] sage: P [0, 1, 2, 3] sage: Q [1, 2, 3, 3]
- sage.matrix.matrix_mod2_dense.to_png(A, filename)#
Saves the matrix
A
to filename as a 1-bit PNG image.INPUT:
A
- a matrix over GF(2)filename
- a string for a file in a writable position
EXAMPLES:
sage: from sage.matrix.matrix_mod2_dense import from_png, to_png sage: A = random_matrix(GF(2),10,10) sage: fn = tmp_filename() sage: to_png(A, fn) sage: B = from_png(fn) sage: A == B True
- sage.matrix.matrix_mod2_dense.unpickle_matrix_mod2_dense_v2(r, c, data, size, immutable=False)#
Deserialize a matrix encoded in the string
s
.INPUT:
r
– number of rows of matrixc
– number of columns of matrixs
– a stringsize
– length of the strings
immutable
– (default:False
) whether the matrix is immutable or not
EXAMPLES:
sage: A = random_matrix(GF(2),100,101) sage: _, (r,c,s,s2,i) = A.__reduce__() sage: from sage.matrix.matrix_mod2_dense import unpickle_matrix_mod2_dense_v2 sage: unpickle_matrix_mod2_dense_v2(r,c,s,s2,i) == A True sage: loads(dumps(A)) == A True