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
|
### =========================================================================
### Read/write blocks of array data
### -------------------------------------------------------------------------
###
### - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
### read_block() and write_block()
###
### Must return an ordinay array and propagate the dimnames.
setGeneric("read_block", signature="x",
function(x, viewport)
{
stopifnot(is(viewport, "ArrayViewport"),
identical(refdim(viewport), dim(x)))
ans <- standardGeneric("read_block")
check_returned_array(ans, dim(viewport), "read_block", class(x))
}
)
### Must return 'x' (possibly modified if it's an in-memory object).
setGeneric("write_block", signature="x",
function(x, viewport, block)
{
stopifnot(is(viewport, "ArrayViewport"),
identical(refdim(viewport), dim(x)),
is.array(block),
identical(dim(block), dim(viewport)))
standardGeneric("write_block")
}
)
### Work on any object 'x' that supports extract_array() e.g. an ordinary
### array or a DelayedArray, HDF5ArraySeed, or DelayedOp object, etc...
### Propagate the dimnames.
setMethod("read_block", "ANY",
function(x, viewport)
{
## We use extract_array_by_Nindex() instead of extract_array() to
## propagate the dimnames.
Nindex <- makeNindexFromArrayViewport(viewport)
extract_array_by_Nindex(x, Nindex)
}
)
setMethod("write_block", "ANY",
function(x, viewport, block)
{
Nindex <- makeNindexFromArrayViewport(viewport)
replace_by_Nindex(x, Nindex, block)
}
)
### - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
### read_sparse_block() and write_sparse_block()
###
### Like extract_sparse_array(), which it is based on, read_sparse_block()
### should be called only on an array-like object 'x' for which 'is_sparse(x)'
### is TRUE. For the sake of efficiency, this is NOT checked and is the
### responsibility of the user. See SparseArraySeed-class.R
### Must return a SparseArraySeed object.
setGeneric("read_sparse_block", signature="x",
function(x, viewport)
{
stopifnot(is(viewport, "ArrayViewport"),
identical(refdim(viewport), dim(x)))
ans <- standardGeneric("read_sparse_block")
## TODO: Display a more user/developper-friendly error by
## doing something like the read_block() generic where
## check_returned_array() is used to display a long and
## detailed error message.
stopifnot(is(ans, "SparseArraySeed"),
identical(dim(ans), dim(viewport)))
ans
}
)
setMethod("read_sparse_block", "ANY",
function(x, viewport)
{
index <- makeNindexFromArrayViewport(viewport, expand.RangeNSBS=TRUE)
extract_sparse_array(x, index)
}
)
### The default "read_sparse_block" method above would work just fine on a
### SparseArraySeed object but we overwrite it with a method that is slightly
### more efficient (can be 2x to 3x faster on big SparseArraySeed objects
### with hundreds of thousands of nonzero elements).
.read_sparse_block_from_SparseArraySeed <- function(x, viewport)
{
taind <- t(x@aind)
keep_idx <- which(colAlls(taind >= start(viewport) &
taind <= end(viewport)))
taind <- taind[ , keep_idx, drop=FALSE]
offsets <- start(viewport) - 1L
x0_aind <- t(taind - offsets)
x0_nzdata <- x@nzdata[keep_idx]
BiocGenerics:::replaceSlots(x, dim=dim(viewport),
aind=x0_aind,
nzdata=x0_nzdata,
check=FALSE)
}
setMethod("read_sparse_block", "SparseArraySeed",
.read_sparse_block_from_SparseArraySeed
)
### 'sparse_block' must be a SparseArraySeed object.
### Must return 'x' (possibly modified if it's an in-memory object).
setGeneric("write_sparse_block", signature="x",
function(x, viewport, sparse_block)
{
stopifnot(is(viewport, "ArrayViewport"),
identical(refdim(viewport), dim(x)),
is(sparse_block, "SparseArraySeed"),
identical(dim(sparse_block), dim(viewport)))
standardGeneric("write_sparse_block")
}
)
setMethod("write_sparse_block", "ANY",
function(x, viewport, sparse_block)
{
block <- sparse2dense(sparse_block)
write_block(x, viewport, block)
}
)
|