Sparse matrices over \(\ZZ/n\ZZ\) for \(n\) small#

This is a compiled implementation of sparse matrices over \(\ZZ/n\ZZ\) for \(n\) small.

Todo

move vectors into a Cython vector class - add _add_ and _mul_ methods.

EXAMPLES:

sage: a = matrix(Integers(37),3,3,range(9),sparse=True); a
[0 1 2]
[3 4 5]
[6 7 8]
sage: type(a)
<class 'sage.matrix.matrix_modn_sparse.Matrix_modn_sparse'>
sage: parent(a)
Full MatrixSpace of 3 by 3 sparse matrices over Ring of integers modulo 37
sage: a^2
[15 18 21]
[ 5 17 29]
[32 16  0]
sage: a+a
[ 0  2  4]
[ 6  8 10]
[12 14 16]
sage: b = a.new_matrix(2,3,range(6)); b
[0 1 2]
[3 4 5]
sage: a*b
Traceback (most recent call last):
...
TypeError: unsupported operand parent(s) for *: 'Full MatrixSpace of 3 by 3 sparse matrices over Ring of integers modulo 37' and 'Full MatrixSpace of 2 by 3 sparse matrices over Ring of integers modulo 37'
sage: b*a
[15 18 21]
[ 5 17 29]
>>> from sage.all import *
>>> a = matrix(Integers(Integer(37)),Integer(3),Integer(3),range(Integer(9)),sparse=True); a
[0 1 2]
[3 4 5]
[6 7 8]
>>> type(a)
<class 'sage.matrix.matrix_modn_sparse.Matrix_modn_sparse'>
>>> parent(a)
Full MatrixSpace of 3 by 3 sparse matrices over Ring of integers modulo 37
>>> a**Integer(2)
[15 18 21]
[ 5 17 29]
[32 16  0]
>>> a+a
[ 0  2  4]
[ 6  8 10]
[12 14 16]
>>> b = a.new_matrix(Integer(2),Integer(3),range(Integer(6))); b
[0 1 2]
[3 4 5]
>>> a*b
Traceback (most recent call last):
...
TypeError: unsupported operand parent(s) for *: 'Full MatrixSpace of 3 by 3 sparse matrices over Ring of integers modulo 37' and 'Full MatrixSpace of 2 by 3 sparse matrices over Ring of integers modulo 37'
>>> b*a
[15 18 21]
[ 5 17 29]
sage: TestSuite(a).run()
sage: TestSuite(b).run()
>>> from sage.all import *
>>> TestSuite(a).run()
>>> TestSuite(b).run()
sage: a.echelonize(); a
[ 1  0 36]
[ 0  1  2]
[ 0  0  0]
sage: b.echelonize(); b
[ 1  0 36]
[ 0  1  2]
sage: a.pivots()
(0, 1)
sage: b.pivots()
(0, 1)
sage: a.rank()
2
sage: b.rank()
2
sage: a[2,2] = 5
sage: a.rank()
3
>>> from sage.all import *
>>> a.echelonize(); a
[ 1  0 36]
[ 0  1  2]
[ 0  0  0]
>>> b.echelonize(); b
[ 1  0 36]
[ 0  1  2]
>>> a.pivots()
(0, 1)
>>> b.pivots()
(0, 1)
>>> a.rank()
2
>>> b.rank()
2
>>> a[Integer(2),Integer(2)] = Integer(5)
>>> a.rank()
3
class sage.matrix.matrix_modn_sparse.Matrix_modn_sparse[source]#

Bases: Matrix_sparse

Create a sparse matrix over the integers modulo n.

INPUT:

  • parent – a matrix space over the integers modulo n

  • entries – see matrix()

  • copy – ignored (for backwards compatibility)

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

density()[source]#

Return the density of self, i.e., the ratio of the number of nonzero entries of self to the total size of self.

EXAMPLES:

sage: A = matrix(QQ,3,3,[0,1,2,3,0,0,6,7,8],sparse=True)
sage: A.density()
2/3
>>> from sage.all import *
>>> A = matrix(QQ,Integer(3),Integer(3),[Integer(0),Integer(1),Integer(2),Integer(3),Integer(0),Integer(0),Integer(6),Integer(7),Integer(8)],sparse=True)
>>> A.density()
2/3

Notice that the density parameter does not ensure the density of a matrix; it is only an upper bound.

sage: A = random_matrix(GF(127), 200, 200, density=0.3, sparse=True)
sage: density_sum = float(A.density())
sage: total = 1
sage: expected_density = 1.0 - (199/200)^60
sage: expected_density
0.2597...
sage: while abs(density_sum/total - expected_density) > 0.001:
....:     A = random_matrix(GF(127), 200, 200, density=0.3, sparse=True)
....:     density_sum += float(A.density())
....:     total += 1
>>> from sage.all import *
>>> A = random_matrix(GF(Integer(127)), Integer(200), Integer(200), density=RealNumber('0.3'), sparse=True)
>>> density_sum = float(A.density())
>>> total = Integer(1)
>>> expected_density = RealNumber('1.0') - (Integer(199)/Integer(200))**Integer(60)
>>> expected_density
0.2597...
>>> while abs(density_sum/total - expected_density) > RealNumber('0.001'):
...     A = random_matrix(GF(Integer(127)), Integer(200), Integer(200), density=RealNumber('0.3'), sparse=True)
...     density_sum += float(A.density())
...     total += Integer(1)
determinant(algorithm=None)[source]#

Return the determinant of this matrix.

INPUT:

  • algorithm – either "linbox" (default) or "generic".

EXAMPLES:

sage: A = matrix(GF(3), 4, range(16), sparse=True)
sage: B = identity_matrix(GF(3), 4, sparse=True)
sage: (A + B).det()
2
sage: (A + B).det(algorithm="linbox")
2
sage: (A + B).det(algorithm="generic")
2
sage: (A + B).det(algorithm="hey")
Traceback (most recent call last):
...
ValueError: no algorithm 'hey'

sage: matrix(GF(11), 1, 2, sparse=True).det()
Traceback (most recent call last):
...
ValueError: self must be a square matrix
>>> from sage.all import *
>>> A = matrix(GF(Integer(3)), Integer(4), range(Integer(16)), sparse=True)
>>> B = identity_matrix(GF(Integer(3)), Integer(4), sparse=True)
>>> (A + B).det()
2
>>> (A + B).det(algorithm="linbox")
2
>>> (A + B).det(algorithm="generic")
2
>>> (A + B).det(algorithm="hey")
Traceback (most recent call last):
...
ValueError: no algorithm 'hey'

>>> matrix(GF(Integer(11)), Integer(1), Integer(2), sparse=True).det()
Traceback (most recent call last):
...
ValueError: self must be a square matrix
matrix_from_columns(cols)[source]#

Return the matrix constructed from self using columns with indices in the columns list.

EXAMPLES:

sage: M = MatrixSpace(GF(127),3,3,sparse=True)
sage: A = M(range(9)); A
[0 1 2]
[3 4 5]
[6 7 8]
sage: A.matrix_from_columns([2,1])
[2 1]
[5 4]
[8 7]
>>> from sage.all import *
>>> M = MatrixSpace(GF(Integer(127)),Integer(3),Integer(3),sparse=True)
>>> A = M(range(Integer(9))); A
[0 1 2]
[3 4 5]
[6 7 8]
>>> A.matrix_from_columns([Integer(2),Integer(1)])
[2 1]
[5 4]
[8 7]
matrix_from_rows(rows)[source]#

Return the matrix constructed from self using rows with indices in the rows list.

INPUT:

  • rows – list or tuple of row indices

EXAMPLES:

sage: M = MatrixSpace(GF(127),3,3,sparse=True)
sage: A = M(range(9)); A
[0 1 2]
[3 4 5]
[6 7 8]
sage: A.matrix_from_rows([2,1])
[6 7 8]
[3 4 5]
>>> from sage.all import *
>>> M = MatrixSpace(GF(Integer(127)),Integer(3),Integer(3),sparse=True)
>>> A = M(range(Integer(9))); A
[0 1 2]
[3 4 5]
[6 7 8]
>>> A.matrix_from_rows([Integer(2),Integer(1)])
[6 7 8]
[3 4 5]
p[source]#
rank(algorithm=None)[source]#

Return the rank of this matrix.

INPUT:

  • algorithm – either "linbox" (only available for matrices over prime fields) or "generic"

EXAMPLES:

sage: A = matrix(GF(127), 2, 2, sparse=True)
sage: A[0,0] = 34
sage: A[0,1] = 102
sage: A[1,0] = 55
sage: A[1,1] = 74
sage: A.rank()
2

sage: A._clear_cache()
sage: A.rank(algorithm="generic")
2
sage: A._clear_cache()
sage: A.rank(algorithm="hey")
Traceback (most recent call last):
...
ValueError: no algorithm 'hey'
>>> from sage.all import *
>>> A = matrix(GF(Integer(127)), Integer(2), Integer(2), sparse=True)
>>> A[Integer(0),Integer(0)] = Integer(34)
>>> A[Integer(0),Integer(1)] = Integer(102)
>>> A[Integer(1),Integer(0)] = Integer(55)
>>> A[Integer(1),Integer(1)] = Integer(74)
>>> A.rank()
2

>>> A._clear_cache()
>>> A.rank(algorithm="generic")
2
>>> A._clear_cache()
>>> A.rank(algorithm="hey")
Traceback (most recent call last):
...
ValueError: no algorithm 'hey'

REFERENCES:

Note

For very sparse matrices Gaussian elimination is faster because it barely has anything to do. If the fill in needs to be considered, ‘Symbolic Reordering’ is usually much faster.

swap_rows(r1, r2)[source]#
transpose()[source]#

Return the transpose of self.

EXAMPLES:

sage: A = matrix(GF(127),3,3,[0,1,0,2,0,0,3,0,0],sparse=True)
sage: A
[0 1 0]
[2 0 0]
[3 0 0]
sage: A.transpose()
[0 2 3]
[1 0 0]
[0 0 0]
>>> from sage.all import *
>>> A = matrix(GF(Integer(127)),Integer(3),Integer(3),[Integer(0),Integer(1),Integer(0),Integer(2),Integer(0),Integer(0),Integer(3),Integer(0),Integer(0)],sparse=True)
>>> A
[0 1 0]
[2 0 0]
[3 0 0]
>>> A.transpose()
[0 2 3]
[1 0 0]
[0 0 0]

.T is a convenient shortcut for the transpose:

sage: A.T
[0 2 3]
[1 0 0]
[0 0 0]
>>> from sage.all import *
>>> A.T
[0 2 3]
[1 0 0]
[0 0 0]