writeH5AD: Write H5AD

View source: R/write.R

writeH5ADR Documentation

Write H5AD

Description

Write a H5AD file from a SingleCellExperiment object.

Usage

writeH5AD(
  sce,
  file,
  X_name = NULL,
  skip_assays = FALSE,
  compression = c("none", "gzip", "lzf"),
  version = NULL,
  verbose = NULL,
  ...
)

Arguments

sce

A SingleCellExperiment object.

file

String containing a path to write the new .h5ad file.

X_name

Name of the assay to use as the primary matrix (X) of the AnnData object. If NULL, the first assay of sce will be used by default.

skip_assays

Logical scalar indicating whether assay matrices should be ignored when writing to file.

compression

Type of compression when writing the new .h5ad file.

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 SCE2AnnData

assays,colData,rowData,reducedDims,metadata,colPairs,rowPairs

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.

Details

Skipping assays

Setting skip_assays = TRUE can occasionally be useful if the matrices in sce are stored in a format that is not amenable for efficient conversion to a numpy-compatible format. In such cases, it can be better to create an empty placeholder dataset in file and fill it in R afterwards.

DelayedArray assays

If sce contains any DelayedArray matrices as assays writeH5AD() will write them to disk using the rhdf5 package directly rather than via Python to avoid instantiating them in memory. However there is currently an issue which prevents this being done for sparse DelayedArray matrices.

Known conversion issues

Coercion to factors

The anndata package automatically converts some character vectors to factors when saving .h5ad files. This can effect columns of rowData(sce) and colData(sce) which may change type when the .h5ad file is read back into R.

Environment

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

Value

A NULL is invisibly returned.

Author(s)

Luke Zappia

Aaron Lun

See Also

readH5AD(), to read a SingleCellExperiment file from a H5AD file.

SCE2AnnData(), for developers to create an AnnData object from a SingleCellExperiment.

Examples

# Using the Zeisel brain dataset
if (requireNamespace("scRNAseq", quietly = TRUE)) {
    library(scRNAseq)
    sce <- ZeiselBrainData()

    # Writing to a H5AD file
    temp <- tempfile(fileext = ".h5ad")
    writeH5AD(sce, temp)
}

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