AnnData-Conversion | R Documentation |
Conversion between Python AnnData objects and SingleCellExperiment::SingleCellExperiment objects.
AnnData2SCE(
adata,
X_name = NULL,
layers = TRUE,
uns = TRUE,
var = TRUE,
obs = TRUE,
varm = TRUE,
obsm = TRUE,
varp = TRUE,
obsp = TRUE,
raw = FALSE,
skip_assays = FALSE,
hdf5_backed = TRUE,
verbose = NULL
)
SCE2AnnData(
sce,
X_name = NULL,
assays = TRUE,
colData = TRUE,
rowData = TRUE,
varm = TRUE,
reducedDims = TRUE,
metadata = TRUE,
colPairs = TRUE,
rowPairs = TRUE,
skip_assays = FALSE,
verbose = NULL
)
adata |
A reticulate reference to a Python AnnData object. |
X_name |
For |
layers , uns , var , obs , varm , obsm , varp , obsp , raw |
Arguments specifying how
these slots are converted. If |
skip_assays |
Logical scalar indicating whether to skip conversion of
any assays in |
hdf5_backed |
Logical scalar indicating whether HDF5-backed matrices
in |
verbose |
Logical scalar indicating whether to print progress messages.
If |
sce |
A SingleCellExperiment::SingleCellExperiment object. |
assays , colData , rowData , reducedDims , metadata , colPairs , rowPairs |
Arguments specifying how these slots are converted. If |
These functions assume that an appropriate Python environment has already been loaded. As such, they are largely intended for developer use, most typically inside a basilisk context.
The conversion is not entirely lossless. The current mapping is shown below (also at https://tinyurl.com/AnnData2SCE):
In SCE2AnnData()
, matrices are converted to a numpy-friendly format.
Sparse matrices are converted to
Matrix::dgCMatrix objects while all
other matrices are converted into ordinary matrices. If skip_assays = TRUE
,
empty sparse matrices are created instead and the user is expected to fill in
the assays on the Python side.
For AnnData2SCE()
, a warning is raised if there is no corresponding R
format for a matrix in the AnnData
object, and an empty sparse matrix is
created instead as a placeholder. If skip_assays = NA
, no warning is
emitted but variables are created in the
int_metadata()
of the output to
specify which assays were skipped.
If skip_assays = TRUE
, empty sparse matrices are created for all assays,
regardless of whether they might be convertible to an R format or not.
In both cases, the user is expected to fill in the assays on the R side.
metadata
/uns
conversionWe attempt to convert between items in the
SingleCellExperiment::SingleCellExperiment
metadata()
slot and the AnnData
uns
slot. If
an item cannot be converted a warning will be raised.
uns
conversionValues stored in the varm
slot of an AnnData
object are stored in a
column of rowData()
in a
SingleCellExperiment::SingleCellExperiment
as a S4Vectors::DataFrame-class of matrices.
If this column is present an attempt is made to transfer this information
when converting from
SingleCellExperiment::SingleCellExperiment
to AnnData
.
SpatialExperiment
conversionIn SCE2AnnData()
, if sce
is a SpatialExperiment::SpatialExperiment
object, the spatial coordinates are added to the reducedDims
slot before
conversion to an AnnData
object.
AnnData2SCE()
will return a
SingleCellExperiment::SingleCellExperiment
containing the equivalent data from adata
.
SCE2AnnData()
will return a reticulate reference to an AnnData object
containing the content of sce
.
Luke Zappia
Aaron Lun
writeH5AD()
and readH5AD()
for dealing directly with H5AD files.
if (requireNamespace("scRNAseq", quietly = TRUE)) {
library(basilisk)
library(scRNAseq)
seger <- SegerstolpePancreasData()
# These functions are designed to be run inside
# a specified Python environment
roundtrip <- basiliskRun(fun = function(sce) {
# Convert SCE to AnnData:
adata <- zellkonverter::SCE2AnnData(sce)
# Maybe do some work in Python on 'adata':
# BLAH BLAH BLAH
# Convert back to an SCE:
zellkonverter::AnnData2SCE(adata)
}, env = zellkonverterAnnDataEnv(), sce = seger)
}
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.