DelayedMatrixStats ports the matrixStats API to work with DelayedMatrix objects from the DelayedArray package. It provides high-performing functions operating on rows and columns of DelayedMatrix objects, including all subclasses such as RleArray (from the DelayedArray package) and HDF5Array (from the HDF5Array) as well as supporting all types of seeds, such as matrix (from the base package) and Matrix (from the Matrix package).
The DelayedArray package allows developers to store array-like data using in-memory or on-disk representations (e.g., in HDF5 files) and provides a common and familiar array-like interface for interacting with these data.
The DelayedMatrixStats package is designed to make life easier for Bioconductor developers wanting to use DelayedArray by providing a rich set of column-wise and row-wise summary functions.
We briefly demonstrate and explain these two features using a simple example. We’ll simulate some (unrealistic) RNA-seq read counts data from 10,000 genes and 20 samples and store it on disk as a HDF5Array:
library(DelayedArray)
x <- do.call(cbind, lapply(1:20, function(j) {
rpois(n = 10000, lambda = sample(20:40, 10000, replace = TRUE))
}))
colnames(x) <- paste0("S", 1:20)
x <- realize(x, "HDF5Array")
x
#> <10000 x 20> HDF5Matrix object of type "integer":
#> S1 S2 S3 S4 ... S17 S18 S19 S20
#> [1,] 28 18 20 31 . 27 35 20 24
#> [2,] 26 42 35 42 . 22 23 34 28
#> [3,] 19 31 36 42 . 25 29 42 34
#> [4,] 39 24 36 20 . 37 41 27 59
#> [5,] 38 42 13 33 . 17 40 42 28
#> ... . . . . . . . . .
#> [9996,] 24 35 28 20 . 27 29 16 36
#> [9997,] 43 32 38 17 . 44 28 29 41
#> [9998,] 28 45 19 25 . 22 42 25 22
#> [9999,] 30 25 36 33 . 24 14 22 28
#> [10000,] 35 34 43 22 . 17 37 19 18Suppose you wish to compute the standard deviation of the read counts for each gene.
You might think to use apply() like in the
following:
system.time(row_sds <- apply(x, 1, sd))
#> user system elapsed
#> 134.979 8.302 143.322
head(row_sds)
#> [1] 8.519946 7.319045 8.041177 10.236544 9.786672 8.708586This works, but takes quite a while.
Or perhaps you already know that the matrixStats
package provides a rowSds() function:
matrixStats::rowSds(x)
#> Error in rowVars(x, rows = rows, cols = cols, na.rm = na.rm, refine = refine, : Argument 'x' must be a matrix or a vectorUnfortunately (and perhaps unsurprisingly) this doesn’t work. matrixStats is designed for use on in-memory matrix objects. Well, why don’t we just first realize our data in-memory and then use matrixStats
system.time(row_sds <- matrixStats::rowSds(as.matrix(x)))
#> user system elapsed
#> 0.008 0.000 0.008
head(row_sds)
#> [1] 8.519946 7.319045 8.041177 10.236544 9.786672 8.708586This works and is many times faster than the
apply()-based approach! However, it rather defeats the
purpose of using a HDF5Array for storing the data since we have
to bring all the data into memory at once to compute the result.
Instead, we can use DelayedMatrixStats::rowSds(), which
has the speed benefits of matrixStats::rowSds()1 but
without having to load the entire data into memory at once2:
library(DelayedMatrixStats)
system.time(row_sds <- rowSds(x))
#> user system elapsed
#> 0.017 0.001 0.017
head(row_sds)
#> [1] 8.519946 7.319045 8.041177 10.236544 9.786672 8.708586Finally, by using DelayedMatrixStats
we can use the same code, (colMedians(x)) regardless of
whether the input is an ordinary matrix or a
DelayedMatrix. This is useful for packages wishing to support
both types of objects, e.g., packages wanting to retain backward
compatibility or during a transition period from matrix-based
to DelayeMatrix-based objects.
The initial release of DelayedMatrixStats supports the complete column-wise and row-wise API matrixStats API3. Please see the matrixStats vignette (available online) for a summary these methods. The following table documents the API coverage and availability of ‘seed-aware’ methods in the current version of DelayedMatrixStats, where:
| Method | Block processing | base::matrix optimized | Matrix::dgCMatrix optimized | Matrix::lgCMatrix optimized | DelayedArray::RleArray (SolidRleArraySeed) optimized | DelayedArray::RleArray (ChunkedRleArraySeed) optimized | HDF5Array::HDF5Matrix optimized | base::data.frame optimized | S4Vectors::DataFrame optimized |
|---|---|---|---|---|---|---|---|---|---|
colAlls() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colAnyMissings() |
❌ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colAnyNAs() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colAnys() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colAvgsPerRowSet() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colCollapse() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colCounts() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colCummaxs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colCummins() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colCumprods() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colCumsums() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colDiffs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colIQRDiffs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colIQRs() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colLogSumExps() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colMadDiffs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colMads() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colMaxs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colMeans2() |
✔ | ✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ |
colMedians() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colMins() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colOrderStats() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colProds() |
✔ | ✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ |
colQuantiles() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colRanges() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colRanks() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colSdDiffs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colSds() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colSums2() |
✔ | ✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ |
colTabulates() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colVarDiffs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colVars() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colWeightedMads() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colWeightedMeans() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colWeightedMedians() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
colWeightedSds() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colWeightedVars() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
colsum() |
☑️ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowAlls() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowAnyMissings() |
❌ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowAnyNAs() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowAnys() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowAvgsPerColSet() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowCollapse() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowCounts() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowCummaxs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowCummins() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowCumprods() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowCumsums() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowDiffs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowIQRDiffs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowIQRs() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowLogSumExps() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowMadDiffs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowMads() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowMaxs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowMeans2() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowMedians() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowMins() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowOrderStats() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowProds() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowQuantiles() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowRanges() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowRanks() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowSdDiffs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowSds() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowSums2() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowTabulates() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowVarDiffs() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowVars() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowWeightedMads() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowWeightedMeans() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowWeightedMedians() |
✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowWeightedSds() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowWeightedVars() |
✔ | ✔ | ✔ | ✔ | ❌ | ❌ | ❌ | ❌ | ❌ |
rowsum() |
☑️ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ | ❌ |
As well as offering a familiar API, DelayedMatrixStats provides ‘seed-aware’ methods that are optimized for specific types of DelayedMatrix objects.
To illustrate this idea, we will compare two ways of computing the column sums of a DelayedMatrix object:
force_block_processing argumentWe will demonstrate this by computing the column sums matrices with 20,000 rows and 600 columns where the data have different structure and are stored in DelayedMatrix objects with different types of seed:
We use the microbenchmark package to measure running time and the profmem package to measure the total memory allocations of each method.
In each case, the ‘seed-aware’ method is many times faster and allocates substantially lower total memory.
library(DelayedMatrixStats)
library(sparseMatrixStats)
library(microbenchmark)
library(profmem)
set.seed(666)
# -----------------------------------------------------------------------------
# Dense with values in (0, 1)
# Fast, memory-efficient column sums of DelayedMatrix with ordinary matrix seed
#
# Generate some data
dense_matrix <- matrix(runif(20000 * 600),
nrow = 20000,
ncol = 600)
# Benchmark
dm_matrix <- DelayedArray(dense_matrix)
class(seed(dm_matrix))
#> [1] "matrix" "array"
dm_matrix
#> <20000 x 600> DelayedMatrix object of type "double":
#> [,1] [,2] [,3] ... [,599] [,600]
#> [1,] 0.7743685 0.6601787 0.4098798 . 0.89118118 0.05776471
#> [2,] 0.1972242 0.8436035 0.9198450 . 0.31799523 0.63099417
#> [3,] 0.9780138 0.2017589 0.4696158 . 0.31783791 0.02830454
#> [4,] 0.2013274 0.8797239 0.6474768 . 0.55217184 0.09678816
#> [5,] 0.3612444 0.8158778 0.5928599 . 0.08530977 0.39224147
#> ... . . . . . .
#> [19996,] 0.19490291 0.07763570 0.56391725 . 0.09703424 0.62659353
#> [19997,] 0.61182993 0.01910121 0.04046034 . 0.59708388 0.88389731
#> [19998,] 0.12932744 0.21155070 0.19344085 . 0.51682032 0.13378223
#> [19999,] 0.18985573 0.41716539 0.35110782 . 0.62939661 0.94601427
#> [20000,] 0.87889047 0.25308041 0.54666920 . 0.81630322 0.73272217
microbenchmark(
block_processing = colSums2(dm_matrix, force_block_processing = TRUE),
seed_aware = colSums2(dm_matrix),
times = 10)
#> Unit: milliseconds
#> expr min lq mean median uq max neval
#> block_processing 38.91601 39.37248 106.19388 41.63752 43.33834 369.23688 10
#> seed_aware 12.01726 12.18392 12.18829 12.19826 12.23082 12.24316 10
total(profmem(colSums2(dm_matrix, force_block_processing = TRUE)))
#> [1] 96101408
total(profmem(colSums2(dm_matrix)))
#> [1] 6624
# -----------------------------------------------------------------------------
# Sparse (60% zero) with values in (0, 1)
# Fast, memory-efficient column sums of DelayedMatrix with ordinary matrix seed
#
# Generate some data
sparse_matrix <- dense_matrix
zero_idx <- sample(length(sparse_matrix), 0.6 * length(sparse_matrix))
sparse_matrix[zero_idx] <- 0
# Benchmark
dm_dgCMatrix <- DelayedArray(Matrix(sparse_matrix, sparse = TRUE))
class(seed(dm_dgCMatrix))
#> [1] "dgCMatrix"
#> attr(,"package")
#> [1] "Matrix"
dm_dgCMatrix
#> <20000 x 600> sparse DelayedMatrix object of type "double":
#> [,1] [,2] [,3] ... [,599] [,600]
#> [1,] 0.7743685 0.0000000 0.0000000 . 0.89118118 0.00000000
#> [2,] 0.1972242 0.0000000 0.9198450 . 0.00000000 0.00000000
#> [3,] 0.9780138 0.0000000 0.4696158 . 0.31783791 0.00000000
#> [4,] 0.0000000 0.8797239 0.6474768 . 0.55217184 0.00000000
#> [5,] 0.3612444 0.0000000 0.0000000 . 0.08530977 0.39224147
#> ... . . . . . .
#> [19996,] 0.1949029 0.0776357 0.0000000 . 0.09703424 0.00000000
#> [19997,] 0.0000000 0.0000000 0.0000000 . 0.00000000 0.88389731
#> [19998,] 0.0000000 0.2115507 0.1934408 . 0.00000000 0.00000000
#> [19999,] 0.1898557 0.0000000 0.3511078 . 0.62939661 0.94601427
#> [20000,] 0.8788905 0.2530804 0.0000000 . 0.00000000 0.73272217
microbenchmark(
block_processing = colSums2(dm_dgCMatrix, force_block_processing = TRUE),
seed_aware = colSums2(dm_dgCMatrix),
times = 10)
#> Unit: milliseconds
#> expr min lq mean median uq max neval
#> block_processing 74.35313 79.42809 113.78854 80.80639 82.68261 416.08231 10
#> seed_aware 19.13179 19.23066 19.28002 19.27613 19.30205 19.55056 10
total(profmem(colSums2(dm_dgCMatrix, force_block_processing = TRUE)))
#> [1] 153762792
total(profmem(colSums2(dm_dgCMatrix)))
#> [1] 5408
# -----------------------------------------------------------------------------
# Dense with values in {0, 100} featuring runs of identical values
# Fast, memory-efficient column sums of DelayedMatrix with Rle-based seed
#
# Generate some data
runs <- rep(sample(100, 500000, replace = TRUE), rpois(500000, 100))
runs <- runs[seq_len(20000 * 600)]
runs_matrix <- matrix(runs,
nrow = 20000,
ncol = 600)
# Benchmark
dm_rle <- RleArray(Rle(runs),
dim = c(20000, 600))
class(seed(dm_rle))
#> [1] "SolidRleArraySeed"
#> attr(,"package")
#> [1] "DelayedArray"
dm_rle
#> <20000 x 600> RleMatrix object of type "integer":
#> [,1] [,2] [,3] [,4] ... [,597] [,598] [,599] [,600]
#> [1,] 20 65 9 50 . 43 74 85 14
#> [2,] 20 65 9 50 . 43 74 85 14
#> [3,] 20 65 9 50 . 43 74 85 14
#> [4,] 20 65 9 50 . 43 74 85 14
#> [5,] 20 65 9 50 . 43 74 85 14
#> ... . . . . . . . . .
#> [19996,] 65 9 50 16 . 74 85 14 55
#> [19997,] 65 9 50 16 . 74 85 14 55
#> [19998,] 65 9 50 16 . 74 85 14 55
#> [19999,] 65 9 50 16 . 74 85 14 55
#> [20000,] 65 9 50 16 . 74 85 14 55
microbenchmark(
block_processing = colSums2(dm_rle, force_block_processing = TRUE),
seed_aware = colSums2(dm_rle),
times = 10)
#> Unit: milliseconds
#> expr min lq mean median uq
#> block_processing 313.542707 314.465297 348.974150 317.334165 318.852255
#> seed_aware 2.774472 3.043012 4.274285 3.123883 3.704185
#> max neval
#> 642.07149 10
#> 13.49501 10
total(profmem(colSums2(dm_rle, force_block_processing = TRUE)))
#> [1] 144257360
total(profmem(colSums2(dm_rle)))
#> [1] 57400The development of ‘seed-aware’ methods is ongoing work (see the Roadmap), and for now only a few methods and seed-types have a ‘seed-aware’ method.
An extensive set of benchmarks is under development at http://peterhickey.org/BenchmarkingDelayedMatrixStats/.
A key feature of a DelayedArray is the ability to register
‘delayed operations’. For example, let’s compute
sin(dm_matrix):
This instantaneous because the operation is not actually performed,
rather it is registered and only performed when the object is
realized. All methods in DelayedMatrixStats
will correctly realise these delayed operations before computing the
final result. For example, let’s compute
colSums2(sin_dm_matrix) and compare check we get the
correct answer:
The initial version of DelayedMatrixStats provides complete coverage of the matrixStats column-wise and row-wise API4, allowing package developers to use these functions with DelayedMatrix objects as well as with ordinary matrix objects. This should simplify package development and assist authors to support to their software for large datasets stored in disk-backed data structures such as HDF5Array. Such large datasets are increasingly common with the rise of single-cell genomics.
Future releases of DelayedMatrixStats
will improve the performance of these methods, specifically by
developing additional ‘seed-aware’ methods. The plan is to prioritise
commonly used methods (e.g.,
colMeans2()/rowMeans2(),
colSums2()/rowSums2(), etc.) and the
development of ‘seed-aware’ methods for the HDF5Matrix class.
To do so, we will leverage the beachmat
package. Proof-of-concept code has shown that this can greatly increase
the performance when analysing such disk-backed data.
Importantly, all package developers using methods from DelayedMatrixStats will immediately gain from performance improvements to these low-level routines. By using DelayedMatrixStats, package developers will be able to focus on higher level programming tasks and address important scientific questions and technological challenges in high-throughput biology.
sessionInfo()
#> R version 4.5.2 (2025-10-31)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.3 LTS
#>
#> Matrix products: default
#> BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
#> LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_US.UTF-8 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
#>
#> time zone: Etc/UTC
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats4 stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] profmem_0.7.0 microbenchmark_1.5.0
#> [3] sparseMatrixStats_1.23.0 DelayedMatrixStats_1.32.0
#> [5] DelayedArray_0.37.0 SparseArray_1.11.1
#> [7] S4Arrays_1.11.0 abind_1.4-8
#> [9] IRanges_2.45.0 S4Vectors_0.49.0
#> [11] MatrixGenerics_1.23.0 matrixStats_1.5.0
#> [13] BiocGenerics_0.56.0 generics_0.1.4
#> [15] Matrix_1.7-4 BiocStyle_2.38.0
#>
#> loaded via a namespace (and not attached):
#> [1] jsonlite_2.0.0 compiler_4.5.2 BiocManager_1.30.26
#> [4] Rcpp_1.1.0 rhdf5filters_1.23.0 jquerylib_0.1.4
#> [7] yaml_2.3.10 fastmap_1.2.0 lattice_0.22-7
#> [10] R6_2.6.1 XVector_0.51.0 knitr_1.50
#> [13] maketools_1.3.2 h5mread_1.3.0 bslib_0.9.0
#> [16] rlang_1.1.6 cachem_1.1.0 HDF5Array_1.39.0
#> [19] xfun_0.54 sass_0.4.10 sys_3.4.3
#> [22] cli_3.6.5 Rhdf5lib_1.33.0 digest_0.6.37
#> [25] grid_4.5.2 rhdf5_2.55.7 lifecycle_1.0.4
#> [28] evaluate_1.0.5 buildtools_1.0.0 rmarkdown_2.30
#> [31] tools_4.5.2 htmltools_0.5.8.1In fact, it currently uses
matrixStats::rowSds() under the hood.↩︎
In this case, it loads blocks of data row-by-row. The
amount of data loaded into memory at any one time is controlled by the
default block size global setting; see
?DelayedArray::getAutoBlockSize for details. Notably, if
the data are small enough (and the default block size is large enough)
then all the data is loaded as a single block, but this approach
generalizes and still works when the data are too large to be loaded
into memory in one block.↩︎
Some of the API is covered via inheritance to functionality in DelayedArray↩︎
Some of the API is covered via inheritance to functionality in DelayedArray↩︎