Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/matrix/matrix_integer_sparse.pyx
4036 views
1
"""
2
Sparse integer matrices.
3
4
AUTHORS:
5
-- William Stein (2007-02-21)
6
-- Soroosh Yazdani (2007-02-21)
7
8
TESTS:
9
sage: a = matrix(ZZ,2,range(4), sparse=True)
10
sage: TestSuite(a).run()
11
sage: Matrix(ZZ,0,0,sparse=True).inverse()
12
[]
13
"""
14
15
##############################################################################
16
# Copyright (C) 2007 William Stein <[email protected]>
17
# Distributed under the terms of the GNU General Public License (GPL)
18
# The full text of the GPL is available at:
19
# http://www.gnu.org/licenses/
20
##############################################################################
21
22
include '../modules/binary_search.pxi'
23
include '../modules/vector_integer_sparse_h.pxi'
24
include '../modules/vector_integer_sparse_c.pxi'
25
include '../modules/vector_modn_sparse_h.pxi'
26
include '../modules/vector_modn_sparse_c.pxi'
27
28
include '../ext/stdsage.pxi'
29
30
from sage.rings.integer cimport Integer
31
from matrix cimport Matrix
32
33
from matrix_modn_sparse cimport Matrix_modn_sparse
34
from sage.structure.element cimport ModuleElement, RingElement, Element, Vector
35
36
import matrix_space
37
38
from sage.rings.integer_ring import ZZ
39
from sage.rings.finite_rings.integer_mod_ring import IntegerModRing
40
41
42
cdef class Matrix_integer_sparse(matrix_sparse.Matrix_sparse):
43
44
########################################################################
45
# LEVEL 1 functionality
46
# * __cinit__
47
# * __dealloc__
48
# * __init__
49
# * set_unsafe
50
# * get_unsafe
51
# * __richcmp__ -- always the same
52
# * __hash__ -- always simple
53
########################################################################
54
def __cinit__(self, parent, entries, copy, coerce):
55
self._initialized = False
56
# set the parent, nrows, ncols, etc.
57
matrix_sparse.Matrix_sparse.__init__(self, parent)
58
59
self._matrix = <mpz_vector*> sage_malloc(parent.nrows()*sizeof(mpz_vector))
60
if self._matrix == NULL:
61
raise MemoryError("error allocating sparse matrix")
62
63
# initialize the rows
64
for i from 0 <= i < parent.nrows():
65
mpz_vector_init(&self._matrix[i], self._ncols, 0)
66
# record that rows have been initialized
67
self._initialized = True
68
69
def __dealloc__(self):
70
cdef Py_ssize_t i
71
if self._initialized:
72
for i from 0 <= i < self._nrows:
73
mpz_vector_clear(&self._matrix[i])
74
sage_free(self._matrix)
75
76
def __init__(self, parent, entries, copy, coerce):
77
"""
78
Create a sparse matrix over the integers.
79
80
INPUT:
81
parent -- a matrix space
82
entries -- * a Python list of triples (i,j,x), where 0 <= i < nrows,
83
0 <= j < ncols, and x is coercible to an int. The i,j
84
entry of self is set to x. The x's can be 0.
85
* Alternatively, entries can be a list of *all* the entries
86
of the sparse matrix (so they would be mostly 0).
87
copy -- ignored
88
coerce -- ignored
89
"""
90
cdef Py_ssize_t i, j, k
91
cdef Integer z
92
cdef void** X
93
94
# fill in entries in the dict case
95
if entries is None: return
96
if isinstance(entries, dict):
97
R = self.base_ring()
98
for ij, x in entries.iteritems():
99
z = R(x)
100
if z != 0:
101
i, j = ij # nothing better to do since this is user input, which may be bogus.
102
if i < 0 or j < 0 or i >= self._nrows or j >= self._ncols:
103
raise IndexError("invalid entries list")
104
mpz_vector_set_entry(&self._matrix[i], j, z.value)
105
106
elif isinstance(entries, list):
107
108
# Dense input format -- fill in entries
109
if len(entries) != self._nrows * self._ncols:
110
raise TypeError("list of entries must be a dictionary of (i,j):x or a dense list of n * m elements")
111
seq = PySequence_Fast(entries,"expected a list")
112
X = PySequence_Fast_ITEMS(seq)
113
k = 0
114
R = self.base_ring()
115
# Get fast access to the entries list.
116
for i from 0 <= i < self._nrows:
117
for j from 0 <= j < self._ncols:
118
z = R(<object>X[k])
119
if z != 0:
120
mpz_vector_set_entry(&self._matrix[i], j, z.value)
121
k = k + 1
122
123
else:
124
125
# fill in entries in the scalar case
126
z = Integer(entries)
127
if z == 0:
128
return
129
if self._nrows != self._ncols:
130
raise TypeError("matrix must be square to initialize with a scalar.")
131
for i from 0 <= i < self._nrows:
132
mpz_vector_set_entry(&self._matrix[i], i, z.value)
133
134
135
cdef set_unsafe(self, Py_ssize_t i, Py_ssize_t j, x):
136
mpz_vector_set_entry(&self._matrix[i], j, (<Integer> x).value)
137
138
cdef get_unsafe(self, Py_ssize_t i, Py_ssize_t j):
139
cdef Integer x
140
x = Integer()
141
mpz_vector_get_entry(&x.value, &self._matrix[i], j)
142
return x
143
144
def __richcmp__(Matrix self, right, int op): # always need for mysterious reasons.
145
return self._richcmp(right, op)
146
147
def __hash__(self):
148
return self._hash()
149
150
151
########################################################################
152
# LEVEL 2 functionality
153
# * def _pickle
154
# * def _unpickle
155
# * cdef _add_
156
# * cdef _sub_
157
# * cdef _mul_
158
# * cdef _cmp_c_impl
159
# * __neg__
160
# * __invert__
161
# * __copy__
162
# * _multiply_classical
163
# * _matrix_times_matrix_
164
# * _list -- list of underlying elements (need not be a copy)
165
# * x _dict -- sparse dictionary of underlying elements (need not be a copy)
166
########################################################################
167
# def _pickle(self):
168
# def _unpickle(self, data, int version): # use version >= 0
169
# cpdef ModuleElement _add_(self, ModuleElement right):
170
# cdef _mul_(self, Matrix right):
171
# cdef int _cmp_c_impl(self, Matrix right) except -2:
172
# def __neg__(self):
173
# def __invert__(self):
174
# def __copy__(self):
175
# def _multiply_classical(left, matrix.Matrix _right):
176
# def _list(self):
177
178
cpdef ModuleElement _lmul_(self, RingElement right):
179
"""
180
EXAMPLES::
181
182
sage: a = matrix(ZZ,2,range(6), sparse=True)
183
sage: 3 * a
184
[ 0 3 6]
185
[ 9 12 15]
186
"""
187
cdef Py_ssize_t i
188
cdef mpz_vector* self_row, *M_row
189
cdef Matrix_integer_sparse M
190
cdef Integer _x
191
_x = Integer(right)
192
M = Matrix_integer_sparse.__new__(Matrix_integer_sparse, self._parent, None, None, None)
193
for i from 0 <= i < self._nrows:
194
self_row = &self._matrix[i]
195
M_row = &M._matrix[i]
196
mpz_vector_scalar_multiply(M_row, self_row, _x.value)
197
return M
198
199
cpdef ModuleElement _add_(self, ModuleElement right):
200
cdef Py_ssize_t i, j
201
cdef mpz_vector* self_row, *M_row
202
cdef Matrix_integer_sparse M
203
204
M = Matrix_integer_sparse.__new__(Matrix_integer_sparse, self._parent, None, None, None)
205
cdef mpz_t mul
206
mpz_init_set_si(mul,1)
207
for i from 0 <= i < self._nrows:
208
mpz_vector_clear(&M._matrix[i])
209
add_mpz_vector_init(&M._matrix[i], &self._matrix[i], &(<Matrix_integer_sparse>right)._matrix[i], mul)
210
mpz_clear(mul)
211
return M
212
213
cpdef ModuleElement _sub_(self, ModuleElement right):
214
cdef Py_ssize_t i, j
215
cdef mpz_vector* self_row, *M_row
216
cdef Matrix_integer_sparse M
217
218
M = Matrix_integer_sparse.__new__(Matrix_integer_sparse, self._parent, None, None, None)
219
cdef mpz_t mul
220
mpz_init_set_si(mul,-1)
221
for i from 0 <= i < self._nrows:
222
mpz_vector_clear(&M._matrix[i])
223
add_mpz_vector_init(&M._matrix[i], &self._matrix[i], &(<Matrix_integer_sparse>right)._matrix[i], mul)
224
mpz_clear(mul)
225
return M
226
227
def _dict(self):
228
"""
229
Unsafe version of the dict method, mainly for internal use.
230
This may return the dict of elements, but as an *unsafe*
231
reference to the underlying dict of the object. It might
232
be dangerous if you change entries of the returned dict.
233
"""
234
d = self.fetch('dict')
235
if not d is None:
236
return d
237
238
cdef Py_ssize_t i, j, k
239
d = {}
240
for i from 0 <= i < self._nrows:
241
for j from 0 <= j < self._matrix[i].num_nonzero:
242
x = Integer()
243
mpz_set((<Integer>x).value, self._matrix[i].entries[j])
244
d[(int(i),int(self._matrix[i].positions[j]))] = x
245
self.cache('dict', d)
246
return d
247
248
########################################################################
249
# LEVEL 3 functionality (Optional)
250
# * cdef _sub_
251
# * __deepcopy__
252
# * __invert__
253
# * Matrix windows -- only if you need strassen for that base
254
# * Other functions (list them here):
255
########################################################################
256
257
def _nonzero_positions_by_row(self, copy=True):
258
"""
259
Returns the list of pairs (i,j) such that self[i,j] != 0.
260
261
It is safe to change the resulting list (unless you give the option copy=False).
262
263
EXAMPLE::
264
sage: M = Matrix(ZZ, [[0,0,0,1,0,0,0,0],[0,1,0,0,0,0,1,0]], sparse=True); M
265
[0 0 0 1 0 0 0 0]
266
[0 1 0 0 0 0 1 0]
267
sage: M._nonzero_positions_by_row()
268
[(0, 3), (1, 1), (1, 6)]
269
270
"""
271
x = self.fetch('nonzero_positions')
272
if not x is None:
273
if copy:
274
return list(x)
275
return x
276
nzp = []
277
cdef Py_ssize_t i, j
278
for i from 0 <= i < self._nrows:
279
for j from 0 <= j < self._matrix[i].num_nonzero:
280
nzp.append((i,self._matrix[i].positions[j]))
281
self.cache('nonzero_positions', nzp)
282
if copy:
283
return list(nzp)
284
return nzp
285
286
def _nonzero_positions_by_column(self, copy=True):
287
"""
288
Returns the list of pairs (i,j) such that self[i,j] != 0, but
289
sorted by columns, i.e., column j=0 entries occur first, then
290
column j=1 entries, etc.
291
292
It is safe to change the resulting list (unless you give the option copy=False).
293
294
EXAMPLE::
295
sage: M = Matrix(ZZ, [[0,0,0,1,0,0,0,0],[0,1,0,0,0,0,1,0]], sparse=True); M
296
[0 0 0 1 0 0 0 0]
297
[0 1 0 0 0 0 1 0]
298
sage: M._nonzero_positions_by_column()
299
[(1, 1), (0, 3), (1, 6)]
300
301
"""
302
x = self.fetch('nonzero_positions_by_column')
303
if not x is None:
304
if copy:
305
return list(x)
306
return x
307
nzc = [[] for _ in xrange(self._ncols)]
308
cdef Py_ssize_t i, j
309
for i from 0 <= i < self._nrows:
310
for j from 0 <= j < self._matrix[i].num_nonzero:
311
p = self._matrix[i].positions[j]
312
nzc[p].append((i,p))
313
nzc = sum(nzc,[])
314
self.cache('nonzero_positions_by_column', nzc)
315
if copy:
316
return list(nzc)
317
return nzc
318
319
def _mod_int(self, modulus):
320
"""
321
Helper function in reducing matrices mod n.
322
323
INPUT:
324
325
- `modulus` - a number
326
327
OUTPUT:
328
329
This matrix, over `ZZ/nZZ`.
330
331
TESTS::
332
333
sage: M = Matrix(ZZ, sparse=True)
334
sage: B = M._mod_int(7)
335
sage: B.parent()
336
Full MatrixSpace of 0 by 0 sparse matrices over Ring of integers modulo 7
337
338
"""
339
return self._mod_int_c(modulus)
340
341
cdef _mod_int_c(self, mod_int p):
342
cdef Py_ssize_t i, j
343
cdef Matrix_modn_sparse res
344
cdef mpz_vector* self_row
345
cdef c_vector_modint* res_row
346
res = Matrix_modn_sparse.__new__(Matrix_modn_sparse, matrix_space.MatrixSpace(
347
IntegerModRing(p), self._nrows, self._ncols, sparse=True), None, None, None)
348
for i from 0 <= i < self._nrows:
349
self_row = &(self._matrix[i])
350
res_row = &(res.rows[i])
351
for j from 0 <= j < self_row.num_nonzero:
352
set_entry(res_row, self_row.positions[j], mpz_fdiv_ui(self_row.entries[j], p))
353
return res
354
355
356
def rational_reconstruction(self, N):
357
"""
358
Use rational reconstruction to lift self to a matrix over the
359
rational numbers (if possible), where we view self as a matrix
360
modulo N.
361
362
EXAMPLES:
363
sage: A = matrix(ZZ, 3, 4, [(1/3)%500, 2, 3, (-4)%500, 7, 2, 2, 3, 4, 3, 4, (5/7)%500], sparse=True)
364
sage: A.rational_reconstruction(500)
365
[1/3 2 3 -4]
366
[ 7 2 2 3]
367
[ 4 3 4 5/7]
368
369
TEST:
370
371
Check that ticket #9345 is fixed::
372
373
sage: A = random_matrix(ZZ, 3, 3, sparse = True)
374
sage: A.rational_reconstruction(0)
375
Traceback (most recent call last):
376
...
377
ZeroDivisionError: The modulus cannot be zero
378
"""
379
import misc
380
return misc.matrix_integer_sparse_rational_reconstruction(self, N)
381
382
def _linbox_sparse(self):
383
cdef Py_ssize_t i, j
384
v = ['%s %s M'%(self._nrows, self._ncols)]
385
d = self._dict()
386
for ij, x in d.iteritems():
387
v.append('%s %s %s'%(ij[0]+1,ij[1]+1,x))
388
v.append('0 0 0\n')
389
return '\n'.join(v)
390
391
def _right_kernel_matrix(self, **kwds):
392
r"""
393
Returns a pair that includes a matrix of basis vectors
394
for the right kernel of ``self``.
395
396
INPUT:
397
398
- ``algorithm`` - determines which algorithm to use, options are:
399
400
- 'pari' - use the ``matkerint()`` function from the PARI library
401
- 'padic' - use the p-adic algorithm from the IML library
402
- 'default' - use a heuristic to decide which of the two above
403
routines is fastest. This is the default value.
404
405
- ``proof`` - this is passed to the p-adic IML algorithm.
406
If not specified, the global flag for linear algebra will be used.
407
408
OUTPUT:
409
410
Returns a pair. First item is the string is either
411
'computed-pari-int' or 'computed-iml-int', which identifies
412
the nature of the basis vectors.
413
414
Second item is a matrix whose rows are a basis for the right kernel,
415
over the integers, as computed by either the IML or PARI libraries.
416
417
EXAMPLES::
418
419
sage: A = matrix(ZZ, [[4, 7, 9, 7, 5, 0],
420
... [1, 0, 5, 8, 9, 1],
421
... [0, 1, 0, 1, 9, 7],
422
... [4, 7, 6, 5, 1, 4]],
423
... sparse = True)
424
425
sage: result = A._right_kernel_matrix(algorithm='pari')
426
sage: result[0]
427
'computed-pari-int'
428
sage: X = result[1]; X
429
[-26 31 -30 21 2 -10]
430
[-47 -13 48 -14 -11 18]
431
sage: A*X.transpose() == zero_matrix(ZZ, 4, 2)
432
True
433
434
sage: result = A._right_kernel_matrix(algorithm='padic')
435
sage: result[0]
436
'computed-iml-int'
437
sage: X = result[1]; X
438
[-469 214 -30 119 -37 0]
439
[ 370 -165 18 -91 30 -2]
440
sage: A*X.transpose() == zero_matrix(ZZ, 4, 2)
441
True
442
443
sage: result = A._right_kernel_matrix(algorithm='default')
444
sage: result[0]
445
'computed-pari-int'
446
sage: result[1]
447
[-26 31 -30 21 2 -10]
448
[-47 -13 48 -14 -11 18]
449
450
sage: result = A._right_kernel_matrix()
451
sage: result[0]
452
'computed-pari-int'
453
sage: result[1]
454
[-26 31 -30 21 2 -10]
455
[-47 -13 48 -14 -11 18]
456
457
With the 'default' given as the algorithm, several heuristics are
458
used to determine if PARI or IML ('padic') is used. The code has
459
exact details, but roughly speaking, relevant factors are: the
460
absolute size of the matrix, or the relative dimensions, or the
461
magnitude of the entries. ::
462
463
sage: A = random_matrix(ZZ, 18, 11, sparse=True)
464
sage: A._right_kernel_matrix(algorithm='default')[0]
465
'computed-pari-int'
466
sage: A = random_matrix(ZZ, 18, 11, x = 10^200, sparse=True)
467
sage: A._right_kernel_matrix(algorithm='default')[0]
468
'computed-iml-int'
469
sage: A = random_matrix(ZZ, 60, 60, sparse=True)
470
sage: A._right_kernel_matrix(algorithm='default')[0]
471
'computed-iml-int'
472
sage: A = random_matrix(ZZ, 60, 55, sparse=True)
473
sage: A._right_kernel_matrix(algorithm='default')[0]
474
'computed-pari-int'
475
476
TESTS:
477
478
We test three trivial cases. PARI is used for small matrices,
479
but we let the heuristic decide that. ::
480
481
sage: A = matrix(ZZ, 0, 2, sparse=True)
482
sage: A._right_kernel_matrix()[1]
483
[1 0]
484
[0 1]
485
sage: A = matrix(ZZ, 2, 0, sparse=True)
486
sage: A._right_kernel_matrix()[1].parent()
487
Full MatrixSpace of 0 by 0 dense matrices over Integer Ring
488
sage: A = zero_matrix(ZZ, 4, 3, sparse=True)
489
sage: A._right_kernel_matrix()[1]
490
[1 0 0]
491
[0 1 0]
492
[0 0 1]
493
"""
494
return self.dense_matrix()._right_kernel_matrix(**kwds)
495
496
hermite_form = Matrix.echelon_form
497
498
def elementary_divisors(self, algorithm='pari'):
499
"""
500
Return the elementary divisors of self, in order.
501
502
The elementary divisors are the invariants of the finite
503
abelian group that is the cokernel of *left* multiplication by
504
this matrix. They are ordered in reverse by divisibility.
505
506
INPUT:
507
self -- matrix
508
algorithm -- (default: 'pari')
509
'pari': works robustly, but is slower.
510
'linbox' -- use linbox (currently off, broken)
511
512
OUTPUT:
513
list of integers
514
515
EXAMPLES:
516
sage: matrix(3, range(9),sparse=True).elementary_divisors()
517
[1, 3, 0]
518
sage: M = matrix(ZZ, 3, [1,5,7, 3,6,9, 0,1,2], sparse=True)
519
sage: M.elementary_divisors()
520
[1, 1, 6]
521
522
This returns a copy, which is safe to change:
523
sage: edivs = M.elementary_divisors()
524
sage: edivs.pop()
525
6
526
sage: M.elementary_divisors()
527
[1, 1, 6]
528
529
SEE ALSO: smith_form
530
"""
531
return self.dense_matrix().elementary_divisors(algorithm=algorithm)
532
533
def smith_form(self):
534
r"""
535
Returns matrices S, U, and V such that S = U\*self\*V, and S is in
536
Smith normal form. Thus S is diagonal with diagonal entries the
537
ordered elementary divisors of S.
538
539
This version is for sparse matrices and simply makes the matrix
540
dense and calls the version for dense integer matrices.
541
542
.. warning::
543
544
The elementary_divisors function, which returns the
545
diagonal entries of S, is VASTLY faster than this function.
546
547
The elementary divisors are the invariants of the finite abelian
548
group that is the cokernel of this matrix. They are ordered in
549
reverse by divisibility.
550
551
EXAMPLES::
552
553
sage: A = MatrixSpace(IntegerRing(), 3, sparse=True)(range(9))
554
sage: D, U, V = A.smith_form()
555
sage: D
556
[1 0 0]
557
[0 3 0]
558
[0 0 0]
559
sage: U
560
[ 0 1 0]
561
[ 0 -1 1]
562
[-1 2 -1]
563
sage: V
564
[-1 4 1]
565
[ 1 -3 -2]
566
[ 0 0 1]
567
sage: U*A*V
568
[1 0 0]
569
[0 3 0]
570
[0 0 0]
571
572
It also makes sense for nonsquare matrices::
573
574
sage: A = Matrix(ZZ,3,2,range(6), sparse=True)
575
sage: D, U, V = A.smith_form()
576
sage: D
577
[1 0]
578
[0 2]
579
[0 0]
580
sage: U
581
[ 0 1 0]
582
[ 0 -1 1]
583
[-1 2 -1]
584
sage: V
585
[-1 3]
586
[ 1 -2]
587
sage: U * A * V
588
[1 0]
589
[0 2]
590
[0 0]
591
592
The examples above show that Trac ticket #10626 has been implemented.
593
594
595
.. seealso::
596
597
:meth:`elementary_divisors`
598
"""
599
return self.dense_matrix().smith_form()
600
601