Write a H5AD file from a SingleCellExperiment object.

writeH5AD(
  sce,
  file,
  X_name = NULL,
  skip_assays = FALSE,
  compression = c("none", "gzip", "lzf"),
  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.

verbose

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

...

Arguments passed on to SCE2AnnData

varm

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.

assays

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.

colData

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.

rowData

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.

reducedDims

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.

metadata

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.

colPairs

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.

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.

Value

A NULL is invisibly returned.

Details

Environment

When first run, this function will instantiate a conda environment containing all of the necessary dependencies. This will not be performed on any subsequent run or if any other zellkonverter function has been run prior to this one.

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.

See also

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

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

Author

Luke Zappia

Aaron Lun

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)
}
#> snapshotDate(): 2021-10-18
#> see ?scRNAseq and browseVignettes('scRNAseq') for documentation
#> downloading 1 resources
#> retrieving 1 resource
#> loading from cache
#> see ?scRNAseq and browseVignettes('scRNAseq') for documentation
#> downloading 1 resources
#> retrieving 1 resource
#> loading from cache
#> see ?scRNAseq and browseVignettes('scRNAseq') for documentation
#> downloading 1 resources
#> retrieving 1 resource
#> loading from cache
#> snapshotDate(): 2021-10-18
#> see ?scRNAseq and browseVignettes('scRNAseq') for documentation
#> loading from cache
#>  Using the 'counts' assay as the X matrix