DelayedRandomArray 1.2.0
The DelayedRandomArray package implements DelayedArray
subclasses containing dynamically sampled random values.
Specifically, the actual values are never fully held in memory but are generated when the relevant part of the array is accessed.
This allows users to create very large arrays of random values that would not otherwise be possible by filling an ordinary matrix.
To install the package, follow the instructions on DelayedRandomArray landing page. Using the package is then as simple as:
library(DelayedRandomArray)
X <- RandomUnifArray(c(1e6, 1e6))
X
## <1000000 x 1000000> matrix of class RandomUnifMatrix and type "double":
## [,1] [,2] [,3] ... [,999999] [,1000000]
## [1,] 0.5486587 0.7946366 0.4656066 . 0.92717022 0.52624366
## [2,] 0.9399642 0.1731777 0.5014395 . 0.04463803 0.96987919
## [3,] 0.6352189 0.3435728 0.5650094 . 0.62403823 0.40307523
## [4,] 0.7279476 0.3939002 0.1544186 . 0.12475696 0.07679103
## [5,] 0.6184666 0.3158507 0.8395698 . 0.43707783 0.39171108
## ... . . . . . .
## [999996,] 0.5255953 0.3782828 0.7988924 . 0.24270297 0.06832743
## [999997,] 0.7784576 0.2602471 0.3683895 . 0.10298012 0.79112390
## [999998,] 0.5669219 0.9443055 0.2206555 . 0.89146668 0.66201774
## [999999,] 0.7244116 0.4783655 0.8127185 . 0.57545293 0.63719682
## [1000000,] 0.3973751 0.2477242 0.4124443 . 0.16440406 0.88812259
The resulting array can be used in any pipeline that is compatible with DelayedArray
objects.
This object occupies only 64 MB in memory,
whereas an ordinary matrix
would require 8 PB instead.
Almost every distribution in stats is available here. To list a few:
RandomNormArray(c(100, 50))
## <100 x 50> matrix of class RandomNormMatrix and type "double":
## [,1] [,2] [,3] ... [,49] [,50]
## [1,] 0.6401411 -0.7366411 -0.9154968 . -0.6251722 -0.8757643
## [2,] -1.1723617 -1.0610589 -0.2041376 . 0.6053334 0.4260002
## [3,] -1.2055987 1.5136896 -2.0441965 . -0.6830776 -0.9729994
## [4,] -1.5669036 0.5797026 0.1806591 . 0.2033038 1.6915929
## [5,] 0.5717801 -0.7158436 -0.4933458 . 0.1308775 0.1372781
## ... . . . . . .
## [96,] -2.11338480 0.56088063 -1.18857215 . 1.42113786 -0.87659859
## [97,] -0.02235439 -1.67029808 0.95417893 . -1.17865682 -2.46154226
## [98,] 0.74278272 -0.76828394 1.59538154 . -0.99879417 -0.38246881
## [99,] 2.41121757 0.52644887 -1.01282662 . -1.18409547 -0.08717828
## [100,] 1.00909005 -0.08316495 1.48806640 . 0.30389343 -0.13052843
RandomPoisArray(c(100, 50), lambda=5)
## <100 x 50> matrix of class RandomPoisMatrix and type "double":
## [,1] [,2] [,3] ... [,49] [,50]
## [1,] 5 4 1 . 4 4
## [2,] 3 3 5 . 7 9
## [3,] 3 4 8 . 3 5
## [4,] 5 4 4 . 4 4
## [5,] 5 7 2 . 6 6
## ... . . . . . .
## [96,] 4 4 7 . 8 2
## [97,] 3 6 3 . 4 5
## [98,] 4 4 5 . 9 1
## [99,] 6 4 1 . 3 8
## [100,] 6 8 4 . 4 4
RandomGammaArray(c(100, 50), shape=2, rate=5)
## <100 x 50> matrix of class RandomGammaMatrix and type "double":
## [,1] [,2] [,3] ... [,49] [,50]
## [1,] 0.01931223 0.57654351 0.33603508 . 0.32832657 0.40863873
## [2,] 0.51068744 0.45074450 0.27984757 . 0.68515184 0.16867340
## [3,] 0.08598661 0.30274657 0.08524634 . 0.80770591 0.87356179
## [4,] 0.75439707 0.25162604 0.09103332 . 0.22328078 0.42061398
## [5,] 0.66749154 0.32985811 0.47924038 . 0.44249425 0.06677582
## ... . . . . . .
## [96,] 0.21356222 0.42673965 0.20241269 . 0.16834888 0.07391251
## [97,] 0.49998389 0.57424574 0.30711896 . 1.10778428 0.18182112
## [98,] 0.13339240 0.56606440 0.46550126 . 0.75404519 0.31050333
## [99,] 0.16921423 0.08051658 1.04282548 . 0.18811712 0.22344999
## [100,] 0.33662015 0.91105495 0.22563113 . 0.42782138 0.92481571
RandomWeibullArray(c(100, 50), shape=5)
## <100 x 50> matrix of class RandomWeibullMatrix and type "double":
## [,1] [,2] [,3] ... [,49] [,50]
## [1,] 0.9267396 1.1430982 1.3022031 . 0.8552672 0.6156038
## [2,] 1.0771462 1.0625598 0.8323211 . 0.9920387 1.2290863
## [3,] 0.7607703 0.8590234 0.9159920 . 0.7412455 1.0485659
## [4,] 0.7645334 1.1998753 0.5505189 . 1.4657181 0.5902842
## [5,] 0.7952678 1.3456427 0.9217310 . 0.9631565 1.1653904
## ... . . . . . .
## [96,] 0.6029694 1.0623164 0.5418988 . 0.7872284 0.8268339
## [97,] 0.9859221 1.0395097 0.8047293 . 1.0054900 0.9482861
## [98,] 0.8662814 0.8509869 1.0654178 . 1.2349706 1.0817435
## [99,] 1.1181221 0.8386192 1.0563748 . 0.6563652 1.0544186
## [100,] 0.8790426 1.2065749 0.8585798 . 1.1955416 1.1093309
Distributional parameters can either be scalars:
RandomNormArray(c(100, 50), mean=1)
## <100 x 50> matrix of class RandomNormMatrix and type "double":
## [,1] [,2] [,3] ... [,49] [,50]
## [1,] -0.75316026 -0.09056193 2.29413495 . -1.4290391 1.6181027
## [2,] 2.82712685 0.55856375 3.05192959 . -1.0576624 0.5648332
## [3,] 1.12325156 -0.15931940 0.97908338 . -0.1129238 1.5325238
## [4,] -1.58436801 0.07534559 1.43058657 . 0.8009343 -0.2716130
## [5,] -0.33740846 -0.91460647 1.40798643 . 1.9963799 1.9680257
## ... . . . . . .
## [96,] 1.81452006 1.22671535 0.71127699 . 0.8479829 0.9859512
## [97,] -1.36690064 0.74984040 0.26047137 . -0.5319456 1.7162226
## [98,] 0.93775413 -0.36857639 0.04192479 . 2.7149064 2.1340845
## [99,] 0.80042780 0.38389342 1.20213512 . 1.0715856 -0.8960479
## [100,] 0.59805583 0.62864113 3.09767299 . 0.5095343 1.6679576
Or vectors, which are recycled along the length of the array:
RandomNormArray(c(100, 50), mean=1:100)
## <100 x 50> matrix of class RandomNormMatrix and type "double":
## [,1] [,2] [,3] ... [,49] [,50]
## [1,] -0.3760225 -0.1295281 -0.3036742 . -1.369089 3.088067
## [2,] 2.4189640 2.5627345 2.7464447 . 3.301008 2.407192
## [3,] 3.8947744 4.6213657 3.1366911 . 4.390639 2.652099
## [4,] 3.3851539 4.7432649 2.2405681 . 4.902984 5.426490
## [5,] 5.2067178 2.3901385 4.4961530 . 6.581792 6.418634
## ... . . . . . .
## [96,] 95.92437 96.65958 97.06213 . 96.40126 96.40516
## [97,] 98.56048 96.71971 97.79066 . 97.84698 98.17902
## [98,] 97.68315 98.59804 98.83767 . 97.38438 98.79016
## [99,] 98.32998 99.22381 98.16929 . 97.11744 101.21516
## [100,] 100.53798 100.31805 100.91001 . 99.66872 99.74978
Or other arrays of the same dimensions, which are used to sample the corresponding parts of the random array:
means <- RandomNormArray(c(100, 50))
RandomPoisArray(c(100, 50), lambda=2^means)
## <100 x 50> matrix of class RandomPoisMatrix and type "double":
## [,1] [,2] [,3] ... [,49] [,50]
## [1,] 1 1 1 . 1 2
## [2,] 3 1 1 . 3 0
## [3,] 3 0 0 . 0 2
## [4,] 1 1 0 . 2 0
## [5,] 0 0 2 . 1 2
## ... . . . . . .
## [96,] 0 0 1 . 0 0
## [97,] 0 0 1 . 1 3
## [98,] 4 0 2 . 3 1
## [99,] 0 2 2 . 4 1
## [100,] 0 1 0 . 2 1
For example, a hypothetical simulation of a million-cell single-cell RNA-seq dataset might look like this:
ngenes <- 20000
log.abundances <- runif(ngenes, -2, 5)
nclusters <- 20 # define 20 clusters and their population means.
cluster.means <- matrix(2^rnorm(ngenes*nclusters, log.abundances, sd=2), ncol=nclusters)
ncells <- 1e6
clusters <- sample(nclusters, ncells, replace=TRUE) # randomly allocate cells
cell.means <- DelayedArray(cluster.means)[,clusters]
dispersions <- 0.05 + 10/cell.means # typical mean variance trend.
y <- RandomNbinomArray(c(ngenes, ncells), mu=cell.means, size=1/dispersions)
y
## <20000 x 1000000> matrix of class RandomNbinomMatrix and type "double":
## [,1] [,2] [,3] ... [,999999] [,1000000]
## [1,] 5 0 11 . 4 12
## [2,] 0 0 0 . 0 0
## [3,] 0 0 0 . 0 0
## [4,] 1 32 0 . 4 10
## [5,] 0 16 0 . 9 0
## ... . . . . . .
## [19996,] 0 0 0 . 3 0
## [19997,] 0 23 6 . 7 6
## [19998,] 202 4 1 . 290 135
## [19999,] 36 0 15 . 25 7
## [20000,] 1 0 0 . 0 0
Each random DelayedArray
s is broken into contiguous rectangular chunks of identical size and shape.
Each chunk is assigned a seed at construction time that is used to initialize a random number stream (using the PCG32 generator from the dqrng package).
When the user accesses any part of the array, we generate the random numbers in the overlapping chunks and return the desired values.
This provides efficient random access to any subarray without the need to use any jump-ahead functionality.
The chunking scheme determines the efficiency of accessing our random DelayedArray
s.
Chunks that are too large require unnecessary number generation when a subarray is requested, while chunks that are too small would increase memory usage and book-keeping overhead.
The “best” choice also depends on the downstream access pattern, if such information is known.
For example, in a matrix where each column is a chunk, retrieval of a column would be very efficient while retrieval of a single row would be very slow.
The default chunk dimensions are set to the square root of the array dimensions (or 100, whichever is larger), providing a reasonable compromise between all of these considerations.
This can also be manually specified with the chunkdim=
argument.
# Row-wise chunks:
RandomUnifArray(c(1000, 500), chunkdim=c(1, 500))
## <1000 x 500> matrix of class RandomUnifMatrix and type "double":
## [,1] [,2] [,3] ... [,499] [,500]
## [1,] 0.12552122 0.12602151 0.87428992 . 0.09904807 0.53380128
## [2,] 0.85929620 0.55253760 0.36919103 . 0.33372818 0.44120995
## [3,] 0.03062967 0.57340039 0.81855662 . 0.09862350 0.77172037
## [4,] 0.07057872 0.92139388 0.80500709 . 0.61162464 0.29782295
## [5,] 0.40548009 0.07988773 0.47774225 . 0.98474657 0.97728792
## ... . . . . . .
## [996,] 0.45958159 0.70677583 0.41605786 . 0.57063688 0.77065296
## [997,] 0.37620219 0.79944334 0.79581162 . 0.80025884 0.27668267
## [998,] 0.79355880 0.18670673 0.24851211 . 0.86731980 0.15844791
## [999,] 0.13865887 0.45470416 0.11304151 . 0.88390126 0.02011783
## [1000,] 0.02060416 0.89029306 0.95057059 . 0.21307731 0.48039245
# Column-wise chunks:
RandomUnifArray(c(1000, 500), chunkdim=c(1000, 1))
## <1000 x 500> matrix of class RandomUnifMatrix and type "double":
## [,1] [,2] [,3] ... [,499] [,500]
## [1,] 0.8497022 0.4271245 0.9977609 . 0.8648180 0.7704341
## [2,] 0.3171040 0.2048571 0.7795784 . 0.2564039 0.7521729
## [3,] 0.0190134 0.6070411 0.1670825 . 0.6454981 0.1366991
## [4,] 0.5540234 0.2317951 0.4837381 . 0.0372647 0.9248257
## [5,] 0.4856773 0.4325407 0.1739521 . 0.9451870 0.5676749
## ... . . . . . .
## [996,] 0.49502304 0.63233814 0.33384929 . 0.30780392 0.51634514
## [997,] 0.41042033 0.72667507 0.31332396 . 0.78063057 0.74293030
## [998,] 0.38939918 0.88260553 0.99172190 . 0.09292509 0.54855504
## [999,] 0.66750675 0.13490573 0.30198703 . 0.73342482 0.11530273
## [1000,] 0.41536984 0.02050215 0.04207997 . 0.01914548 0.87434696
Unlike other chunk-based DelayedArray
s, the actual values of the random DelayedArray
are dependent on the chunk parameters.
This is because the sampling is done within each chunk and any alteration to the chunk shape or size will rearrange the stream of random numbers within the array.
Thus, even when the seed is set, a different chunkdim
will yield different results:
set.seed(199)
RandomUnifArray(c(10, 5), chunkdim=c(1, 5))
## <10 x 5> matrix of class RandomUnifMatrix and type "double":
## [,1] [,2] [,3] [,4] [,5]
## [1,] 0.4432813900 0.7119991907 0.7257159729 0.9557158216 0.2104000037
## [2,] 0.3134128384 0.0670862596 0.3807734565 0.1526811279 0.1971345709
## [3,] 0.5525409468 0.7807079460 0.9905007351 0.0898367104 0.8541555854
## [4,] 0.8914577304 0.3764661429 0.6028462166 0.2576166289 0.0145520803
## [5,] 0.5042253651 0.2618340398 0.2405915416 0.6800763716 0.5603335327
## [6,] 0.3953960796 0.1314800435 0.9008538304 0.1165704445 0.9035755624
## [7,] 0.6270407308 0.0001234491 0.3542078000 0.3546805161 0.9760325255
## [8,] 0.2882098067 0.9133890264 0.8018615395 0.7615402588 0.6458653482
## [9,] 0.2233627134 0.2957882064 0.7371763836 0.3001928469 0.8730950402
## [10,] 0.6481108814 0.5491673953 0.6821353873 0.4434014931 0.1461723484
set.seed(199)
RandomUnifArray(c(10, 5), chunkdim=c(10, 1))
## <10 x 5> matrix of class RandomUnifMatrix and type "double":
## [,1] [,2] [,3] [,4] [,5]
## [1,] 0.44328139 0.31341284 0.55254095 0.89145773 0.50422537
## [2,] 0.71199919 0.06708626 0.78070795 0.37646614 0.26183404
## [3,] 0.72571597 0.38077346 0.99050074 0.60284622 0.24059154
## [4,] 0.95571582 0.15268113 0.08983671 0.25761663 0.68007637
## [5,] 0.21040000 0.19713457 0.85415559 0.01455208 0.56033353
## [6,] 0.63126383 0.10118984 0.76830283 0.17583353 0.62183470
## [7,] 0.87478047 0.26076972 0.58018989 0.92146566 0.90680388
## [8,] 0.55180537 0.20464839 0.38008429 0.56606812 0.75101891
## [9,] 0.40323090 0.20083487 0.28559824 0.67769322 0.11943279
## [10,] 0.13404760 0.94959186 0.61432837 0.67675354 0.64178865
Like any other random process, the seed should be set to achieve reproducible results.
We stress that the R-level seed only needs to be set before construction of the random DelayedArray
; it is not necessary to set the seed during its use.
This is because the class itself will define further seeds (one per chunk) and store these in the object for use in per-chunk sampling.
set.seed(999)
RandomNormArray(c(10, 5))
## <10 x 5> matrix of class RandomNormMatrix and type "double":
## [,1] [,2] [,3] [,4] [,5]
## [1,] 0.07299897 -0.62967402 0.48686043 -1.00406271 0.09839211
## [2,] -1.87549362 0.11653501 0.98496315 -0.85515988 -1.44932474
## [3,] 0.62150210 -0.24704984 -0.69766033 -1.60089511 0.94752723
## [4,] -1.68035676 0.20939894 -3.29663662 0.57270346 0.94413669
## [5,] 0.22933664 1.12333362 -0.77725398 1.59886411 1.07042538
## [6,] 0.43082191 -0.22740387 1.25381398 0.32842239 0.20448699
## [7,] -0.71424097 0.27773134 1.53748664 0.01761569 1.15470498
## [8,] 0.35908991 0.35618136 -0.49375892 -0.38652760 0.95591577
## [9,] 0.06472627 0.79333115 -0.73803086 -0.56434493 0.22473943
## [10,] -2.41389122 0.76273286 0.47567035 -0.08524682 -1.17849232
set.seed(999)
RandomNormArray(c(10, 5))
## <10 x 5> matrix of class RandomNormMatrix and type "double":
## [,1] [,2] [,3] [,4] [,5]
## [1,] 0.07299897 -0.62967402 0.48686043 -1.00406271 0.09839211
## [2,] -1.87549362 0.11653501 0.98496315 -0.85515988 -1.44932474
## [3,] 0.62150210 -0.24704984 -0.69766033 -1.60089511 0.94752723
## [4,] -1.68035676 0.20939894 -3.29663662 0.57270346 0.94413669
## [5,] 0.22933664 1.12333362 -0.77725398 1.59886411 1.07042538
## [6,] 0.43082191 -0.22740387 1.25381398 0.32842239 0.20448699
## [7,] -0.71424097 0.27773134 1.53748664 0.01761569 1.15470498
## [8,] 0.35908991 0.35618136 -0.49375892 -0.38652760 0.95591577
## [9,] 0.06472627 0.79333115 -0.73803086 -0.56434493 0.22473943
## [10,] -2.41389122 0.76273286 0.47567035 -0.08524682 -1.17849232
For certain distributions, it is possible to indicate that the array is sparse. This does not change the result or efficiency of the sampling process, but can still be useful as it allows downstream functions to use more efficient sparse algorithms. Of course, this is only relevant if the distributional parameters are such that sparsity is actually observed.
RandomPoisArray(c(1e6, 1e6), lambda=0.5) # dense by default
## <1000000 x 1000000> matrix of class RandomPoisMatrix and type "double":
## [,1] [,2] [,3] ... [,999999] [,1000000]
## [1,] 0 0 1 . 1 0
## [2,] 0 1 1 . 0 1
## [3,] 0 1 1 . 0 0
## [4,] 0 1 0 . 1 0
## [5,] 0 0 0 . 1 0
## ... . . . . . .
## [999996,] 1 0 0 . 1 1
## [999997,] 0 1 0 . 1 1
## [999998,] 2 0 1 . 0 0
## [999999,] 1 0 1 . 0 1
## [1000000,] 0 0 0 . 1 0
RandomPoisArray(c(1e6, 1e6), lambda=0.5, sparse=TRUE) # treat as sparse
## <1000000 x 1000000> sparse matrix of class RandomPoisMatrix and type "double":
## [,1] [,2] [,3] ... [,999999] [,1000000]
## [1,] 0 1 1 . 1 0
## [2,] 0 1 1 . 1 0
## [3,] 0 1 1 . 0 0
## [4,] 1 2 0 . 0 0
## [5,] 0 1 0 . 0 2
## ... . . . . . .
## [999996,] 0 0 1 . 0 0
## [999997,] 0 0 0 . 1 0
## [999998,] 3 0 0 . 1 1
## [999999,] 1 0 0 . 1 0
## [1000000,] 0 2 1 . 1 0
sessionInfo()
## R version 4.1.1 (2021-08-10)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.14-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] DelayedRandomArray_1.2.0 DelayedArray_0.20.0 IRanges_2.28.0
## [4] S4Vectors_0.32.0 MatrixGenerics_1.6.0 matrixStats_0.61.0
## [7] BiocGenerics_0.40.0 Matrix_1.3-4 BiocStyle_2.22.0
##
## loaded via a namespace (and not attached):
## [1] Rcpp_1.0.7 knitr_1.36 magrittr_2.0.1
## [4] lattice_0.20-45 R6_2.5.1 rlang_0.4.12
## [7] fastmap_1.1.0 stringr_1.4.0 tools_4.1.1
## [10] grid_4.1.1 xfun_0.27 dqrng_0.3.0
## [13] jquerylib_0.1.4 htmltools_0.5.2 yaml_2.2.1
## [16] digest_0.6.28 bookdown_0.24 BiocManager_1.30.16
## [19] sass_0.4.0 evaluate_0.14 rmarkdown_2.11
## [22] stringi_1.7.5 compiler_4.1.1 bslib_0.3.1
## [25] jsonlite_1.7.2