Contents

1 Overview

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).

2 How can DelayedMatrixStats help me?

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> matrix of class HDF5Matrix and type "integer":
#>           S1  S2  S3  S4 ... S17 S18 S19 S20
#>     [1,]  32  27  32  24   .  18  17  21  32
#>     [2,]  27  20  29  38   .  30  30  22  13
#>     [3,]  40  28  22  17   .  28  31  27  38
#>     [4,]  32  24  30  43   .  24  29  32  31
#>     [5,]  23  16  49  40   .  29  21  50  27
#>      ...   .   .   .   .   .   .   .   .   .
#>  [9996,]  56  28  29  39   .  26  14  38  46
#>  [9997,]  26  31  24  27   .  26  41  21  27
#>  [9998,]  15  37  36  25   .  37  21  15  22
#>  [9999,]  18  22  27  32   .  22  38  29  38
#> [10000,]  20  26  26  39   .  28  40  38  30

Suppose 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 
#> 125.430   2.552 127.987
head(row_sds)
#> [1]  5.353553  8.192423  7.192833  7.103743 10.205649  7.677033

This 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, center = center, : Argument 'x' must be a matrix or a vector.

Unfortunately (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.015   0.000   0.015
head(row_sds)
#> [1]  5.353553  8.192423  7.192833  7.103743 10.205649  7.677033

This 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 In fact, it currently uses matrixStats::rowSds() under the hood. but without having to load the entire data into memory at once2 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.:

library(DelayedMatrixStats)

system.time(row_sds <- rowSds(x))
#>    user  system elapsed 
#>   0.032   0.000   0.032
head(row_sds)
#> [1]  5.353553  8.192423  7.192833  7.103743 10.205649  7.677033

Finally, 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.

3 Supported methods

The initial release of DelayedMatrixStats supports the complete column-wise and row-wise API matrixStats API3 Some of the API is covered via inheritance to functionality in DelayedArray. 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() ☑️

4 ‘Seed-aware’ methods

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:

  1. The ‘block-processing’ strategy. This was developed in the DelayedArray package and is available for all methods in the DelayedMatrixStats through the force_block_processing argument
  2. The ‘seed-aware’ strategy. This is implemented in the DelayedMatrixStats and is optimized for both speed and memory but only for DelayedMatrix objects with certain types of seed.

We 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> matrix of class DelayedMatrix and 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
#>  block_processing 110.81405 113.85501 214.55605 117.25571 127.35435 624.82345
#>        seed_aware  29.63172  30.89894  31.68722  31.71754  32.68202  33.03397
#>  neval cld
#>     10   b
#>     10  a
total(profmem(colSums2(dm_matrix, force_block_processing = TRUE)))
#> Error in profmem_begin(threshold = threshold): Profiling of memory allocations is not supported on this R system (capabilities('profmem') reports FALSE). See help('tracemem'). To enable memory profiling for R on Linux, R needs to be configured and built using './configure --enable-memory-profiling'.
total(profmem(colSums2(dm_matrix)))
#> Error in profmem_begin(threshold = threshold): Profiling of memory allocations is not supported on this R system (capabilities('profmem') reports FALSE). See help('tracemem'). To enable memory profiling for R on Linux, R needs to be configured and built using './configure --enable-memory-profiling'.

# -----------------------------------------------------------------------------
# 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 matrix of class DelayedMatrix and 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
#>  block_processing 277.71071 289.95044 487.07898 308.74965 787.52181 814.22108
#>        seed_aware  17.03122  17.94872  18.34263  18.38639  18.67344  19.78566
#>  neval cld
#>     10   b
#>     10  a
total(profmem(colSums2(dm_dgCMatrix, force_block_processing = TRUE)))
#> Error in profmem_begin(threshold = threshold): Profiling of memory allocations is not supported on this R system (capabilities('profmem') reports FALSE). See help('tracemem'). To enable memory profiling for R on Linux, R needs to be configured and built using './configure --enable-memory-profiling'.
total(profmem(colSums2(dm_dgCMatrix)))
#> Error in profmem_begin(threshold = threshold): Profiling of memory allocations is not supported on this R system (capabilities('profmem') reports FALSE). See help('tracemem'). To enable memory profiling for R on Linux, R needs to be configured and built using './configure --enable-memory-profiling'.

# -----------------------------------------------------------------------------
# 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> matrix of class RleMatrix and 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 770.084952 785.36034 792.697903 793.669717 797.334121
#>        seed_aware   5.665917   5.75912   7.821419   6.199941   7.109136
#>        max neval cld
#>  813.43786    10   b
#>   21.45908    10  a
total(profmem(colSums2(dm_rle, force_block_processing = TRUE)))
#> Error in profmem_begin(threshold = threshold): Profiling of memory allocations is not supported on this R system (capabilities('profmem') reports FALSE). See help('tracemem'). To enable memory profiling for R on Linux, R needs to be configured and built using './configure --enable-memory-profiling'.
total(profmem(colSums2(dm_rle)))
#> Error in profmem_begin(threshold = threshold): Profiling of memory allocations is not supported on this R system (capabilities('profmem') reports FALSE). See help('tracemem'). To enable memory profiling for R on Linux, R needs to be configured and built using './configure --enable-memory-profiling'.

The 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/.

5 Delayed operations

A key feature of a DelayedArray is the ability to register ‘delayed operations’. For example, let’s compute sin(dm_matrix):

system.time(sin_dm_matrix <- sin(dm_matrix))
#>    user  system elapsed 
#>   0.006   0.000   0.006

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:

all.equal(colSums2(sin_dm_matrix), colSums(sin(as.matrix(dm_matrix))))
#> [1] TRUE

6 Roadmap

The initial version of DelayedMatrixStats provides complete coverage of the matrixStats column-wise and row-wise API4 Some of the API is covered via inheritance to functionality in DelayedArray, 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.