This class is an abstract representation of an AnnData object. It is
intended to be used as a base class for concrete implementations of
AnnData objects, such as InMemoryAnnData
or HDF5AnnData
.
The following functions can be used to create an object that inherits
from AbstractAnnData
:
AnnData()
: Create an in-memory AnnData object.read_h5ad()
: Create an HDF5-backed AnnData object.from_Seurat()
: Create an in-memory AnnData object from a Seurat object.from_SingleCellExperiment()
: Create an in-memory AnnData object from a SingleCellExperiment object.
Active bindings
X
NULL or an observation x variable matrix (without dimnames) consistent with the number of rows of
obs
andvar
.layers
The layers slot. Must be NULL or a named list with with all elements having the dimensions consistent with
obs
andvar
.obs
A
data.frame
with columns containing information about observations. The number of rows ofobs
defines the observation dimension of the AnnData object.var
A
data.frame
with columns containing information about variables. The number of rows ofvar
defines the variable dimension of the AnnData object.obs_names
Either NULL or a vector of unique identifiers used to identify each row of
obs
and to act as an index into the observation dimension of the AnnData object. For compatibility with R representations,obs_names
should be a character vector.var_names
Either NULL or a vector of unique identifiers used to identify each row of
var
and to act as an index into the variable dimension of the AnnData object. For compatibility with R representations,var_names
should be a character vector.obsm
The obsm slot. Must be
NULL
or a named list with with all elements having the same number of rows asobs
.varm
The varm slot. Must be
NULL
or a named list with with all elements having the same number of rows asvar
.obsp
The obsp slot. Must be
NULL
or a named list with with all elements having the same number of rows and columns asobs
.varp
The varp slot. Must be
NULL
or a named list with with all elements having the same number of rows and columns asvar
.uns
The uns slot. Must be
NULL
or a named list.
Methods
Method print()
Print a summary of the AnnData object. print()
methods should be implemented so that they are not
computationally expensive.
Method to_Seurat()
Convert to Seurat
See to_Seurat()
for more details on the conversion and each of the parameters.
Usage
AbstractAnnData$to_Seurat(
assay_name = "RNA",
layers_mapping = NULL,
reduction_mapping = NULL,
graph_mapping = NULL,
misc_mapping = NULL
)
Arguments
assay_name
The name of the assay to use as the main data
layers_mapping
A named list mapping Seurat layers to AnnData layers
reduction_mapping
A named list mapping Seurat reductions to AnnData obsm/varm
graph_mapping
A named list mapping Seurat graphs to AnnData obsp/varp
misc_mapping
A named list mapping Seurat misc to AnnData uns
Method to_HDF5AnnData()
Convert to an HDF5 Backed AnnData
Arguments
file
The path to the HDF5 file
compression
The compression algorithm to use when writing the HDF5 file. Can be one of
"none"
,"gzip"
or"lzf"
. Defaults to"none"
.mode
The mode to open the HDF5 file.
a
creates a new file or opens an existing one for read/write.r+
opens an existing file for read/write.w
creates a file, truncating any existing onesw-
/x
are synonyms creating a file and failing if it already exists.
Method write_h5ad()
Write the AnnData object to an H5AD file.
Arguments
path
The path to the H5AD file
compression
The compression algorithm to use when writing the HDF5 file. Can be one of
"none"
,"gzip"
or"lzf"
. Defaults to"none"
.mode
The mode to open the HDF5 file.
a
creates a new file or opens an existing one for read/write.r+
opens an existing file for read/write.w
creates a file, truncating any existing onesw-
/x
are synonyms creating a file and failing if it already exists.
Examples
adata <- AnnData(
X = matrix(1:5, 3L, 5L),
layers = list(
A = matrix(5:1, 3L, 5L),
B = matrix(letters[1:5], 3L, 5L)
),
obs = data.frame(row.names = LETTERS[1:3], cell = 1:3),
var = data.frame(row.names = letters[1:5], gene = 1:5)
)
h5ad_file <- tempfile(fileext = ".h5ad")
adata$write_h5ad(h5ad_file)
Examples
## ------------------------------------------------
## Method `AbstractAnnData$write_h5ad`
## ------------------------------------------------
adata <- AnnData(
X = matrix(1:5, 3L, 5L),
layers = list(
A = matrix(5:1, 3L, 5L),
B = matrix(letters[1:5], 3L, 5L)
),
obs = data.frame(row.names = LETTERS[1:3], cell = 1:3),
var = data.frame(row.names = letters[1:5], gene = 1:5)
)
h5ad_file <- tempfile(fileext = ".h5ad")
adata$write_h5ad(h5ad_file)
#> [1] "/tmp/Rtmp3F3EQq/file1d144b02a93c.h5ad"