Base class for sparse matrices#

class sage.matrix.matrix_sparse.Matrix_sparse[source]#

Bases: Matrix

antitranspose()[source]#

Return the antitranspose of self, without changing self.

This is the mirror image along the other diagonal.

EXAMPLES:

sage: M = MatrixSpace(QQ, 2, sparse=True)
sage: A = M([1,2,3,4]); A
[1 2]
[3 4]
sage: A.antitranspose()
[4 2]
[3 1]
>>> from sage.all import *
>>> M = MatrixSpace(QQ, Integer(2), sparse=True)
>>> A = M([Integer(1),Integer(2),Integer(3),Integer(4)]); A
[1 2]
[3 4]
>>> A.antitranspose()
[4 2]
[3 1]

See also

transpose()

apply_map(phi, R=None, sparse=True)[source]#

Apply the given map phi (an arbitrary Python function or callable object) to this matrix.

If R is not given, automatically determine the base ring of the resulting matrix.

INPUT:

  • phi – arbitrary Python function or callable object

  • R – (optional) ring

  • sparse – (default: True) whether to return a sparse or a dense matrix

OUTPUT: a matrix over R

EXAMPLES:

sage: m = matrix(ZZ, 10000, {(1,2): 17}, sparse=True)

sage: # needs sage.rings.finite_rings
sage: k.<a> = GF(9)
sage: f = lambda x: k(x)
sage: n = m.apply_map(f)
sage: n.parent()
Full MatrixSpace of 10000 by 10000 sparse matrices
 over Finite Field in a of size 3^2
sage: n[1, 2]
2
>>> from sage.all import *
>>> m = matrix(ZZ, Integer(10000), {(Integer(1),Integer(2)): Integer(17)}, sparse=True)

>>> # needs sage.rings.finite_rings
>>> k = GF(Integer(9), names=('a',)); (a,) = k._first_ngens(1)
>>> f = lambda x: k(x)
>>> n = m.apply_map(f)
>>> n.parent()
Full MatrixSpace of 10000 by 10000 sparse matrices
 over Finite Field in a of size 3^2
>>> n[Integer(1), Integer(2)]
2

An example where the codomain is explicitly specified.

sage: n = m.apply_map(lambda x: x%3, GF(3))
sage: n.parent()
Full MatrixSpace of 10000 by 10000 sparse matrices
 over Finite Field of size 3
sage: n[1, 2]
2
>>> from sage.all import *
>>> n = m.apply_map(lambda x: x%Integer(3), GF(Integer(3)))
>>> n.parent()
Full MatrixSpace of 10000 by 10000 sparse matrices
 over Finite Field of size 3
>>> n[Integer(1), Integer(2)]
2

If we did not specify the codomain, the resulting matrix in the above case ends up over \(\ZZ\) again:

sage: n = m.apply_map(lambda x: x%3)
sage: n.parent()
Full MatrixSpace of 10000 by 10000 sparse matrices over Integer Ring
sage: n[1, 2]
2
>>> from sage.all import *
>>> n = m.apply_map(lambda x: x%Integer(3))
>>> n.parent()
Full MatrixSpace of 10000 by 10000 sparse matrices over Integer Ring
>>> n[Integer(1), Integer(2)]
2

If self is subdivided, the result will be as well:

sage: m = matrix(2, 2, [0, 0, 3, 0])
sage: m.subdivide(None, 1); m
[0|0]
[3|0]
sage: m.apply_map(lambda x: x*x)
[0|0]
[9|0]
>>> from sage.all import *
>>> m = matrix(Integer(2), Integer(2), [Integer(0), Integer(0), Integer(3), Integer(0)])
>>> m.subdivide(None, Integer(1)); m
[0|0]
[3|0]
>>> m.apply_map(lambda x: x*x)
[0|0]
[9|0]

If the map sends zero to a non-zero value, then it may be useful to get the result as a dense matrix.

sage: m = matrix(ZZ, 3, 3, [0] * 7 + [1,2], sparse=True); m
[0 0 0]
[0 0 0]
[0 1 2]
sage: parent(m)
Full MatrixSpace of 3 by 3 sparse matrices over Integer Ring
sage: n = m.apply_map(lambda x: x+polygen(QQ), sparse=False); n
[    x     x     x]
[    x     x     x]
[    x x + 1 x + 2]
sage: parent(n)
Full MatrixSpace of 3 by 3 dense matrices over Univariate Polynomial Ring in x over Rational Field
>>> from sage.all import *
>>> m = matrix(ZZ, Integer(3), Integer(3), [Integer(0)] * Integer(7) + [Integer(1),Integer(2)], sparse=True); m
[0 0 0]
[0 0 0]
[0 1 2]
>>> parent(m)
Full MatrixSpace of 3 by 3 sparse matrices over Integer Ring
>>> n = m.apply_map(lambda x: x+polygen(QQ), sparse=False); n
[    x     x     x]
[    x     x     x]
[    x x + 1 x + 2]
>>> parent(n)
Full MatrixSpace of 3 by 3 dense matrices over Univariate Polynomial Ring in x over Rational Field
apply_morphism(phi)[source]#

Apply the morphism phi to the coefficients of this sparse matrix.

The resulting matrix is over the codomain of phi.

INPUT:

  • phi – a morphism, so phi is callable and

    phi.domain() and phi.codomain() are defined. The codomain must be a ring.

OUTPUT: a matrix over the codomain of phi

EXAMPLES:

sage: m = matrix(ZZ, 3, range(9), sparse=True)
sage: phi = ZZ.hom(GF(5))
sage: m.apply_morphism(phi)
[0 1 2]
[3 4 0]
[1 2 3]
sage: m.apply_morphism(phi).parent()
Full MatrixSpace of 3 by 3 sparse matrices
 over Finite Field of size 5
>>> from sage.all import *
>>> m = matrix(ZZ, Integer(3), range(Integer(9)), sparse=True)
>>> phi = ZZ.hom(GF(Integer(5)))
>>> m.apply_morphism(phi)
[0 1 2]
[3 4 0]
[1 2 3]
>>> m.apply_morphism(phi).parent()
Full MatrixSpace of 3 by 3 sparse matrices
 over Finite Field of size 5
augment(right, subdivide=False)[source]#

Return the augmented matrix of the form:

[self | right].

EXAMPLES:

sage: M = MatrixSpace(QQ, 2, 2, sparse=True)
sage: A = M([1,2, 3,4])
sage: A
[1 2]
[3 4]
sage: N = MatrixSpace(QQ, 2, 1, sparse=True)
sage: B = N([9,8])
sage: B
[9]
[8]
sage: A.augment(B)
[1 2 9]
[3 4 8]
sage: B.augment(A)
[9 1 2]
[8 3 4]
>>> from sage.all import *
>>> M = MatrixSpace(QQ, Integer(2), Integer(2), sparse=True)
>>> A = M([Integer(1),Integer(2), Integer(3),Integer(4)])
>>> A
[1 2]
[3 4]
>>> N = MatrixSpace(QQ, Integer(2), Integer(1), sparse=True)
>>> B = N([Integer(9),Integer(8)])
>>> B
[9]
[8]
>>> A.augment(B)
[1 2 9]
[3 4 8]
>>> B.augment(A)
[9 1 2]
[8 3 4]

A vector may be augmented to a matrix.

sage: A = matrix(QQ, 3, 4, range(12), sparse=True)
sage: v = vector(QQ, 3, range(3), sparse=True)
sage: A.augment(v)
[ 0  1  2  3  0]
[ 4  5  6  7  1]
[ 8  9 10 11  2]
>>> from sage.all import *
>>> A = matrix(QQ, Integer(3), Integer(4), range(Integer(12)), sparse=True)
>>> v = vector(QQ, Integer(3), range(Integer(3)), sparse=True)
>>> A.augment(v)
[ 0  1  2  3  0]
[ 4  5  6  7  1]
[ 8  9 10 11  2]

The subdivide option will add a natural subdivision between self and right. For more details about how subdivisions are managed when augmenting, see sage.matrix.matrix1.Matrix.augment().

sage: A = matrix(QQ, 3, 5, range(15), sparse=True)
sage: B = matrix(QQ, 3, 3, range(9), sparse=True)
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]
>>> from sage.all import *
>>> A = matrix(QQ, Integer(3), Integer(5), range(Integer(15)), sparse=True)
>>> B = matrix(QQ, Integer(3), Integer(3), range(Integer(9)), sparse=True)
>>> 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]
change_ring(ring)[source]#

Return the matrix obtained by coercing the entries of this matrix into the given ring.

Always returns a copy (unless self is immutable, in which case returns self).

EXAMPLES:

sage: x = polygen(ZZ, 'x')
sage: A = matrix(QQ['x,y'], 2, [0,-1,2*x,-2], sparse=True); A
[  0  -1]
[2*x  -2]
sage: A.change_ring(QQ['x,y,z'])
[  0  -1]
[2*x  -2]
>>> from sage.all import *
>>> x = polygen(ZZ, 'x')
>>> A = matrix(QQ['x,y'], Integer(2), [Integer(0),-Integer(1),Integer(2)*x,-Integer(2)], sparse=True); A
[  0  -1]
[2*x  -2]
>>> A.change_ring(QQ['x,y,z'])
[  0  -1]
[2*x  -2]

Subdivisions are preserved when changing rings:

sage: A.subdivide([2],[]); A
[  0  -1]
[2*x  -2]
[-------]
sage: A.change_ring(RR['x,y'])
[                 0  -1.00000000000000]
[2.00000000000000*x  -2.00000000000000]
[-------------------------------------]
>>> from sage.all import *
>>> A.subdivide([Integer(2)],[]); A
[  0  -1]
[2*x  -2]
[-------]
>>> A.change_ring(RR['x,y'])
[                 0  -1.00000000000000]
[2.00000000000000*x  -2.00000000000000]
[-------------------------------------]
charpoly(var='x', **kwds)[source]#

Return the characteristic polynomial of this matrix.

Note

the generic sparse charpoly implementation in Sage is to just compute the charpoly of the corresponding dense matrix, so this could use a lot of memory. In particular, for this matrix, the charpoly will be computed using a dense algorithm.

EXAMPLES:

sage: A = matrix(ZZ, 4, range(16), sparse=True)
sage: A.charpoly()
x^4 - 30*x^3 - 80*x^2
sage: A.charpoly('y')
y^4 - 30*y^3 - 80*y^2
sage: A.charpoly()
x^4 - 30*x^3 - 80*x^2
>>> from sage.all import *
>>> A = matrix(ZZ, Integer(4), range(Integer(16)), sparse=True)
>>> A.charpoly()
x^4 - 30*x^3 - 80*x^2
>>> A.charpoly('y')
y^4 - 30*y^3 - 80*y^2
>>> A.charpoly()
x^4 - 30*x^3 - 80*x^2
density()[source]#

Return the density of the matrix.

By density we understand the ratio of the number of nonzero positions and the number self.nrows() * self.ncols(), i.e. the number of possible nonzero positions.

EXAMPLES:

sage: a = matrix([[],[],[],[]], sparse=True); a.density()
0
sage: a = matrix(5000,5000,{(1,2): 1}); a.density()
1/25000000
>>> from sage.all import *
>>> a = matrix([[],[],[],[]], sparse=True); a.density()
0
>>> a = matrix(Integer(5000),Integer(5000),{(Integer(1),Integer(2)): Integer(1)}); a.density()
1/25000000
determinant(**kwds)[source]#

Return the determinant of this matrix.

Note

the generic sparse determinant implementation in Sage is to just compute the determinant of the corresponding dense matrix, so this could use a lot of memory. In particular, for this matrix, the determinant will be computed using a dense algorithm.

EXAMPLES:

sage: A = matrix(ZZ, 4, range(16), sparse=True)
sage: B = A + identity_matrix(ZZ, 4, sparse=True)
sage: B.det()
-49
>>> from sage.all import *
>>> A = matrix(ZZ, Integer(4), range(Integer(16)), sparse=True)
>>> B = A + identity_matrix(ZZ, Integer(4), sparse=True)
>>> B.det()
-49
matrix_from_rows_and_columns(rows, columns)[source]#

Return the matrix constructed from self from the given rows and columns.

EXAMPLES:

sage: M = MatrixSpace(Integers(8),3,3, sparse=True)
sage: A = M(range(9)); A
[0 1 2]
[3 4 5]
[6 7 0]
sage: A.matrix_from_rows_and_columns([1], [0,2])
[3 5]
sage: A.matrix_from_rows_and_columns([1,2], [1,2])
[4 5]
[7 0]
>>> from sage.all import *
>>> M = MatrixSpace(Integers(Integer(8)),Integer(3),Integer(3), sparse=True)
>>> A = M(range(Integer(9))); A
[0 1 2]
[3 4 5]
[6 7 0]
>>> A.matrix_from_rows_and_columns([Integer(1)], [Integer(0),Integer(2)])
[3 5]
>>> A.matrix_from_rows_and_columns([Integer(1),Integer(2)], [Integer(1),Integer(2)])
[4 5]
[7 0]

Note that row and column indices can be reordered or repeated:

sage: A.matrix_from_rows_and_columns([2,1], [2,1])
[0 7]
[5 4]
>>> from sage.all import *
>>> A.matrix_from_rows_and_columns([Integer(2),Integer(1)], [Integer(2),Integer(1)])
[0 7]
[5 4]

For example here we take from row 1 columns 2 then 0 twice, and do this 3 times.

sage: A.matrix_from_rows_and_columns([1,1,1],[2,0,0])
[5 3 3]
[5 3 3]
[5 3 3]
>>> from sage.all import *
>>> A.matrix_from_rows_and_columns([Integer(1),Integer(1),Integer(1)],[Integer(2),Integer(0),Integer(0)])
[5 3 3]
[5 3 3]
[5 3 3]

We can efficiently extract large submatrices:

sage: A = random_matrix(ZZ, 100000, density=.00005, sparse=True)  # long time (4s on sage.math, 2012)
sage: B = A[50000:,:50000]        # long time
sage: count = 0
sage: for i, j in A.nonzero_positions():  # long time
....:     if i >= 50000 and j < 50000:
....:         assert B[i-50000, j] == A[i, j]
....:         count += 1
sage: count == sum(1 for _ in B.nonzero_positions())  # long time
True
>>> from sage.all import *
>>> A = random_matrix(ZZ, Integer(100000), density=RealNumber('.00005'), sparse=True)  # long time (4s on sage.math, 2012)
>>> B = A[Integer(50000):,:Integer(50000)]        # long time
>>> count = Integer(0)
>>> for i, j in A.nonzero_positions():  # long time
...     if i >= Integer(50000) and j < Integer(50000):
...         assert B[i-Integer(50000), j] == A[i, j]
...         count += Integer(1)
>>> count == sum(Integer(1) for _ in B.nonzero_positions())  # long time
True

We must pass in a list of indices:

sage: A = random_matrix(ZZ,100,density=.02,sparse=True)
sage: A.matrix_from_rows_and_columns(1,[2,3])
Traceback (most recent call last):
...
TypeError: 'sage.rings.integer.Integer' object is not iterable

sage: A.matrix_from_rows_and_columns([1,2],3)
Traceback (most recent call last):
...
TypeError: 'sage.rings.integer.Integer' object is not iterable
>>> from sage.all import *
>>> A = random_matrix(ZZ,Integer(100),density=RealNumber('.02'),sparse=True)
>>> A.matrix_from_rows_and_columns(Integer(1),[Integer(2),Integer(3)])
Traceback (most recent call last):
...
TypeError: 'sage.rings.integer.Integer' object is not iterable

>>> A.matrix_from_rows_and_columns([Integer(1),Integer(2)],Integer(3))
Traceback (most recent call last):
...
TypeError: 'sage.rings.integer.Integer' object is not iterable

AUTHORS:

  • Jaap Spies (2006-02-18)

  • Didier Deshommes: some Pyrex speedups implemented

  • Jason Grout: sparse matrix optimizations

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: M = MatrixSpace(QQ, 2, sparse=True)
sage: A = M([1,2,3,4]); A
[1 2]
[3 4]
sage: B = A.transpose(); B
[1 3]
[2 4]
>>> from sage.all import *
>>> M = MatrixSpace(QQ, Integer(2), sparse=True)
>>> A = M([Integer(1),Integer(2),Integer(3),Integer(4)]); A
[1 2]
[3 4]
>>> B = A.transpose(); B
[1 3]
[2 4]

.T is a convenient shortcut for the transpose:

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

See also

antitranspose()