H5SparseMatrix-class: HDF5 sparse matrices as DelayedMatrix objects

H5SparseMatrix-classR Documentation

HDF5 sparse matrices as DelayedMatrix objects

Description

The H5SparseMatrix class is a DelayedMatrix subclass for representing and operating on an HDF5 sparse matrix stored in CSR/CSC/Yale format.

All the operations available for DelayedMatrix objects work on H5SparseMatrix objects.

Usage

## Constructor function:
H5SparseMatrix(filepath, group)

Arguments

filepath

The path (as a single string) to the HDF5 file (.h5 or .h5ad) where the sparse matrix is located.

group

The name of the group in the HDF5 file where the sparse matrix is stored.

Value

An H5SparseMatrix object.

See Also

  • HDF5Array objects for representing conventional (a.k.a. dense) HDF5 datasets as DelayedArray 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.

  • DelayedMatrix objects in the DelayedArray package.

  • The H5SparseMatrixSeed helper class.

  • h5ls to list the content of an HDF5 file (.h5 or .h5ad).

  • SparseArray objects in the SparseArray package.

Examples

library(zellkonverter)
h5ad_file <- system.file("extdata", "example_anndata.h5ad",
                         package="zellkonverter")
rhdf5::h5ls(h5ad_file)

M <- H5SparseMatrix(h5ad_file, "/obsp/connectivities")
M

class(M)  # H5SparseMatrix
is(M, "DelayedMatrix")  # TRUE

seed(M)
class(seed(M))  # CSC_H5SparseMatrixSeed

dim(M)
path(M)
is_sparse(M)  # TRUE

## Use coercion to load the full dataset into memory:
as.matrix(M)          # as ordinary array (usually not recommended)
as(M, "dgCMatrix")    # as dgCMatrix
as(M, "SparseArray")  # as SparseArray object (most efficient)
SparseArray(M)        # equivalent to 'as(M, "SparseArray")'

Bioconductor/HDF5Array documentation built on Oct. 31, 2024, 9:16 a.m.