_construct.py 30 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947
  1. """Functions to construct sparse matrices
  2. """
  3. __docformat__ = "restructuredtext en"
  4. __all__ = ['spdiags', 'eye', 'identity', 'kron', 'kronsum',
  5. 'hstack', 'vstack', 'bmat', 'rand', 'random', 'diags', 'block_diag']
  6. import numbers
  7. from functools import partial
  8. import numpy as np
  9. from scipy._lib._util import check_random_state, rng_integers
  10. from ._sputils import upcast, get_index_dtype, isscalarlike
  11. from ._sparsetools import csr_hstack
  12. from ._csr import csr_matrix
  13. from ._csc import csc_matrix
  14. from ._bsr import bsr_matrix
  15. from ._coo import coo_matrix
  16. from ._dia import dia_matrix
  17. from ._base import issparse
  18. def spdiags(data, diags, m=None, n=None, format=None):
  19. """
  20. Return a sparse matrix from diagonals.
  21. Parameters
  22. ----------
  23. data : array_like
  24. Matrix diagonals stored row-wise
  25. diags : sequence of int or an int
  26. Diagonals to set:
  27. * k = 0 the main diagonal
  28. * k > 0 the kth upper diagonal
  29. * k < 0 the kth lower diagonal
  30. m, n : int, tuple, optional
  31. Shape of the result. If `n` is None and `m` is a given tuple,
  32. the shape is this tuple. If omitted, the matrix is square and
  33. its shape is len(data[0]).
  34. format : str, optional
  35. Format of the result. By default (format=None) an appropriate sparse
  36. matrix format is returned. This choice is subject to change.
  37. See Also
  38. --------
  39. diags : more convenient form of this function
  40. dia_matrix : the sparse DIAgonal format.
  41. Examples
  42. --------
  43. >>> import numpy as np
  44. >>> from scipy.sparse import spdiags
  45. >>> data = np.array([[1, 2, 3, 4], [1, 2, 3, 4], [1, 2, 3, 4]])
  46. >>> diags = np.array([0, -1, 2])
  47. >>> spdiags(data, diags, 4, 4).toarray()
  48. array([[1, 0, 3, 0],
  49. [1, 2, 0, 4],
  50. [0, 2, 3, 0],
  51. [0, 0, 3, 4]])
  52. """
  53. if m is None and n is None:
  54. m = n = len(data[0])
  55. elif n is None:
  56. m, n = m
  57. return dia_matrix((data, diags), shape=(m, n)).asformat(format)
  58. def diags(diagonals, offsets=0, shape=None, format=None, dtype=None):
  59. """
  60. Construct a sparse matrix from diagonals.
  61. Parameters
  62. ----------
  63. diagonals : sequence of array_like
  64. Sequence of arrays containing the matrix diagonals,
  65. corresponding to `offsets`.
  66. offsets : sequence of int or an int, optional
  67. Diagonals to set:
  68. - k = 0 the main diagonal (default)
  69. - k > 0 the kth upper diagonal
  70. - k < 0 the kth lower diagonal
  71. shape : tuple of int, optional
  72. Shape of the result. If omitted, a square matrix large enough
  73. to contain the diagonals is returned.
  74. format : {"dia", "csr", "csc", "lil", ...}, optional
  75. Matrix format of the result. By default (format=None) an
  76. appropriate sparse matrix format is returned. This choice is
  77. subject to change.
  78. dtype : dtype, optional
  79. Data type of the matrix.
  80. See Also
  81. --------
  82. spdiags : construct matrix from diagonals
  83. Notes
  84. -----
  85. This function differs from `spdiags` in the way it handles
  86. off-diagonals.
  87. The result from `diags` is the sparse equivalent of::
  88. np.diag(diagonals[0], offsets[0])
  89. + ...
  90. + np.diag(diagonals[k], offsets[k])
  91. Repeated diagonal offsets are disallowed.
  92. .. versionadded:: 0.11
  93. Examples
  94. --------
  95. >>> from scipy.sparse import diags
  96. >>> diagonals = [[1, 2, 3, 4], [1, 2, 3], [1, 2]]
  97. >>> diags(diagonals, [0, -1, 2]).toarray()
  98. array([[1, 0, 1, 0],
  99. [1, 2, 0, 2],
  100. [0, 2, 3, 0],
  101. [0, 0, 3, 4]])
  102. Broadcasting of scalars is supported (but shape needs to be
  103. specified):
  104. >>> diags([1, -2, 1], [-1, 0, 1], shape=(4, 4)).toarray()
  105. array([[-2., 1., 0., 0.],
  106. [ 1., -2., 1., 0.],
  107. [ 0., 1., -2., 1.],
  108. [ 0., 0., 1., -2.]])
  109. If only one diagonal is wanted (as in `numpy.diag`), the following
  110. works as well:
  111. >>> diags([1, 2, 3], 1).toarray()
  112. array([[ 0., 1., 0., 0.],
  113. [ 0., 0., 2., 0.],
  114. [ 0., 0., 0., 3.],
  115. [ 0., 0., 0., 0.]])
  116. """
  117. # if offsets is not a sequence, assume that there's only one diagonal
  118. if isscalarlike(offsets):
  119. # now check that there's actually only one diagonal
  120. if len(diagonals) == 0 or isscalarlike(diagonals[0]):
  121. diagonals = [np.atleast_1d(diagonals)]
  122. else:
  123. raise ValueError("Different number of diagonals and offsets.")
  124. else:
  125. diagonals = list(map(np.atleast_1d, diagonals))
  126. offsets = np.atleast_1d(offsets)
  127. # Basic check
  128. if len(diagonals) != len(offsets):
  129. raise ValueError("Different number of diagonals and offsets.")
  130. # Determine shape, if omitted
  131. if shape is None:
  132. m = len(diagonals[0]) + abs(int(offsets[0]))
  133. shape = (m, m)
  134. # Determine data type, if omitted
  135. if dtype is None:
  136. dtype = np.common_type(*diagonals)
  137. # Construct data array
  138. m, n = shape
  139. M = max([min(m + offset, n - offset) + max(0, offset)
  140. for offset in offsets])
  141. M = max(0, M)
  142. data_arr = np.zeros((len(offsets), M), dtype=dtype)
  143. K = min(m, n)
  144. for j, diagonal in enumerate(diagonals):
  145. offset = offsets[j]
  146. k = max(0, offset)
  147. length = min(m + offset, n - offset, K)
  148. if length < 0:
  149. raise ValueError("Offset %d (index %d) out of bounds" % (offset, j))
  150. try:
  151. data_arr[j, k:k+length] = diagonal[...,:length]
  152. except ValueError as e:
  153. if len(diagonal) != length and len(diagonal) != 1:
  154. raise ValueError(
  155. "Diagonal length (index %d: %d at offset %d) does not "
  156. "agree with matrix size (%d, %d)." % (
  157. j, len(diagonal), offset, m, n)) from e
  158. raise
  159. return dia_matrix((data_arr, offsets), shape=(m, n)).asformat(format)
  160. def identity(n, dtype='d', format=None):
  161. """Identity matrix in sparse format
  162. Returns an identity matrix with shape (n,n) using a given
  163. sparse format and dtype.
  164. Parameters
  165. ----------
  166. n : int
  167. Shape of the identity matrix.
  168. dtype : dtype, optional
  169. Data type of the matrix
  170. format : str, optional
  171. Sparse format of the result, e.g., format="csr", etc.
  172. Examples
  173. --------
  174. >>> from scipy.sparse import identity
  175. >>> identity(3).toarray()
  176. array([[ 1., 0., 0.],
  177. [ 0., 1., 0.],
  178. [ 0., 0., 1.]])
  179. >>> identity(3, dtype='int8', format='dia')
  180. <3x3 sparse matrix of type '<class 'numpy.int8'>'
  181. with 3 stored elements (1 diagonals) in DIAgonal format>
  182. """
  183. return eye(n, n, dtype=dtype, format=format)
  184. def eye(m, n=None, k=0, dtype=float, format=None):
  185. """Sparse matrix with ones on diagonal
  186. Returns a sparse (m x n) matrix where the kth diagonal
  187. is all ones and everything else is zeros.
  188. Parameters
  189. ----------
  190. m : int
  191. Number of rows in the matrix.
  192. n : int, optional
  193. Number of columns. Default: `m`.
  194. k : int, optional
  195. Diagonal to place ones on. Default: 0 (main diagonal).
  196. dtype : dtype, optional
  197. Data type of the matrix.
  198. format : str, optional
  199. Sparse format of the result, e.g., format="csr", etc.
  200. Examples
  201. --------
  202. >>> import numpy as np
  203. >>> from scipy import sparse
  204. >>> sparse.eye(3).toarray()
  205. array([[ 1., 0., 0.],
  206. [ 0., 1., 0.],
  207. [ 0., 0., 1.]])
  208. >>> sparse.eye(3, dtype=np.int8)
  209. <3x3 sparse matrix of type '<class 'numpy.int8'>'
  210. with 3 stored elements (1 diagonals) in DIAgonal format>
  211. """
  212. if n is None:
  213. n = m
  214. m,n = int(m),int(n)
  215. if m == n and k == 0:
  216. # fast branch for special formats
  217. if format in ['csr', 'csc']:
  218. idx_dtype = get_index_dtype(maxval=n)
  219. indptr = np.arange(n+1, dtype=idx_dtype)
  220. indices = np.arange(n, dtype=idx_dtype)
  221. data = np.ones(n, dtype=dtype)
  222. cls = {'csr': csr_matrix, 'csc': csc_matrix}[format]
  223. return cls((data,indices,indptr),(n,n))
  224. elif format == 'coo':
  225. idx_dtype = get_index_dtype(maxval=n)
  226. row = np.arange(n, dtype=idx_dtype)
  227. col = np.arange(n, dtype=idx_dtype)
  228. data = np.ones(n, dtype=dtype)
  229. return coo_matrix((data, (row, col)), (n, n))
  230. diags = np.ones((1, max(0, min(m + k, n))), dtype=dtype)
  231. return spdiags(diags, k, m, n).asformat(format)
  232. def kron(A, B, format=None):
  233. """kronecker product of sparse matrices A and B
  234. Parameters
  235. ----------
  236. A : sparse or dense matrix
  237. first matrix of the product
  238. B : sparse or dense matrix
  239. second matrix of the product
  240. format : str, optional
  241. format of the result (e.g. "csr")
  242. Returns
  243. -------
  244. kronecker product in a sparse matrix format
  245. Examples
  246. --------
  247. >>> import numpy as np
  248. >>> from scipy import sparse
  249. >>> A = sparse.csr_matrix(np.array([[0, 2], [5, 0]]))
  250. >>> B = sparse.csr_matrix(np.array([[1, 2], [3, 4]]))
  251. >>> sparse.kron(A, B).toarray()
  252. array([[ 0, 0, 2, 4],
  253. [ 0, 0, 6, 8],
  254. [ 5, 10, 0, 0],
  255. [15, 20, 0, 0]])
  256. >>> sparse.kron(A, [[1, 2], [3, 4]]).toarray()
  257. array([[ 0, 0, 2, 4],
  258. [ 0, 0, 6, 8],
  259. [ 5, 10, 0, 0],
  260. [15, 20, 0, 0]])
  261. """
  262. B = coo_matrix(B)
  263. if (format is None or format == "bsr") and 2*B.nnz >= B.shape[0] * B.shape[1]:
  264. # B is fairly dense, use BSR
  265. A = csr_matrix(A,copy=True)
  266. output_shape = (A.shape[0]*B.shape[0], A.shape[1]*B.shape[1])
  267. if A.nnz == 0 or B.nnz == 0:
  268. # kronecker product is the zero matrix
  269. return coo_matrix(output_shape).asformat(format)
  270. B = B.toarray()
  271. data = A.data.repeat(B.size).reshape(-1,B.shape[0],B.shape[1])
  272. data = data * B
  273. return bsr_matrix((data,A.indices,A.indptr), shape=output_shape)
  274. else:
  275. # use COO
  276. A = coo_matrix(A)
  277. output_shape = (A.shape[0]*B.shape[0], A.shape[1]*B.shape[1])
  278. if A.nnz == 0 or B.nnz == 0:
  279. # kronecker product is the zero matrix
  280. return coo_matrix(output_shape).asformat(format)
  281. # expand entries of a into blocks
  282. row = A.row.repeat(B.nnz)
  283. col = A.col.repeat(B.nnz)
  284. data = A.data.repeat(B.nnz)
  285. if max(A.shape[0]*B.shape[0], A.shape[1]*B.shape[1]) > np.iinfo('int32').max:
  286. row = row.astype(np.int64)
  287. col = col.astype(np.int64)
  288. row *= B.shape[0]
  289. col *= B.shape[1]
  290. # increment block indices
  291. row,col = row.reshape(-1,B.nnz),col.reshape(-1,B.nnz)
  292. row += B.row
  293. col += B.col
  294. row,col = row.reshape(-1),col.reshape(-1)
  295. # compute block entries
  296. data = data.reshape(-1,B.nnz) * B.data
  297. data = data.reshape(-1)
  298. return coo_matrix((data,(row,col)), shape=output_shape).asformat(format)
  299. def kronsum(A, B, format=None):
  300. """kronecker sum of sparse matrices A and B
  301. Kronecker sum of two sparse matrices is a sum of two Kronecker
  302. products kron(I_n,A) + kron(B,I_m) where A has shape (m,m)
  303. and B has shape (n,n) and I_m and I_n are identity matrices
  304. of shape (m,m) and (n,n), respectively.
  305. Parameters
  306. ----------
  307. A
  308. square matrix
  309. B
  310. square matrix
  311. format : str
  312. format of the result (e.g. "csr")
  313. Returns
  314. -------
  315. kronecker sum in a sparse matrix format
  316. Examples
  317. --------
  318. """
  319. A = coo_matrix(A)
  320. B = coo_matrix(B)
  321. if A.shape[0] != A.shape[1]:
  322. raise ValueError('A is not square')
  323. if B.shape[0] != B.shape[1]:
  324. raise ValueError('B is not square')
  325. dtype = upcast(A.dtype, B.dtype)
  326. L = kron(eye(B.shape[0],dtype=dtype), A, format=format)
  327. R = kron(B, eye(A.shape[0],dtype=dtype), format=format)
  328. return (L+R).asformat(format) # since L + R is not always same format
  329. def _compressed_sparse_stack(blocks, axis):
  330. """
  331. Stacking fast path for CSR/CSC matrices
  332. (i) vstack for CSR, (ii) hstack for CSC.
  333. """
  334. other_axis = 1 if axis == 0 else 0
  335. data = np.concatenate([b.data for b in blocks])
  336. constant_dim = blocks[0].shape[other_axis]
  337. idx_dtype = get_index_dtype(arrays=[b.indptr for b in blocks],
  338. maxval=max(data.size, constant_dim))
  339. indices = np.empty(data.size, dtype=idx_dtype)
  340. indptr = np.empty(sum(b.shape[axis] for b in blocks) + 1, dtype=idx_dtype)
  341. last_indptr = idx_dtype(0)
  342. sum_dim = 0
  343. sum_indices = 0
  344. for b in blocks:
  345. if b.shape[other_axis] != constant_dim:
  346. raise ValueError(f'incompatible dimensions for axis {other_axis}')
  347. indices[sum_indices:sum_indices+b.indices.size] = b.indices
  348. sum_indices += b.indices.size
  349. idxs = slice(sum_dim, sum_dim + b.shape[axis])
  350. indptr[idxs] = b.indptr[:-1]
  351. indptr[idxs] += last_indptr
  352. sum_dim += b.shape[axis]
  353. last_indptr += b.indptr[-1]
  354. indptr[-1] = last_indptr
  355. if axis == 0:
  356. return csr_matrix((data, indices, indptr),
  357. shape=(sum_dim, constant_dim))
  358. else:
  359. return csc_matrix((data, indices, indptr),
  360. shape=(constant_dim, sum_dim))
  361. def _stack_along_minor_axis(blocks, axis):
  362. """
  363. Stacking fast path for CSR/CSC matrices along the minor axis
  364. (i) hstack for CSR, (ii) vstack for CSC.
  365. """
  366. n_blocks = len(blocks)
  367. if n_blocks == 0:
  368. raise ValueError('Missing block matrices')
  369. if n_blocks == 1:
  370. return blocks[0]
  371. # check for incompatible dimensions
  372. other_axis = 1 if axis == 0 else 0
  373. other_axis_dims = set(b.shape[other_axis] for b in blocks)
  374. if len(other_axis_dims) > 1:
  375. raise ValueError(f'Mismatching dimensions along axis {other_axis}: '
  376. f'{other_axis_dims}')
  377. constant_dim, = other_axis_dims
  378. # Do the stacking
  379. indptr_list = [b.indptr for b in blocks]
  380. data_cat = np.concatenate([b.data for b in blocks])
  381. # Need to check if any indices/indptr, would be too large post-
  382. # concatenation for np.int32:
  383. # - The max value of indices is the output array's stacking-axis length - 1
  384. # - The max value in indptr is the number of non-zero entries. This is
  385. # exceedingly unlikely to require int64, but is checked out of an
  386. # abundance of caution.
  387. sum_dim = sum(b.shape[axis] for b in blocks)
  388. nnz = sum(len(b.indices) for b in blocks)
  389. idx_dtype = get_index_dtype(maxval=max(sum_dim - 1, nnz))
  390. stack_dim_cat = np.array([b.shape[axis] for b in blocks], dtype=idx_dtype)
  391. if data_cat.size > 0:
  392. indptr_cat = np.concatenate(indptr_list).astype(idx_dtype)
  393. indices_cat = (np.concatenate([b.indices for b in blocks])
  394. .astype(idx_dtype))
  395. indptr = np.empty(constant_dim + 1, dtype=idx_dtype)
  396. indices = np.empty_like(indices_cat)
  397. data = np.empty_like(data_cat)
  398. csr_hstack(n_blocks, constant_dim, stack_dim_cat,
  399. indptr_cat, indices_cat, data_cat,
  400. indptr, indices, data)
  401. else:
  402. indptr = np.zeros(constant_dim + 1, dtype=idx_dtype)
  403. indices = np.empty(0, dtype=idx_dtype)
  404. data = np.empty(0, dtype=data_cat.dtype)
  405. if axis == 0:
  406. return csc_matrix((data, indices, indptr),
  407. shape=(sum_dim, constant_dim))
  408. else:
  409. return csr_matrix((data, indices, indptr),
  410. shape=(constant_dim, sum_dim))
  411. def hstack(blocks, format=None, dtype=None):
  412. """
  413. Stack sparse matrices horizontally (column wise)
  414. Parameters
  415. ----------
  416. blocks
  417. sequence of sparse matrices with compatible shapes
  418. format : str
  419. sparse format of the result (e.g., "csr")
  420. by default an appropriate sparse matrix format is returned.
  421. This choice is subject to change.
  422. dtype : dtype, optional
  423. The data-type of the output matrix. If not given, the dtype is
  424. determined from that of `blocks`.
  425. See Also
  426. --------
  427. vstack : stack sparse matrices vertically (row wise)
  428. Examples
  429. --------
  430. >>> from scipy.sparse import coo_matrix, hstack
  431. >>> A = coo_matrix([[1, 2], [3, 4]])
  432. >>> B = coo_matrix([[5], [6]])
  433. >>> hstack([A,B]).toarray()
  434. array([[1, 2, 5],
  435. [3, 4, 6]])
  436. """
  437. return bmat([blocks], format=format, dtype=dtype)
  438. def vstack(blocks, format=None, dtype=None):
  439. """
  440. Stack sparse matrices vertically (row wise)
  441. Parameters
  442. ----------
  443. blocks
  444. sequence of sparse matrices with compatible shapes
  445. format : str, optional
  446. sparse format of the result (e.g., "csr")
  447. by default an appropriate sparse matrix format is returned.
  448. This choice is subject to change.
  449. dtype : dtype, optional
  450. The data-type of the output matrix. If not given, the dtype is
  451. determined from that of `blocks`.
  452. See Also
  453. --------
  454. hstack : stack sparse matrices horizontally (column wise)
  455. Examples
  456. --------
  457. >>> from scipy.sparse import coo_matrix, vstack
  458. >>> A = coo_matrix([[1, 2], [3, 4]])
  459. >>> B = coo_matrix([[5, 6]])
  460. >>> vstack([A, B]).toarray()
  461. array([[1, 2],
  462. [3, 4],
  463. [5, 6]])
  464. """
  465. return bmat([[b] for b in blocks], format=format, dtype=dtype)
  466. def bmat(blocks, format=None, dtype=None):
  467. """
  468. Build a sparse matrix from sparse sub-blocks
  469. Parameters
  470. ----------
  471. blocks : array_like
  472. Grid of sparse matrices with compatible shapes.
  473. An entry of None implies an all-zero matrix.
  474. format : {'bsr', 'coo', 'csc', 'csr', 'dia', 'dok', 'lil'}, optional
  475. The sparse format of the result (e.g. "csr"). By default an
  476. appropriate sparse matrix format is returned.
  477. This choice is subject to change.
  478. dtype : dtype, optional
  479. The data-type of the output matrix. If not given, the dtype is
  480. determined from that of `blocks`.
  481. Returns
  482. -------
  483. bmat : sparse matrix
  484. See Also
  485. --------
  486. block_diag, diags
  487. Examples
  488. --------
  489. >>> from scipy.sparse import coo_matrix, bmat
  490. >>> A = coo_matrix([[1, 2], [3, 4]])
  491. >>> B = coo_matrix([[5], [6]])
  492. >>> C = coo_matrix([[7]])
  493. >>> bmat([[A, B], [None, C]]).toarray()
  494. array([[1, 2, 5],
  495. [3, 4, 6],
  496. [0, 0, 7]])
  497. >>> bmat([[A, None], [None, C]]).toarray()
  498. array([[1, 2, 0],
  499. [3, 4, 0],
  500. [0, 0, 7]])
  501. """
  502. blocks = np.asarray(blocks, dtype='object')
  503. if blocks.ndim != 2:
  504. raise ValueError('blocks must be 2-D')
  505. M,N = blocks.shape
  506. # check for fast path cases
  507. if (format in (None, 'csr') and all(isinstance(b, csr_matrix)
  508. for b in blocks.flat)):
  509. if N > 1:
  510. # stack along columns (axis 1):
  511. blocks = [[_stack_along_minor_axis(blocks[b, :], 1)]
  512. for b in range(M)] # must have shape: (M, 1)
  513. blocks = np.asarray(blocks, dtype='object')
  514. # stack along rows (axis 0):
  515. A = _compressed_sparse_stack(blocks[:, 0], 0)
  516. if dtype is not None:
  517. A = A.astype(dtype)
  518. return A
  519. elif (format in (None, 'csc') and all(isinstance(b, csc_matrix)
  520. for b in blocks.flat)):
  521. if M > 1:
  522. # stack along rows (axis 0):
  523. blocks = [[_stack_along_minor_axis(blocks[:, b], 0)
  524. for b in range(N)]] # must have shape: (1, N)
  525. blocks = np.asarray(blocks, dtype='object')
  526. # stack along columns (axis 1):
  527. A = _compressed_sparse_stack(blocks[0, :], 1)
  528. if dtype is not None:
  529. A = A.astype(dtype)
  530. return A
  531. block_mask = np.zeros(blocks.shape, dtype=bool)
  532. brow_lengths = np.zeros(M, dtype=np.int64)
  533. bcol_lengths = np.zeros(N, dtype=np.int64)
  534. # convert everything to COO format
  535. for i in range(M):
  536. for j in range(N):
  537. if blocks[i,j] is not None:
  538. A = coo_matrix(blocks[i,j])
  539. blocks[i,j] = A
  540. block_mask[i,j] = True
  541. if brow_lengths[i] == 0:
  542. brow_lengths[i] = A.shape[0]
  543. elif brow_lengths[i] != A.shape[0]:
  544. msg = (f'blocks[{i},:] has incompatible row dimensions. '
  545. f'Got blocks[{i},{j}].shape[0] == {A.shape[0]}, '
  546. f'expected {brow_lengths[i]}.')
  547. raise ValueError(msg)
  548. if bcol_lengths[j] == 0:
  549. bcol_lengths[j] = A.shape[1]
  550. elif bcol_lengths[j] != A.shape[1]:
  551. msg = (f'blocks[:,{j}] has incompatible column '
  552. f'dimensions. '
  553. f'Got blocks[{i},{j}].shape[1] == {A.shape[1]}, '
  554. f'expected {bcol_lengths[j]}.')
  555. raise ValueError(msg)
  556. nnz = sum(block.nnz for block in blocks[block_mask])
  557. if dtype is None:
  558. all_dtypes = [blk.dtype for blk in blocks[block_mask]]
  559. dtype = upcast(*all_dtypes) if all_dtypes else None
  560. row_offsets = np.append(0, np.cumsum(brow_lengths))
  561. col_offsets = np.append(0, np.cumsum(bcol_lengths))
  562. shape = (row_offsets[-1], col_offsets[-1])
  563. data = np.empty(nnz, dtype=dtype)
  564. idx_dtype = get_index_dtype(maxval=max(shape))
  565. row = np.empty(nnz, dtype=idx_dtype)
  566. col = np.empty(nnz, dtype=idx_dtype)
  567. nnz = 0
  568. ii, jj = np.nonzero(block_mask)
  569. for i, j in zip(ii, jj):
  570. B = blocks[i, j]
  571. idx = slice(nnz, nnz + B.nnz)
  572. data[idx] = B.data
  573. np.add(B.row, row_offsets[i], out=row[idx], dtype=idx_dtype)
  574. np.add(B.col, col_offsets[j], out=col[idx], dtype=idx_dtype)
  575. nnz += B.nnz
  576. return coo_matrix((data, (row, col)), shape=shape).asformat(format)
  577. def block_diag(mats, format=None, dtype=None):
  578. """
  579. Build a block diagonal sparse matrix from provided matrices.
  580. Parameters
  581. ----------
  582. mats : sequence of matrices
  583. Input matrices.
  584. format : str, optional
  585. The sparse format of the result (e.g., "csr"). If not given, the matrix
  586. is returned in "coo" format.
  587. dtype : dtype specifier, optional
  588. The data-type of the output matrix. If not given, the dtype is
  589. determined from that of `blocks`.
  590. Returns
  591. -------
  592. res : sparse matrix
  593. Notes
  594. -----
  595. .. versionadded:: 0.11.0
  596. See Also
  597. --------
  598. bmat, diags
  599. Examples
  600. --------
  601. >>> from scipy.sparse import coo_matrix, block_diag
  602. >>> A = coo_matrix([[1, 2], [3, 4]])
  603. >>> B = coo_matrix([[5], [6]])
  604. >>> C = coo_matrix([[7]])
  605. >>> block_diag((A, B, C)).toarray()
  606. array([[1, 2, 0, 0],
  607. [3, 4, 0, 0],
  608. [0, 0, 5, 0],
  609. [0, 0, 6, 0],
  610. [0, 0, 0, 7]])
  611. """
  612. row = []
  613. col = []
  614. data = []
  615. r_idx = 0
  616. c_idx = 0
  617. for a in mats:
  618. if isinstance(a, (list, numbers.Number)):
  619. a = coo_matrix(a)
  620. nrows, ncols = a.shape
  621. if issparse(a):
  622. a = a.tocoo()
  623. row.append(a.row + r_idx)
  624. col.append(a.col + c_idx)
  625. data.append(a.data)
  626. else:
  627. a_row, a_col = np.divmod(np.arange(nrows*ncols), ncols)
  628. row.append(a_row + r_idx)
  629. col.append(a_col + c_idx)
  630. data.append(a.ravel())
  631. r_idx += nrows
  632. c_idx += ncols
  633. row = np.concatenate(row)
  634. col = np.concatenate(col)
  635. data = np.concatenate(data)
  636. return coo_matrix((data, (row, col)),
  637. shape=(r_idx, c_idx),
  638. dtype=dtype).asformat(format)
  639. def random(m, n, density=0.01, format='coo', dtype=None,
  640. random_state=None, data_rvs=None):
  641. """Generate a sparse matrix of the given shape and density with randomly
  642. distributed values.
  643. Parameters
  644. ----------
  645. m, n : int
  646. shape of the matrix
  647. density : real, optional
  648. density of the generated matrix: density equal to one means a full
  649. matrix, density of 0 means a matrix with no non-zero items.
  650. format : str, optional
  651. sparse matrix format.
  652. dtype : dtype, optional
  653. type of the returned matrix values.
  654. random_state : {None, int, `numpy.random.Generator`,
  655. `numpy.random.RandomState`}, optional
  656. If `seed` is None (or `np.random`), the `numpy.random.RandomState`
  657. singleton is used.
  658. If `seed` is an int, a new ``RandomState`` instance is used,
  659. seeded with `seed`.
  660. If `seed` is already a ``Generator`` or ``RandomState`` instance then
  661. that instance is used.
  662. This random state will be used
  663. for sampling the sparsity structure, but not necessarily for sampling
  664. the values of the structurally nonzero entries of the matrix.
  665. data_rvs : callable, optional
  666. Samples a requested number of random values.
  667. This function should take a single argument specifying the length
  668. of the ndarray that it will return. The structurally nonzero entries
  669. of the sparse random matrix will be taken from the array sampled
  670. by this function. By default, uniform [0, 1) random values will be
  671. sampled using the same random state as is used for sampling
  672. the sparsity structure.
  673. Returns
  674. -------
  675. res : sparse matrix
  676. Notes
  677. -----
  678. Only float types are supported for now.
  679. Examples
  680. --------
  681. >>> from scipy.sparse import random
  682. >>> from scipy import stats
  683. >>> from numpy.random import default_rng
  684. >>> rng = default_rng()
  685. >>> rvs = stats.poisson(25, loc=10).rvs
  686. >>> S = random(3, 4, density=0.25, random_state=rng, data_rvs=rvs)
  687. >>> S.A
  688. array([[ 36., 0., 33., 0.], # random
  689. [ 0., 0., 0., 0.],
  690. [ 0., 0., 36., 0.]])
  691. >>> from scipy.sparse import random
  692. >>> from scipy.stats import rv_continuous
  693. >>> class CustomDistribution(rv_continuous):
  694. ... def _rvs(self, size=None, random_state=None):
  695. ... return random_state.standard_normal(size)
  696. >>> X = CustomDistribution(seed=rng)
  697. >>> Y = X() # get a frozen version of the distribution
  698. >>> S = random(3, 4, density=0.25, random_state=rng, data_rvs=Y.rvs)
  699. >>> S.A
  700. array([[ 0. , 0. , 0. , 0. ], # random
  701. [ 0.13569738, 1.9467163 , -0.81205367, 0. ],
  702. [ 0. , 0. , 0. , 0. ]])
  703. """
  704. if density < 0 or density > 1:
  705. raise ValueError("density expected to be 0 <= density <= 1")
  706. dtype = np.dtype(dtype)
  707. mn = m * n
  708. tp = np.intc
  709. if mn > np.iinfo(tp).max:
  710. tp = np.int64
  711. if mn > np.iinfo(tp).max:
  712. msg = """\
  713. Trying to generate a random sparse matrix such as the product of dimensions is
  714. greater than %d - this is not supported on this machine
  715. """
  716. raise ValueError(msg % np.iinfo(tp).max)
  717. # Number of non zero values
  718. k = int(round(density * m * n))
  719. random_state = check_random_state(random_state)
  720. if data_rvs is None:
  721. if np.issubdtype(dtype, np.integer):
  722. def data_rvs(n):
  723. return rng_integers(random_state,
  724. np.iinfo(dtype).min,
  725. np.iinfo(dtype).max,
  726. n,
  727. dtype=dtype)
  728. elif np.issubdtype(dtype, np.complexfloating):
  729. def data_rvs(n):
  730. return (random_state.uniform(size=n) +
  731. random_state.uniform(size=n) * 1j)
  732. else:
  733. data_rvs = partial(random_state.uniform, 0., 1.)
  734. ind = random_state.choice(mn, size=k, replace=False)
  735. j = np.floor(ind * 1. / m).astype(tp, copy=False)
  736. i = (ind - j * m).astype(tp, copy=False)
  737. vals = data_rvs(k).astype(dtype, copy=False)
  738. return coo_matrix((vals, (i, j)), shape=(m, n)).asformat(format,
  739. copy=False)
  740. def rand(m, n, density=0.01, format="coo", dtype=None, random_state=None):
  741. """Generate a sparse matrix of the given shape and density with uniformly
  742. distributed values.
  743. Parameters
  744. ----------
  745. m, n : int
  746. shape of the matrix
  747. density : real, optional
  748. density of the generated matrix: density equal to one means a full
  749. matrix, density of 0 means a matrix with no non-zero items.
  750. format : str, optional
  751. sparse matrix format.
  752. dtype : dtype, optional
  753. type of the returned matrix values.
  754. random_state : {None, int, `numpy.random.Generator`,
  755. `numpy.random.RandomState`}, optional
  756. If `seed` is None (or `np.random`), the `numpy.random.RandomState`
  757. singleton is used.
  758. If `seed` is an int, a new ``RandomState`` instance is used,
  759. seeded with `seed`.
  760. If `seed` is already a ``Generator`` or ``RandomState`` instance then
  761. that instance is used.
  762. Returns
  763. -------
  764. res : sparse matrix
  765. Notes
  766. -----
  767. Only float types are supported for now.
  768. See Also
  769. --------
  770. scipy.sparse.random : Similar function that allows a user-specified random
  771. data source.
  772. Examples
  773. --------
  774. >>> from scipy.sparse import rand
  775. >>> matrix = rand(3, 4, density=0.25, format="csr", random_state=42)
  776. >>> matrix
  777. <3x4 sparse matrix of type '<class 'numpy.float64'>'
  778. with 3 stored elements in Compressed Sparse Row format>
  779. >>> matrix.toarray()
  780. array([[0.05641158, 0. , 0. , 0.65088847],
  781. [0. , 0. , 0. , 0.14286682],
  782. [0. , 0. , 0. , 0. ]])
  783. """
  784. return random(m, n, density, format, dtype, random_state)