read_block: Read array blocks

View source: R/read_block.R

read_blockR Documentation

Read array blocks

Description

Use read_block to read a block of data from an array-like object.

Note that this function is typically used in the context of block processing of on-disk objects (e.g. DelayedArray objects), often in combination with write_block.

Usage

read_block(x, viewport, as.sparse=NA)

## Internal generic function used by read_block() when is_sparse(x)
## is FALSE:
read_block_as_dense(x, viewport)

Arguments

x

An array-like object.

This can be an ordinary array, a SparseArray object from the SparseArray package, a dgCMatrix object from the Matrix package, a DelayedArray object from the DelayedArray package, or any object with an array semantic (i.e. an object for which dim(x) is not NULL).

viewport

An ArrayViewport object compatible with x, that is, such that refdim(viewport) is identical to dim(x).

as.sparse

Can be FALSE, TRUE, or NA.

If FALSE, the block is returned as an ordinary array (a.k.a. dense array).

If TRUE, it's returned as a SparseArray object.

If NA (the default), the block is returned as an ordinary array if is_sparse(x) is FALSE and as a SparseArray object otherwise. In other words, using as.sparse=NA is equivalent to using as.sparse=is_sparse(x). This preserves sparsity and is the most efficient way to read a block.

Note that when returned as a 2D SparseArray object with numeric or logical data, a block can easily and efficiently be coerced to a sparseMatrix derivative from the Matrix package with as(block, "sparseMatrix"). This will return a dgCMatrix object if type(block) is "double" or "integer", and a lgCMatrix object if it's "logical".

Details

read_block() delegates to 2 internal generic functions for reading a block:

  • read_block_as_dense: used when is_sparse(x) is FALSE.

  • read_block_as_sparse (defined in the SparseArray package): used when is_sparse(x) is TRUE.

Note that these 2 internal generic functions are not meant to be called directly by the end user. The end user should always call the higher-level user-facing read_block() function instead.

Value

A block of data. More precisely, the data from x that belongs to the block delimited by the specified viewport.

The block of data is returned either as an ordinary (dense) array or as a SparseArray object from the SparseArray package.

Note that the returned block of data is guaranteed to have the same type as x and the same dimensions as the viewport. More formally, if block is the value returned by read_block(x, viewport), then:

    identical(type(block), type(x))

and

    identical(dim(block), dim(viewport))

are always TRUE.

See Also

  • ArrayGrid for ArrayGrid and ArrayViewport objects.

  • is_sparse to check whether an object uses a sparse representation of the data or not.

  • SparseArray objects implemented in the SparseArray package.

  • S4Arrays::type to get the type of the elements of an array-like object.

  • The read_block_as_sparse internal generic function defined in the SparseArray package and used by read_block() when is_sparse(x) is TRUE.

  • write_block to write a block of data to an array-like object.

  • blockApply and family, in the DelayedArray package, for convenient block processing of an array-like object.

  • dgCMatrix and lgCMatrix objects implemented in the Matrix package.

  • DelayedArray objects implemented in the DelayedArray package.

  • array and matrix objects in base R.

Examples

## Please note that, although educative, the examples below are somewhat
## artificial and do not illustrate real-world usage of read_block().
## See '?RealizationSink' in the DelayedArray package for more realistic
## read_block/write_block examples.

## ---------------------------------------------------------------------
## BASIC EXAMPLE 1: READ A BLOCK FROM AN ORDINARY MATRIX (DENSE)
## ---------------------------------------------------------------------
m1 <- matrix(1:30, ncol=5)
m1

## Define the viewport on 'm1' to read the data from:
block1_dim <- c(4, 3)
viewport1 <- ArrayViewport(dim(m1), IRanges(c(3, 2), width=block1_dim))
viewport1

## Read the block:
block1 <- read_block(m1, viewport1)  # same as m1[3:6, 2:4, drop=FALSE]
block1

## Use 'as.sparse=TRUE' to read the block as sparse object:
block1b <- read_block(m1, viewport1, as.sparse=TRUE)
block1b
is_sparse(block1b)  # TRUE
class(block1b)      # an SVT_SparseArray object

## Sanity checks:
stopifnot(identical(type(m1), type(block1)))
stopifnot(identical(dim(viewport1), dim(block1)))
stopifnot(identical(m1[3:6, 2:4, drop=FALSE], block1))
stopifnot(is(block1b, "SparseArray"))
stopifnot(identical(type(m1), type(block1b)))
stopifnot(identical(dim(viewport1), dim(block1b)))
stopifnot(identical(block1, as.array(block1b)))

## ---------------------------------------------------------------------
## BASIC EXAMPLE 2: READ A BLOCK FROM A SPARSE MATRIX
## ---------------------------------------------------------------------
m2 <- rsparsematrix(12, 20, density=0.2,
                    rand.x=function(n) sample(25, n, replace=TRUE))
m2
is_sparse(m2)  # TRUE

## Define the viewport on 'm2' to read the data from:
block2_dim <- c(2, 20)
viewport2 <- ArrayViewport(dim(m2), IRanges(c(1, 1), width=block2_dim))
viewport2

## By default, read_block() preserves sparsity:
block2 <- read_block(m2, viewport2)
block2
is_sparse(block2)  # TRUE
class(block2)      # an SVT_SparseArray object

## Use 'as.sparse=FALSE' to force read_block() to return an ordinary
## matrix or array:
block2b <- read_block(m2, viewport2, as.sparse=FALSE)
block2b
as(block2b, "sparseMatrix")

## Sanity checks:
stopifnot(is(block2, "SparseArray"))
stopifnot(identical(type(m2), type(block2)))
stopifnot(identical(dim(viewport2), dim(block2)))
stopifnot(identical(type(m2), type(block2b)))
stopifnot(identical(dim(viewport2), dim(block2b)))
stopifnot(identical(block2b, as.array(block2)))

## ---------------------------------------------------------------------
## BASIC EXAMPLE 3: READ A BLOCK FROM A 3D ARRAY
## ---------------------------------------------------------------------
a3 <- array(1:60, dim=5:3)

## Define the viewport on 'a3' to read the data from:
block3_dim <- c(2, 4, 1)
viewport3 <- ArrayViewport(dim(a3), IRanges(c(1, 1, 3), width=block3_dim))
viewport3

## Read the block:
block3 <- read_block(a3, viewport3)  # same as a3[1:2, 1:4, 3, drop=FALSE]
block3

## Note that unlike [, read_block() never drops dimensions.

## Sanity checks:
stopifnot(identical(type(a3), type(block3)))
stopifnot(identical(dim(viewport3), dim(block3)))
stopifnot(identical(a3[1:2, 1:4, 3, drop=FALSE], block3))

## ---------------------------------------------------------------------
## BASIC EXAMPLE 4: READ AND PROCESS BLOCKS DEFINED BY A GRID
## ---------------------------------------------------------------------
a4 <- array(runif(120), dim=6:4)

## Define a grid of 2x3x2 blocks on 'a4':
grid4 <- RegularArrayGrid(dim(a4), spacings=c(2,3,2))
grid4
nblock <- length(grid4)  # number of blocks

## Walk on the grid and print the corresponding blocks:
for (bid in seq_len(nblock)) {
    viewport <- grid4[[bid]]
    block <- read_block(a4, viewport)
    cat("====== Block ", bid, "/", nblock, " ======\n", sep="")
    print(block)
}

## Walk on the grid and compute the sum of each block:
block_sums <- sapply(grid4,
    function(viewport) sum(read_block(a4, viewport))
)
block_sums

## Sanity checks:
stopifnot(identical(length(block_sums), nblock))
stopifnot(all.equal(sum(block_sums), sum(a4)))

## ---------------------------------------------------------------------
## THE read_block/write_block COMBO
## ---------------------------------------------------------------------
## See '?write_block' for examples that use the read_block/write_block
## combo.

Bioconductor/S4Arrays documentation built on Oct. 30, 2024, 12:13 p.m.