RealizationSink objects
RealizationSink-class.RdUse a RealizationSink object in combination with
write_block() to write blocks of array
data to disk.
RealizationSink is a virtual class with various concrete subclasses that support writing data into specific formats.
sinkApply() is a convenience function for walking on a
RealizationSink object, typically for the purpose of filling it
with blocks of data.
Note that write_block() is typically used
inside the callback function passed to sinkApply().
Usage
## Walk on a RealizationSink derivative:
sinkApply(sink, FUN, ..., grid=NULL, verbose=NA)
## Backend-agnostic RealizationSink constructor:
AutoRealizationSink(dim, dimnames=NULL, type="double", as.sparse=FALSE)
## Get/set the "automatic realization backend":
getAutoRealizationBackend()
setAutoRealizationBackend(BACKEND=NULL)
registeredRealizationBackends()Arguments
- sink
A **writable** array-like object, typically a RealizationSink derivative. Some important notes:
DelayedArray objects are NEVER writable, even when they don't carry delayed operations (e.g. HDF5Array objects from the HDF5Array package), and even when they don't carry delayed operations **and** have all their data in memory (e.g. RleArray objects). In other words, there are NO exceptions.
RealizationSink is a **virtual** class so
sinkwill always be a RealizationSink **derivative**, that is, an object that belongs to a **concrete** subclass of the RealizationSink class (e.g. an HDF5RealizationSink object from the HDF5Array package).RealizationSink derivatives are considered array-like objects i.e. they have dimensions and possibly dimnames.
Although
write_block()andsinkApply()will typically be used on a RealizationSink derivative, they can also be used on an ordinary array or other writable in-memory array-like objects like dgCMatrix objects from the Matrix package.- FUN
The callback function to apply to each **viewport** of the grid used to walk on
sink.sinkApply()will performsink <- FUN(sink, viewport, ...)on each viewport, soFUNmust take at least two arguments, typicallysinkandviewport(but the exact names can differ).The function is expected to return its 1st argument (
sink) possibly modified (e.g. whenFUNcontains a call towrite_block(), which is typically the case).- ...
Additional arguments passed to
FUN.- grid
The grid used for the walk, that is, an ArrayGrid object that defines the viewports to walk on. It must be compatible with the geometry of
sink. If not specified, an automatic grid is created by callingdefaultSinkAutoGrid(sink), and used. See?defaultSinkAutoGridfor more information.- verbose
Whether block processing progress should be displayed or not. If set to
NA(the default), verbosity is controlled byDelayedArray:::get_verbose_block_processing(). SettingverbosetoTRUEorFALSEoverrides this.- dim
The dimensions (specified as an integer vector) of the RealizationSink derivative to create.
- dimnames
The dimnames (specified as a list of character vectors or NULLs) of the RealizationSink derivative to create.
- type
The type of the data that will be written to the RealizationSink derivative to create.
- as.sparse
Whether the data should be written as sparse or not to the RealizationSink derivative to create. Not all realization backends support this.
- BACKEND
NULL(the default), or a single string specifying the name of a realization backend e.g."HDF5Array"or"RleArray"etc...
Details
*** The RealizationSink API ***
The DelayedArray package provides a simple API for writing blocks of array data to disk (or to memory): the "RealizationSink API". This API allows the developper to write code that is agnostic about the particular on-disk (or in-memory) format being used to store the data.
Here is how to use it:
Create a realization sink.
Write blocks of array data to the realization sink with one or several calls to
write_block().Close the realization sink with
close().Coerce the realization sink to DelayedArray.
A realization sink is formally represented by a RealizationSink derivative. Note that RealizationSink is a virtual class with various concrete subclasses like HDF5RealizationSink from the HDF5Array package, or RleRealizationSink. Each subclass implements the "RealizationSink API" for a specific realization backend.
To create a realization sink, use the specific constructor function.
This function should be named as the class itself e.g.
HDF5RealizationSink().
To create a realization sink in a backend-agnostic way, use
AutoRealizationSink(). It will create a RealizationSink derivative
for the current automatic realization backend (see below).
Once writing to the realization sink is completed, the RealizationSink
derivative must be closed (with close(sink)), then coerced to
DelayedArray (with as(sink, "DelayedArray"). What
specific DelayedArray derivative this coercion will return
depends on the specific class of the RealizationSink derivative. For
example, if sink is an HDF5RealizationSink
object from the HDF5Array package, then as(sink, "DelayedArray")
will return an HDF5Array instance (the
HDF5Array class is a DelayedArray subclass).
*** The automatic realization backend ***
The automatic realization backend is a user-controlled global
setting that indicates what specific RealizationSink derivative
AutoRealizationSink() should return.
In the context of block processing of a DelayedArray object,
this controls where/how realization happens e.g. as an ordinary array
if not set (i.e. set to NULL), or as an HDF5Array
object if set to "HDF5Array", or as an RleArray object
if set to "RleArray", etc...
Use getAutoRealizationBackend() or setAutoRealizationBackend()
to get or set the automatic realization backend.
Use registeredRealizationBackends() to get the list of realization
backends that are currently registered.
*** Cross realization backend compatibility ***
Two important things to keep in mind for developers aiming at writing code that is compatible across realization backends:
Realization backends don't necessarily support concurrent writing.
More precisely: Even though it is safe to assume that any DelayedArray object will support concurrent
read_block()calls, it is not so safe to assume that any RealizationSink derivative will support concurrent calls towrite_block(). For example, at the moment, HDF5RealizationSink objects do not support concurrent writing.This means that in order to remain compatible across realization backends, code that contains calls to
write_block()should NOT be parallelized.Some realization backends are "linear write only", that is, they don't support random write access, only linear write access.
Such backends will provide a relization sink where the blocks of data must be written in linear order (i.e. by ascending rank). Furthermore, the geometry of the blocks must also be compatible with linear write access, that is, they must have a "first-dim-grows-first" shape. Concretely this means that the grid used to walk on the relization sink must be created with something like:
colAutoGrid(sink)for a two-dimensional sink, or with something like:
defaultSinkAutoGrid(sink)for a sink with an arbitrary number of dimensions.
See
?defaultSinkAutoGridfor more information.For obvious reasons, "linear write only" realization backends do not support concurrent writing.
Value
For sinkApply(), its 1st argument (sink) possibly
modified (e.g. when callback function FUN contains a call to
write_block(), which is typically the case).
For AutoRealizationSink(), a RealizationSink derivative with the
class associated with the current automatic realization backend.
For getAutoRealizationBackend, NULL (no backend set yet)
or a single string specifying the name of the automatic realization
backend currently in use.
For registeredRealizationBackends, a data frame with 1 row
per registered realization backend.
See also
read_blockandwrite_blockin the S4Arrays package.ArrayGrid in the S4Arrays package for the formal representation of grids and viewports.
defaultSinkAutoGridto create an automatic grid on a RealizationSink derivative.blockApplyand family for convenient block processing of an array-like object.HDF5RealizationSink objects in the HDF5Array package.
HDF5-dump-management in the HDF5Array package to control the location and physical properties of automatically created HDF5 datasets.
RleArray objects.
DelayedArray objects.
array objects in base R.
Examples
## ---------------------------------------------------------------------
## USING THE "RealizationSink API": EXAMPLE 1
## ---------------------------------------------------------------------
## -- STEP 1 --
## Create a realization sink. Note that instead of creating a
## realization sink by calling a backend-specific sink constructor
## (e.g. HDF5Array::HDF5RealizationSink), we set the "automatic
## realization backend" to "HDF5Array" and use backend-agnostic
## constructor AutoRealizationSink():
setAutoRealizationBackend("HDF5Array")
sink <- AutoRealizationSink(c(35L, 50L, 8L))
dim(sink)
#> [1] 35 50 8
## -- STEP 2 --
## Define the grid of viewports to walk on. Here we define a grid made
## of very small viewports on 'sink'. Note that, for real-world use cases,
## block processing will typically use grids made of much bigger
## viewports, usually obtained with defaultSinkAutoGrid().
## Also please note that this grid would not be compatible with "linear
## write only" realization backends. See "Cross realization backend
## compatibility" above in this man page for more information.
sink_grid <- RegularArrayGrid(dim(sink), spacings=c(20, 20, 4))
## -- STEP 3 --
## Walk on the grid, and, for each viewport, write random data to it.
for (bid in seq_along(sink_grid)) {
viewport <- sink_grid[[bid]]
block <- array(runif(length(viewport)), dim=dim(viewport))
sink <- write_block(sink, viewport, block)
}
## -- An alternative to STEP 3 --
FUN <- function(sink, viewport) {
block <- array(runif(length(viewport)), dim=dim(viewport))
write_block(sink, viewport, block)
}
sink <- sinkApply(sink, FUN, grid=sink_grid, verbose=TRUE)
#> \ processing viewport 1/12 ...
#> ok
#> \ processing viewport 2/12 ...
#> ok
#> \ processing viewport 3/12 ...
#> ok
#> \ processing viewport 4/12 ...
#> ok
#> \ processing viewport 5/12 ...
#> ok
#> \ processing viewport 6/12 ...
#> ok
#> \ processing viewport 7/12 ...
#> ok
#> \ processing viewport 8/12 ...
#> ok
#> \ processing viewport 9/12 ...
#> ok
#> \ processing viewport 10/12 ...
#> ok
#> \ processing viewport 11/12 ...
#> ok
#> \ processing viewport 12/12 ...
#> ok
## -- STEP 4 --
## Close the sink and turn it into a DelayedArray object:
close(sink)
A <- as(sink, "DelayedArray")
A
#> <35 x 50 x 8> HDF5Array object of type "double":
#> ,,1
#> [,1] [,2] [,3] ... [,49] [,50]
#> [1,] 0.58525713 0.32895605 0.12428716 . 0.8776101 0.1236370
#> [2,] 0.52524119 0.76037904 0.06078881 . 0.8527295 0.4081144
#> ... . . . . . .
#> [34,] 0.88538694 0.91002592 0.10621641 . 0.5801200 0.9171059
#> [35,] 0.05126328 0.48823863 0.82587009 . 0.9908749 0.6969871
#>
#> ...
#>
#> ,,8
#> [,1] [,2] [,3] ... [,49] [,50]
#> [1,] 0.98410378 0.93937357 0.56783511 . 0.64869086 0.47334866
#> [2,] 0.44225391 0.21997686 0.08355569 . 0.20312493 0.07917907
#> ... . . . . . .
#> [34,] 0.6065857 0.4595153 0.4679179 . 0.4527369 0.6612981
#> [35,] 0.9599744 0.4690177 0.1221533 . 0.8188423 0.7592227
#>
setAutoRealizationBackend() # restore default (NULL)
## ---------------------------------------------------------------------
## USING THE "RealizationSink API": EXAMPLE 2
## ---------------------------------------------------------------------
## Say we have a 3D array and want to collapse its 3rd dimension by
## summing the array elements that are stacked vertically, that is, we
## want to compute the matrix M obtained by doing sum(A[i, j, ]) for all
## valid i and j. This is very easy to do with an ordinary array:
collapse_3rd_dim <- function(a) apply(a, MARGIN=1:2, sum)
## or, in a slightly more efficient way:
collapse_3rd_dim <- function(a) {
m <- matrix(0, nrow=nrow(a), ncol=ncol(a))
for (z in seq_len(dim(a)[[3]]))
m <- m + a[ , , z]
m
}
## With a toy 3D array:
a <- array(runif(8000), dim=c(25, 40, 8))
dim(collapse_3rd_dim(a))
#> [1] 25 40
stopifnot(identical(sum(a), sum(collapse_3rd_dim(a)))) # sanity check
## Now say that A is so big that even M wouldn't fit in memory. This is
## a situation where we'd want to compute M block by block:
## -- STEP 1 --
## Create the 2D realization sink:
setAutoRealizationBackend("HDF5Array")
sink <- AutoRealizationSink(dim(a)[1:2])
dim(sink)
#> [1] 25 40
## -- STEP 2 --
## Define two grids: one for 'sink' and one for 'a'. Since we're going
## to walk on the two grids simultaneously, read a block from 'a' and
## write it to 'sink', we need to make sure that we define grids that
## are "aligned". More precisely, the two grids must have the same number
## of viewports, and the viewports in one must correspond to the viewports
## in the other one:
sink_grid <- colAutoGrid(sink, ncol=10)
a_spacings <- c(dim(sink_grid[[1L]]), dim(a)[[3]])
a_grid <- RegularArrayGrid(dim(a), spacings=a_spacings)
dims(sink_grid) # dimensions of the individual viewports
#> [,1] [,2]
#> [1,] 25 10
#> [2,] 25 10
#> [3,] 25 10
#> [4,] 25 10
dims(a_grid) # dimensions of the individual viewports
#> [,1] [,2] [,3]
#> [1,] 25 10 8
#> [2,] 25 10 8
#> [3,] 25 10 8
#> [4,] 25 10 8
## Let's check that our two grids are actually "aligned":
stopifnot(identical(length(sink_grid), length(a_grid)))
stopifnot(identical(dims(sink_grid), dims(a_grid)[ , 1:2, drop=FALSE]))
## -- STEP 3 --
## Walk on the two grids simultaneously:
for (bid in seq_along(sink_grid)) {
## Read block from 'a'.
a_viewport <- a_grid[[bid]]
block <- read_block(a, a_viewport)
## Collapse it.
block <- collapse_3rd_dim(block)
## Write the collapsed block to 'sink'.
sink_viewport <- sink_grid[[bid]]
sink <- write_block(sink, sink_viewport, block)
}
## -- An alternative to STEP 3 --
FUN <- function(sink, sink_viewport) {
## Read block from 'a'.
bid <- currentBlockId()
a_viewport <- a_grid[[bid]]
block <- read_block(a, a_viewport)
## Collapse it.
block <- collapse_3rd_dim(block)
## Write the collapsed block to 'sink'.
write_block(sink, sink_viewport, block)
}
sink <- sinkApply(sink, FUN, grid=sink_grid, verbose=TRUE)
#> \ processing viewport 1/4 ...
#> ok
#> \ processing viewport 2/4 ...
#> ok
#> \ processing viewport 3/4 ...
#> ok
#> \ processing viewport 4/4 ...
#> ok
## -- STEP 4 --
## Close the sink and turn it into a DelayedArray object:
close(sink)
M <- as(sink, "DelayedArray")
M
#> <25 x 40> HDF5Matrix object of type "double":
#> [,1] [,2] [,3] ... [,39] [,40]
#> [1,] 4.891130 2.691342 3.844899 . 3.801618 4.728601
#> [2,] 5.223214 4.866491 4.198213 . 5.106603 4.619065
#> [3,] 4.387821 4.074536 4.972230 . 3.137642 1.622522
#> [4,] 3.797170 3.338020 3.103377 . 4.336885 4.184482
#> [5,] 4.727747 4.038604 4.316469 . 3.526290 2.712939
#> ... . . . . . .
#> [21,] 4.227927 3.142827 3.992997 . 2.336606 4.843255
#> [22,] 3.788819 3.796723 3.844115 . 4.083134 3.226500
#> [23,] 4.225157 5.019150 4.490189 . 5.669049 4.756903
#> [24,] 2.989142 3.395270 3.981459 . 2.689478 3.644247
#> [25,] 3.783744 3.847474 4.051030 . 4.473004 3.612669
## Sanity check:
stopifnot(identical(collapse_3rd_dim(a), as.array(M)))
setAutoRealizationBackend() # restore default (NULL)
## ---------------------------------------------------------------------
## USING THE "RealizationSink API": AN ADVANCED EXAMPLE
## ---------------------------------------------------------------------
## Say we have 2 matrices with the same number of columns. Each column
## represents a biological sample:
library(HDF5Array)
R <- as(matrix(runif(75000), ncol=1000), "HDF5Array") # 75 rows
G <- as(matrix(runif(250000), ncol=1000), "HDF5Array") # 250 rows
## Say we want to compute the matrix U obtained by applying the same
## binary functions FUN() to all samples i.e. U is defined as:
##
## U[ , j] <- FUN(R[ , j], G[ , j]) for 1 <= j <= 1000
##
## Note that FUN() should return a vector of constant length, say 200,
## so U will be a 200x1000 matrix. A naive implementation would be:
##
## pFUN <- function(r, g) {
## stopifnot(ncol(r) == ncol(g)) # sanity check
## sapply(seq_len(ncol(r)), function(j) FUN(r[ , j], g[ , j]))
## }
##
## But because U is going to be too big to fit in memory, we can't
## just do pFUN(R, G). So we want to compute U block by block and
## write the blocks to disk as we go. The blocks will be made of full
## columns. Also since we need to walk on 2 matrices at the same time
## (R and G), we can't use blockApply() or blockReduce() so we'll use
## a "for" loop.
## Before we get to the "for" loop, we need 4 things:
## 1. Two grids of blocks, one on R and one on G. The blocks in the
## two grids must contain the same number of columns. We arbitrarily
## choose to use blocks of 150 columns:
R_grid <- colAutoGrid(R, ncol=150)
G_grid <- colAutoGrid(G, ncol=150)
## 2. The function pFUN(). It will take 2 blocks as input, 1 from R
## and 1 from G, apply FUN() to all the samples in the blocks,
## and return a matrix with one columns per sample:
pFUN <- function(r, g) {
stopifnot(ncol(r) == ncol(g)) # sanity check
## Return a matrix with 200 rows with random values. Completely
## artificial sorry. A realistic example would actually need to
## apply the same binary function to r[ ,j] and g[ , j] for
## 1 <= j <= ncol(r).
matrix(runif(200 * ncol(r)), nrow=200)
}
## 3. A RealizationSink derivative where to write the matrices returned
## by pFUN() as we go:
setAutoRealizationBackend("HDF5Array")
U_sink <- AutoRealizationSink(c(200L, 1000L))
## 4. Finally, we create a grid on U_sink with viewports that contain
## the same number of columns as the corresponding blocks in R and G:
U_grid <- colAutoGrid(U_sink, ncol=150)
## Note that the three grids should have the same number of viewports:
stopifnot(length(U_grid) == length(R_grid))
stopifnot(length(U_grid) == length(G_grid))
## 5. Now we can proceed. We use a "for" loop to walk on R and G
## simultaneously, block by block, apply pFUN(), and write the
## output of pFUN() to U_sink:
for (bid in seq_along(U_grid)) {
R_block <- read_block(R, R_grid[[bid]])
G_block <- read_block(G, G_grid[[bid]])
U_block <- pFUN(R_block, G_block)
U_sink <- write_block(U_sink, U_grid[[bid]], U_block)
}
## An alternative to the "for" loop is to use sinkApply():
FUN <- function(U_sink, U_viewport) {
bid <- currentBlockId()
R_block <- read_block(R, R_grid[[bid]])
G_block <- read_block(G, G_grid[[bid]])
U_block <- pFUN(R_block, G_block)
write_block(U_sink, U_viewport, U_block)
}
U_sink <- sinkApply(U_sink, FUN, grid=U_grid, verbose=TRUE)
#> \ processing viewport 1/7 ...
#> ok
#> \ processing viewport 2/7 ...
#> ok
#> \ processing viewport 3/7 ...
#> ok
#> \ processing viewport 4/7 ...
#> ok
#> \ processing viewport 5/7 ...
#> ok
#> \ processing viewport 6/7 ...
#> ok
#> \ processing viewport 7/7 ...
#> ok
close(U_sink)
U <- as(U_sink, "DelayedArray")
U
#> <200 x 1000> HDF5Matrix object of type "double":
#> [,1] [,2] [,3] ... [,999] [,1000]
#> [1,] 0.844809629 0.536645272 0.964087917 . 0.78816070 0.26438504
#> [2,] 0.955293005 0.174254668 0.136979119 . 0.92319025 0.32704078
#> [3,] 0.002817264 0.144103834 0.526218672 . 0.96381486 0.66496762
#> [4,] 0.143393669 0.453105410 0.074382507 . 0.27267989 0.36118448
#> [5,] 0.119990439 0.301102219 0.893128434 . 0.05651871 0.88406591
#> ... . . . . . .
#> [196,] 0.6704552 0.3816505 0.6479495 . 0.06030346 0.36984943
#> [197,] 0.6216413 0.7311405 0.1014469 . 0.59756641 0.40552255
#> [198,] 0.8637928 0.9678953 0.4578831 . 0.33012191 0.60925925
#> [199,] 0.7480010 0.5937706 0.9522801 . 0.56477554 0.12757410
#> [200,] 0.1811271 0.3787599 0.8476047 . 0.63253921 0.54073128
setAutoRealizationBackend() # restore default (NULL)
## ---------------------------------------------------------------------
## VERY BASIC (BUT ALSO VERY ARTIFICIAL) USAGE OF THE
## read_block()/write_block() COMBO
## ---------------------------------------------------------------------
###### On an ordinary matrix ######
m1 <- matrix(1:30, ncol=5)
## Define a viewport on 'm1':
block1_dim <- c(4, 3)
viewport1 <- ArrayViewport(dim(m1), IRanges(c(3, 2), width=block1_dim))
## Read/tranform/write:
block1 <- read_block(m1, viewport1)
write_block(m1, viewport1, block1 + 1000L)
#> [,1] [,2] [,3] [,4] [,5]
#> [1,] 1 7 13 19 25
#> [2,] 2 8 14 20 26
#> [3,] 3 1009 1015 1021 27
#> [4,] 4 1010 1016 1022 28
#> [5,] 5 1011 1017 1023 29
#> [6,] 6 1012 1018 1024 30
## Define another viewport on 'm1':
viewport1b <- ArrayViewport(dim(m1), IRanges(c(1, 3), width=block1_dim))
## Read/tranform/write:
write_block(m1, viewport1b, block1 + 1000L)
#> [,1] [,2] [,3] [,4] [,5]
#> [1,] 1 7 1009 1015 1021
#> [2,] 2 8 1010 1016 1022
#> [3,] 3 9 1011 1017 1023
#> [4,] 4 10 1012 1018 1024
#> [5,] 5 11 17 23 29
#> [6,] 6 12 18 24 30
## No-op:
m <- write_block(m1, viewport1, read_block(m1, viewport1))
stopifnot(identical(m1, m))
########## On a 3D array ##########
a3 <- array(1:60, 5:3)
## Define a viewport on 'a3':
block3_dim <- c(2, 4, 1)
viewport3 <- ArrayViewport(dim(a3), IRanges(c(1, 1, 3), width=block3_dim))
## Read/tranform/write:
block3 <- read_block(a3, viewport3)
write_block(a3, viewport3, block3 + 1000L)
#> , , 1
#>
#> [,1] [,2] [,3] [,4]
#> [1,] 1 6 11 16
#> [2,] 2 7 12 17
#> [3,] 3 8 13 18
#> [4,] 4 9 14 19
#> [5,] 5 10 15 20
#>
#> , , 2
#>
#> [,1] [,2] [,3] [,4]
#> [1,] 21 26 31 36
#> [2,] 22 27 32 37
#> [3,] 23 28 33 38
#> [4,] 24 29 34 39
#> [5,] 25 30 35 40
#>
#> , , 3
#>
#> [,1] [,2] [,3] [,4]
#> [1,] 1041 1046 1051 1056
#> [2,] 1042 1047 1052 1057
#> [3,] 43 48 53 58
#> [4,] 44 49 54 59
#> [5,] 45 50 55 60
#>
## Define another viewport on 'a3':
viewport3b <- ArrayViewport(dim(a3), IRanges(c(3, 1, 3), width=block3_dim))
## Read/tranform/write:
write_block(a3, viewport3b, block3 + 1000L)
#> , , 1
#>
#> [,1] [,2] [,3] [,4]
#> [1,] 1 6 11 16
#> [2,] 2 7 12 17
#> [3,] 3 8 13 18
#> [4,] 4 9 14 19
#> [5,] 5 10 15 20
#>
#> , , 2
#>
#> [,1] [,2] [,3] [,4]
#> [1,] 21 26 31 36
#> [2,] 22 27 32 37
#> [3,] 23 28 33 38
#> [4,] 24 29 34 39
#> [5,] 25 30 35 40
#>
#> , , 3
#>
#> [,1] [,2] [,3] [,4]
#> [1,] 41 46 51 56
#> [2,] 42 47 52 57
#> [3,] 1041 1046 1051 1056
#> [4,] 1042 1047 1052 1057
#> [5,] 45 50 55 60
#>
## No-op:
a <- write_block(a3, viewport3, read_block(a3, viewport3))
stopifnot(identical(a3, a))
## ---------------------------------------------------------------------
## LESS BASIC (BUT STILL VERY ARTIFICIAL) USAGE OF THE
## read_block()/write_block() COMBO
## ---------------------------------------------------------------------
grid1 <- RegularArrayGrid(dim(m1), spacings=c(3L, 2L))
grid1
#> 2 x 3 RegularArrayGrid object on a 6 x 5 array:
#> [,1] [,2] [,3]
#> [1,] [1-3,1-2] [1-3,3-4] [1-3,5]
#> [2,] [4-6,1-2] [4-6,3-4] [4-6,5]
length(grid1) # number of blocks defined by the grid
#> [1] 6
read_block(m1, grid1[[3L]]) # read 3rd block
#> [,1] [,2]
#> [1,] 13 19
#> [2,] 14 20
#> [3,] 15 21
read_block(m1, grid1[[1L, 3L]])
#> [,1]
#> [1,] 25
#> [2,] 26
#> [3,] 27
## Walk on the grid, colum by column:
m1a <- m1
for (bid in seq_along(grid1)) {
viewport <- grid1[[bid]]
block <- read_block(m1a, viewport)
block <- bid * 1000L + block
m1a <- write_block(m1a, viewport, block)
}
m1a
#> [,1] [,2] [,3] [,4] [,5]
#> [1,] 1001 1007 3013 3019 5025
#> [2,] 1002 1008 3014 3020 5026
#> [3,] 1003 1009 3015 3021 5027
#> [4,] 2004 2010 4016 4022 6028
#> [5,] 2005 2011 4017 4023 6029
#> [6,] 2006 2012 4018 4024 6030
## Walk on the grid, row by row:
m1b <- m1
for (i in seq_len(dim(grid1)[[1]])) {
for (j in seq_len(dim(grid1)[[2]])) {
viewport <- grid1[[i, j]]
block <- read_block(m1b, viewport)
block <- (i * 10L + j) * 1000L + block
m1b <- write_block(m1b, viewport, block)
}
}
m1b
#> [,1] [,2] [,3] [,4] [,5]
#> [1,] 11001 11007 12013 12019 13025
#> [2,] 11002 11008 12014 12020 13026
#> [3,] 11003 11009 12015 12021 13027
#> [4,] 21004 21010 22016 22022 23028
#> [5,] 21005 21011 22017 22023 23029
#> [6,] 21006 21012 22018 22024 23030
## ---------------------------------------------------------------------
## registeredRealizationBackends() AND FAMILY
## ---------------------------------------------------------------------
getAutoRealizationBackend() # no backend set yet
#> NULL
registeredRealizationBackends()
#> BACKEND package
#> 1 RleArray DelayedArray
#> 2 HDF5Array HDF5Array
#> 3 TENxMatrix HDF5Array
setAutoRealizationBackend("HDF5Array")
getAutoRealizationBackend() # backend is set to "HDF5Array"
#> [1] "HDF5Array"
registeredRealizationBackends()
#> BACKEND package
#> 1 RleArray DelayedArray
#> 2 -> HDF5Array HDF5Array <-
#> 3 TENxMatrix HDF5Array
getHDF5DumpChunkLength()
#> [1] 1000000
setHDF5DumpChunkLength(500L)
getHDF5DumpChunkShape()
#> [1] "scale"
sink <- AutoRealizationSink(c(120L, 50L))
class(sink) # HDF5-specific realization sink
#> [1] "HDF5RealizationSink"
#> attr(,"package")
#> [1] "HDF5Array"
dim(sink)
#> [1] 120 50
chunkdim(sink)
#> [1] 34 14
grid <- defaultSinkAutoGrid(sink, block.length=600)
for (bid in seq_along(grid)) {
viewport <- grid[[bid]]
block <- 101 * bid + runif(length(viewport))
dim(block) <- dim(viewport)
sink <- write_block(sink, viewport, block)
}
close(sink)
A <- as(sink, "DelayedArray")
A
#> <120 x 50> HDF5Matrix object of type "double":
#> [,1] [,2] [,3] ... [,49] [,50]
#> [1,] 101.0416 101.7265 101.6960 . 1313.964 1313.440
#> [2,] 101.9350 101.0209 101.4228 . 1313.785 1313.632
#> [3,] 101.9617 101.8247 101.0887 . 1313.049 1313.394
#> [4,] 101.2810 101.2333 101.8391 . 1313.177 1313.691
#> [5,] 101.8044 101.1290 101.1236 . 1313.680 1313.959
#> ... . . . . . .
#> [116,] 404.3529 404.1235 404.8089 . 1616.797 1616.137
#> [117,] 404.3167 404.2574 404.3650 . 1616.177 1616.372
#> [118,] 404.2191 404.2059 404.2975 . 1616.011 1616.577
#> [119,] 404.0483 404.6489 404.1236 . 1616.923 1616.870
#> [120,] 404.5980 404.8793 404.8491 . 1616.448 1616.730
setAutoRealizationBackend() # restore default (NULL)