rsparsematrix.Rd
Generate a random sparse matrix efficiently. The default has rounded
gaussian non-zero entries, and rand.x = NULL
generates random
pattern matrices, i.e. inheriting from nsparseMatrix
.
number of rows and columns, i.e., the matrix
dimension (dim
).
optional number in \([0,1]\), the density is the
proportion of non-zero entries among all matrix entries. If
specified it determines the default for nnz
, otherwise
nnz
needs to be specified.
number of non-zero entries, for a sparse matrix typically
considerably smaller than nrow*ncol
. Must be specified if
density
is not.
logical indicating if result should be a matrix of
class symmetricMatrix
. Note that in the symmetric
case, nnz
denotes the number of non zero entries of the upper
(or lower) part of the matrix, including the diagonal.
NULL
or the random number generator for the x
slot, a
function
such that rand.x(n)
generates a
numeric vector of length n
. Typical examples are
rand.x = rnorm
, or rand.x = runif
; the default is nice
for didactical purposes.
optionally further arguments passed to
sparseMatrix()
, notably repr
.
The algorithm first samples “encoded” \((i,j)\)s without
replacement, via one dimensional indices, if not symmetric
sample.int(nrow*ncol, nnz)
, then—if rand.x
is
not NULL
—gets x <- rand.x(nnz)
and calls
sparseMatrix(i=i, j=j, x=x, ..)
. When
rand.x=NULL
, sparseMatrix(i=i, j=j, ..)
will
return a pattern matrix (i.e., inheriting from
nsparseMatrix
).
a sparseMatrix
, say M
of dimension (nrow,
ncol), i.e., with dim(M) == c(nrow, ncol)
, if symmetric
is not true, with nzM <- nnzero(M)
fulfilling
nzM <= nnz
and typically, nzM == nnz
.
set.seed(17)# to be reproducible
M <- rsparsematrix(8, 12, nnz = 30) # small example, not very sparse
M
#> 8 x 12 sparse Matrix of class "dgCMatrix"
#>
#> [1,] -0.710 . . . 0.54 . . -0.62 . . -1.2 .
#> [2,] . . . . . 0.12 -0.720 . . -0.96 . .
#> [3,] . . 1.20 -0.83 -0.23 . . . . . . .
#> [4,] 0.052 . -0.45 . . . . 1.80 . . . .
#> [5,] -0.780 . 0.25 . 2.40 -0.60 . . . . . .
#> [6,] -0.880 . 0.26 -3.60 . . 0.074 . . . . 1.7
#> [7,] . . . . 1.20 0.38 . . . -0.33 . .
#> [8,] -0.990 -0.21 . . . . -1.400 0.32 1.9 . . .
M1 <- rsparsematrix(1000, 20, nnz = 123, rand.x = runif)
summary(M1)
#> 1000 x 20 sparse Matrix of class "dgCMatrix", with 123 entries
#> i j x
#> 1 30 1 0.121030300
#> 2 256 1 0.171285824
#> 3 471 1 0.658082934
#> 4 490 1 0.672053193
#> 5 515 1 0.563469009
#> 6 799 1 0.754721688
#> 7 914 1 0.214298567
#> 8 112 2 0.230193013
#> 9 146 2 0.520770952
#> 10 238 2 0.895551263
#> 11 460 2 0.747560311
#> 12 494 2 0.338654988
#> 13 669 2 0.633221050
#> 14 753 2 0.029495869
#> 15 26 3 0.857548872
#> 16 103 3 0.877821950
#> 17 123 3 0.373133136
#> 18 163 3 0.365224795
#> 19 395 3 0.244816811
#> 20 425 3 0.474136304
#> 21 697 3 0.463969983
#> 22 858 3 0.272504956
#> 23 875 3 0.840738221
#> 24 157 4 0.277547142
#> 25 317 4 0.394847750
#> 26 400 4 0.991649291
#> 27 510 4 0.362585672
#> 28 514 4 0.520276636
#> 29 634 4 0.808294520
#> 30 850 4 0.007048685
#> 31 934 4 0.536796092
#> 32 30 5 0.166392594
#> 33 263 5 0.306517981
#> 34 301 5 0.834864402
#> 35 884 5 0.772908976
#> 36 967 5 0.792207911
#> 37 268 6 0.697779252
#> 38 290 6 0.080128319
#> 39 434 6 0.307286057
#> 40 535 6 0.539794671
#> 41 759 6 0.651017611
#> 42 866 6 0.069991217
#> 43 167 7 0.720624894
#> 44 394 7 0.798488620
#> 45 502 7 0.837568737
#> 46 580 7 0.182004452
#> 47 705 7 0.289417504
#> 48 723 7 0.114946217
#> 49 216 8 0.780055751
#> 50 409 8 0.961610656
#> 51 631 8 0.188212461
#> 52 361 9 0.849970878
#> 53 411 9 0.110761594
#> 54 482 9 0.918999826
#> 55 564 9 0.106341774
#> 56 634 9 0.966994960
#> 57 817 9 0.080746639
#> 58 842 9 0.936882774
#> 59 940 9 0.723470934
#> 60 989 9 0.964498797
#> 61 47 10 0.359321831
#> 62 405 10 0.023504461
#> 63 507 10 0.112909265
#> 64 717 10 0.748889274
#> 65 734 10 0.642421125
#> 66 329 11 0.111511966
#> 67 438 11 0.498936020
#> 68 492 11 0.691674638
#> 69 570 11 0.246009281
#> 70 748 11 0.161394011
#> 71 753 11 0.872274340
#> 72 783 11 0.100716683
#> 73 977 11 0.462361495
#> 74 355 12 0.807811439
#> 75 645 12 0.652533316
#> 76 672 12 0.437335808
#> 77 697 12 0.308880390
#> 78 864 12 0.614445118
#> 79 949 12 0.469374262
#> 80 167 13 0.640432993
#> 81 318 13 0.598486273
#> 82 382 13 0.373182097
#> 83 400 13 0.973641726
#> 84 547 13 0.338046568
#> 85 601 13 0.117673368
#> 86 729 13 0.159017836
#> 87 996 13 0.350850822
#> 88 307 14 0.998227660
#> 89 385 14 0.323179815
#> 90 453 14 0.351955581
#> 91 526 14 0.706791886
#> 92 680 14 0.538912676
#> 93 916 14 0.176190115
#> 94 926 14 0.276029720
#> 95 11 15 0.118834062
#> 96 161 15 0.343942903
#> 97 211 15 0.850774512
#> 98 323 15 0.133445217
#> 99 478 15 0.788471288
#> 100 597 15 0.145946067
#> 101 608 15 0.197623800
#> 102 886 15 0.862795090
#> 103 275 16 0.714235556
#> 104 510 16 0.019910095
#> 105 527 16 0.639534228
#> 106 695 16 0.428969391
#> 107 765 16 0.905435610
#> 108 103 17 0.899292917
#> 109 736 17 0.150202421
#> 110 101 18 0.505359340
#> 111 513 18 0.784202072
#> 112 739 18 0.770542073
#> 113 941 18 0.728378107
#> 114 225 19 0.216264430
#> 115 558 19 0.311900054
#> 116 680 19 0.807296937
#> 117 948 19 0.314206922
#> 118 154 20 0.226125864
#> 119 341 20 0.942557140
#> 120 637 20 0.635889334
#> 121 650 20 0.914174104
#> 122 767 20 0.673403587
#> 123 952 20 0.392104291
## a random *symmetric* Matrix
(S9 <- rsparsematrix(9, 9, nnz = 10, symmetric=TRUE)) # dsCMatrix
#> 9 x 9 sparse Matrix of class "dsCMatrix"
#>
#> [1,] -0.038 . . 1.50 0.17 . 0.96 . .
#> [2,] . . . . . . -0.48 . .
#> [3,] . . . . . . . 1.3 .
#> [4,] 1.500 . . . . -0.48 . . .
#> [5,] 0.170 . . . . . . 2.4 .
#> [6,] . . . -0.48 . . -0.48 . .
#> [7,] 0.960 -0.48 . . . -0.48 -0.33 . .
#> [8,] . . 1.3 . 2.40 . . . .
#> [9,] . . . . . . . . .
nnzero(S9)# ~ 20: as 'nnz' only counts one "triangle"
#> [1] 18
## a random patter*n* aka boolean Matrix (no 'x' slot):
(n7 <- rsparsematrix(5, 12, nnz = 10, rand.x = NULL))
#> 5 x 12 sparse Matrix of class "ngCMatrix"
#>
#> [1,] . . . . . . . . . . . .
#> [2,] . . . . . . . . | . . .
#> [3,] . . | . . . | . . . . .
#> [4,] | . . . . . . | . | . .
#> [5,] . | . . . | . . . | . |
## a [T]riplet representation sparseMatrix:
T2 <- rsparsematrix(40, 12, nnz = 99, repr = "T")
head(T2)
#> 6 x 12 sparse Matrix of class "dgTMatrix"
#>
#> [1,] 1.8 . -2.8 . . -0.77 . -1.20 . . . 0.27
#> [2,] . -0.71 . . . . . . . 0.96 . -0.32
#> [3,] . . . . . . . . . . -0.051 .
#> [4,] . . . . . . . 0.17 . . -0.027 .
#> [5,] . . . . . . . . . . -1.400 -0.94
#> [6,] . . . -0.26 . . . . . . . .