TENxMatrixSeed objects
TENxMatrixSeed-class.Rd
TENxMatrixSeed is a low-level helper class that is a direct extension of the H5SparseMatrixSeed class. It is used to represent a pointer to an HDF5 sparse matrix that is stored in the CSR/CSC/Yale format ("Compressed Sparse Row") and follows the 10x Genomics convention for storing the dimensions of the matrix.
Note that a TENxMatrixSeed object is not intended to be used directly.
Most end users will typically create and manipulate a higher-level
TENxMatrix object instead. See ?TENxMatrix
for
more information.
Arguments
- filepath, group
See
?TENxMatrix
for a description of these arguments.
Details
A TENxMatrixSeed object supports the same limited set of methods as
an H5SparseMatrixSeed object. See ?H5SparseMatrixSeed
for the details.
TENxMatrixSeed vs TENxMatrix objects
In order to have access to the full set of operations that are available
for DelayedMatrix objects, a TENxMatrixSeed object
first needs to be wrapped in a DelayedMatrix object,
typically by calling the DelayedArray()
constructor on it.
This is what the TENxMatrix()
constructor function does.
Note that the result of this wrapping is a TENxMatrix object, which is just a TENxMatrixSeed object wrapped in a DelayedMatrix object.
See also
TENxMatrix objects.
H5SparseMatrixSeed objects.
The
TENxBrainData
dataset (in the TENxBrainData package).h5ls
to list the content of an HDF5 file.
Examples
## The 1.3 Million Brain Cell Dataset from 10x Genomics is available
## via ExperimentHub:
library(ExperimentHub)
hub <- ExperimentHub()
query(hub, "TENxBrainData")
#> ExperimentHub with 8 records
#> # snapshotDate(): 2025-04-12
#> # $dataprovider: 10X Genomics
#> # $species: Mus musculus
#> # $rdataclass: character
#> # additional mcols(): taxonomyid, genome, description,
#> # coordinate_1_based, maintainer, rdatadateadded, preparerclass, tags,
#> # rdatapath, sourceurl, sourcetype
#> # retrieve records with, e.g., 'object[["EH1039"]]'
#>
#> title
#> EH1039 | Brain scRNA-seq data, 'HDF5-based 10X Genomics' format
#> EH1040 | Brain scRNA-seq data, 'dense matrix' format
#> EH1041 | Brain scRNA-seq data, sample (column) annotation
#> EH1042 | Brain scRNA-seq data, gene (row) annotation
#> EH1689 | Brain scRNA-seq data 20k subset, 'HDF5-based 10x Genomics' format
#> EH1690 | Brain scRNA-seq data 20k subset, 'dense matrix' format
#> EH1691 | Brain scRNA-seq data 20k subset, sample (column) annotation
#> EH1692 | Brain scRNA-seq data 20k subset, gene (row) annotation
fname <- hub[["EH1039"]]
#> see ?TENxBrainData and browseVignettes('TENxBrainData') for documentation
#> loading from cache
## 'fname' is an HDF5 file. Use h5ls() to list its content:
h5ls(fname)
#> group name otype dclass dim
#> 0 / mm10 H5I_GROUP
#> 1 /mm10 barcodes H5I_DATASET STRING 1306127
#> 2 /mm10 data H5I_DATASET INTEGER 2624828308
#> 3 /mm10 gene_names H5I_DATASET STRING 27998
#> 4 /mm10 genes H5I_DATASET STRING 27998
#> 5 /mm10 indices H5I_DATASET INTEGER 2624828308
#> 6 /mm10 indptr H5I_DATASET INTEGER 1306128
#> 7 /mm10 shape H5I_DATASET INTEGER 2
## The 1.3 Million Brain Cell Dataset is represented by the "mm10"
## group. We point the TENxMatrixSeed() constructor to this group
## to create a TENxMatrixSeed object representing the dataset:
seed <- TENxMatrixSeed(fname, group="mm10")
seed
#> <27998 x 1306127> sparse TENxMatrixSeed object of type "integer":
#> # dirname: /github/home/.cache/R/ExperimentHub
#> # basename: 7b946c52a4b_1039
#> # group: /mm10
path(seed)
#> [1] "/github/home/.cache/R/ExperimentHub/7b946c52a4b_1039"
dim(seed)
#> [1] 27998 1306127
is_sparse(seed)
#> [1] TRUE
sparsity(seed)
#> [1] 0.9282225
DelayedArray(seed)
#> <27998 x 1306127> sparse TENxMatrix object of type "integer":
#> AAACCTGAGATAGGAG-1 ... TTTGTCATCTGAAAGA-133
#> ENSMUSG00000051951 0 . 0
#> ENSMUSG00000089699 0 . 0
#> ENSMUSG00000102343 0 . 0
#> ENSMUSG00000025900 0 . 0
#> ENSMUSG00000109048 0 . 0
#> ... . . .
#> ENSMUSG00000079808 0 . 0
#> ENSMUSG00000095041 1 . 0
#> ENSMUSG00000063897 0 . 0
#> ENSMUSG00000096730 0 . 0
#> ENSMUSG00000095742 0 . 0
stopifnot(class(DelayedArray(seed)) == "TENxMatrix")