HDF5 datasets as DelayedArray objects
HDF5Array-class.Rd
The HDF5Array class is a DelayedArray subclass for representing and operating on a conventional (a.k.a. dense) HDF5 dataset.
All the operations available for DelayedArray objects work on HDF5Array objects.
Arguments
- filepath
The path (as a single string or H5File object) to the HDF5 file (
.h5
or.h5ad
) where the dataset is located.Note that you must create and use an H5File object if the HDF5 file to access is stored in an Amazon S3 bucket. See
?H5File
for how to do this.Also please note that H5File objects must NOT be used in the context of parallel evaluation at the moment.
- name
The name of the dataset in the HDF5 file.
- as.sparse
Whether the HDF5 dataset should be flagged as sparse or not, that is, whether it should be considered sparse (and treated as such) or not. Note that HDF5 doesn't natively support sparse storage at the moment so HDF5 datasets cannot be stored in a sparse format, only in a dense one. However a dataset stored in a dense format can still contain a lot of zeros. Using
as.sparse=TRUE
on such dataset will enable some optimizations that can lead to a lower memory footprint (and possibly better performance) when operating on the HDF5Array.IMPORTANT NOTE: If the dataset is in the 10x Genomics format (i.e. if it uses the HDF5-based sparse matrix representation from 10x Genomics), you should use the
TENxMatrix()
constructor instead of theHDF5Array()
constructor.- type
By default the
type
of the returned object is inferred from the H5 datatype of the HDF5 dataset. This can be overridden by specifying thetype
argument. The specified type must be an R atomic type (e.g."integer"
) or"list"
.
Note
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 (a.k.a. dense) HDF5 representation.
If your dataset uses the conventional (a.k.a. dense) HDF5 representation,
use the HDF5Array()
constructor documented here.
But if your dataset uses the HDF5 sparse matrix representation from
10x Genomics, use the TENxMatrix()
constructor instead.
See also
H5File objects.
H5SparseMatrix objects for representing HDF5 sparse matrices as DelayedMatrix objects.
H5ADMatrix objects for representing h5ad central matrices (or matrices in the
/layers
group) as DelayedMatrix objects.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
andloadHDF5SummarizedExperiment
in this package (the HDF5Array package) for saving/loading an HDF5-based SummarizedExperiment object to/from disk.The HDF5ArraySeed helper class.
h5ls
to list the content of an HDF5 file (.h5
or.h5ad
).
Examples
## ---------------------------------------------------------------------
## A. CONSTRUCTION
## ---------------------------------------------------------------------
## With a local file:
toy_h5 <- system.file("extdata", "toy.h5", package="HDF5Array")
h5ls(toy_h5)
#> group name otype dclass dim
#> 0 / M1 H5I_DATASET FLOAT 10000 x 150
#> 1 / M2 H5I_DATASET FLOAT 150 x 200
HDF5Array(toy_h5, "M2")
#> <150 x 200> HDF5Matrix object of type "double":
#> [,1] [,2] [,3] ... [,199] [,200]
#> [1,] 9.0325451 -0.8372894 -3.0655634 . 4.202503 7.212745
#> [2,] 0.5603696 1.4590217 1.6268596 . 14.854123 3.481012
#> [3,] 8.1094489 2.9110584 12.7574103 . 0.520843 8.876323
#> [4,] 5.8518510 6.4912588 13.8403073 . 8.795605 12.342113
#> [5,] -0.8206218 10.1040618 10.6219819 . 8.756145 3.130079
#> ... . . . . . .
#> [146,] 11.2827878 -3.7808660 5.4538427 . -0.762771429 6.070777178
#> [147,] -2.3159837 10.5172435 -2.9781124 . 0.824398901 -3.154179584
#> [148,] 0.9635204 13.8644578 11.6860867 . 5.761504602 1.688992344
#> [149,] 12.6951503 11.6835890 13.4213585 . 13.590451800 2.134752800
#> [150,] -1.8871983 9.9152485 8.6534266 . -0.007439521 7.164514768
HDF5Array(toy_h5, "M2", type="integer")
#> <150 x 200> HDF5Matrix object of type "integer":
#> [,1] [,2] [,3] [,4] ... [,197] [,198] [,199] [,200]
#> [1,] 9 0 -3 14 . -1 11 4 7
#> [2,] 0 1 1 14 . 12 6 14 3
#> [3,] 8 2 12 6 . 10 0 0 8
#> [4,] 5 6 13 12 . -1 6 8 12
#> [5,] 0 10 10 -4 . 9 14 8 3
#> ... . . . . . . . . .
#> [146,] 11 -3 5 0 . 1 6 0 6
#> [147,] -2 10 -2 -1 . 3 9 0 -3
#> [148,] 0 13 11 -2 . 2 -4 5 1
#> [149,] 12 11 13 8 . 13 3 13 2
#> [150,] -1 9 8 9 . 10 -2 0 7
HDF5Array(toy_h5, "M2", type="complex")
#> <150 x 200> HDF5Matrix object of type "complex":
#> [,1] [,2] [,3] ... [,199]
#> [1,] 9.0325451+0i -0.8372894+0i -3.0655634+0i . 4.202503+0i
#> [2,] 0.5603696+0i 1.4590217+0i 1.6268596+0i . 14.854123+0i
#> [3,] 8.1094489+0i 2.9110584+0i 12.7574103+0i . 0.520843+0i
#> [4,] 5.8518510+0i 6.4912588+0i 13.8403073+0i . 8.795605+0i
#> [5,] -0.8206218+0i 10.1040618+0i 10.6219819+0i . 8.756145+0i
#> ... . . . . .
#> [146,] 11.2827878+0i -3.7808660+0i 5.4538427+0i . -0.762771429+0i
#> [147,] -2.3159837+0i 10.5172435+0i -2.9781124+0i . 0.824398901+0i
#> [148,] 0.9635204+0i 13.8644578+0i 11.6860867+0i . 5.761504602+0i
#> [149,] 12.6951503+0i 11.6835890+0i 13.4213585+0i . 13.590451800+0i
#> [150,] -1.8871983+0i 9.9152485+0i 8.6534266+0i . -0.007439521+0i
#> [,200]
#> [1,] 7.212745+0i
#> [2,] 3.481012+0i
#> [3,] 8.876323+0i
#> [4,] 12.342113+0i
#> [5,] 3.130079+0i
#> ... .
#> [146,] 6.070777178+0i
#> [147,] -3.154179584+0i
#> [148,] 1.688992344+0i
#> [149,] 2.134752800+0i
#> [150,] 7.164514768+0i
## With a file stored in an Amazon S3 bucket:
if (Sys.info()[["sysname"]] != "Darwin") {
public_S3_url <-
"https://rhdf5-public.s3.eu-central-1.amazonaws.com/rhdf5ex_t_float_3d.h5"
h5file <- H5File(public_S3_url, s3=TRUE)
h5ls(h5file)
HDF5Array(h5file, "a1")
}
#> <5 x 10 x 2> HDF5Array object of type "double":
#> ,,1
#> [,1] [,2] [,3] ... [,9] [,10]
#> [1,] 0.122377600 0.473829743 0.244448476 . 0.06294312 0.25262216
#> [2,] 0.257908699 0.009508783 0.387372283 . 0.52225184 0.94390880
#> [3,] 0.886223679 0.211665150 0.404533111 . 0.95703000 0.18594012
#> [4,] 0.549336179 0.071420332 0.617852765 . 0.20320586 0.48317974
#> [5,] 0.179772913 0.557769188 0.553050760 . 0.72190719 0.18467657
#>
#> ,,2
#> [,1] [,2] [,3] ... [,9] [,10]
#> [1,] 0.6848327 0.1282442 0.7906603 . 0.40603368 0.01085798
#> [2,] 0.6579706 0.4897645 0.3274960 . 0.47773518 0.38283946
#> [3,] 0.7180033 0.4929566 0.2159982 . 0.48220542 0.56122233
#> [4,] 0.8445332 0.8016226 0.8083879 . 0.26790274 0.17247477
#> [5,] 0.4086140 0.3476558 0.7067179 . 0.18221174 0.99788310
#>
## ---------------------------------------------------------------------
## B. BASIC MANIPULATION
## ---------------------------------------------------------------------
library(h5vcData)
tally_file <- system.file("extdata", "example.tally.hfs5",
package="h5vcData")
h5ls(tally_file)
#> group name otype dclass dim
#> 0 / ExampleStudy H5I_GROUP
#> 1 /ExampleStudy 16 H5I_GROUP
#> 2 /ExampleStudy/16 Counts H5I_DATASET INTEGER 12 x 6 x 2 x 90354753
#> 3 /ExampleStudy/16 Coverages H5I_DATASET INTEGER 6 x 2 x 90354753
#> 4 /ExampleStudy/16 Deletions H5I_DATASET INTEGER 6 x 2 x 90354753
#> 5 /ExampleStudy/16 Reference H5I_DATASET INTEGER 90354753
#> 6 /ExampleStudy 22 H5I_GROUP
#> 7 /ExampleStudy/22 Counts H5I_DATASET INTEGER 12 x 6 x 2 x 51304566
#> 8 /ExampleStudy/22 Coverages H5I_DATASET INTEGER 6 x 2 x 51304566
#> 9 /ExampleStudy/22 Deletions H5I_DATASET INTEGER 6 x 2 x 51304566
#> 10 /ExampleStudy/22 Reference H5I_DATASET INTEGER 51304566
## Pick up "Coverages" dataset for Human chromosome 16:
name <- "/ExampleStudy/16/Coverages"
cvg <- HDF5Array(tally_file, name)
cvg
#> <6 x 2 x 90354753> HDF5Array object of type "integer":
#> ,,1
#> [,1] [,2]
#> [1,] 0 0
#> [2,] 0 0
#> ... . .
#> [5,] 0 0
#> [6,] 0 0
#>
#> ...
#>
#> ,,90354753
#> [,1] [,2]
#> [1,] 0 0
#> [2,] 0 0
#> ... . .
#> [5,] 0 0
#> [6,] 0 0
#>
is(cvg, "DelayedArray") # TRUE
#> [1] TRUE
seed(cvg)
#> An object of class "HDF5ArraySeed"
#> Slot "filepath":
#> [1] "/usr/local/lib/R/site-library/h5vcData/extdata/example.tally.hfs5"
#>
#> Slot "name":
#> [1] "/ExampleStudy/16/Coverages"
#>
#> Slot "as_sparse":
#> [1] FALSE
#>
#> Slot "type":
#> [1] NA
#>
#> Slot "dim":
#> [1] 6 2 90354753
#>
#> Slot "chunkdim":
#> [1] 1 1 88238
#>
#> Slot "first_val":
#> [1] 0
#>
path(cvg)
#> [1] "/usr/local/lib/R/site-library/h5vcData/extdata/example.tally.hfs5"
chunkdim(cvg)
#> [1] 1 1 88238
## The data in the dataset looks sparse. In this case it is recommended
## to set 'as.sparse' to TRUE when constructing the HDF5Array object.
## This will make block processing (used in operations like sum()) more
## memory efficient and likely faster:
cvg0 <- HDF5Array(tally_file, name, as.sparse=TRUE)
is_sparse(cvg0) # TRUE
#> [1] TRUE
## Note that we can also flag the HDF5Array object as sparse after
## creation:
is_sparse(cvg) <- TRUE
cvg # same as 'cvg0'
#> <6 x 2 x 90354753> sparse HDF5Array object of type "integer":
#> ,,1
#> [,1] [,2]
#> [1,] 0 0
#> [2,] 0 0
#> ... . .
#> [5,] 0 0
#> [6,] 0 0
#>
#> ...
#>
#> ,,90354753
#> [,1] [,2]
#> [1,] 0 0
#> [2,] 0 0
#> ... . .
#> [5,] 0 0
#> [6,] 0 0
#>
## dim/dimnames:
dim(cvg0)
#> [1] 6 2 90354753
dimnames(cvg0)
#> NULL
dimnames(cvg0) <- list(paste0("s", 1:6), c("+", "-"), NULL)
dimnames(cvg0)
#> [[1]]
#> [1] "s1" "s2" "s3" "s4" "s5" "s6"
#>
#> [[2]]
#> [1] "+" "-"
#>
#> [[3]]
#> NULL
#>
## ---------------------------------------------------------------------
## C. SLICING (A.K.A. SUBSETTING)
## ---------------------------------------------------------------------
cvg1 <- cvg0[ , , 29000001:29000007]
cvg1
#> <6 x 2 x 7> sparse DelayedArray object of type "integer":
#> ,,1
#> + -
#> s1 19 40
#> s2 24 49
#> .. . .
#> s5 15 24
#> s6 22 30
#>
#> ...
#>
#> ,,7
#> + -
#> s1 18 30
#> s2 24 41
#> .. . .
#> s5 15 17
#> s6 21 18
#>
dim(cvg1)
#> [1] 6 2 7
as.array(cvg1)
#> , , 1
#>
#> + -
#> s1 19 40
#> s2 24 49
#> s3 7 20
#> s4 22 49
#> s5 15 24
#> s6 22 30
#>
#> , , 2
#>
#> + -
#> s1 19 36
#> s2 23 47
#> s3 7 20
#> s4 19 46
#> s5 15 24
#> s6 22 28
#>
#> , , 3
#>
#> + -
#> s1 19 36
#> s2 23 47
#> s3 7 20
#> s4 19 46
#> s5 15 24
#> s6 22 27
#>
#> , , 4
#>
#> + -
#> s1 18 36
#> s2 23 47
#> s3 7 19
#> s4 19 45
#> s5 15 23
#> s6 23 27
#>
#> , , 5
#>
#> + -
#> s1 19 31
#> s2 23 42
#> s3 7 14
#> s4 19 42
#> s5 15 17
#> s6 23 19
#>
#> , , 6
#>
#> + -
#> s1 19 30
#> s2 23 42
#> s3 7 14
#> s4 18 42
#> s5 15 18
#> s6 22 19
#>
#> , , 7
#>
#> + -
#> s1 18 30
#> s2 24 41
#> s3 7 14
#> s4 18 41
#> s5 15 17
#> s6 21 18
#>
stopifnot(identical(dim(as.array(cvg1)), dim(cvg1)))
stopifnot(identical(dimnames(as.array(cvg1)), dimnames(cvg1)))
cvg2 <- cvg0[ , "+", 29000001:29000007]
cvg2
#> <6 x 7> sparse DelayedMatrix object of type "integer":
#> [,1] [,2] [,3] [,4] [,5] [,6] [,7]
#> s1 19 19 19 18 19 19 18
#> s2 24 23 23 23 23 23 24
#> s3 7 7 7 7 7 7 7
#> s4 22 19 19 19 19 18 18
#> s5 15 15 15 15 15 15 15
#> s6 22 22 22 23 23 22 21
as.matrix(cvg2)
#> [,1] [,2] [,3] [,4] [,5] [,6] [,7]
#> s1 19 19 19 18 19 19 18
#> s2 24 23 23 23 23 23 24
#> s3 7 7 7 7 7 7 7
#> s4 22 19 19 19 19 18 18
#> s5 15 15 15 15 15 15 15
#> s6 22 22 22 23 23 22 21
## ---------------------------------------------------------------------
## D. 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)
#> Loading required package: GenomicRanges
#> Loading required package: GenomeInfoDb
#> Loading required package: Biobase
#> Welcome to Bioconductor
#>
#> Vignettes contain introductory material; view with
#> 'browseVignettes()'. To cite Bioconductor, see
#> 'citation("Biobase")', and for packages 'citation("pkgname")'.
#>
#> Attaching package: ‘Biobase’
#> The following object is masked from ‘package:SparseArray’:
#>
#> rowMedians
#> The following object is masked from ‘package:MatrixGenerics’:
#>
#> rowMedians
#> The following objects are masked from ‘package:matrixStats’:
#>
#> anyMissing, rowMedians
pcvg <- cvg0[ , 1, ] # coverage on plus strand
mcvg <- cvg0[ , 2, ] # coverage on minus strand
nrow(pcvg) # nb of samples
#> [1] 6
ncol(pcvg) # length of Human chromosome 16
#> [1] 90354753
## 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
#> class: SummarizedExperiment
#> dim: 90354753 6
#> metadata(0):
#> assays(2): pcvg mcvg
#> rownames: NULL
#> rowData names(0):
#> colnames(6): s1 s2 ... s5 s6
#> colData names(0):
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
#> StitchedGPos object with 90354753 positions and 0 metadata columns:
#> seqnames pos strand
#> <Rle> <integer> <Rle>
#> [1] 16 1 *
#> [2] 16 2 *
#> [3] 16 3 *
#> [4] 16 4 *
#> [5] 16 5 *
#> ... ... ... ...
#> [90354749] 16 90354749 *
#> [90354750] 16 90354750 *
#> [90354751] 16 90354751 *
#> [90354752] 16 90354752 *
#> [90354753] 16 90354753 *
#> -------
#> seqinfo: 1 sequence from an unspecified genome; no seqlengths
rowRanges(se) <- gpos
se
#> class: RangedSummarizedExperiment
#> dim: 90354753 6
#> metadata(0):
#> assays(2): pcvg mcvg
#> rownames: NULL
#> rowData names(0):
#> colnames(6): s1 s2 ... s5 s6
#> colData names(0):
stopifnot(validObject(se))
assays(se)$pcvg
#> <90354753 x 6> sparse DelayedMatrix object of type "integer":
#> s1 s2 s3 s4 s5 s6
#> [1,] 0 0 0 0 0 0
#> [2,] 0 0 0 0 0 0
#> [3,] 0 0 0 0 0 0
#> [4,] 0 0 0 0 0 0
#> [5,] 0 0 0 0 0 0
#> ... . . . . . .
#> [90354749,] 0 0 0 0 0 0
#> [90354750,] 0 0 0 0 0 0
#> [90354751,] 0 0 0 0 0 0
#> [90354752,] 0 0 0 0 0 0
#> [90354753,] 0 0 0 0 0 0
assays(se)$mcvg
#> <90354753 x 6> sparse DelayedMatrix object of type "integer":
#> s1 s2 s3 s4 s5 s6
#> [1,] 0 0 0 0 0 0
#> [2,] 0 0 0 0 0 0
#> [3,] 0 0 0 0 0 0
#> [4,] 0 0 0 0 0 0
#> [5,] 0 0 0 0 0 0
#> ... . . . . . .
#> [90354749,] 0 0 0 0 0 0
#> [90354750,] 0 0 0 0 0 0
#> [90354751,] 0 0 0 0 0 0
#> [90354752,] 0 0 0 0 0 0
#> [90354753,] 0 0 0 0 0 0