File: DelayedArray-class.Rd

package info (click to toggle)
r-bioc-delayedarray 0.8.0%2Bdfsg-2
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 980 kB
  • sloc: ansic: 93; makefile: 2; sh: 1
file content (412 lines) | stat: -rw-r--r-- 14,207 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
\name{DelayedArray-class}
\docType{class}

\alias{class:DelayedArray}
\alias{DelayedArray-class}

\alias{class:DelayedMatrix}
\alias{DelayedMatrix-class}
\alias{DelayedMatrix}

\alias{coerce,DelayedArray,DelayedMatrix-method}
\alias{coerce,DelayedMatrix,DelayedArray-method}

\alias{dim,DelayedArray-method}
\alias{dimnames,DelayedArray-method}
\alias{extract_array,DelayedArray-method}

\alias{new_DelayedArray}
\alias{DelayedArray}
\alias{DelayedArray,ANY-method}
\alias{DelayedArray,DelayedArray-method}
\alias{DelayedArray,DelayedOp-method}

\alias{class:DelayedArray1}
\alias{DelayedArray1-class}
\alias{DelayedArray1}
\alias{updateObject,DelayedArray-method}

\alias{type}

\alias{nseed}
\alias{nseed,ANY-method}

\alias{seed}
\alias{seed,DelayedOp-method}
\alias{seed<-}
\alias{seed<-,DelayedOp-method}

\alias{path}
\alias{path,DelayedArray-method}
\alias{path<-,DelayedArray-method}

\alias{aperm}
\alias{aperm.DelayedArray}
\alias{aperm,DelayedArray-method}

\alias{dim<-,DelayedArray-method}
\alias{dimnames<-,DelayedArray-method}
\alias{names,DelayedArray-method}
\alias{names<-,DelayedArray-method}

\alias{drop,DelayedArray-method}
\alias{[,DelayedArray-method}
\alias{[<-,DelayedArray-method}

\alias{coerce,DelayedArray,SparseArraySeed-method}
\alias{coerce,DelayedMatrix,dgCMatrix-method}
\alias{coerce,DelayedMatrix,sparseMatrix-method}

\alias{[[,DelayedArray-method}

\alias{show,DelayedArray-method}

\alias{c,DelayedArray-method}
\alias{splitAsList,DelayedArray-method}
\alias{split.DelayedArray}
\alias{split,DelayedArray,ANY-method}

% Internal stuff
\alias{matrixClass}
\alias{matrixClass,DelayedArray-method}

\title{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
nseed(x)            # seed counter
path(object, ...)   # path getter
type(x)
}

\arguments{
  \item{seed}{
    An array-like object.
  }
  \item{x, object}{
    A DelayedArray object. For \code{type()}, \code{x} can also be any
    array-like object, that is, any object for which \code{dim(x)} is not
    NULL.
  }
  \item{...}{
    Additional arguments passed to methods.
  }
}

\section{In-memory versus on-disk realization}{
  To \emph{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 \code{as.array} on it. However this realizes the
  full object at once \emph{in memory} which could require too much memory
  if the object is big. A big DelayedArray object is preferrably realized
  \emph{on disk} e.g. by calling \code{\link[HDF5Array]{writeHDF5Array}} on
  it (this function is defined in the \pkg{HDF5Array} package) or coercing it
  to an \link[HDF5Array]{HDF5Array} object with \code{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 \code{?\link[HDF5Array]{writeHDF5Array}} in the \pkg{HDF5Array} package
  for more information about this.
}

\section{Accessors}{
  DelayedArray objects support the same set of getters as ordinary arrays
  i.e. \code{dim()}, \code{length()}, and \code{dimnames()}.
  In addition, they support \code{seed()}, \code{nseed()}, \code{path()},
  and \code{type()}.

  \code{type()} is the DelayedArray equivalent of \code{typeof()} (or
  \code{storage.mode()}) for ordinary arrays and vectors. Note that, for
  convenience and consistency, \code{type()} also supports ordinary arrays
  and vectors. It also supports any array-like object, that is, any object
  \code{x} for which \code{dim(x)} is not NULL.

  \code{dimnames()}, \code{seed()}, and \code{path()} also work as setters.
}

\section{Subsetting}{
  A DelayedArray object can be subsetted with \code{[} like an ordinary array,
  but with the following differences:
  \itemize{
    \item \emph{Multi-dimensional single bracket subsetting} (i.e. subsetting
          of the form \code{x[i_1, i_2, ..., i_n]} with one (possibly missing)
          subscript per dimension) returns a DelayedArray object where the
          subsetting is actually delayed. So it's a very light operation.
          One notable exception to this is when \code{drop=TRUE} and the
          result has only one dimension, in which case it is returned as an
          ordinary vector (atomic or list).
          Note that NAs in the subscripts are not supported.

    \item \emph{Linear single bracket subsetting} (a.k.a. 1D-style subsetting,
          that is, subsetting of the form \code{x[i]}) only works if the
          subscript \code{i} is a numeric vector at the moment. Furthermore,
          \code{i} cannot contain NAs and all the indices in it must be >= 1
          and <= \code{length(x)} for now. It returns an atomic vector of the
          same length as \code{i}. This is NOT a delayed operation (block
          processing is triggered).
  }

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

  Subassignment to a DelayedArray object with \code{[<-} is also supported
  like with an ordinary array, but with the following restrictions:
  \itemize{
    \item \emph{Multi-dimensional subassignment} (i.e. subassignment of the
          form \code{x[i_1, i_2, ..., i_n] <- value} with one (possibly
          missing) subscript per dimension) only accepts a replacement
          value (a.k.a. right value) that is an array-like object (e.g.
          ordinary array, dgCMatrix object, DelayedArray object, etc...)
          or an ordinary vector (atomic or list) of length 1.

    \item \emph{Linear subassignment} (a.k.a. 1D-style subassignment, that
          is, subassignment of the form \code{x[i] <- value}) only works if
          the subscript \code{i} is a logical DelayedArray object of the same
          dimensions as \code{x} and if the replacement value is an ordinary
          vector (atomic or list) of length 1.

    \item \emph{Filling with a vector}, that is, subassignment of the form
          \code{x[] <- v} where \code{v} is an ordinary vector (atomic or
          list), is only supported if the length of the vector is a divisor
          of \code{nrow(x)}.
  }
  These 3 forms of subassignment are implemented as delayed operations so
  are very light.

  Single value replacement (\code{x[[...]] <- value}) is not supported yet.
}

\seealso{
  \itemize{
    \item \code{\link{realize}} for realizing a DelayedArray object in memory
          or on disk.

    \item \link{block_processing} for more information about block processing
          of an array-like object.

    \item \link{DelayedArray-utils} for common operations on
          \link{DelayedArray} objects.

    \item \link{DelayedArray-stats} for statistical functions on
          DelayedArray objects.

    \item \link{DelayedMatrix-stats} for \link{DelayedMatrix} row/col
          summarization.

    \item \link{RleArray} objects.

    \item \link[HDF5Array]{HDF5Array} objects in the \pkg{HDF5Array} package.

    \item \link[S4Vectors]{DataFrame} objects in the \pkg{S4Vectors} package.

    \item \link[base]{array} objects in base R.
  }
}

\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 DelayedArray object is **always** treated as a
## "read-only" object so will never be modified by the operations
## we perform on A:
stopifnot(identical(a, seed(A)))
type(A)

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

## Linear single bracket subsetting:
A[11:20]
A[A <= 1e-5]
stopifnot(identical(a[a <= 1e-5], A[A <= 1e-5]))

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

A[2:5, 1:2, ] <- array(1:40, c(4, 2, 5))
a[2:5, 1:2, ] <- array(1:40, c(4, 2, 5))
stopifnot(identical(a, as.array(A)))

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

B <- crazy(A)  # very fast! (all operations are delayed)
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 memory-efficient
colSums(A2)

## ---------------------------------------------------------------------
## C. AN 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)                   # an HDF5ArraySeed object

B3 <- crazy(A3)            # very fast! (all operations are delayed)
B3                         # not an HDF5Array object anymore because
                           # now it carries delayed operations
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)
path(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. MODIFY THE PATH OF A DelayedArray OBJECT
## ---------------------------------------------------------------------
## This can be useful if the file containing the array data is on a
## shared partition but the exact path to the partition depends on the
## machine from which the data is being accessed.
## For example:

\dontrun{
library(HDF5Array)
A <- HDF5Array("/path/to/lab_data/my_precious_data.h5")
path(A)

## Operate on A...
## Now A carries delayed operations.
## Make sure path(A) still works:
path(A)

## Save A:
save(A, file="A.rda")

## A.rda should be small (it doesn't contain the array data).
## Send it to a co-worker that has access to my_precious_data.h5.

## Co-worker loads it:
load("A.rda")
path(A)

## A is broken because path(A) is incorrect for co-worker:
A  # error!

## Co-worker fixes the path (in this case this is better done using the
## dirname() setter rather than the path() setter):
dirname(A) <- "E:/other/path/to/lab_data"

## A "works" again:
A
}

## ---------------------------------------------------------------------
## F. WRAP A SPARSE MATRIX IN A DelayedArray OBJECT
## ---------------------------------------------------------------------
\dontrun{
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))) %*% p

## 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 check the automatic block size:
getAutoBlockSize()

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

## Increasing the block size increases the speed but also memory usage:
setAutoBlockSize(2e8)
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))

setAutoBlockSize()
}
}
\keyword{classes}
\keyword{methods}