readH5AD: Read H5AD

View source: R/read.R

readH5ADR Documentation

Read H5AD

Description

Reads a H5AD file and returns a SingleCellExperiment object.

Usage

readH5AD(
  file,
  X_name = NULL,
  use_hdf5 = FALSE,
  reader = c("python", "R"),
  version = NULL,
  verbose = NULL,
  ...
)

Arguments

file

String containing a path to a .h5ad file.

X_name

Name used when saving X as an assay. If NULL looks for an X_name value in uns, otherwise uses "X".

use_hdf5

Logical scalar indicating whether assays should be loaded as HDF5-based matrices from the HDF5Array package.

reader

Which HDF5 reader to use. Either "python" for reading with the anndata Python package via reticulate or "R" for zellkonverter's native R reader.

version

A string giving the version of the anndata Python library to use. Allowed values are available in .AnnDataVersions. By default the latest version is used.

verbose

Logical scalar indicating whether to print progress messages. If NULL uses getOption("zellkonverter.verbose").

...

Arguments passed on to AnnData2SCE

layers,uns,var,obs,varm,obsm,varp,obsp,raw

Arguments specifying how these slots are converted. If TRUE everything in that slot is converted, if FALSE nothing is converted and if a character vector only those items or columns are converted.

skip_assays

Logical scalar indicating whether to skip conversion of any assays in sce or adata, replacing them with empty sparse matrices instead.

Details

Setting use_hdf5 = TRUE allows for very large datasets to be efficiently represented on machines with little memory. However, this comes at the cost of access speed as data needs to be fetched from the HDF5 file upon request.

Setting reader = "R" will use an experimental native R reader instead of reading the file into Python and converting the result. This avoids the need for a Python environment and some of the issues with conversion but is still under development and is likely to return slightly different output.

See AnnData-Environment for more details on zellkonverter Python environments.

Value

A SingleCellExperiment object is returned.

Author(s)

Luke Zappia

Aaron Lun

See Also

writeH5AD(), to write a SingleCellExperiment object to a H5AD file.

AnnData2SCE(), for developers to convert existing AnnData instances to a SingleCellExperiment.

Examples

library(SummarizedExperiment)

file <- system.file("extdata", "krumsiek11.h5ad", package = "zellkonverter")
sce <- readH5AD(file)
class(assay(sce))

sce2 <- readH5AD(file, use_hdf5 = TRUE)
class(assay(sce2))

sce3 <- readH5AD(file, reader = "R")

theislab/zellkonverter documentation built on Oct. 24, 2024, 7:43 a.m.