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
|
### =========================================================================
### H5SparseMatrix objects
### -------------------------------------------------------------------------
###
setClass("H5SparseMatrix",
contains="DelayedMatrix",
representation(seed="H5SparseMatrixSeed")
)
### - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
### Constructor
###
setMethod("DelayedArray", "H5SparseMatrixSeed",
function(seed) new_DelayedArray(seed, Class="H5SparseMatrix")
)
### Works directly on an H5SparseMatrixSeed derivative, in which case it must
### be called with a single argument.
H5SparseMatrix <- function(filepath, group)
{
if (is(filepath, "H5SparseMatrixSeed")) {
if (!missing(group))
stop(wmsg("H5SparseMatrix() must be called with a single argument ",
"when passed an H5SparseMatrixSeed object"))
seed <- filepath
} else {
seed <- H5SparseMatrixSeed(filepath, group)
}
DelayedArray(seed)
}
### - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
### Taking advantage of sparsity
###
setMethod("nzcount", "H5SparseMatrix", function(x) nzcount(x@seed))
setMethod("read_sparse_block", "H5SparseMatrix",
function(x, viewport) read_sparse_block(x@seed, viewport)
)
setMethod("extractNonzeroDataByCol", "H5SparseMatrix",
function(x, j) extractNonzeroDataByCol(x@seed, j)
)
setMethod("extractNonzeroDataByRow", "H5SparseMatrix",
function(x, i) extractNonzeroDataByCol(x@seed, i)
)
|