DelayedArray-class {DelayedArray}R Documentation

DelayedArray objects

Description

Wrapping an array-like object (typically an on-disk object) in a DelayedArray object allows one to perform common array operations on it without loading the object in memory. In order to reduce memory usage and optimize performance, operations on the object are either delayed or executed using a block processing mechanism.

Usage

DelayedArray(seed)  # constructor function
seed(x)             # seed getter
type(x)

Arguments

seed

An array-like object.

x

A DelayedArray object. (Can also be an ordinary array in case of type.)

In-memory versus on-disk realization

To realize a DelayedArray object (i.e. to trigger execution of the delayed operations carried by the object and return the result as an ordinary array), call as.array on it. However this realizes the full object at once in memory which could require too much memory if the object is big. A big DelayedArray object is preferrably realized on disk e.g. by calling writeHDF5Array on it (this function is defined in the HDF5Array package) or coercing it to an HDF5Array object with as(x, "HDF5Array"). Other on-disk backends can be supported. This uses a block-processing strategy so that the full object is not realized at once in memory. Instead the object is processed block by block i.e. the blocks are realized in memory and written to disk one at a time. See ?writeHDF5Array in the HDF5Array package for more information about this.

Accessors

DelayedArray objects support the same set of getters as ordinary arrays i.e. dim(), length(), and dimnames(). In addition, they support type(), which is the DelayedArray equivalent of typeof() or storage.mode() for ordinary arrays. Note that, for convenience and consistency, type() also works on ordinary arrays.

Only dimnames() is supported as a setter.

Subsetting

A DelayedArray object can be subsetted with [ like an ordinary array but with the following differences:

Subsetting with [[ is supported but only the linear form of it at the moment i.e. the x[[i]] form where i is a single numeric value >= 1 and <= length(x). It is equivalent to x[i].

DelayedArray objects support only 2 forms of subassignment at the moment: x[i] <- value and x[] <- value. The former is supported only when the subscript i is a logical DelayedArray object with the same dimensions as x and when value is a scalar (i.e. an atomic vector of length 1). The latter is supported only when value is an atomic vector and length(value) is a divisor of nrow(x). Both are delayed operations so are very light.

Single value replacement (x[[...]] <- value) is not supported.

Binding

Binding DelayedArray objects along the rows (or columns) is supported via the rbind and arbind (or cbind and acbind) methods for DelayedArray objects. All these operations are delayed.

See Also

Examples

## ---------------------------------------------------------------------
## A. WRAP AN ORDINARY ARRAY IN A DelayedArray OBJECT
## ---------------------------------------------------------------------
a <- array(runif(1500000), dim=c(10000, 30, 5))
A <- DelayedArray(a)
A
## The seed of A is treated as a "read-only" object so won't change when
## we start operating on A:
stopifnot(identical(a, seed(A)))
type(A)

## Multi-dimensional single bracket subsetting:
m <- a[11:20 , 5, ]        # a matrix
A[11:20 , 5, ]             # not a DelayedMatrix (still 3 dimensions)
M <- drop(A[11:20 , 5, ])  # a DelayedMatrix object
stopifnot(identical(m, as.array(M)))
stopifnot(identical(a, seed(M)))

## Linear single bracket subsetting:
A[11:20]
A[which(A <= 1e-5)]

## Subassignment:
A[A < 0.2] <- NA
a[a < 0.2] <- NA
stopifnot(identical(a, as.array(A)))

## Other operations:
toto <- function(x) (5 * x[ , , 1] ^ 3 + 1L) * log(x[, , 2])
b <- toto(a)
head(b)

B <- toto(A)  # very fast! (operations are delayed)
B  # still 3 dimensions (subsetting a DelayedArray object never drops
   # dimensions)
B <- drop(B)
B

cs <- colSums(b)
CS <- colSums(B)
stopifnot(identical(cs, CS))

## ---------------------------------------------------------------------
## B. WRAP A DataFrame OBJECT IN A DelayedArray OBJECT
## ---------------------------------------------------------------------

## Generate random coverage and score along an imaginary chromosome:
cov <- Rle(sample(20, 5000, replace=TRUE), sample(6, 5000, replace=TRUE))
score <- Rle(sample(100, nrun(cov), replace=TRUE), runLength(cov))

DF <- DataFrame(cov, score)
A2 <- DelayedArray(DF)
A2
seed(A2)  # 'DF'

## Coercion of a DelayedMatrix object to DataFrame produces a DataFrame
## object with Rle columns:
as(A2, "DataFrame")
stopifnot(identical(DF, as(A2, "DataFrame")))

t(A2)  # transposition is delayed so is very fast and very memory
       # efficient
stopifnot(identical(DF, seed(t(A2))))  # the "seed" is still the same
colSums(A2)

## ---------------------------------------------------------------------
## C. A HDF5Array OBJECT IS A (PARTICULAR KIND OF) DelayedArray OBJECT
## ---------------------------------------------------------------------
library(HDF5Array)
A3 <- as(a, "HDF5Array")   # write 'a' to an HDF5 file
A3
is(A3, "DelayedArray")     # TRUE
seed(A3)                   # a HDF5ArraySeed object

B3 <- toto(A3)             # very fast! (operations are delayed)

B3                         # not a HDF5Array object because now it
                           # carries delayed operations
B3 <- drop(B3)

CS3 <- colSums(B3)
stopifnot(identical(cs, CS3))

## ---------------------------------------------------------------------
## D. PERFORM THE DELAYED OPERATIONS
## ---------------------------------------------------------------------
as(B3, "HDF5Array")        # "realize" 'B3' on disk

## If this is just an intermediate result, you can either keep going
## with B3 or replace it with its "realized" version:
B3 <- as(B3, "HDF5Array")  # no more delayed operations on new 'B3'
seed(B3)

## For convenience, realize() can be used instead of explicit coercion.
## The current "realization backend" controls where realization
## happens e.g. in memory if set to NULL or in an HDF5 file if set
## to "HDF5Array":
D <- cbind(B3, exp(B3))
D
setRealizationBackend("HDF5Array")
D <- realize(D)
D
## See '?realize' for more information about "realization backends".

## ---------------------------------------------------------------------
## E. BIND DelayedArray OBJECTS
## ---------------------------------------------------------------------

## rbind/cbind

library(HDF5Array)
toy_h5 <- system.file("extdata", "toy.h5", package="HDF5Array")
h5ls(toy_h5)

M1 <- HDF5Array(toy_h5, "M1")
M2 <- HDF5Array(toy_h5, "M2")

M12 <- rbind(M1, t(M2))
M12
colMeans(M12)

## arbind/acbind

example(acbind)  # to create arrays a1, a2, a3

A1 <- DelayedArray(a1)
A2 <- DelayedArray(a2)
A3 <- DelayedArray(a3)
A <- arbind(A1, A2, A3)
A

## Sanity check:
stopifnot(identical(arbind(a1, a2, a3), as.array(A)))

## ---------------------------------------------------------------------
## F. WRAP A SPARSE MATRIX IN A DelayedArray OBJECT
## ---------------------------------------------------------------------
## Not run: 
library(Matrix)
M <- 75000L
N <- 1800L
p <- sparseMatrix(sample(M, 9000000, replace=TRUE),
                  sample(N, 9000000, replace=TRUE),
                  x=runif(9000000), dims=c(M, N))
P <- DelayedArray(p)
P
p2 <- as(P, "sparseMatrix")
stopifnot(identical(p, p2))

## The following is based on the following post by Murat Tasan on the
## R-help mailing list:
##   https://stat.ethz.ch/pipermail/r-help/2017-May/446702.html

## As pointed out by Murat, the straight-forward row normalization
## directly on sparse matrix 'p' would consume too much memory:
row_normalized_p <- p / rowSums(p^2)  # consumes too much memory
## because the rowSums() result is being recycled (appropriately) into a
## *dense* matrix with dimensions equal to dim(p).

## Murat came up with the following solution that is very fast and memory
## efficient:
row_normalized_p1 <- Diagonal(x=1/sqrt(Matrix::rowSums(p^2))) 

## With a DelayedArray object, the straight-forward approach uses a
## block processing strategy behind the scene so it doesn't consume
## too much memory.

## First, let's see the block processing in action:
DelayedArray:::set_verbose_block_processing(TRUE)
## and set block size to a bigger value than the default:
getOption("DelayedArray.block.size")
options(DelayedArray.block.size=80e6)

row_normalized_P <- P / sqrt(DelayedArray::rowSums(P^2))

## Increasing the block size increases the speed but also memory usage:
options(DelayedArray.block.size=200e6)
row_normalized_P2 <- P / sqrt(DelayedArray::rowSums(P^2))
stopifnot(all.equal(row_normalized_P, row_normalized_P2))

## Back to sparse representation:
DelayedArray:::set_verbose_block_processing(FALSE)
row_normalized_p2 <- as(row_normalized_P, "sparseMatrix")
stopifnot(all.equal(row_normalized_p1, row_normalized_p2))

options(DelayedArray.block.size=10e6)

## End(Not run)

[Package DelayedArray version 0.4.1 Index]