Path: blob/master/sage/matrix/matrix_integer_sparse.pyx
4036 views
"""1Sparse integer matrices.23AUTHORS:4-- William Stein (2007-02-21)5-- Soroosh Yazdani (2007-02-21)67TESTS:8sage: a = matrix(ZZ,2,range(4), sparse=True)9sage: TestSuite(a).run()10sage: Matrix(ZZ,0,0,sparse=True).inverse()11[]12"""1314##############################################################################15# Copyright (C) 2007 William Stein <[email protected]>16# Distributed under the terms of the GNU General Public License (GPL)17# The full text of the GPL is available at:18# http://www.gnu.org/licenses/19##############################################################################2021include '../modules/binary_search.pxi'22include '../modules/vector_integer_sparse_h.pxi'23include '../modules/vector_integer_sparse_c.pxi'24include '../modules/vector_modn_sparse_h.pxi'25include '../modules/vector_modn_sparse_c.pxi'2627include '../ext/stdsage.pxi'2829from sage.rings.integer cimport Integer30from matrix cimport Matrix3132from matrix_modn_sparse cimport Matrix_modn_sparse33from sage.structure.element cimport ModuleElement, RingElement, Element, Vector3435import matrix_space3637from sage.rings.integer_ring import ZZ38from sage.rings.finite_rings.integer_mod_ring import IntegerModRing394041cdef class Matrix_integer_sparse(matrix_sparse.Matrix_sparse):4243########################################################################44# LEVEL 1 functionality45# * __cinit__46# * __dealloc__47# * __init__48# * set_unsafe49# * get_unsafe50# * __richcmp__ -- always the same51# * __hash__ -- always simple52########################################################################53def __cinit__(self, parent, entries, copy, coerce):54self._initialized = False55# set the parent, nrows, ncols, etc.56matrix_sparse.Matrix_sparse.__init__(self, parent)5758self._matrix = <mpz_vector*> sage_malloc(parent.nrows()*sizeof(mpz_vector))59if self._matrix == NULL:60raise MemoryError("error allocating sparse matrix")6162# initialize the rows63for i from 0 <= i < parent.nrows():64mpz_vector_init(&self._matrix[i], self._ncols, 0)65# record that rows have been initialized66self._initialized = True6768def __dealloc__(self):69cdef Py_ssize_t i70if self._initialized:71for i from 0 <= i < self._nrows:72mpz_vector_clear(&self._matrix[i])73sage_free(self._matrix)7475def __init__(self, parent, entries, copy, coerce):76"""77Create a sparse matrix over the integers.7879INPUT:80parent -- a matrix space81entries -- * a Python list of triples (i,j,x), where 0 <= i < nrows,820 <= j < ncols, and x is coercible to an int. The i,j83entry of self is set to x. The x's can be 0.84* Alternatively, entries can be a list of *all* the entries85of the sparse matrix (so they would be mostly 0).86copy -- ignored87coerce -- ignored88"""89cdef Py_ssize_t i, j, k90cdef Integer z91cdef void** X9293# fill in entries in the dict case94if entries is None: return95if isinstance(entries, dict):96R = self.base_ring()97for ij, x in entries.iteritems():98z = R(x)99if z != 0:100i, j = ij # nothing better to do since this is user input, which may be bogus.101if i < 0 or j < 0 or i >= self._nrows or j >= self._ncols:102raise IndexError("invalid entries list")103mpz_vector_set_entry(&self._matrix[i], j, z.value)104105elif isinstance(entries, list):106107# Dense input format -- fill in entries108if len(entries) != self._nrows * self._ncols:109raise TypeError("list of entries must be a dictionary of (i,j):x or a dense list of n * m elements")110seq = PySequence_Fast(entries,"expected a list")111X = PySequence_Fast_ITEMS(seq)112k = 0113R = self.base_ring()114# Get fast access to the entries list.115for i from 0 <= i < self._nrows:116for j from 0 <= j < self._ncols:117z = R(<object>X[k])118if z != 0:119mpz_vector_set_entry(&self._matrix[i], j, z.value)120k = k + 1121122else:123124# fill in entries in the scalar case125z = Integer(entries)126if z == 0:127return128if self._nrows != self._ncols:129raise TypeError("matrix must be square to initialize with a scalar.")130for i from 0 <= i < self._nrows:131mpz_vector_set_entry(&self._matrix[i], i, z.value)132133134cdef set_unsafe(self, Py_ssize_t i, Py_ssize_t j, x):135mpz_vector_set_entry(&self._matrix[i], j, (<Integer> x).value)136137cdef get_unsafe(self, Py_ssize_t i, Py_ssize_t j):138cdef Integer x139x = Integer()140mpz_vector_get_entry(&x.value, &self._matrix[i], j)141return x142143def __richcmp__(Matrix self, right, int op): # always need for mysterious reasons.144return self._richcmp(right, op)145146def __hash__(self):147return self._hash()148149150########################################################################151# LEVEL 2 functionality152# * def _pickle153# * def _unpickle154# * cdef _add_155# * cdef _sub_156# * cdef _mul_157# * cdef _cmp_c_impl158# * __neg__159# * __invert__160# * __copy__161# * _multiply_classical162# * _matrix_times_matrix_163# * _list -- list of underlying elements (need not be a copy)164# * x _dict -- sparse dictionary of underlying elements (need not be a copy)165########################################################################166# def _pickle(self):167# def _unpickle(self, data, int version): # use version >= 0168# cpdef ModuleElement _add_(self, ModuleElement right):169# cdef _mul_(self, Matrix right):170# cdef int _cmp_c_impl(self, Matrix right) except -2:171# def __neg__(self):172# def __invert__(self):173# def __copy__(self):174# def _multiply_classical(left, matrix.Matrix _right):175# def _list(self):176177cpdef ModuleElement _lmul_(self, RingElement right):178"""179EXAMPLES::180181sage: a = matrix(ZZ,2,range(6), sparse=True)182sage: 3 * a183[ 0 3 6]184[ 9 12 15]185"""186cdef Py_ssize_t i187cdef mpz_vector* self_row, *M_row188cdef Matrix_integer_sparse M189cdef Integer _x190_x = Integer(right)191M = Matrix_integer_sparse.__new__(Matrix_integer_sparse, self._parent, None, None, None)192for i from 0 <= i < self._nrows:193self_row = &self._matrix[i]194M_row = &M._matrix[i]195mpz_vector_scalar_multiply(M_row, self_row, _x.value)196return M197198cpdef ModuleElement _add_(self, ModuleElement right):199cdef Py_ssize_t i, j200cdef mpz_vector* self_row, *M_row201cdef Matrix_integer_sparse M202203M = Matrix_integer_sparse.__new__(Matrix_integer_sparse, self._parent, None, None, None)204cdef mpz_t mul205mpz_init_set_si(mul,1)206for i from 0 <= i < self._nrows:207mpz_vector_clear(&M._matrix[i])208add_mpz_vector_init(&M._matrix[i], &self._matrix[i], &(<Matrix_integer_sparse>right)._matrix[i], mul)209mpz_clear(mul)210return M211212cpdef ModuleElement _sub_(self, ModuleElement right):213cdef Py_ssize_t i, j214cdef mpz_vector* self_row, *M_row215cdef Matrix_integer_sparse M216217M = Matrix_integer_sparse.__new__(Matrix_integer_sparse, self._parent, None, None, None)218cdef mpz_t mul219mpz_init_set_si(mul,-1)220for i from 0 <= i < self._nrows:221mpz_vector_clear(&M._matrix[i])222add_mpz_vector_init(&M._matrix[i], &self._matrix[i], &(<Matrix_integer_sparse>right)._matrix[i], mul)223mpz_clear(mul)224return M225226def _dict(self):227"""228Unsafe version of the dict method, mainly for internal use.229This may return the dict of elements, but as an *unsafe*230reference to the underlying dict of the object. It might231be dangerous if you change entries of the returned dict.232"""233d = self.fetch('dict')234if not d is None:235return d236237cdef Py_ssize_t i, j, k238d = {}239for i from 0 <= i < self._nrows:240for j from 0 <= j < self._matrix[i].num_nonzero:241x = Integer()242mpz_set((<Integer>x).value, self._matrix[i].entries[j])243d[(int(i),int(self._matrix[i].positions[j]))] = x244self.cache('dict', d)245return d246247########################################################################248# LEVEL 3 functionality (Optional)249# * cdef _sub_250# * __deepcopy__251# * __invert__252# * Matrix windows -- only if you need strassen for that base253# * Other functions (list them here):254########################################################################255256def _nonzero_positions_by_row(self, copy=True):257"""258Returns the list of pairs (i,j) such that self[i,j] != 0.259260It is safe to change the resulting list (unless you give the option copy=False).261262EXAMPLE::263sage: M = Matrix(ZZ, [[0,0,0,1,0,0,0,0],[0,1,0,0,0,0,1,0]], sparse=True); M264[0 0 0 1 0 0 0 0]265[0 1 0 0 0 0 1 0]266sage: M._nonzero_positions_by_row()267[(0, 3), (1, 1), (1, 6)]268269"""270x = self.fetch('nonzero_positions')271if not x is None:272if copy:273return list(x)274return x275nzp = []276cdef Py_ssize_t i, j277for i from 0 <= i < self._nrows:278for j from 0 <= j < self._matrix[i].num_nonzero:279nzp.append((i,self._matrix[i].positions[j]))280self.cache('nonzero_positions', nzp)281if copy:282return list(nzp)283return nzp284285def _nonzero_positions_by_column(self, copy=True):286"""287Returns the list of pairs (i,j) such that self[i,j] != 0, but288sorted by columns, i.e., column j=0 entries occur first, then289column j=1 entries, etc.290291It is safe to change the resulting list (unless you give the option copy=False).292293EXAMPLE::294sage: M = Matrix(ZZ, [[0,0,0,1,0,0,0,0],[0,1,0,0,0,0,1,0]], sparse=True); M295[0 0 0 1 0 0 0 0]296[0 1 0 0 0 0 1 0]297sage: M._nonzero_positions_by_column()298[(1, 1), (0, 3), (1, 6)]299300"""301x = self.fetch('nonzero_positions_by_column')302if not x is None:303if copy:304return list(x)305return x306nzc = [[] for _ in xrange(self._ncols)]307cdef Py_ssize_t i, j308for i from 0 <= i < self._nrows:309for j from 0 <= j < self._matrix[i].num_nonzero:310p = self._matrix[i].positions[j]311nzc[p].append((i,p))312nzc = sum(nzc,[])313self.cache('nonzero_positions_by_column', nzc)314if copy:315return list(nzc)316return nzc317318def _mod_int(self, modulus):319"""320Helper function in reducing matrices mod n.321322INPUT:323324- `modulus` - a number325326OUTPUT:327328This matrix, over `ZZ/nZZ`.329330TESTS::331332sage: M = Matrix(ZZ, sparse=True)333sage: B = M._mod_int(7)334sage: B.parent()335Full MatrixSpace of 0 by 0 sparse matrices over Ring of integers modulo 7336337"""338return self._mod_int_c(modulus)339340cdef _mod_int_c(self, mod_int p):341cdef Py_ssize_t i, j342cdef Matrix_modn_sparse res343cdef mpz_vector* self_row344cdef c_vector_modint* res_row345res = Matrix_modn_sparse.__new__(Matrix_modn_sparse, matrix_space.MatrixSpace(346IntegerModRing(p), self._nrows, self._ncols, sparse=True), None, None, None)347for i from 0 <= i < self._nrows:348self_row = &(self._matrix[i])349res_row = &(res.rows[i])350for j from 0 <= j < self_row.num_nonzero:351set_entry(res_row, self_row.positions[j], mpz_fdiv_ui(self_row.entries[j], p))352return res353354355def rational_reconstruction(self, N):356"""357Use rational reconstruction to lift self to a matrix over the358rational numbers (if possible), where we view self as a matrix359modulo N.360361EXAMPLES:362sage: A = matrix(ZZ, 3, 4, [(1/3)%500, 2, 3, (-4)%500, 7, 2, 2, 3, 4, 3, 4, (5/7)%500], sparse=True)363sage: A.rational_reconstruction(500)364[1/3 2 3 -4]365[ 7 2 2 3]366[ 4 3 4 5/7]367368TEST:369370Check that ticket #9345 is fixed::371372sage: A = random_matrix(ZZ, 3, 3, sparse = True)373sage: A.rational_reconstruction(0)374Traceback (most recent call last):375...376ZeroDivisionError: The modulus cannot be zero377"""378import misc379return misc.matrix_integer_sparse_rational_reconstruction(self, N)380381def _linbox_sparse(self):382cdef Py_ssize_t i, j383v = ['%s %s M'%(self._nrows, self._ncols)]384d = self._dict()385for ij, x in d.iteritems():386v.append('%s %s %s'%(ij[0]+1,ij[1]+1,x))387v.append('0 0 0\n')388return '\n'.join(v)389390def _right_kernel_matrix(self, **kwds):391r"""392Returns a pair that includes a matrix of basis vectors393for the right kernel of ``self``.394395INPUT:396397- ``algorithm`` - determines which algorithm to use, options are:398399- 'pari' - use the ``matkerint()`` function from the PARI library400- 'padic' - use the p-adic algorithm from the IML library401- 'default' - use a heuristic to decide which of the two above402routines is fastest. This is the default value.403404- ``proof`` - this is passed to the p-adic IML algorithm.405If not specified, the global flag for linear algebra will be used.406407OUTPUT:408409Returns a pair. First item is the string is either410'computed-pari-int' or 'computed-iml-int', which identifies411the nature of the basis vectors.412413Second item is a matrix whose rows are a basis for the right kernel,414over the integers, as computed by either the IML or PARI libraries.415416EXAMPLES::417418sage: A = matrix(ZZ, [[4, 7, 9, 7, 5, 0],419... [1, 0, 5, 8, 9, 1],420... [0, 1, 0, 1, 9, 7],421... [4, 7, 6, 5, 1, 4]],422... sparse = True)423424sage: result = A._right_kernel_matrix(algorithm='pari')425sage: result[0]426'computed-pari-int'427sage: X = result[1]; X428[-26 31 -30 21 2 -10]429[-47 -13 48 -14 -11 18]430sage: A*X.transpose() == zero_matrix(ZZ, 4, 2)431True432433sage: result = A._right_kernel_matrix(algorithm='padic')434sage: result[0]435'computed-iml-int'436sage: X = result[1]; X437[-469 214 -30 119 -37 0]438[ 370 -165 18 -91 30 -2]439sage: A*X.transpose() == zero_matrix(ZZ, 4, 2)440True441442sage: result = A._right_kernel_matrix(algorithm='default')443sage: result[0]444'computed-pari-int'445sage: result[1]446[-26 31 -30 21 2 -10]447[-47 -13 48 -14 -11 18]448449sage: result = A._right_kernel_matrix()450sage: result[0]451'computed-pari-int'452sage: result[1]453[-26 31 -30 21 2 -10]454[-47 -13 48 -14 -11 18]455456With the 'default' given as the algorithm, several heuristics are457used to determine if PARI or IML ('padic') is used. The code has458exact details, but roughly speaking, relevant factors are: the459absolute size of the matrix, or the relative dimensions, or the460magnitude of the entries. ::461462sage: A = random_matrix(ZZ, 18, 11, sparse=True)463sage: A._right_kernel_matrix(algorithm='default')[0]464'computed-pari-int'465sage: A = random_matrix(ZZ, 18, 11, x = 10^200, sparse=True)466sage: A._right_kernel_matrix(algorithm='default')[0]467'computed-iml-int'468sage: A = random_matrix(ZZ, 60, 60, sparse=True)469sage: A._right_kernel_matrix(algorithm='default')[0]470'computed-iml-int'471sage: A = random_matrix(ZZ, 60, 55, sparse=True)472sage: A._right_kernel_matrix(algorithm='default')[0]473'computed-pari-int'474475TESTS:476477We test three trivial cases. PARI is used for small matrices,478but we let the heuristic decide that. ::479480sage: A = matrix(ZZ, 0, 2, sparse=True)481sage: A._right_kernel_matrix()[1]482[1 0]483[0 1]484sage: A = matrix(ZZ, 2, 0, sparse=True)485sage: A._right_kernel_matrix()[1].parent()486Full MatrixSpace of 0 by 0 dense matrices over Integer Ring487sage: A = zero_matrix(ZZ, 4, 3, sparse=True)488sage: A._right_kernel_matrix()[1]489[1 0 0]490[0 1 0]491[0 0 1]492"""493return self.dense_matrix()._right_kernel_matrix(**kwds)494495hermite_form = Matrix.echelon_form496497def elementary_divisors(self, algorithm='pari'):498"""499Return the elementary divisors of self, in order.500501The elementary divisors are the invariants of the finite502abelian group that is the cokernel of *left* multiplication by503this matrix. They are ordered in reverse by divisibility.504505INPUT:506self -- matrix507algorithm -- (default: 'pari')508'pari': works robustly, but is slower.509'linbox' -- use linbox (currently off, broken)510511OUTPUT:512list of integers513514EXAMPLES:515sage: matrix(3, range(9),sparse=True).elementary_divisors()516[1, 3, 0]517sage: M = matrix(ZZ, 3, [1,5,7, 3,6,9, 0,1,2], sparse=True)518sage: M.elementary_divisors()519[1, 1, 6]520521This returns a copy, which is safe to change:522sage: edivs = M.elementary_divisors()523sage: edivs.pop()5246525sage: M.elementary_divisors()526[1, 1, 6]527528SEE ALSO: smith_form529"""530return self.dense_matrix().elementary_divisors(algorithm=algorithm)531532def smith_form(self):533r"""534Returns matrices S, U, and V such that S = U\*self\*V, and S is in535Smith normal form. Thus S is diagonal with diagonal entries the536ordered elementary divisors of S.537538This version is for sparse matrices and simply makes the matrix539dense and calls the version for dense integer matrices.540541.. warning::542543The elementary_divisors function, which returns the544diagonal entries of S, is VASTLY faster than this function.545546The elementary divisors are the invariants of the finite abelian547group that is the cokernel of this matrix. They are ordered in548reverse by divisibility.549550EXAMPLES::551552sage: A = MatrixSpace(IntegerRing(), 3, sparse=True)(range(9))553sage: D, U, V = A.smith_form()554sage: D555[1 0 0]556[0 3 0]557[0 0 0]558sage: U559[ 0 1 0]560[ 0 -1 1]561[-1 2 -1]562sage: V563[-1 4 1]564[ 1 -3 -2]565[ 0 0 1]566sage: U*A*V567[1 0 0]568[0 3 0]569[0 0 0]570571It also makes sense for nonsquare matrices::572573sage: A = Matrix(ZZ,3,2,range(6), sparse=True)574sage: D, U, V = A.smith_form()575sage: D576[1 0]577[0 2]578[0 0]579sage: U580[ 0 1 0]581[ 0 -1 1]582[-1 2 -1]583sage: V584[-1 3]585[ 1 -2]586sage: U * A * V587[1 0]588[0 2]589[0 0]590591The examples above show that Trac ticket #10626 has been implemented.592593594.. seealso::595596:meth:`elementary_divisors`597"""598return self.dense_matrix().smith_form()599600601