_base.py 44 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331
  1. """Base class for sparse matrices"""
  2. from warnings import warn
  3. import numpy as np
  4. from ._sputils import (asmatrix, check_reshape_kwargs, check_shape,
  5. get_sum_dtype, isdense, isintlike, isscalarlike,
  6. matrix, validateaxis)
  7. __all__ = ['spmatrix', 'isspmatrix', 'issparse',
  8. 'SparseWarning', 'SparseEfficiencyWarning']
  9. class SparseWarning(Warning):
  10. pass
  11. class SparseFormatWarning(SparseWarning):
  12. pass
  13. class SparseEfficiencyWarning(SparseWarning):
  14. pass
  15. # The formats that we might potentially understand.
  16. _formats = {'csc': [0, "Compressed Sparse Column"],
  17. 'csr': [1, "Compressed Sparse Row"],
  18. 'dok': [2, "Dictionary Of Keys"],
  19. 'lil': [3, "List of Lists"],
  20. 'dod': [4, "Dictionary of Dictionaries"],
  21. 'sss': [5, "Symmetric Sparse Skyline"],
  22. 'coo': [6, "COOrdinate"],
  23. 'lba': [7, "Linpack BAnded"],
  24. 'egd': [8, "Ellpack-itpack Generalized Diagonal"],
  25. 'dia': [9, "DIAgonal"],
  26. 'bsr': [10, "Block Sparse Row"],
  27. 'msr': [11, "Modified compressed Sparse Row"],
  28. 'bsc': [12, "Block Sparse Column"],
  29. 'msc': [13, "Modified compressed Sparse Column"],
  30. 'ssk': [14, "Symmetric SKyline"],
  31. 'nsk': [15, "Nonsymmetric SKyline"],
  32. 'jad': [16, "JAgged Diagonal"],
  33. 'uss': [17, "Unsymmetric Sparse Skyline"],
  34. 'vbr': [18, "Variable Block Row"],
  35. 'und': [19, "Undefined"]
  36. }
  37. # These univariate ufuncs preserve zeros.
  38. _ufuncs_with_fixed_point_at_zero = frozenset([
  39. np.sin, np.tan, np.arcsin, np.arctan, np.sinh, np.tanh, np.arcsinh,
  40. np.arctanh, np.rint, np.sign, np.expm1, np.log1p, np.deg2rad,
  41. np.rad2deg, np.floor, np.ceil, np.trunc, np.sqrt])
  42. MAXPRINT = 50
  43. class spmatrix:
  44. """ This class provides a base class for all sparse matrices. It
  45. cannot be instantiated. Most of the work is provided by subclasses.
  46. """
  47. __array_priority__ = 10.1
  48. ndim = 2
  49. @property
  50. def _bsr_container(self):
  51. from ._bsr import bsr_matrix
  52. return bsr_matrix
  53. @property
  54. def _coo_container(self):
  55. from ._coo import coo_matrix
  56. return coo_matrix
  57. @property
  58. def _csc_container(self):
  59. from ._csc import csc_matrix
  60. return csc_matrix
  61. @property
  62. def _csr_container(self):
  63. from ._csr import csr_matrix
  64. return csr_matrix
  65. @property
  66. def _dia_container(self):
  67. from ._dia import dia_matrix
  68. return dia_matrix
  69. @property
  70. def _dok_container(self):
  71. from ._dok import dok_matrix
  72. return dok_matrix
  73. @property
  74. def _lil_container(self):
  75. from ._lil import lil_matrix
  76. return lil_matrix
  77. _is_array = False
  78. def __init__(self, maxprint=MAXPRINT):
  79. self._shape = None
  80. if self.__class__.__name__ == 'spmatrix':
  81. raise ValueError("This class is not intended"
  82. " to be instantiated directly.")
  83. self.maxprint = maxprint
  84. def set_shape(self, shape):
  85. """See `reshape`."""
  86. # Make sure copy is False since this is in place
  87. # Make sure format is unchanged because we are doing a __dict__ swap
  88. new_matrix = self.reshape(shape, copy=False).asformat(self.format)
  89. self.__dict__ = new_matrix.__dict__
  90. def get_shape(self):
  91. """Get shape of a matrix."""
  92. return self._shape
  93. shape = property(fget=get_shape, fset=set_shape)
  94. def reshape(self, *args, **kwargs):
  95. """reshape(self, shape, order='C', copy=False)
  96. Gives a new shape to a sparse matrix without changing its data.
  97. Parameters
  98. ----------
  99. shape : length-2 tuple of ints
  100. The new shape should be compatible with the original shape.
  101. order : {'C', 'F'}, optional
  102. Read the elements using this index order. 'C' means to read and
  103. write the elements using C-like index order; e.g., read entire first
  104. row, then second row, etc. 'F' means to read and write the elements
  105. using Fortran-like index order; e.g., read entire first column, then
  106. second column, etc.
  107. copy : bool, optional
  108. Indicates whether or not attributes of self should be copied
  109. whenever possible. The degree to which attributes are copied varies
  110. depending on the type of sparse matrix being used.
  111. Returns
  112. -------
  113. reshaped_matrix : sparse matrix
  114. A sparse matrix with the given `shape`, not necessarily of the same
  115. format as the current object.
  116. See Also
  117. --------
  118. numpy.matrix.reshape : NumPy's implementation of 'reshape' for
  119. matrices
  120. """
  121. # If the shape already matches, don't bother doing an actual reshape
  122. # Otherwise, the default is to convert to COO and use its reshape
  123. shape = check_shape(args, self.shape)
  124. order, copy = check_reshape_kwargs(kwargs)
  125. if shape == self.shape:
  126. if copy:
  127. return self.copy()
  128. else:
  129. return self
  130. return self.tocoo(copy=copy).reshape(shape, order=order, copy=False)
  131. def resize(self, shape):
  132. """Resize the matrix in-place to dimensions given by ``shape``
  133. Any elements that lie within the new shape will remain at the same
  134. indices, while non-zero elements lying outside the new shape are
  135. removed.
  136. Parameters
  137. ----------
  138. shape : (int, int)
  139. number of rows and columns in the new matrix
  140. Notes
  141. -----
  142. The semantics are not identical to `numpy.ndarray.resize` or
  143. `numpy.resize`. Here, the same data will be maintained at each index
  144. before and after reshape, if that index is within the new bounds. In
  145. numpy, resizing maintains contiguity of the array, moving elements
  146. around in the logical matrix but not within a flattened representation.
  147. We give no guarantees about whether the underlying data attributes
  148. (arrays, etc.) will be modified in place or replaced with new objects.
  149. """
  150. # As an inplace operation, this requires implementation in each format.
  151. raise NotImplementedError(
  152. '{}.resize is not implemented'.format(type(self).__name__))
  153. def astype(self, dtype, casting='unsafe', copy=True):
  154. """Cast the matrix elements to a specified type.
  155. Parameters
  156. ----------
  157. dtype : string or numpy dtype
  158. Typecode or data-type to which to cast the data.
  159. casting : {'no', 'equiv', 'safe', 'same_kind', 'unsafe'}, optional
  160. Controls what kind of data casting may occur.
  161. Defaults to 'unsafe' for backwards compatibility.
  162. 'no' means the data types should not be cast at all.
  163. 'equiv' means only byte-order changes are allowed.
  164. 'safe' means only casts which can preserve values are allowed.
  165. 'same_kind' means only safe casts or casts within a kind,
  166. like float64 to float32, are allowed.
  167. 'unsafe' means any data conversions may be done.
  168. copy : bool, optional
  169. If `copy` is `False`, the result might share some memory with this
  170. matrix. If `copy` is `True`, it is guaranteed that the result and
  171. this matrix do not share any memory.
  172. """
  173. dtype = np.dtype(dtype)
  174. if self.dtype != dtype:
  175. return self.tocsr().astype(
  176. dtype, casting=casting, copy=copy).asformat(self.format)
  177. elif copy:
  178. return self.copy()
  179. else:
  180. return self
  181. @classmethod
  182. def _ascontainer(cls, X, **kwargs):
  183. if cls._is_array:
  184. return np.asarray(X, **kwargs)
  185. else:
  186. return asmatrix(X, **kwargs)
  187. @classmethod
  188. def _container(cls, X, **kwargs):
  189. if cls._is_array:
  190. return np.array(X, **kwargs)
  191. else:
  192. return matrix(X, **kwargs)
  193. def asfptype(self):
  194. """Upcast matrix to a floating point format (if necessary)"""
  195. fp_types = ['f', 'd', 'F', 'D']
  196. if self.dtype.char in fp_types:
  197. return self
  198. else:
  199. for fp_type in fp_types:
  200. if self.dtype <= np.dtype(fp_type):
  201. return self.astype(fp_type)
  202. raise TypeError('cannot upcast [%s] to a floating '
  203. 'point format' % self.dtype.name)
  204. def __iter__(self):
  205. for r in range(self.shape[0]):
  206. yield self[r, :]
  207. def getmaxprint(self):
  208. """Maximum number of elements to display when printed."""
  209. return self.maxprint
  210. def count_nonzero(self):
  211. """Number of non-zero entries, equivalent to
  212. np.count_nonzero(a.toarray())
  213. Unlike getnnz() and the nnz property, which return the number of stored
  214. entries (the length of the data attribute), this method counts the
  215. actual number of non-zero entries in data.
  216. """
  217. raise NotImplementedError("count_nonzero not implemented for %s." %
  218. self.__class__.__name__)
  219. def getnnz(self, axis=None):
  220. """Number of stored values, including explicit zeros.
  221. Parameters
  222. ----------
  223. axis : None, 0, or 1
  224. Select between the number of values across the whole matrix, in
  225. each column, or in each row.
  226. See also
  227. --------
  228. count_nonzero : Number of non-zero entries
  229. """
  230. raise NotImplementedError("getnnz not implemented for %s." %
  231. self.__class__.__name__)
  232. @property
  233. def nnz(self):
  234. """Number of stored values, including explicit zeros.
  235. See also
  236. --------
  237. count_nonzero : Number of non-zero entries
  238. """
  239. return self.getnnz()
  240. def getformat(self):
  241. """Format of a matrix representation as a string."""
  242. return getattr(self, 'format', 'und')
  243. def __repr__(self):
  244. _, format_name = _formats[self.getformat()]
  245. sparse_cls = 'array' if self._is_array else 'matrix'
  246. return f"<%dx%d sparse {sparse_cls} of type '%s'\n" \
  247. "\twith %d stored elements in %s format>" % \
  248. (self.shape + (self.dtype.type, self.nnz, format_name))
  249. def __str__(self):
  250. maxprint = self.getmaxprint()
  251. A = self.tocoo()
  252. # helper function, outputs "(i,j) v"
  253. def tostr(row, col, data):
  254. triples = zip(list(zip(row, col)), data)
  255. return '\n'.join([(' %s\t%s' % t) for t in triples])
  256. if self.nnz > maxprint:
  257. half = maxprint // 2
  258. out = tostr(A.row[:half], A.col[:half], A.data[:half])
  259. out += "\n :\t:\n"
  260. half = maxprint - maxprint//2
  261. out += tostr(A.row[-half:], A.col[-half:], A.data[-half:])
  262. else:
  263. out = tostr(A.row, A.col, A.data)
  264. return out
  265. def __bool__(self): # Simple -- other ideas?
  266. if self.shape == (1, 1):
  267. return self.nnz != 0
  268. else:
  269. raise ValueError("The truth value of an array with more than one "
  270. "element is ambiguous. Use a.any() or a.all().")
  271. __nonzero__ = __bool__
  272. # What should len(sparse) return? For consistency with dense matrices,
  273. # perhaps it should be the number of rows? But for some uses the number of
  274. # non-zeros is more important. For now, raise an exception!
  275. def __len__(self):
  276. raise TypeError("sparse matrix length is ambiguous; use getnnz()"
  277. " or shape[0]")
  278. def asformat(self, format, copy=False):
  279. """Return this matrix in the passed format.
  280. Parameters
  281. ----------
  282. format : {str, None}
  283. The desired matrix format ("csr", "csc", "lil", "dok", "array", ...)
  284. or None for no conversion.
  285. copy : bool, optional
  286. If True, the result is guaranteed to not share data with self.
  287. Returns
  288. -------
  289. A : This matrix in the passed format.
  290. """
  291. if format is None or format == self.format:
  292. if copy:
  293. return self.copy()
  294. else:
  295. return self
  296. else:
  297. try:
  298. convert_method = getattr(self, 'to' + format)
  299. except AttributeError as e:
  300. raise ValueError('Format {} is unknown.'.format(format)) from e
  301. # Forward the copy kwarg, if it's accepted.
  302. try:
  303. return convert_method(copy=copy)
  304. except TypeError:
  305. return convert_method()
  306. ###################################################################
  307. # NOTE: All arithmetic operations use csr_matrix by default.
  308. # Therefore a new sparse matrix format just needs to define a
  309. # .tocsr() method to provide arithmetic support. Any of these
  310. # methods can be overridden for efficiency.
  311. ####################################################################
  312. def multiply(self, other):
  313. """Point-wise multiplication by another matrix
  314. """
  315. return self.tocsr().multiply(other)
  316. def maximum(self, other):
  317. """Element-wise maximum between this and another matrix."""
  318. return self.tocsr().maximum(other)
  319. def minimum(self, other):
  320. """Element-wise minimum between this and another matrix."""
  321. return self.tocsr().minimum(other)
  322. def dot(self, other):
  323. """Ordinary dot product
  324. Examples
  325. --------
  326. >>> import numpy as np
  327. >>> from scipy.sparse import csr_matrix
  328. >>> A = csr_matrix([[1, 2, 0], [0, 0, 3], [4, 0, 5]])
  329. >>> v = np.array([1, 0, -1])
  330. >>> A.dot(v)
  331. array([ 1, -3, -1], dtype=int64)
  332. """
  333. if np.isscalar(other):
  334. return self * other
  335. else:
  336. return self @ other
  337. def power(self, n, dtype=None):
  338. """Element-wise power."""
  339. return self.tocsr().power(n, dtype=dtype)
  340. def __eq__(self, other):
  341. return self.tocsr().__eq__(other)
  342. def __ne__(self, other):
  343. return self.tocsr().__ne__(other)
  344. def __lt__(self, other):
  345. return self.tocsr().__lt__(other)
  346. def __gt__(self, other):
  347. return self.tocsr().__gt__(other)
  348. def __le__(self, other):
  349. return self.tocsr().__le__(other)
  350. def __ge__(self, other):
  351. return self.tocsr().__ge__(other)
  352. def __abs__(self):
  353. return abs(self.tocsr())
  354. def __round__(self, ndigits=0):
  355. return round(self.tocsr(), ndigits=ndigits)
  356. def _add_sparse(self, other):
  357. return self.tocsr()._add_sparse(other)
  358. def _add_dense(self, other):
  359. return self.tocoo()._add_dense(other)
  360. def _sub_sparse(self, other):
  361. return self.tocsr()._sub_sparse(other)
  362. def _sub_dense(self, other):
  363. return self.todense() - other
  364. def _rsub_dense(self, other):
  365. # note: this can't be replaced by other + (-self) for unsigned types
  366. return other - self.todense()
  367. def __add__(self, other): # self + other
  368. if isscalarlike(other):
  369. if other == 0:
  370. return self.copy()
  371. # Now we would add this scalar to every element.
  372. raise NotImplementedError('adding a nonzero scalar to a '
  373. 'sparse matrix is not supported')
  374. elif isspmatrix(other):
  375. if other.shape != self.shape:
  376. raise ValueError("inconsistent shapes")
  377. return self._add_sparse(other)
  378. elif isdense(other):
  379. other = np.broadcast_to(other, self.shape)
  380. return self._add_dense(other)
  381. else:
  382. return NotImplemented
  383. def __radd__(self,other): # other + self
  384. return self.__add__(other)
  385. def __sub__(self, other): # self - other
  386. if isscalarlike(other):
  387. if other == 0:
  388. return self.copy()
  389. raise NotImplementedError('subtracting a nonzero scalar from a '
  390. 'sparse matrix is not supported')
  391. elif isspmatrix(other):
  392. if other.shape != self.shape:
  393. raise ValueError("inconsistent shapes")
  394. return self._sub_sparse(other)
  395. elif isdense(other):
  396. other = np.broadcast_to(other, self.shape)
  397. return self._sub_dense(other)
  398. else:
  399. return NotImplemented
  400. def __rsub__(self,other): # other - self
  401. if isscalarlike(other):
  402. if other == 0:
  403. return -self.copy()
  404. raise NotImplementedError('subtracting a sparse matrix from a '
  405. 'nonzero scalar is not supported')
  406. elif isdense(other):
  407. other = np.broadcast_to(other, self.shape)
  408. return self._rsub_dense(other)
  409. else:
  410. return NotImplemented
  411. def _mul_dispatch(self, other):
  412. """`np.matrix`-compatible mul, i.e. `dot` or `NotImplemented`
  413. interpret other and call one of the following
  414. self._mul_scalar()
  415. self._mul_vector()
  416. self._mul_multivector()
  417. self._mul_sparse_matrix()
  418. """
  419. # This method has to be different from `__mul__` because it is also
  420. # called by sparse array classes via matmul, while their mul is
  421. # elementwise.
  422. M, N = self.shape
  423. if other.__class__ is np.ndarray:
  424. # Fast path for the most common case
  425. if other.shape == (N,):
  426. return self._mul_vector(other)
  427. elif other.shape == (N, 1):
  428. return self._mul_vector(other.ravel()).reshape(M, 1)
  429. elif other.ndim == 2 and other.shape[0] == N:
  430. return self._mul_multivector(other)
  431. if isscalarlike(other):
  432. # scalar value
  433. return self._mul_scalar(other)
  434. if issparse(other):
  435. if self.shape[1] != other.shape[0]:
  436. raise ValueError('dimension mismatch')
  437. return self._mul_sparse_matrix(other)
  438. # If it's a list or whatever, treat it like a matrix
  439. other_a = np.asanyarray(other)
  440. if other_a.ndim == 0 and other_a.dtype == np.object_:
  441. # Not interpretable as an array; return NotImplemented so that
  442. # other's __rmul__ can kick in if that's implemented.
  443. return NotImplemented
  444. try:
  445. other.shape
  446. except AttributeError:
  447. other = other_a
  448. if other.ndim == 1 or other.ndim == 2 and other.shape[1] == 1:
  449. # dense row or column vector
  450. if other.shape != (N,) and other.shape != (N, 1):
  451. raise ValueError('dimension mismatch')
  452. result = self._mul_vector(np.ravel(other))
  453. if isinstance(other, np.matrix):
  454. result = self._ascontainer(result)
  455. if other.ndim == 2 and other.shape[1] == 1:
  456. # If 'other' was an (nx1) column vector, reshape the result
  457. result = result.reshape(-1, 1)
  458. return result
  459. elif other.ndim == 2:
  460. ##
  461. # dense 2D array or matrix ("multivector")
  462. if other.shape[0] != self.shape[1]:
  463. raise ValueError('dimension mismatch')
  464. result = self._mul_multivector(np.asarray(other))
  465. if isinstance(other, np.matrix):
  466. result = self._ascontainer(result)
  467. return result
  468. else:
  469. raise ValueError('could not interpret dimensions')
  470. def __mul__(self, other):
  471. return self._mul_dispatch(other)
  472. # by default, use CSR for __mul__ handlers
  473. def _mul_scalar(self, other):
  474. return self.tocsr()._mul_scalar(other)
  475. def _mul_vector(self, other):
  476. return self.tocsr()._mul_vector(other)
  477. def _mul_multivector(self, other):
  478. return self.tocsr()._mul_multivector(other)
  479. def _mul_sparse_matrix(self, other):
  480. return self.tocsr()._mul_sparse_matrix(other)
  481. def _rmul_dispatch(self, other):
  482. if isscalarlike(other):
  483. return self._mul_scalar(other)
  484. else:
  485. # Don't use asarray unless we have to
  486. try:
  487. tr = other.transpose()
  488. except AttributeError:
  489. tr = np.asarray(other).transpose()
  490. ret = self.transpose()._mul_dispatch(tr)
  491. if ret is NotImplemented:
  492. return NotImplemented
  493. return ret.transpose()
  494. def __rmul__(self, other): # other * self
  495. return self._rmul_dispatch(other)
  496. #######################
  497. # matmul (@) operator #
  498. #######################
  499. def __matmul__(self, other):
  500. if isscalarlike(other):
  501. raise ValueError("Scalar operands are not allowed, "
  502. "use '*' instead")
  503. return self._mul_dispatch(other)
  504. def __rmatmul__(self, other):
  505. if isscalarlike(other):
  506. raise ValueError("Scalar operands are not allowed, "
  507. "use '*' instead")
  508. return self._rmul_dispatch(other)
  509. ####################
  510. # Other Arithmetic #
  511. ####################
  512. def _divide(self, other, true_divide=False, rdivide=False):
  513. if isscalarlike(other):
  514. if rdivide:
  515. if true_divide:
  516. return np.true_divide(other, self.todense())
  517. else:
  518. return np.divide(other, self.todense())
  519. if true_divide and np.can_cast(self.dtype, np.float_):
  520. return self.astype(np.float_)._mul_scalar(1./other)
  521. else:
  522. r = self._mul_scalar(1./other)
  523. scalar_dtype = np.asarray(other).dtype
  524. if (np.issubdtype(self.dtype, np.integer) and
  525. np.issubdtype(scalar_dtype, np.integer)):
  526. return r.astype(self.dtype)
  527. else:
  528. return r
  529. elif isdense(other):
  530. if not rdivide:
  531. if true_divide:
  532. return np.true_divide(self.todense(), other)
  533. else:
  534. return np.divide(self.todense(), other)
  535. else:
  536. if true_divide:
  537. return np.true_divide(other, self.todense())
  538. else:
  539. return np.divide(other, self.todense())
  540. elif isspmatrix(other):
  541. if rdivide:
  542. return other._divide(self, true_divide, rdivide=False)
  543. self_csr = self.tocsr()
  544. if true_divide and np.can_cast(self.dtype, np.float_):
  545. return self_csr.astype(np.float_)._divide_sparse(other)
  546. else:
  547. return self_csr._divide_sparse(other)
  548. else:
  549. return NotImplemented
  550. def __truediv__(self, other):
  551. return self._divide(other, true_divide=True)
  552. def __div__(self, other):
  553. # Always do true division
  554. return self._divide(other, true_divide=True)
  555. def __rtruediv__(self, other):
  556. # Implementing this as the inverse would be too magical -- bail out
  557. return NotImplemented
  558. def __rdiv__(self, other):
  559. # Implementing this as the inverse would be too magical -- bail out
  560. return NotImplemented
  561. def __neg__(self):
  562. return -self.tocsr()
  563. def __iadd__(self, other):
  564. return NotImplemented
  565. def __isub__(self, other):
  566. return NotImplemented
  567. def __imul__(self, other):
  568. return NotImplemented
  569. def __idiv__(self, other):
  570. return self.__itruediv__(other)
  571. def __itruediv__(self, other):
  572. return NotImplemented
  573. def __pow__(self, other):
  574. M, N = self.shape[0], self.shape[1]
  575. if M != N:
  576. raise TypeError('matrix is not square')
  577. if isintlike(other):
  578. other = int(other)
  579. if other < 0:
  580. raise ValueError('exponent must be >= 0')
  581. if other == 0:
  582. from ._construct import eye
  583. E = eye(M, dtype=self.dtype)
  584. if self._is_array:
  585. from ._arrays import dia_array
  586. E = dia_array(E)
  587. return E
  588. elif other == 1:
  589. return self.copy()
  590. else:
  591. tmp = self.__pow__(other//2)
  592. if (other % 2):
  593. return self @ tmp @ tmp
  594. else:
  595. return tmp @ tmp
  596. elif isscalarlike(other):
  597. raise ValueError('exponent must be an integer')
  598. else:
  599. return NotImplemented
  600. def __getattr__(self, attr):
  601. if attr == 'A':
  602. if self._is_array:
  603. warn(np.VisibleDeprecationWarning(
  604. "Please use `.todense()` instead"
  605. ))
  606. return self.toarray()
  607. elif attr == 'T':
  608. return self.transpose()
  609. elif attr == 'H':
  610. if self._is_array:
  611. warn(np.VisibleDeprecationWarning(
  612. "Please use `.conj().T` instead"
  613. ))
  614. return self.getH()
  615. elif attr == 'real':
  616. return self._real()
  617. elif attr == 'imag':
  618. return self._imag()
  619. elif attr == 'size':
  620. return self.getnnz()
  621. else:
  622. raise AttributeError(attr + " not found")
  623. def transpose(self, axes=None, copy=False):
  624. """
  625. Reverses the dimensions of the sparse matrix.
  626. Parameters
  627. ----------
  628. axes : None, optional
  629. This argument is in the signature *solely* for NumPy
  630. compatibility reasons. Do not pass in anything except
  631. for the default value.
  632. copy : bool, optional
  633. Indicates whether or not attributes of `self` should be
  634. copied whenever possible. The degree to which attributes
  635. are copied varies depending on the type of sparse matrix
  636. being used.
  637. Returns
  638. -------
  639. p : `self` with the dimensions reversed.
  640. See Also
  641. --------
  642. numpy.matrix.transpose : NumPy's implementation of 'transpose'
  643. for matrices
  644. """
  645. return self.tocsr(copy=copy).transpose(axes=axes, copy=False)
  646. def conj(self, copy=True):
  647. """Element-wise complex conjugation.
  648. If the matrix is of non-complex data type and `copy` is False,
  649. this method does nothing and the data is not copied.
  650. Parameters
  651. ----------
  652. copy : bool, optional
  653. If True, the result is guaranteed to not share data with self.
  654. Returns
  655. -------
  656. A : The element-wise complex conjugate.
  657. """
  658. if np.issubdtype(self.dtype, np.complexfloating):
  659. return self.tocsr(copy=copy).conj(copy=False)
  660. elif copy:
  661. return self.copy()
  662. else:
  663. return self
  664. def conjugate(self, copy=True):
  665. return self.conj(copy=copy)
  666. conjugate.__doc__ = conj.__doc__
  667. # Renamed conjtranspose() -> getH() for compatibility with dense matrices
  668. def getH(self):
  669. """Return the Hermitian transpose of this matrix.
  670. See Also
  671. --------
  672. numpy.matrix.getH : NumPy's implementation of `getH` for matrices
  673. """
  674. return self.transpose().conj()
  675. def _real(self):
  676. return self.tocsr()._real()
  677. def _imag(self):
  678. return self.tocsr()._imag()
  679. def nonzero(self):
  680. """nonzero indices
  681. Returns a tuple of arrays (row,col) containing the indices
  682. of the non-zero elements of the matrix.
  683. Examples
  684. --------
  685. >>> from scipy.sparse import csr_matrix
  686. >>> A = csr_matrix([[1,2,0],[0,0,3],[4,0,5]])
  687. >>> A.nonzero()
  688. (array([0, 0, 1, 2, 2]), array([0, 1, 2, 0, 2]))
  689. """
  690. # convert to COOrdinate format
  691. A = self.tocoo()
  692. nz_mask = A.data != 0
  693. return (A.row[nz_mask], A.col[nz_mask])
  694. def getcol(self, j):
  695. """Returns a copy of column j of the matrix, as an (m x 1) sparse
  696. matrix (column vector).
  697. """
  698. # Spmatrix subclasses should override this method for efficiency.
  699. # Post-multiply by a (n x 1) column vector 'a' containing all zeros
  700. # except for a_j = 1
  701. n = self.shape[1]
  702. if j < 0:
  703. j += n
  704. if j < 0 or j >= n:
  705. raise IndexError("index out of bounds")
  706. col_selector = self._csc_container(([1], [[j], [0]]),
  707. shape=(n, 1), dtype=self.dtype)
  708. return self @ col_selector
  709. def getrow(self, i):
  710. """Returns a copy of row i of the matrix, as a (1 x n) sparse
  711. matrix (row vector).
  712. """
  713. # Spmatrix subclasses should override this method for efficiency.
  714. # Pre-multiply by a (1 x m) row vector 'a' containing all zeros
  715. # except for a_i = 1
  716. m = self.shape[0]
  717. if i < 0:
  718. i += m
  719. if i < 0 or i >= m:
  720. raise IndexError("index out of bounds")
  721. row_selector = self._csr_container(([1], [[0], [i]]),
  722. shape=(1, m), dtype=self.dtype)
  723. return row_selector @ self
  724. # The following dunder methods cannot be implemented.
  725. #
  726. # def __array__(self):
  727. # # Sparse matrices rely on NumPy wrapping them in object arrays under
  728. # # the hood to make unary ufuncs work on them. So we cannot raise
  729. # # TypeError here - which would be handy to not give users object
  730. # # arrays they probably don't want (they're looking for `.toarray()`).
  731. # #
  732. # # Conversion with `toarray()` would also break things because of the
  733. # # behavior discussed above, plus we want to avoid densification by
  734. # # accident because that can too easily blow up memory.
  735. #
  736. # def __array_ufunc__(self):
  737. # # We cannot implement __array_ufunc__ due to mismatching semantics.
  738. # # See gh-7707 and gh-7349 for details.
  739. #
  740. # def __array_function__(self):
  741. # # We cannot implement __array_function__ due to mismatching semantics.
  742. # # See gh-10362 for details.
  743. def todense(self, order=None, out=None):
  744. """
  745. Return a dense matrix representation of this matrix.
  746. Parameters
  747. ----------
  748. order : {'C', 'F'}, optional
  749. Whether to store multi-dimensional data in C (row-major)
  750. or Fortran (column-major) order in memory. The default
  751. is 'None', which provides no ordering guarantees.
  752. Cannot be specified in conjunction with the `out`
  753. argument.
  754. out : ndarray, 2-D, optional
  755. If specified, uses this array (or `numpy.matrix`) as the
  756. output buffer instead of allocating a new array to
  757. return. The provided array must have the same shape and
  758. dtype as the sparse matrix on which you are calling the
  759. method.
  760. Returns
  761. -------
  762. arr : numpy.matrix, 2-D
  763. A NumPy matrix object with the same shape and containing
  764. the same data represented by the sparse matrix, with the
  765. requested memory order. If `out` was passed and was an
  766. array (rather than a `numpy.matrix`), it will be filled
  767. with the appropriate values and returned wrapped in a
  768. `numpy.matrix` object that shares the same memory.
  769. """
  770. return self._ascontainer(self.toarray(order=order, out=out))
  771. def toarray(self, order=None, out=None):
  772. """
  773. Return a dense ndarray representation of this matrix.
  774. Parameters
  775. ----------
  776. order : {'C', 'F'}, optional
  777. Whether to store multidimensional data in C (row-major)
  778. or Fortran (column-major) order in memory. The default
  779. is 'None', which provides no ordering guarantees.
  780. Cannot be specified in conjunction with the `out`
  781. argument.
  782. out : ndarray, 2-D, optional
  783. If specified, uses this array as the output buffer
  784. instead of allocating a new array to return. The provided
  785. array must have the same shape and dtype as the sparse
  786. matrix on which you are calling the method. For most
  787. sparse types, `out` is required to be memory contiguous
  788. (either C or Fortran ordered).
  789. Returns
  790. -------
  791. arr : ndarray, 2-D
  792. An array with the same shape and containing the same
  793. data represented by the sparse matrix, with the requested
  794. memory order. If `out` was passed, the same object is
  795. returned after being modified in-place to contain the
  796. appropriate values.
  797. """
  798. return self.tocoo(copy=False).toarray(order=order, out=out)
  799. # Any sparse matrix format deriving from spmatrix must define one of
  800. # tocsr or tocoo. The other conversion methods may be implemented for
  801. # efficiency, but are not required.
  802. def tocsr(self, copy=False):
  803. """Convert this matrix to Compressed Sparse Row format.
  804. With copy=False, the data/indices may be shared between this matrix and
  805. the resultant csr_matrix.
  806. """
  807. return self.tocoo(copy=copy).tocsr(copy=False)
  808. def todok(self, copy=False):
  809. """Convert this matrix to Dictionary Of Keys format.
  810. With copy=False, the data/indices may be shared between this matrix and
  811. the resultant dok_matrix.
  812. """
  813. return self.tocoo(copy=copy).todok(copy=False)
  814. def tocoo(self, copy=False):
  815. """Convert this matrix to COOrdinate format.
  816. With copy=False, the data/indices may be shared between this matrix and
  817. the resultant coo_matrix.
  818. """
  819. return self.tocsr(copy=False).tocoo(copy=copy)
  820. def tolil(self, copy=False):
  821. """Convert this matrix to List of Lists format.
  822. With copy=False, the data/indices may be shared between this matrix and
  823. the resultant lil_matrix.
  824. """
  825. return self.tocsr(copy=False).tolil(copy=copy)
  826. def todia(self, copy=False):
  827. """Convert this matrix to sparse DIAgonal format.
  828. With copy=False, the data/indices may be shared between this matrix and
  829. the resultant dia_matrix.
  830. """
  831. return self.tocoo(copy=copy).todia(copy=False)
  832. def tobsr(self, blocksize=None, copy=False):
  833. """Convert this matrix to Block Sparse Row format.
  834. With copy=False, the data/indices may be shared between this matrix and
  835. the resultant bsr_matrix.
  836. When blocksize=(R, C) is provided, it will be used for construction of
  837. the bsr_matrix.
  838. """
  839. return self.tocsr(copy=False).tobsr(blocksize=blocksize, copy=copy)
  840. def tocsc(self, copy=False):
  841. """Convert this matrix to Compressed Sparse Column format.
  842. With copy=False, the data/indices may be shared between this matrix and
  843. the resultant csc_matrix.
  844. """
  845. return self.tocsr(copy=copy).tocsc(copy=False)
  846. def copy(self):
  847. """Returns a copy of this matrix.
  848. No data/indices will be shared between the returned value and current
  849. matrix.
  850. """
  851. return self.__class__(self, copy=True)
  852. def sum(self, axis=None, dtype=None, out=None):
  853. """
  854. Sum the matrix elements over a given axis.
  855. Parameters
  856. ----------
  857. axis : {-2, -1, 0, 1, None} optional
  858. Axis along which the sum is computed. The default is to
  859. compute the sum of all the matrix elements, returning a scalar
  860. (i.e., `axis` = `None`).
  861. dtype : dtype, optional
  862. The type of the returned matrix and of the accumulator in which
  863. the elements are summed. The dtype of `a` is used by default
  864. unless `a` has an integer dtype of less precision than the default
  865. platform integer. In that case, if `a` is signed then the platform
  866. integer is used while if `a` is unsigned then an unsigned integer
  867. of the same precision as the platform integer is used.
  868. .. versionadded:: 0.18.0
  869. out : np.matrix, optional
  870. Alternative output matrix in which to place the result. It must
  871. have the same shape as the expected output, but the type of the
  872. output values will be cast if necessary.
  873. .. versionadded:: 0.18.0
  874. Returns
  875. -------
  876. sum_along_axis : np.matrix
  877. A matrix with the same shape as `self`, with the specified
  878. axis removed.
  879. See Also
  880. --------
  881. numpy.matrix.sum : NumPy's implementation of 'sum' for matrices
  882. """
  883. validateaxis(axis)
  884. # We use multiplication by a matrix of ones to achieve this.
  885. # For some sparse matrix formats more efficient methods are
  886. # possible -- these should override this function.
  887. m, n = self.shape
  888. # Mimic numpy's casting.
  889. res_dtype = get_sum_dtype(self.dtype)
  890. if axis is None:
  891. # sum over rows and columns
  892. return (
  893. self @ self._ascontainer(np.ones((n, 1), dtype=res_dtype))
  894. ).sum(dtype=dtype, out=out)
  895. if axis < 0:
  896. axis += 2
  897. # axis = 0 or 1 now
  898. if axis == 0:
  899. # sum over columns
  900. ret = self._ascontainer(
  901. np.ones((1, m), dtype=res_dtype)
  902. ) @ self
  903. else:
  904. # sum over rows
  905. ret = self @ self._ascontainer(
  906. np.ones((n, 1), dtype=res_dtype)
  907. )
  908. if out is not None and out.shape != ret.shape:
  909. raise ValueError("dimensions do not match")
  910. return ret.sum(axis=axis, dtype=dtype, out=out)
  911. def mean(self, axis=None, dtype=None, out=None):
  912. """
  913. Compute the arithmetic mean along the specified axis.
  914. Returns the average of the matrix elements. The average is taken
  915. over all elements in the matrix by default, otherwise over the
  916. specified axis. `float64` intermediate and return values are used
  917. for integer inputs.
  918. Parameters
  919. ----------
  920. axis : {-2, -1, 0, 1, None} optional
  921. Axis along which the mean is computed. The default is to compute
  922. the mean of all elements in the matrix (i.e., `axis` = `None`).
  923. dtype : data-type, optional
  924. Type to use in computing the mean. For integer inputs, the default
  925. is `float64`; for floating point inputs, it is the same as the
  926. input dtype.
  927. .. versionadded:: 0.18.0
  928. out : np.matrix, optional
  929. Alternative output matrix in which to place the result. It must
  930. have the same shape as the expected output, but the type of the
  931. output values will be cast if necessary.
  932. .. versionadded:: 0.18.0
  933. Returns
  934. -------
  935. m : np.matrix
  936. See Also
  937. --------
  938. numpy.matrix.mean : NumPy's implementation of 'mean' for matrices
  939. """
  940. def _is_integral(dtype):
  941. return (np.issubdtype(dtype, np.integer) or
  942. np.issubdtype(dtype, np.bool_))
  943. validateaxis(axis)
  944. res_dtype = self.dtype.type
  945. integral = _is_integral(self.dtype)
  946. # output dtype
  947. if dtype is None:
  948. if integral:
  949. res_dtype = np.float64
  950. else:
  951. res_dtype = np.dtype(dtype).type
  952. # intermediate dtype for summation
  953. inter_dtype = np.float64 if integral else res_dtype
  954. inter_self = self.astype(inter_dtype)
  955. if axis is None:
  956. return (inter_self / np.array(
  957. self.shape[0] * self.shape[1]))\
  958. .sum(dtype=res_dtype, out=out)
  959. if axis < 0:
  960. axis += 2
  961. # axis = 0 or 1 now
  962. if axis == 0:
  963. return (inter_self * (1.0 / self.shape[0])).sum(
  964. axis=0, dtype=res_dtype, out=out)
  965. else:
  966. return (inter_self * (1.0 / self.shape[1])).sum(
  967. axis=1, dtype=res_dtype, out=out)
  968. def diagonal(self, k=0):
  969. """Returns the kth diagonal of the matrix.
  970. Parameters
  971. ----------
  972. k : int, optional
  973. Which diagonal to get, corresponding to elements a[i, i+k].
  974. Default: 0 (the main diagonal).
  975. .. versionadded:: 1.0
  976. See also
  977. --------
  978. numpy.diagonal : Equivalent numpy function.
  979. Examples
  980. --------
  981. >>> from scipy.sparse import csr_matrix
  982. >>> A = csr_matrix([[1, 2, 0], [0, 0, 3], [4, 0, 5]])
  983. >>> A.diagonal()
  984. array([1, 0, 5])
  985. >>> A.diagonal(k=1)
  986. array([2, 3])
  987. """
  988. return self.tocsr().diagonal(k=k)
  989. def trace(self, offset=0):
  990. """Returns the sum along diagonals of the sparse matrix.
  991. Parameters
  992. ----------
  993. offset : int, optional
  994. Which diagonal to get, corresponding to elements a[i, i+offset].
  995. Default: 0 (the main diagonal).
  996. """
  997. return self.diagonal(k=offset).sum()
  998. def setdiag(self, values, k=0):
  999. """
  1000. Set diagonal or off-diagonal elements of the array.
  1001. Parameters
  1002. ----------
  1003. values : array_like
  1004. New values of the diagonal elements.
  1005. Values may have any length. If the diagonal is longer than values,
  1006. then the remaining diagonal entries will not be set. If values are
  1007. longer than the diagonal, then the remaining values are ignored.
  1008. If a scalar value is given, all of the diagonal is set to it.
  1009. k : int, optional
  1010. Which off-diagonal to set, corresponding to elements a[i,i+k].
  1011. Default: 0 (the main diagonal).
  1012. """
  1013. M, N = self.shape
  1014. if (k > 0 and k >= N) or (k < 0 and -k >= M):
  1015. raise ValueError("k exceeds matrix dimensions")
  1016. self._setdiag(np.asarray(values), k)
  1017. def _setdiag(self, values, k):
  1018. M, N = self.shape
  1019. if k < 0:
  1020. if values.ndim == 0:
  1021. # broadcast
  1022. max_index = min(M+k, N)
  1023. for i in range(max_index):
  1024. self[i - k, i] = values
  1025. else:
  1026. max_index = min(M+k, N, len(values))
  1027. if max_index <= 0:
  1028. return
  1029. for i, v in enumerate(values[:max_index]):
  1030. self[i - k, i] = v
  1031. else:
  1032. if values.ndim == 0:
  1033. # broadcast
  1034. max_index = min(M, N-k)
  1035. for i in range(max_index):
  1036. self[i, i + k] = values
  1037. else:
  1038. max_index = min(M, N-k, len(values))
  1039. if max_index <= 0:
  1040. return
  1041. for i, v in enumerate(values[:max_index]):
  1042. self[i, i + k] = v
  1043. def _process_toarray_args(self, order, out):
  1044. if out is not None:
  1045. if order is not None:
  1046. raise ValueError('order cannot be specified if out '
  1047. 'is not None')
  1048. if out.shape != self.shape or out.dtype != self.dtype:
  1049. raise ValueError('out array must be same dtype and shape as '
  1050. 'sparse matrix')
  1051. out[...] = 0.
  1052. return out
  1053. else:
  1054. return np.zeros(self.shape, dtype=self.dtype, order=order)
  1055. def isspmatrix(x):
  1056. """Is x of a sparse matrix type?
  1057. Parameters
  1058. ----------
  1059. x
  1060. object to check for being a sparse matrix
  1061. Returns
  1062. -------
  1063. bool
  1064. True if x is a sparse matrix, False otherwise
  1065. Notes
  1066. -----
  1067. issparse and isspmatrix are aliases for the same function.
  1068. Examples
  1069. --------
  1070. >>> from scipy.sparse import csr_matrix, isspmatrix
  1071. >>> isspmatrix(csr_matrix([[5]]))
  1072. True
  1073. >>> from scipy.sparse import isspmatrix
  1074. >>> isspmatrix(5)
  1075. False
  1076. """
  1077. return isinstance(x, spmatrix)
  1078. issparse = isspmatrix