Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1"""Functions to construct sparse matrices 

2""" 

3 

4__docformat__ = "restructuredtext en" 

5 

6__all__ = ['spdiags', 'eye', 'identity', 'kron', 'kronsum', 

7 'hstack', 'vstack', 'bmat', 'rand', 'random', 'diags', 'block_diag'] 

8 

9 

10from functools import partial 

11import numpy as np 

12 

13from scipy._lib._util import check_random_state, rng_integers 

14from .sputils import upcast, get_index_dtype, isscalarlike 

15 

16from .csr import csr_matrix 

17from .csc import csc_matrix 

18from .bsr import bsr_matrix 

19from .coo import coo_matrix 

20from .dia import dia_matrix 

21 

22from .base import issparse 

23 

24 

25def spdiags(data, diags, m, n, format=None): 

26 """ 

27 Return a sparse matrix from diagonals. 

28 

29 Parameters 

30 ---------- 

31 data : array_like 

32 matrix diagonals stored row-wise 

33 diags : diagonals to set 

34 - k = 0 the main diagonal 

35 - k > 0 the k-th upper diagonal 

36 - k < 0 the k-th lower diagonal 

37 m, n : int 

38 shape of the result 

39 format : str, optional 

40 Format of the result. By default (format=None) an appropriate sparse 

41 matrix format is returned. This choice is subject to change. 

42 

43 See Also 

44 -------- 

45 diags : more convenient form of this function 

46 dia_matrix : the sparse DIAgonal format. 

47 

48 Examples 

49 -------- 

50 >>> from scipy.sparse import spdiags 

51 >>> data = np.array([[1, 2, 3, 4], [1, 2, 3, 4], [1, 2, 3, 4]]) 

52 >>> diags = np.array([0, -1, 2]) 

53 >>> spdiags(data, diags, 4, 4).toarray() 

54 array([[1, 0, 3, 0], 

55 [1, 2, 0, 4], 

56 [0, 2, 3, 0], 

57 [0, 0, 3, 4]]) 

58 

59 """ 

60 return dia_matrix((data, diags), shape=(m,n)).asformat(format) 

61 

62 

63def diags(diagonals, offsets=0, shape=None, format=None, dtype=None): 

64 """ 

65 Construct a sparse matrix from diagonals. 

66 

67 Parameters 

68 ---------- 

69 diagonals : sequence of array_like 

70 Sequence of arrays containing the matrix diagonals, 

71 corresponding to `offsets`. 

72 offsets : sequence of int or an int, optional 

73 Diagonals to set: 

74 - k = 0 the main diagonal (default) 

75 - k > 0 the kth upper diagonal 

76 - k < 0 the kth lower diagonal 

77 shape : tuple of int, optional 

78 Shape of the result. If omitted, a square matrix large enough 

79 to contain the diagonals is returned. 

80 format : {"dia", "csr", "csc", "lil", ...}, optional 

81 Matrix format of the result. By default (format=None) an 

82 appropriate sparse matrix format is returned. This choice is 

83 subject to change. 

84 dtype : dtype, optional 

85 Data type of the matrix. 

86 

87 See Also 

88 -------- 

89 spdiags : construct matrix from diagonals 

90 

91 Notes 

92 ----- 

93 This function differs from `spdiags` in the way it handles 

94 off-diagonals. 

95 

96 The result from `diags` is the sparse equivalent of:: 

97 

98 np.diag(diagonals[0], offsets[0]) 

99 + ... 

100 + np.diag(diagonals[k], offsets[k]) 

101 

102 Repeated diagonal offsets are disallowed. 

103 

104 .. versionadded:: 0.11 

105 

106 Examples 

107 -------- 

108 >>> from scipy.sparse import diags 

109 >>> diagonals = [[1, 2, 3, 4], [1, 2, 3], [1, 2]] 

110 >>> diags(diagonals, [0, -1, 2]).toarray() 

111 array([[1, 0, 1, 0], 

112 [1, 2, 0, 2], 

113 [0, 2, 3, 0], 

114 [0, 0, 3, 4]]) 

115 

116 Broadcasting of scalars is supported (but shape needs to be 

117 specified): 

118 

119 >>> diags([1, -2, 1], [-1, 0, 1], shape=(4, 4)).toarray() 

120 array([[-2., 1., 0., 0.], 

121 [ 1., -2., 1., 0.], 

122 [ 0., 1., -2., 1.], 

123 [ 0., 0., 1., -2.]]) 

124 

125 

126 If only one diagonal is wanted (as in `numpy.diag`), the following 

127 works as well: 

128 

129 >>> diags([1, 2, 3], 1).toarray() 

130 array([[ 0., 1., 0., 0.], 

131 [ 0., 0., 2., 0.], 

132 [ 0., 0., 0., 3.], 

133 [ 0., 0., 0., 0.]]) 

134 """ 

135 # if offsets is not a sequence, assume that there's only one diagonal 

136 if isscalarlike(offsets): 

137 # now check that there's actually only one diagonal 

138 if len(diagonals) == 0 or isscalarlike(diagonals[0]): 

139 diagonals = [np.atleast_1d(diagonals)] 

140 else: 

141 raise ValueError("Different number of diagonals and offsets.") 

142 else: 

143 diagonals = list(map(np.atleast_1d, diagonals)) 

144 

145 offsets = np.atleast_1d(offsets) 

146 

147 # Basic check 

148 if len(diagonals) != len(offsets): 

149 raise ValueError("Different number of diagonals and offsets.") 

150 

151 # Determine shape, if omitted 

152 if shape is None: 

153 m = len(diagonals[0]) + abs(int(offsets[0])) 

154 shape = (m, m) 

155 

156 # Determine data type, if omitted 

157 if dtype is None: 

158 dtype = np.common_type(*diagonals) 

159 

160 # Construct data array 

161 m, n = shape 

162 

163 M = max([min(m + offset, n - offset) + max(0, offset) 

164 for offset in offsets]) 

165 M = max(0, M) 

166 data_arr = np.zeros((len(offsets), M), dtype=dtype) 

167 

168 K = min(m, n) 

169 

170 for j, diagonal in enumerate(diagonals): 

171 offset = offsets[j] 

172 k = max(0, offset) 

173 length = min(m + offset, n - offset, K) 

174 if length < 0: 

175 raise ValueError("Offset %d (index %d) out of bounds" % (offset, j)) 

176 try: 

177 data_arr[j, k:k+length] = diagonal[...,:length] 

178 except ValueError: 

179 if len(diagonal) != length and len(diagonal) != 1: 

180 raise ValueError( 

181 "Diagonal length (index %d: %d at offset %d) does not " 

182 "agree with matrix size (%d, %d)." % ( 

183 j, len(diagonal), offset, m, n)) 

184 raise 

185 

186 return dia_matrix((data_arr, offsets), shape=(m, n)).asformat(format) 

187 

188 

189def identity(n, dtype='d', format=None): 

190 """Identity matrix in sparse format 

191 

192 Returns an identity matrix with shape (n,n) using a given 

193 sparse format and dtype. 

194 

195 Parameters 

196 ---------- 

197 n : int 

198 Shape of the identity matrix. 

199 dtype : dtype, optional 

200 Data type of the matrix 

201 format : str, optional 

202 Sparse format of the result, e.g., format="csr", etc. 

203 

204 Examples 

205 -------- 

206 >>> from scipy.sparse import identity 

207 >>> identity(3).toarray() 

208 array([[ 1., 0., 0.], 

209 [ 0., 1., 0.], 

210 [ 0., 0., 1.]]) 

211 >>> identity(3, dtype='int8', format='dia') 

212 <3x3 sparse matrix of type '<class 'numpy.int8'>' 

213 with 3 stored elements (1 diagonals) in DIAgonal format> 

214 

215 """ 

216 return eye(n, n, dtype=dtype, format=format) 

217 

218 

219def eye(m, n=None, k=0, dtype=float, format=None): 

220 """Sparse matrix with ones on diagonal 

221 

222 Returns a sparse (m x n) matrix where the kth diagonal 

223 is all ones and everything else is zeros. 

224 

225 Parameters 

226 ---------- 

227 m : int 

228 Number of rows in the matrix. 

229 n : int, optional 

230 Number of columns. Default: `m`. 

231 k : int, optional 

232 Diagonal to place ones on. Default: 0 (main diagonal). 

233 dtype : dtype, optional 

234 Data type of the matrix. 

235 format : str, optional 

236 Sparse format of the result, e.g., format="csr", etc. 

237 

238 Examples 

239 -------- 

240 >>> from scipy import sparse 

241 >>> sparse.eye(3).toarray() 

242 array([[ 1., 0., 0.], 

243 [ 0., 1., 0.], 

244 [ 0., 0., 1.]]) 

245 >>> sparse.eye(3, dtype=np.int8) 

246 <3x3 sparse matrix of type '<class 'numpy.int8'>' 

247 with 3 stored elements (1 diagonals) in DIAgonal format> 

248 

249 """ 

250 if n is None: 

251 n = m 

252 m,n = int(m),int(n) 

253 

254 if m == n and k == 0: 

255 # fast branch for special formats 

256 if format in ['csr', 'csc']: 

257 idx_dtype = get_index_dtype(maxval=n) 

258 indptr = np.arange(n+1, dtype=idx_dtype) 

259 indices = np.arange(n, dtype=idx_dtype) 

260 data = np.ones(n, dtype=dtype) 

261 cls = {'csr': csr_matrix, 'csc': csc_matrix}[format] 

262 return cls((data,indices,indptr),(n,n)) 

263 elif format == 'coo': 

264 idx_dtype = get_index_dtype(maxval=n) 

265 row = np.arange(n, dtype=idx_dtype) 

266 col = np.arange(n, dtype=idx_dtype) 

267 data = np.ones(n, dtype=dtype) 

268 return coo_matrix((data,(row,col)),(n,n)) 

269 

270 diags = np.ones((1, max(0, min(m + k, n))), dtype=dtype) 

271 return spdiags(diags, k, m, n).asformat(format) 

272 

273 

274def kron(A, B, format=None): 

275 """kronecker product of sparse matrices A and B 

276 

277 Parameters 

278 ---------- 

279 A : sparse or dense matrix 

280 first matrix of the product 

281 B : sparse or dense matrix 

282 second matrix of the product 

283 format : str, optional 

284 format of the result (e.g. "csr") 

285 

286 Returns 

287 ------- 

288 kronecker product in a sparse matrix format 

289 

290 

291 Examples 

292 -------- 

293 >>> from scipy import sparse 

294 >>> A = sparse.csr_matrix(np.array([[0, 2], [5, 0]])) 

295 >>> B = sparse.csr_matrix(np.array([[1, 2], [3, 4]])) 

296 >>> sparse.kron(A, B).toarray() 

297 array([[ 0, 0, 2, 4], 

298 [ 0, 0, 6, 8], 

299 [ 5, 10, 0, 0], 

300 [15, 20, 0, 0]]) 

301 

302 >>> sparse.kron(A, [[1, 2], [3, 4]]).toarray() 

303 array([[ 0, 0, 2, 4], 

304 [ 0, 0, 6, 8], 

305 [ 5, 10, 0, 0], 

306 [15, 20, 0, 0]]) 

307 

308 """ 

309 B = coo_matrix(B) 

310 

311 if (format is None or format == "bsr") and 2*B.nnz >= B.shape[0] * B.shape[1]: 

312 # B is fairly dense, use BSR 

313 A = csr_matrix(A,copy=True) 

314 

315 output_shape = (A.shape[0]*B.shape[0], A.shape[1]*B.shape[1]) 

316 

317 if A.nnz == 0 or B.nnz == 0: 

318 # kronecker product is the zero matrix 

319 return coo_matrix(output_shape) 

320 

321 B = B.toarray() 

322 data = A.data.repeat(B.size).reshape(-1,B.shape[0],B.shape[1]) 

323 data = data * B 

324 

325 return bsr_matrix((data,A.indices,A.indptr), shape=output_shape) 

326 else: 

327 # use COO 

328 A = coo_matrix(A) 

329 output_shape = (A.shape[0]*B.shape[0], A.shape[1]*B.shape[1]) 

330 

331 if A.nnz == 0 or B.nnz == 0: 

332 # kronecker product is the zero matrix 

333 return coo_matrix(output_shape) 

334 

335 # expand entries of a into blocks 

336 row = A.row.repeat(B.nnz) 

337 col = A.col.repeat(B.nnz) 

338 data = A.data.repeat(B.nnz) 

339 

340 if max(A.shape[0]*B.shape[0], A.shape[1]*B.shape[1]) > np.iinfo('int32').max: 

341 row = row.astype(np.int64) 

342 col = col.astype(np.int64) 

343 

344 row *= B.shape[0] 

345 col *= B.shape[1] 

346 

347 # increment block indices 

348 row,col = row.reshape(-1,B.nnz),col.reshape(-1,B.nnz) 

349 row += B.row 

350 col += B.col 

351 row,col = row.reshape(-1),col.reshape(-1) 

352 

353 # compute block entries 

354 data = data.reshape(-1,B.nnz) * B.data 

355 data = data.reshape(-1) 

356 

357 return coo_matrix((data,(row,col)), shape=output_shape).asformat(format) 

358 

359 

360def kronsum(A, B, format=None): 

361 """kronecker sum of sparse matrices A and B 

362 

363 Kronecker sum of two sparse matrices is a sum of two Kronecker 

364 products kron(I_n,A) + kron(B,I_m) where A has shape (m,m) 

365 and B has shape (n,n) and I_m and I_n are identity matrices 

366 of shape (m,m) and (n,n), respectively. 

367 

368 Parameters 

369 ---------- 

370 A 

371 square matrix 

372 B 

373 square matrix 

374 format : str 

375 format of the result (e.g. "csr") 

376 

377 Returns 

378 ------- 

379 kronecker sum in a sparse matrix format 

380 

381 Examples 

382 -------- 

383 

384 

385 """ 

386 A = coo_matrix(A) 

387 B = coo_matrix(B) 

388 

389 if A.shape[0] != A.shape[1]: 

390 raise ValueError('A is not square') 

391 

392 if B.shape[0] != B.shape[1]: 

393 raise ValueError('B is not square') 

394 

395 dtype = upcast(A.dtype, B.dtype) 

396 

397 L = kron(eye(B.shape[0],dtype=dtype), A, format=format) 

398 R = kron(B, eye(A.shape[0],dtype=dtype), format=format) 

399 

400 return (L+R).asformat(format) # since L + R is not always same format 

401 

402 

403def _compressed_sparse_stack(blocks, axis): 

404 """ 

405 Stacking fast path for CSR/CSC matrices 

406 (i) vstack for CSR, (ii) hstack for CSC. 

407 """ 

408 other_axis = 1 if axis == 0 else 0 

409 data = np.concatenate([b.data for b in blocks]) 

410 constant_dim = blocks[0].shape[other_axis] 

411 idx_dtype = get_index_dtype(arrays=[b.indptr for b in blocks], 

412 maxval=max(data.size, constant_dim)) 

413 indices = np.empty(data.size, dtype=idx_dtype) 

414 indptr = np.empty(sum(b.shape[axis] for b in blocks) + 1, dtype=idx_dtype) 

415 last_indptr = idx_dtype(0) 

416 sum_dim = 0 

417 sum_indices = 0 

418 for b in blocks: 

419 if b.shape[other_axis] != constant_dim: 

420 raise ValueError('incompatible dimensions for axis %d' % other_axis) 

421 indices[sum_indices:sum_indices+b.indices.size] = b.indices 

422 sum_indices += b.indices.size 

423 idxs = slice(sum_dim, sum_dim + b.shape[axis]) 

424 indptr[idxs] = b.indptr[:-1] 

425 indptr[idxs] += last_indptr 

426 sum_dim += b.shape[axis] 

427 last_indptr += b.indptr[-1] 

428 indptr[-1] = last_indptr 

429 if axis == 0: 

430 return csr_matrix((data, indices, indptr), 

431 shape=(sum_dim, constant_dim)) 

432 else: 

433 return csc_matrix((data, indices, indptr), 

434 shape=(constant_dim, sum_dim)) 

435 

436 

437def hstack(blocks, format=None, dtype=None): 

438 """ 

439 Stack sparse matrices horizontally (column wise) 

440 

441 Parameters 

442 ---------- 

443 blocks 

444 sequence of sparse matrices with compatible shapes 

445 format : str 

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 

453 See Also 

454 -------- 

455 vstack : stack sparse matrices vertically (row wise) 

456 

457 Examples 

458 -------- 

459 >>> from scipy.sparse import coo_matrix, hstack 

460 >>> A = coo_matrix([[1, 2], [3, 4]]) 

461 >>> B = coo_matrix([[5], [6]]) 

462 >>> hstack([A,B]).toarray() 

463 array([[1, 2, 5], 

464 [3, 4, 6]]) 

465 

466 """ 

467 return bmat([blocks], format=format, dtype=dtype) 

468 

469 

470def vstack(blocks, format=None, dtype=None): 

471 """ 

472 Stack sparse matrices vertically (row wise) 

473 

474 Parameters 

475 ---------- 

476 blocks 

477 sequence of sparse matrices with compatible shapes 

478 format : str, optional 

479 sparse format of the result (e.g., "csr") 

480 by default an appropriate sparse matrix format is returned. 

481 This choice is subject to change. 

482 dtype : dtype, optional 

483 The data-type of the output matrix. If not given, the dtype is 

484 determined from that of `blocks`. 

485 

486 See Also 

487 -------- 

488 hstack : stack sparse matrices horizontally (column wise) 

489 

490 Examples 

491 -------- 

492 >>> from scipy.sparse import coo_matrix, vstack 

493 >>> A = coo_matrix([[1, 2], [3, 4]]) 

494 >>> B = coo_matrix([[5, 6]]) 

495 >>> vstack([A, B]).toarray() 

496 array([[1, 2], 

497 [3, 4], 

498 [5, 6]]) 

499 

500 """ 

501 return bmat([[b] for b in blocks], format=format, dtype=dtype) 

502 

503 

504def bmat(blocks, format=None, dtype=None): 

505 """ 

506 Build a sparse matrix from sparse sub-blocks 

507 

508 Parameters 

509 ---------- 

510 blocks : array_like 

511 Grid of sparse matrices with compatible shapes. 

512 An entry of None implies an all-zero matrix. 

513 format : {'bsr', 'coo', 'csc', 'csr', 'dia', 'dok', 'lil'}, optional 

514 The sparse format of the result (e.g. "csr"). By default an 

515 appropriate sparse matrix format is returned. 

516 This choice is subject to change. 

517 dtype : dtype, optional 

518 The data-type of the output matrix. If not given, the dtype is 

519 determined from that of `blocks`. 

520 

521 Returns 

522 ------- 

523 bmat : sparse matrix 

524 

525 See Also 

526 -------- 

527 block_diag, diags 

528 

529 Examples 

530 -------- 

531 >>> from scipy.sparse import coo_matrix, bmat 

532 >>> A = coo_matrix([[1, 2], [3, 4]]) 

533 >>> B = coo_matrix([[5], [6]]) 

534 >>> C = coo_matrix([[7]]) 

535 >>> bmat([[A, B], [None, C]]).toarray() 

536 array([[1, 2, 5], 

537 [3, 4, 6], 

538 [0, 0, 7]]) 

539 

540 >>> bmat([[A, None], [None, C]]).toarray() 

541 array([[1, 2, 0], 

542 [3, 4, 0], 

543 [0, 0, 7]]) 

544 

545 """ 

546 

547 blocks = np.asarray(blocks, dtype='object') 

548 

549 if blocks.ndim != 2: 

550 raise ValueError('blocks must be 2-D') 

551 

552 M,N = blocks.shape 

553 

554 # check for fast path cases 

555 if (N == 1 and format in (None, 'csr') and all(isinstance(b, csr_matrix) 

556 for b in blocks.flat)): 

557 A = _compressed_sparse_stack(blocks[:,0], 0) 

558 if dtype is not None: 

559 A = A.astype(dtype) 

560 return A 

561 elif (M == 1 and format in (None, 'csc') 

562 and all(isinstance(b, csc_matrix) for b in blocks.flat)): 

563 A = _compressed_sparse_stack(blocks[0,:], 1) 

564 if dtype is not None: 

565 A = A.astype(dtype) 

566 return A 

567 

568 block_mask = np.zeros(blocks.shape, dtype=bool) 

569 brow_lengths = np.zeros(M, dtype=np.int64) 

570 bcol_lengths = np.zeros(N, dtype=np.int64) 

571 

572 # convert everything to COO format 

573 for i in range(M): 

574 for j in range(N): 

575 if blocks[i,j] is not None: 

576 A = coo_matrix(blocks[i,j]) 

577 blocks[i,j] = A 

578 block_mask[i,j] = True 

579 

580 if brow_lengths[i] == 0: 

581 brow_lengths[i] = A.shape[0] 

582 elif brow_lengths[i] != A.shape[0]: 

583 msg = ('blocks[{i},:] has incompatible row dimensions. ' 

584 'Got blocks[{i},{j}].shape[0] == {got}, ' 

585 'expected {exp}.'.format(i=i, j=j, 

586 exp=brow_lengths[i], 

587 got=A.shape[0])) 

588 raise ValueError(msg) 

589 

590 if bcol_lengths[j] == 0: 

591 bcol_lengths[j] = A.shape[1] 

592 elif bcol_lengths[j] != A.shape[1]: 

593 msg = ('blocks[:,{j}] has incompatible row dimensions. ' 

594 'Got blocks[{i},{j}].shape[1] == {got}, ' 

595 'expected {exp}.'.format(i=i, j=j, 

596 exp=bcol_lengths[j], 

597 got=A.shape[1])) 

598 raise ValueError(msg) 

599 

600 nnz = sum(block.nnz for block in blocks[block_mask]) 

601 if dtype is None: 

602 all_dtypes = [blk.dtype for blk in blocks[block_mask]] 

603 dtype = upcast(*all_dtypes) if all_dtypes else None 

604 

605 row_offsets = np.append(0, np.cumsum(brow_lengths)) 

606 col_offsets = np.append(0, np.cumsum(bcol_lengths)) 

607 

608 shape = (row_offsets[-1], col_offsets[-1]) 

609 

610 data = np.empty(nnz, dtype=dtype) 

611 idx_dtype = get_index_dtype(maxval=max(shape)) 

612 row = np.empty(nnz, dtype=idx_dtype) 

613 col = np.empty(nnz, dtype=idx_dtype) 

614 

615 nnz = 0 

616 ii, jj = np.nonzero(block_mask) 

617 for i, j in zip(ii, jj): 

618 B = blocks[i, j] 

619 idx = slice(nnz, nnz + B.nnz) 

620 data[idx] = B.data 

621 row[idx] = B.row + row_offsets[i] 

622 col[idx] = B.col + col_offsets[j] 

623 nnz += B.nnz 

624 

625 return coo_matrix((data, (row, col)), shape=shape).asformat(format) 

626 

627 

628def block_diag(mats, format=None, dtype=None): 

629 """ 

630 Build a block diagonal sparse matrix from provided matrices. 

631 

632 Parameters 

633 ---------- 

634 mats : sequence of matrices 

635 Input matrices. 

636 format : str, optional 

637 The sparse format of the result (e.g., "csr"). If not given, the matrix 

638 is returned in "coo" format. 

639 dtype : dtype specifier, optional 

640 The data-type of the output matrix. If not given, the dtype is 

641 determined from that of `blocks`. 

642 

643 Returns 

644 ------- 

645 res : sparse matrix 

646 

647 Notes 

648 ----- 

649 

650 .. versionadded:: 0.11.0 

651 

652 See Also 

653 -------- 

654 bmat, diags 

655 

656 Examples 

657 -------- 

658 >>> from scipy.sparse import coo_matrix, block_diag 

659 >>> A = coo_matrix([[1, 2], [3, 4]]) 

660 >>> B = coo_matrix([[5], [6]]) 

661 >>> C = coo_matrix([[7]]) 

662 >>> block_diag((A, B, C)).toarray() 

663 array([[1, 2, 0, 0], 

664 [3, 4, 0, 0], 

665 [0, 0, 5, 0], 

666 [0, 0, 6, 0], 

667 [0, 0, 0, 7]]) 

668 

669 """ 

670 nmat = len(mats) 

671 rows = [] 

672 for ia, a in enumerate(mats): 

673 row = [None]*nmat 

674 if issparse(a): 

675 row[ia] = a 

676 else: 

677 row[ia] = coo_matrix(a) 

678 rows.append(row) 

679 return bmat(rows, format=format, dtype=dtype) 

680 

681 

682def random(m, n, density=0.01, format='coo', dtype=None, 

683 random_state=None, data_rvs=None): 

684 """Generate a sparse matrix of the given shape and density with randomly 

685 distributed values. 

686 

687 Parameters 

688 ---------- 

689 m, n : int 

690 shape of the matrix 

691 density : real, optional 

692 density of the generated matrix: density equal to one means a full 

693 matrix, density of 0 means a matrix with no non-zero items. 

694 format : str, optional 

695 sparse matrix format. 

696 dtype : dtype, optional 

697 type of the returned matrix values. 

698 random_state : {numpy.random.RandomState, int}, optional 

699 Random number generator or random seed. If not given, the singleton 

700 numpy.random will be used. This random state will be used 

701 for sampling the sparsity structure, but not necessarily for sampling 

702 the values of the structurally nonzero entries of the matrix. 

703 data_rvs : callable, optional 

704 Samples a requested number of random values. 

705 This function should take a single argument specifying the length 

706 of the ndarray that it will return. The structurally nonzero entries 

707 of the sparse random matrix will be taken from the array sampled 

708 by this function. By default, uniform [0, 1) random values will be 

709 sampled using the same random state as is used for sampling 

710 the sparsity structure. 

711 

712 Returns 

713 ------- 

714 res : sparse matrix 

715 

716 Notes 

717 ----- 

718 Only float types are supported for now. 

719 

720 Examples 

721 -------- 

722 >>> from scipy.sparse import random 

723 >>> from scipy import stats 

724 

725 >>> class CustomRandomState(np.random.RandomState): 

726 ... def randint(self, k): 

727 ... i = np.random.randint(k) 

728 ... return i - i % 2 

729 >>> np.random.seed(12345) 

730 >>> rs = CustomRandomState() 

731 >>> rvs = stats.poisson(25, loc=10).rvs 

732 >>> S = random(3, 4, density=0.25, random_state=rs, data_rvs=rvs) 

733 >>> S.A 

734 array([[ 36., 0., 33., 0.], # random 

735 [ 0., 0., 0., 0.], 

736 [ 0., 0., 36., 0.]]) 

737 

738 >>> from scipy.sparse import random 

739 >>> from scipy.stats import rv_continuous 

740 >>> class CustomDistribution(rv_continuous): 

741 ... def _rvs(self, size=None, random_state=None): 

742 ... return random_state.randn(*size) 

743 >>> X = CustomDistribution(seed=2906) 

744 >>> Y = X() # get a frozen version of the distribution 

745 >>> S = random(3, 4, density=0.25, random_state=2906, data_rvs=Y.rvs) 

746 >>> S.A 

747 array([[ 0. , 0. , 0. , 0. ], 

748 [ 0.13569738, 1.9467163 , -0.81205367, 0. ], 

749 [ 0. , 0. , 0. , 0. ]]) 

750 

751 """ 

752 if density < 0 or density > 1: 

753 raise ValueError("density expected to be 0 <= density <= 1") 

754 dtype = np.dtype(dtype) 

755 

756 mn = m * n 

757 

758 tp = np.intc 

759 if mn > np.iinfo(tp).max: 

760 tp = np.int64 

761 

762 if mn > np.iinfo(tp).max: 

763 msg = """\ 

764Trying to generate a random sparse matrix such as the product of dimensions is 

765greater than %d - this is not supported on this machine 

766""" 

767 raise ValueError(msg % np.iinfo(tp).max) 

768 

769 # Number of non zero values 

770 k = int(round(density * m * n)) 

771 

772 random_state = check_random_state(random_state) 

773 

774 if data_rvs is None: 

775 if np.issubdtype(dtype, np.integer): 

776 def data_rvs(n): 

777 return rng_integers(random_state, 

778 np.iinfo(dtype).min, 

779 np.iinfo(dtype).max, 

780 n, 

781 dtype=dtype) 

782 elif np.issubdtype(dtype, np.complexfloating): 

783 def data_rvs(n): 

784 return (random_state.uniform(size=n) + 

785 random_state.uniform(size=n) * 1j) 

786 else: 

787 data_rvs = partial(random_state.uniform, 0., 1.) 

788 

789 ind = random_state.choice(mn, size=k, replace=False) 

790 

791 j = np.floor(ind * 1. / m).astype(tp, copy=False) 

792 i = (ind - j * m).astype(tp, copy=False) 

793 vals = data_rvs(k).astype(dtype, copy=False) 

794 return coo_matrix((vals, (i, j)), shape=(m, n)).asformat(format, 

795 copy=False) 

796 

797 

798def rand(m, n, density=0.01, format="coo", dtype=None, random_state=None): 

799 """Generate a sparse matrix of the given shape and density with uniformly 

800 distributed values. 

801 

802 Parameters 

803 ---------- 

804 m, n : int 

805 shape of the matrix 

806 density : real, optional 

807 density of the generated matrix: density equal to one means a full 

808 matrix, density of 0 means a matrix with no non-zero items. 

809 format : str, optional 

810 sparse matrix format. 

811 dtype : dtype, optional 

812 type of the returned matrix values. 

813 random_state : {numpy.random.RandomState, int, np.random.Generator}, optional 

814 Random number generator or random seed. If not given, the singleton 

815 numpy.random will be used. 

816 

817 Returns 

818 ------- 

819 res : sparse matrix 

820 

821 Notes 

822 ----- 

823 Only float types are supported for now. 

824 

825 See Also 

826 -------- 

827 scipy.sparse.random : Similar function that allows a user-specified random 

828 data source. 

829 

830 Examples 

831 -------- 

832 >>> from scipy.sparse import rand 

833 >>> matrix = rand(3, 4, density=0.25, format="csr", random_state=42) 

834 >>> matrix 

835 <3x4 sparse matrix of type '<class 'numpy.float64'>' 

836 with 3 stored elements in Compressed Sparse Row format> 

837 >>> matrix.todense() 

838 matrix([[0.05641158, 0. , 0. , 0.65088847], 

839 [0. , 0. , 0. , 0.14286682], 

840 [0. , 0. , 0. , 0. ]]) 

841 

842 """ 

843 return random(m, n, density, format, dtype, random_state)