| HDF5Array-class {HDF5Array} | R Documentation |
The HDF5Array class is a DelayedArray subclass for representing a conventional (i.e. dense) HDF5 dataset.
All the operations available for DelayedArray objects work on HDF5Array objects.
## Constructor function: HDF5Array(filepath, name, type=NA)
filepath |
The path (as a single character string) to the HDF5 file where the dataset is located. |
name |
The name of the dataset in the HDF5 file. |
type |
|
An HDF5Array object.
The 1.3 Million Brain Cell Dataset and other datasets published by 10x Genomics use an HDF5-based sparse matrix representation instead of the conventional (i.e. dense) HDF5 representation.
If your dataset uses the conventional (i.e. dense) HDF5 representation,
use the HDF5Array() constructor.
If your dataset uses the HDF5-based sparse matrix representation from
10x Genomics, use the TENxMatrix() constructor.
TENxMatrix objects for representing 10x Genomics datasets as DelayedMatrix objects.
ReshapedHDF5Array objects for representing HDF5 datasets as DelayedArray objects with a user-supplied upfront virtual reshaping.
DelayedArray objects in the DelayedArray package.
writeHDF5Array for writing an array-like object
to an HDF5 file.
HDF5-dump-management for controlling the location and physical properties of automatically created HDF5 datasets.
saveHDF5SummarizedExperiment and
loadHDF5SummarizedExperiment in this
package (the HDF5Array package) for saving/loading
an HDF5-based SummarizedExperiment
object to/from disk.
The HDF5ArraySeed helper class.
h5ls in the rhdf5 package.
## ---------------------------------------------------------------------
## CONSTRUCTION
## ---------------------------------------------------------------------
library(h5vcData)
tally_file <- system.file("extdata", "example.tally.hfs5",
package="h5vcData")
library(rhdf5) # for h5ls()
h5ls(tally_file)
## Pick up "Coverages" dataset for Human chromosome 16:
cvg0 <- HDF5Array(tally_file, "/ExampleStudy/16/Coverages")
cvg0
is(cvg0, "DelayedArray") # TRUE
seed(cvg0)
path(cvg0)
chunkdim(cvg0)
## ---------------------------------------------------------------------
## dim/dimnames
## ---------------------------------------------------------------------
dim(cvg0)
dimnames(cvg0)
dimnames(cvg0) <- list(paste0("s", 1:6), c("+", "-"), NULL)
dimnames(cvg0)
## ---------------------------------------------------------------------
## SLICING (A.K.A. SUBSETTING)
## ---------------------------------------------------------------------
cvg1 <- cvg0[ , , 29000001:29000007]
cvg1
dim(cvg1)
as.array(cvg1)
stopifnot(identical(dim(as.array(cvg1)), dim(cvg1)))
stopifnot(identical(dimnames(as.array(cvg1)), dimnames(cvg1)))
cvg2 <- cvg0[ , "+", 29000001:29000007]
cvg2
as.matrix(cvg2)
## ---------------------------------------------------------------------
## SummarizedExperiment OBJECTS WITH DELAYED ASSAYS
## ---------------------------------------------------------------------
## DelayedArray objects can be used inside a SummarizedExperiment object
## to hold the assay data and to delay operations on them.
library(SummarizedExperiment)
pcvg <- cvg0[ , 1, ] # coverage on plus strand
mcvg <- cvg0[ , 2, ] # coverage on minus strand
nrow(pcvg) # nb of samples
ncol(pcvg) # length of Human chromosome 16
## The convention for a SummarizedExperiment object is to have 1 column
## per sample so first we need to transpose 'pcvg' and 'mcvg':
pcvg <- t(pcvg)
mcvg <- t(mcvg)
se <- SummarizedExperiment(list(pcvg=pcvg, mcvg=mcvg))
se
stopifnot(validObject(se, complete=TRUE))
## A GPos object can be used to represent the genomic positions along
## the dataset:
gpos <- GPos(GRanges("16", IRanges(1, nrow(se))))
gpos
rowRanges(se) <- gpos
se
stopifnot(validObject(se))
assays(se)$pcvg
assays(se)$mcvg