Dense matrices over F2e for 2e16 using the M4RIE library

The M4RIE library offers two matrix representations:

  1. mzed_t

m x n matrices over F2e are internally represented roughly as m x (en) matrices over F2. Several elements are packed into words such that each element is filled with zeroes until the next power of two. Thus, for example, elements of F23 are represented as [0xxx|0xxx|0xxx|0xxx|...]. This representation is wrapped as Matrix_gf2e_dense in Sage.

Multiplication and elimination both use “Newton-John” tables. These tables are simply all possible multiples of a given row in a matrix such that a scale+add operation is reduced to a table lookup + add. On top of Newton-John multiplication M4RIE implements asymptotically fast Strassen-Winograd multiplication. Elimination uses simple Gaussian elimination which requires O(n3) additions but only O(n22e) multiplications.

  1. mzd_slice_t

m x n matrices over F2e are internally represented as slices of m x n matrices over F2. This representation allows for very fast matrix times matrix products using Karatsuba’s polynomial multiplication for polynomials over matrices. However, it is not feature complete yet and hence not wrapped in Sage for now.

See http://m4ri.sagemath.org for more details on the M4RIE library.

EXAMPLES:

sage: K.<a> = GF(2^8)
sage: A = random_matrix(K, 3,4)
sage: A
[          a^6 + a^5 + a^4 + a^2               a^6 + a^3 + a + 1         a^5 + a^3 + a^2 + a + 1               a^7 + a^6 + a + 1]
[                a^7 + a^6 + a^3             a^7 + a^6 + a^5 + 1         a^5 + a^4 + a^3 + a + 1 a^6 + a^5 + a^4 + a^3 + a^2 + 1]
[              a^6 + a^5 + a + 1                   a^7 + a^3 + 1               a^7 + a^3 + a + 1   a^7 + a^6 + a^3 + a^2 + a + 1]

sage: A.echelon_form()
[                        1                         0                         0     a^6 + a^5 + a^4 + a^2]
[                        0                         1                         0   a^7 + a^5 + a^3 + a + 1]
[                        0                         0                         1 a^6 + a^4 + a^3 + a^2 + 1]

AUTHOR:

Todo

Wrap mzd_slice_t.

REFERENCES:

class sage.matrix.matrix_gf2e_dense.M4RIE_finite_field

Bases: object

A thin wrapper around the M4RIE finite field class such that we can put it in a hash table. This class is not meant for public consumption.

class sage.matrix.matrix_gf2e_dense.Matrix_gf2e_dense

Bases: sage.matrix.matrix_dense.Matrix_dense

Create new matrix over GF(2e) for 2e16.

INPUT:

  • parent – a matrix space over GF(2^e)

  • entries – see matrix()

  • copy – ignored (for backwards compatibility)

  • coerce – if False, assume without checking that the entries lie in the base ring

EXAMPLES:

sage: K.<a> = GF(2^4)
sage: l = [K.random_element() for _ in range(3*4)]; l
[a^2 + 1, a^3 + 1, 0, 0, a, a^3 + a + 1, a + 1, a + 1, a^2, a^3 + a + 1, a^3 + a, a^3 + a]

sage: A = Matrix(K, 3, 4, l); A
[    a^2 + 1     a^3 + 1           0           0]
[          a a^3 + a + 1       a + 1       a + 1]
[        a^2 a^3 + a + 1     a^3 + a     a^3 + a]

sage: A.list()
[a^2 + 1, a^3 + 1, 0, 0, a, a^3 + a + 1, a + 1, a + 1, a^2, a^3 + a + 1, a^3 + a, a^3 + a]

sage: l[0], A[0,0]
(a^2 + 1, a^2 + 1)

sage: A = Matrix(K, 3, 3, a); A
[a 0 0]
[0 a 0]
[0 0 a]
augment(right)

Augments self with right.

INPUT:

  • right - a matrix

EXAMPLES:

sage: K.<a> = GF(2^4)
sage: MS = MatrixSpace(K,3,3)
sage: A = random_matrix(K,3,3)
sage: B = A.augment(MS(1)); B
[              a^2       a^3 + a + 1 a^3 + a^2 + a + 1                 1                 0                 0]
[            a + 1               a^3                 1                 0                 1                 0]
[      a^3 + a + 1     a^3 + a^2 + 1             a + 1                 0                 0                 1]

sage: B.echelonize(); B
[                1                 0                 0           a^2 + a           a^3 + 1           a^3 + a]
[                0                 1                 0     a^3 + a^2 + a a^3 + a^2 + a + 1           a^2 + a]
[                0                 0                 1             a + 1               a^3               a^3]

sage: C = B.matrix_from_columns([3,4,5]); C
[          a^2 + a           a^3 + 1           a^3 + a]
[    a^3 + a^2 + a a^3 + a^2 + a + 1           a^2 + a]
[            a + 1               a^3               a^3]

sage: C == ~A
True

sage: C*A == MS(1)
True
cling(*C)

Pack the matrices over F2 into this matrix over F2e.

Elements in F2e can be represented as ciai where a is a root the minimal polynomial. If this matrix is A then this function writes ciai to the entry A[x,y] where ci is the entry Ci[x,y].

INPUT:

  • C - a list of matrices over GF(2)

EXAMPLES:

sage: K.<a> = GF(2^2)
sage: A = matrix(K, 5, 5)
sage: A0 = random_matrix(GF(2), 5, 5); A0
[0 1 0 1 1]
[0 1 1 1 0]
[0 0 0 1 0]
[0 1 1 0 0]
[0 0 0 1 1]

sage: A1 = random_matrix(GF(2), 5, 5); A1
[0 0 1 1 1]
[1 1 1 1 0]
[0 0 0 1 1]
[1 0 0 0 1]
[1 0 0 1 1]

sage: A.cling(A1, A0); A
[    0     a     1 a + 1 a + 1]
[    1 a + 1 a + 1 a + 1     0]
[    0     0     0 a + 1     1]
[    1     a     a     0     1]
[    1     0     0 a + 1 a + 1]

sage: A0[0,3]*a + A1[0,3], A[0,3]
(a + 1, a + 1)

Slicing and clinging are inverse operations:

sage: B1, B0 = A.slice()
sage: B0 == A0 and B1 == A1
True
echelonize(algorithm='heuristic', reduced=True, **kwds)

Compute the row echelon form of self in place.

INPUT:

  • algorithm - one of the following - heuristic - let M4RIE decide (default) - newton_john - use newton_john table based algorithm - ple - use PLE decomposition - naive - use naive cubic Gaussian elimination (M4RIE implementation) - builtin - use naive cubic Gaussian elimination (Sage implementation)

  • reduced - if True return reduced echelon form. No guarantee is given that the matrix is not reduced if False (default: True)

EXAMPLES:

sage: K.<a> = GF(2^4)
sage: m,n  = 3, 5
sage: A = random_matrix(K, 3, 5); A
[              a^2       a^3 + a + 1 a^3 + a^2 + a + 1             a + 1               a^3]
[                1       a^3 + a + 1     a^3 + a^2 + 1             a + 1           a^3 + 1]
[      a^3 + a + 1 a^3 + a^2 + a + 1           a^2 + a           a^2 + 1           a^2 + a]

sage: A.echelonize(); A
[            1             0             0         a + 1       a^2 + 1]
[            0             1             0           a^2         a + 1]
[            0             0             1 a^3 + a^2 + a           a^3]

sage: K.<a> = GF(2^3)
sage: m,n  = 3, 5
sage: MS = MatrixSpace(K,m,n)
sage: A = random_matrix(K, 3, 5)

sage: copy(A).echelon_form('newton_john')
[          1           0       a + 1           0     a^2 + 1]
[          0           1 a^2 + a + 1           0           a]
[          0           0           0           1 a^2 + a + 1]

sage: copy(A).echelon_form('naive')
[          1           0       a + 1           0     a^2 + 1]
[          0           1 a^2 + a + 1           0           a]
[          0           0           0           1 a^2 + a + 1]

sage: copy(A).echelon_form('builtin')
[          1           0       a + 1           0     a^2 + 1]
[          0           1 a^2 + a + 1           0           a]
[          0           0           0           1 a^2 + a + 1]
randomize(density=1, nonzero=False, *args, **kwds)

Randomize density proportion of the entries of this matrix, leaving the rest unchanged.

INPUT:

  • density - float; proportion (roughly) to be considered for changes

  • nonzero - Bool (default: False); whether the new entries are forced to be non-zero

OUTPUT:

  • None, the matrix is modified in-place

EXAMPLES:

sage: K.<a> = GF(2^4)
sage: A = Matrix(K,3,3)

sage: A.randomize(); A
[              a^2       a^3 + a + 1 a^3 + a^2 + a + 1]
[            a + 1               a^3                 1]
[      a^3 + a + 1     a^3 + a^2 + 1             a + 1]

sage: K.<a> = GF(2^4)
sage: A = random_matrix(K,1000,1000,density=0.1)
sage: float(A.density())
0.0999...

sage: A = random_matrix(K,1000,1000,density=1.0)
sage: float(A.density())
1.0

sage: A = random_matrix(K,1000,1000,density=0.5)
sage: float(A.density())
0.4996...

Note, that the matrix is updated and not zero-ed out before being randomized:

sage: A = matrix(K,1000,1000)
sage: A.randomize(nonzero=False, density=0.1)
sage: float(A.density())
0.0936...

sage: A.randomize(nonzero=False, density=0.05)
sage: float(A.density())
0.135854
rank()

Return the rank of this matrix (cached).

EXAMPLES:

sage: K.<a> = GF(2^4)
sage: A = random_matrix(K, 1000, 1000)
sage: A.rank()
1000

sage: A = matrix(K, 10, 0)
sage: A.rank()
0
slice()

Unpack this matrix into matrices over F2.

Elements in F2e can be represented as ciai where a is a root the minimal polynomial. This function returns a tuple of matrices C whose entry Ci[x,y] is the coefficient of ci in A[x,y] if this matrix is A.

EXAMPLES:

sage: K.<a> = GF(2^2)
sage: A = random_matrix(K, 5, 5); A
[    0 a + 1 a + 1 a + 1     0]
[    1 a + 1     1 a + 1     1]
[a + 1 a + 1     a     1     a]
[    a     1 a + 1     1     0]
[    a     1 a + 1 a + 1     0]

sage: A1,A0 = A.slice()
sage: A0
[0 1 1 1 0]
[0 1 0 1 0]
[1 1 1 0 1]
[1 0 1 0 0]
[1 0 1 1 0]

sage: A1
[0 1 1 1 0]
[1 1 1 1 1]
[1 1 0 1 0]
[0 1 1 1 0]
[0 1 1 1 0]

sage: A0[2,4]*a + A1[2,4], A[2,4]
(a, a)

sage: K.<a> = GF(2^3)
sage: A = random_matrix(K, 5, 5); A
[      a + 1     a^2 + a           1           a     a^2 + a]
[      a + 1     a^2 + a         a^2         a^2     a^2 + 1]
[a^2 + a + 1 a^2 + a + 1           0 a^2 + a + 1     a^2 + 1]
[    a^2 + a           0 a^2 + a + 1           a           a]
[        a^2       a + 1           a     a^2 + 1 a^2 + a + 1]

sage: A0,A1,A2 = A.slice()
sage: A0
[1 0 1 0 0]
[1 0 0 0 1]
[1 1 0 1 1]
[0 0 1 0 0]
[0 1 0 1 1]

Slicing and clinging are inverse operations:

sage: B = matrix(K, 5, 5)
sage: B.cling(A0,A1,A2)
sage: B == A
True
stack(other)

Stack self on top of other.

INPUT:

  • other - a matrix

EXAMPLES:

sage: K.<a> = GF(2^4)
sage: A = random_matrix(K,2,2); A
[              a^2       a^3 + a + 1]
[a^3 + a^2 + a + 1             a + 1]

sage: B = random_matrix(K,2,2); B
[          a^3             1]
[  a^3 + a + 1 a^3 + a^2 + 1]

sage: A.stack(B)
[              a^2       a^3 + a + 1]
[a^3 + a^2 + a + 1             a + 1]
[              a^3                 1]
[      a^3 + a + 1     a^3 + a^2 + 1]

sage: B.stack(A)
[              a^3                 1]
[      a^3 + a + 1     a^3 + a^2 + 1]
[              a^2       a^3 + a + 1]
[a^3 + a^2 + a + 1             a + 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: K.<a> = GF(2^10)
sage: A = random_matrix(K,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
sage.matrix.matrix_gf2e_dense.unpickle_matrix_gf2e_dense_v0(a, base_ring, nrows, ncols)

EXAMPLES:

sage: K.<a> = GF(2^2)
sage: A = random_matrix(K,10,10)
sage: f, s= A.__reduce__()
sage: from sage.matrix.matrix_gf2e_dense import unpickle_matrix_gf2e_dense_v0
sage: f == unpickle_matrix_gf2e_dense_v0
True
sage: f(*s) == A
True

We can still unpickle pickles from before trac ticket #19240:

sage: old_pickle = b'x\x9c\x85RKo\xd3@\x10\xae\xdd$$\xdb&\xe5U\x1e-\x8f\xc2\xc9\x12RD#$\xce\xa0\xb4\x80\x07\xa2\xca\xc2\x07\x0e\xd5\xe2:\x1b\xdb\x8acg\x1c\xa7J\x85*!\xa4\x90\xe6\x07p\xe0\xc4\x01q\xe5\xc4\x19\xf5\xd0?\xc1\x81\xdf\x80\xb8q\x0b\xb3\x8eMS\xa1\x82V;;\xb3\xdf\xce\xf7\xcd\x8e\xe6\xb5j\xf7,GT;V\x1cy\x83\xf4\xe0\x9d\xb0Y\x13\xbc)\x82\x9e`\xfd\xa0\xeb\xd9m_\xf0\xbf1\xbe{\x97\xa1\xa2\x9d\xc6\xf0\x0f\x82,\x7f\x9d\xa1\xaa\x81\n\xb9m\x9c\xd7\xf4\xf1d2\x81-h\xc0#(\x03\x83\x15\xdas\xc9*\xc3\x13x\x0cu0\xd28\x97\x9e*(0\x9f\xfa\x1b\xd0\xd2\x7fH\x82\xb5\xf4\xa2@TO\xe19\x01I\xac\x136\x991\x9f\xa4\xf9&\xcd\x07i\xbe\xcb\xd4ib\t\xba\xa4\xf6\x02zIT\xd1\x8f2(u\x15\xfd\x9d<\xee@\x05V\xd3\x94E*\xb0\x0e\x0fH\xad\xa8\xbf\x97\xa0\r\x03\xfd\xf0\xb8\x1aU\xff\x92\x90\xe8?\xa5\xd6\x814_\xa5\xf9(\xcd\xafc\xe99\xe2\xd9\xa0\x06\xd4\xf5\xcf\xf2\xf2!\xbc\xd4\xdf\x90#\xc0\x8f\r\xccM\x1b\xdd\x8b\xa3\xbe\x1d\xf7#QmYv\x1cF{\xcc\x11\x81\x88<\x9b\xa71\xcf:\xce0\xaf\x9d\x96\xe3\x87a\xbb\xdf\xe5\x8e\x1f\xeeX>\xc3\x82\xb9\xb0\xe9\x05^,6=\xe17\xf1\xcc\xd0\xc0"u\xb0d\xe6wDl\xdd\x1fa)e\x8a\xbc\xc0\xe9U\xbd \x16\x8e\x88X\xc7j\x0b\x9e\x05\xc8L\xe5\x1e%.\x98\x8a5\xc4\xc5\xd9\xf7\xdd\xd0\xdf\x0b\xc2\x8eg\xf93.wZ\xb5\xc1\x94B\xf8\xa2#\x82\x98a\xf9\xffY\x12\xe3v\x18L\xff\x14Fl\xeb\x0ff\x10\xc4\xb0\xa2\xb9y\xcd-\xba%\xcd\xa5\x8ajT\xd1\x92\xa9\x0c\x86x\xb6a\xe6h\xf8\x02<g\xaa\xaf\xf6\xdd%\x89\xae\x13z\xfe \xc6\x0b\xfb1^4p\x99\x1e6\xc6\xd4\xebK\xdbx\xf9\xc4\x8f[Iw\xf8\x89\xef\xcbQf\xcfh\xe3\x95\x8c\xebj&\xb9\xe2.\x8f\x0c\\ui\x89\xf1x\xf4\xd6\xc0kf\xc1\xf1v\xad(\xc4\xeb\x89~\xfa\xf0\x06\xa8\xa4\x7f\x93\xf4\xd7\x0c\xbcE#\xad\x92\xfc\xed\xeao\xefX\\\x03'
sage: loads(old_pickle)
[    0     a]
[a + 1     1]