Title: | AnnData interoperability in R |
---|---|
Description: | Bring the power and flexibility of AnnData to the R ecosystem, allowing you to effortlessly manipulate and analyze your single-cell data. This package lets you work with backed h5ad and zarr files, directly access various slots (e.g. X, obs, var), or convert the data into SingleCellExperiment and Seurat objects. |
Authors: | Robrecht Cannoodt [aut, cre] (<https://orcid.org/0000-0003-3641-729X>, rcannood), Luke Zappia [aut] (<https://orcid.org/0000-0001-7744-8565>, lazappi), Martin Morgan [aut] (<https://orcid.org/0000-0002-5874-8148>, mtmorgan), Louise Deconinck [aut] (<https://orcid.org/0000-0001-8100-6823>, LouiseDck), Danila Bredikhin [ctb] (<https://orcid.org/0000-0001-8089-6983>, gtca), Isaac Virshup [ctb] (<https://orcid.org/0000-0002-1710-8945>, ivirshup), Brian Schilder [ctb] (<https://orcid.org/0000-0001-5949-2191>, bschilder), Chananchida Sang-aram [ctb] (<https://orcid.org/0000-0002-0922-0822>, csangara), Data Intuitive [fnd, cph], Chan Zuckerberg Initiative [fnd] |
Maintainer: | Robrecht Cannoodt <[email protected]> |
License: | MIT + file LICENSE |
Version: | 0.99.0 |
Built: | 2025-02-06 21:38:16 UTC |
Source: | https://github.com/scverse/anndataR |
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.
X
NULL or an observation x variable matrix (without
dimnames) consistent with the number of rows of obs
and var
.
layers
The layers slot. Must be NULL or a named list
with with all elements having the dimensions consistent with
obs
and var
.
obs
A data.frame
with columns containing information
about observations. The number of rows of obs
defines the
observation dimension of the AnnData object.
var
A data.frame
with columns containing information
about variables. The number of rows of var
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 as obs
.
varm
The varm slot. Must be NULL
or a named list with
with all elements having the same number of rows as var
.
obsp
The obsp slot. Must be NULL
or a named list with
with all elements having the same number of rows and columns as obs
.
varp
The varp slot. Must be NULL
or a named list with
with all elements having the same number of rows and columns as var
.
uns
The uns slot. Must be NULL
or a named list.
print()
Print a summary of the AnnData object. print()
methods should be implemented so that they are not
computationally expensive.
AbstractAnnData$print(...)
...
Optional arguments to print method.
shape()
Dimensions (observations x variables) of the AnnData object.
AbstractAnnData$shape()
n_obs()
Number of observations in the AnnData object.
AbstractAnnData$n_obs()
n_vars()
Number of variables in the AnnData object.
AbstractAnnData$n_vars()
obs_keys()
Keys ('column names') of obs
.
AbstractAnnData$obs_keys()
var_keys()
Keys ('column names') of var
.
AbstractAnnData$var_keys()
layers_keys()
Keys (element names) of layers
.
AbstractAnnData$layers_keys()
obsm_keys()
Keys (element names) of obsm
.
AbstractAnnData$obsm_keys()
varm_keys()
Keys (element names) of varm
.
AbstractAnnData$varm_keys()
obsp_keys()
Keys (element names) of obsp
.
AbstractAnnData$obsp_keys()
varp_keys()
Keys (element names) of varp
.
AbstractAnnData$varp_keys()
uns_keys()
Keys (element names) of uns
.
AbstractAnnData$uns_keys()
to_SingleCellExperiment()
Convert to SingleCellExperiment
See to_SingleCellExperiment()
for more details on the conversion.
AbstractAnnData$to_SingleCellExperiment( assays_mapping = NULL, colData_mapping = NULL, rowData_mapping = NULL, reduction_mapping = NULL, colPairs_mapping = NULL, rowPairs_mapping = NULL, metadata_mapping = NULL )
assays_mapping
A named list mapping SingleCellExperiment assays to AnnData layers
colData_mapping
A named list mapping SingleCellExperiment colData to AnnData obs
rowData_mapping
A named list mapping SingleCellExperiment rowData to AnnData var
reduction_mapping
A named list mapping SingleCellExperiment reductions to AnnData obsm/varm
colPairs_mapping
A named list mapping SingleCellExperiment colPairs to AnnData obsp/varp
rowPairs_mapping
A named list mapping SingleCellExperiment rowPairs to AnnData obsp/varp
metadata_mapping
A named list mapping SingleCellExperiment metadata to AnnData uns
A SingleCellExperiment object
to_Seurat()
Convert to Seurat
See to_Seurat()
for more details on the conversion and each of the parameters.
AbstractAnnData$to_Seurat( assay_name = "RNA", layers_mapping = NULL, reduction_mapping = NULL, graph_mapping = NULL, misc_mapping = NULL )
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
A Seurat object
to_InMemoryAnnData()
Convert to an InMemory AnnData
AbstractAnnData$to_InMemoryAnnData()
to_HDF5AnnData()
Convert to an HDF5 Backed AnnData
AbstractAnnData$to_HDF5AnnData( file, compression = c("none", "gzip", "lzf"), mode = c("w-", "r", "r+", "a", "w", "x") )
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 ones
w-
/x
are synonyms creating a file and failing if it already exists.
An HDF5AnnData object
write_h5ad()
Write the AnnData object to an H5AD file.
AbstractAnnData$write_h5ad(path, compression = "gzip", mode = "w")
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 ones
w-
/x
are synonyms creating a file and failing if it already exists.
path
invisibly
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)
clone()
The objects of this class are cloneable with this method.
AbstractAnnData$clone(deep = FALSE)
deep
Whether to make a deep clone.
## ------------------------------------------------ ## 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)
## ------------------------------------------------ ## 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)
For more information on the functionality of an AnnData object, see anndataR-package.
AnnData( X = NULL, obs = NULL, var = NULL, layers = NULL, obsm = NULL, varm = NULL, obsp = NULL, varp = NULL, uns = NULL, shape = NULL )
AnnData( X = NULL, obs = NULL, var = NULL, layers = NULL, obsm = NULL, varm = NULL, obsp = NULL, varp = NULL, uns = NULL, shape = NULL )
X |
Either |
obs |
Either |
var |
Either |
layers |
Either |
obsm |
The obsm slot is used to store multi-dimensional annotation
arrays. It must be either |
varm |
The varm slot is used to store multi-dimensional annotation
arrays. It must be either |
obsp |
The obsp slot is used to store sparse multi-dimensional
annotation arrays. It must be either |
varp |
The varp slot is used to store sparse multi-dimensional
annotation arrays. It must be either |
uns |
The uns slot is used to store unstructured annotation. It must
be either |
shape |
Shape tuple (#observations, #variables). Can be provided
if |
An in-memory AnnData object.
adata <- AnnData( X = matrix(1:12, nrow = 3, ncol = 4), obs = data.frame( row.names = paste0("obs", 1:3), n_counts = c(1, 2, 3), n_cells = c(1, 2, 3) ), var = data.frame( row.names = paste0("var", 1:4), n_cells = c(1, 2, 3, 4) ) ) adata
adata <- AnnData( X = matrix(1:12, nrow = 3, ncol = 4), obs = data.frame( row.names = paste0("obs", 1:3), n_counts = c(1, 2, 3), n_cells = c(1, 2, 3) ), var = data.frame( row.names = paste0("var", 1:4), n_cells = c(1, 2, 3, 4) ) ) adata
from_Seurat()
converts a Seurat object to an AnnData object.
Only one assay can be converted at a time. Arguments are used to configure the conversion.
If NULL
, the functions from_Seurat_guess_*
will be used to guess the mapping.
from_Seurat( seurat_obj, output_class = c("InMemoryAnnData", "HDF5AnnData"), assay_name = NULL, x_mapping = NULL, layers_mapping = NULL, obsm_mapping = NULL, varm_mapping = NULL, obsp_mapping = NULL, varp_mapping = NULL, uns_mapping = NULL, ... )
from_Seurat( seurat_obj, output_class = c("InMemoryAnnData", "HDF5AnnData"), assay_name = NULL, x_mapping = NULL, layers_mapping = NULL, obsm_mapping = NULL, varm_mapping = NULL, obsp_mapping = NULL, varp_mapping = NULL, uns_mapping = NULL, ... )
seurat_obj |
A Seurat object to be converted. |
output_class |
Name of the AnnData class. Must be one of |
assay_name |
The name of the assay to be converted. If |
x_mapping |
A mapping of a Seurat layer to the AnnData |
layers_mapping |
A named list mapping layer names to the names of the layers in the Seurat object. Each item in
the list must be a character vector of length 1. See section " |
obsm_mapping |
A named list mapping reductions to the names of the reductions in the Seurat object. Each item in
the list must be a vector of length 2. See section " |
varm_mapping |
A named list mapping PCA loadings to the names of the PCA loadings in the Seurat object.
Each item in the list must be a character vector of length 1. See section " |
obsp_mapping |
A named list mapping graph names to the names of the graphs in the Seurat object.
Each item in the list must be a character vector of length 1. See section " |
varp_mapping |
A named list mapping miscellaneous data to the names of the data in the Seurat object.
Each item in the list must be a named list with one or two elements. See section " |
uns_mapping |
A named list mapping miscellaneous data to the names of the data in the Seurat object.
Each item in the list must be a named list with one or two elements. See section " |
... |
Additional arguments passed to the generator function. |
For more information on the functionality of an AnnData object, see anndataR-package.
An AnnData object
$X
mappingA mapping of a Seurat layer to the AnnData X
slot. Its value must be NULL
or a character vector of the Seurat
layer name to copy. If NULL
, no data will be copied to the X
slot.
$layers
mappingA named list to map AnnData layers to Seurat layers. Each item in the list must be a character vector of length 1.
The $X
key maps to the X
slot.
Example: layers_mapping = list(counts = "counts", foo = "bar")
.
If NULL
, the internal function from_Seurat_guess_layers
will be used to guess the layer mapping as follows:
All Seurat layers are copied to AnnData layers
by name.
This means that the AnnData X
slot will be NULL
(empty). If you want to copy data to the X
slot,
you must define the layer mapping explicitly.
$obsm
mappingA named list to map Seurat reductions to AnnData $obsm
.
Each item in the list must be a vector of length 2,
where the name corresponds to the name of the resulting $obsm
slot, and the values correspond to the
the location of the data in the Seurat object.
Example: obsm_mapping = list(pca = c("reductions", "pca"), umap = c("reductions", "umap"))
.
If NULL
, the internal function from_Seurat_guess_obsms
will be used to guess the obsm mapping as follows:
All Seurat reductions are prefixed with X_
and copied to AnnData $obsm
.
$varm
mappingA named list to map Seurat reduction loadings to AnnData $varm
.
Each item in the list must be a character vector of length 2, where the name corresponds to the name of the
resulting $varm
slot, and the value corresponds to the location of the data in the Seurat object.
Example: varm_mapping = list(PCs = c("reductions", "pca")
.
If NULL
, the internal function from_Seurat_guess_varms
will be used to guess the varm mapping as follows:
The name of the PCA loadings is copied by name.
$obsp
mappingA named list to map Seurat graphs to AnnData $obsp
.
Example: obsp_mapping = list(nn = "connectivities")
.
If NULL
, the internal function from_Seurat_guess_obsps
will be used to guess the obsp mapping as follows:
All Seurat graphs are copied to $obsp
by name.
$varp
mappingA named list to map Seurat miscellaneous data to AnnData $varp
. The name of each item corresponds to the
resulting $varp
slot, while the value of each item must be a fector which corresponds to the location of the data
in the Seurat object.
Example: varp_mapping = list(foo = c("misc", "foo"))
.
If NULL
, the internal function from_Seurat_guess_varps
will be used to guess the varp mapping as follows:
No data is mapped to $varp
.
$uns
mappingA named list to map Seurat miscellaneous data to AnnData uns
. Each item in the list must be a character of
length 2. The first element must be "misc"
. The second element is the name of the data in the corresponding slot.
Example: uns_mapping = list(foo = c("misc", "foo"))
.
If NULL
, the internal function from_Seurat_guess_uns
will be used to guess the uns mapping as follows:
All Seurat miscellaneous data is copied to uns
by name.
library(Seurat) counts <- matrix(rbinom(20000, 1000, .001), nrow = 100) obj <- CreateSeuratObject(counts = counts) obj <- NormalizeData(obj) obj <- FindVariableFeatures(obj) obj <- ScaleData(obj) obj <- RunPCA(obj, npcs = 10L) obj <- FindNeighbors(obj) obj <- RunUMAP(obj, dims = 1:10) from_Seurat(obj)
library(Seurat) counts <- matrix(rbinom(20000, 1000, .001), nrow = 100) obj <- CreateSeuratObject(counts = counts) obj <- NormalizeData(obj) obj <- FindVariableFeatures(obj) obj <- ScaleData(obj) obj <- RunPCA(obj, npcs = 10L) obj <- FindNeighbors(obj) obj <- RunUMAP(obj, dims = 1:10) from_Seurat(obj)
from_SingleCellExperiment()
converts a
SingleCellExperiment to an AnnData object.
from_SingleCellExperiment( sce, output_class = c("InMemory", "HDF5AnnData"), x_mapping = NULL, layers_mapping = NULL, obs_mapping = NULL, var_mapping = NULL, obsm_mapping = NULL, varm_mapping = NULL, obsp_mapping = NULL, varp_mapping = NULL, uns_mapping = NULL, ... )
from_SingleCellExperiment( sce, output_class = c("InMemory", "HDF5AnnData"), x_mapping = NULL, layers_mapping = NULL, obs_mapping = NULL, var_mapping = NULL, obsm_mapping = NULL, varm_mapping = NULL, obsp_mapping = NULL, varp_mapping = NULL, uns_mapping = NULL, ... )
sce |
An object inheriting from SingleCellExperiment. |
output_class |
Name of the AnnData class. Must be one of |
x_mapping |
Name of the assay in |
layers_mapping |
A named list mapping |
obs_mapping |
A named list mapping |
var_mapping |
A named list mapping |
obsm_mapping |
A named list mapping
|
varm_mapping |
A named list mapping
|
obsp_mapping |
A named list mapping |
varp_mapping |
A named list mapping |
uns_mapping |
A named list mapping |
... |
Additional arguments to pass to the generator function. |
from_SingleCellExperiment()
returns an AnnData object
(e.g., InMemoryAnnData) representing the content of sce
.
## construct an AnnData object from a SingleCellExperiment library(SingleCellExperiment) sce <- SingleCellExperiment( assays = list(counts = matrix(1:5, 5L, 3L)), colData = DataFrame(cell = 1:3), rowData = DataFrame(gene = 1:5) ) from_SingleCellExperiment(sce, "InMemory")
## construct an AnnData object from a SingleCellExperiment library(SingleCellExperiment) sce <- SingleCellExperiment( assays = list(counts = matrix(1:5, 5L, 3L)), colData = DataFrame(cell = 1:3), rowData = DataFrame(gene = 1:5) ) from_SingleCellExperiment(sce, "InMemory")
Generate a dataset with different types of columns and layers
generate_dataset( n_obs = 10L, n_vars = 20L, x_type = "numeric_matrix", layer_types = c("numeric_matrix", "numeric_dense", "numeric_csparse", "numeric_rsparse", "numeric_matrix_with_nas", "numeric_dense_with_nas", "numeric_csparse_with_nas", "numeric_rsparse_with_nas", "integer_matrix", "integer_dense", "integer_csparse", "integer_rsparse", "integer_matrix_with_nas", "integer_dense_with_nas", "integer_csparse_with_nas", "integer_rsparse_with_nas"), obs_types = c("character", "integer", "factor", "factor_ordered", "logical", "numeric", "character_with_nas", "integer_with_nas", "factor_with_nas", "factor_ordered_with_nas", "logical_with_nas", "numeric_with_nas"), var_types = c("character", "integer", "factor", "factor_ordered", "logical", "numeric", "character_with_nas", "integer_with_nas", "factor_with_nas", "factor_ordered_with_nas", "logical_with_nas", "numeric_with_nas"), obsm_types = c("numeric_matrix", "numeric_dense", "numeric_csparse", "numeric_rsparse", "numeric_matrix_with_nas", "numeric_dense_with_nas", "numeric_csparse_with_nas", "numeric_rsparse_with_nas", "integer_matrix", "integer_dense", "integer_csparse", "integer_rsparse", "integer_matrix_with_nas", "integer_dense_with_nas", "integer_csparse_with_nas", "integer_rsparse_with_nas", "character", "integer", "factor", "factor_ordered", "logical", "numeric", "character_with_nas", "integer_with_nas", "factor_with_nas", "factor_ordered_with_nas", "logical_with_nas", "numeric_with_nas"), varm_types = c("numeric_matrix", "numeric_dense", "numeric_csparse", "numeric_rsparse", "numeric_matrix_with_nas", "numeric_dense_with_nas", "numeric_csparse_with_nas", "numeric_rsparse_with_nas", "integer_matrix", "integer_dense", "integer_csparse", "integer_rsparse", "integer_matrix_with_nas", "integer_dense_with_nas", "integer_csparse_with_nas", "integer_rsparse_with_nas", "character", "integer", "factor", "factor_ordered", "logical", "numeric", "character_with_nas", "integer_with_nas", "factor_with_nas", "factor_ordered_with_nas", "logical_with_nas", "numeric_with_nas"), obsp_types = c("numeric_matrix", "numeric_dense", "numeric_csparse", "numeric_rsparse", "numeric_matrix_with_nas", "numeric_dense_with_nas", "numeric_csparse_with_nas", "numeric_rsparse_with_nas", "integer_matrix", "integer_dense", "integer_csparse", "integer_rsparse", "integer_matrix_with_nas", "integer_dense_with_nas", "integer_csparse_with_nas", "integer_rsparse_with_nas"), varp_types = c("numeric_matrix", "numeric_dense", "numeric_csparse", "numeric_rsparse", "numeric_matrix_with_nas", "numeric_dense_with_nas", "numeric_csparse_with_nas", "numeric_rsparse_with_nas", "integer_matrix", "integer_dense", "integer_csparse", "integer_rsparse", "integer_matrix_with_nas", "integer_dense_with_nas", "integer_csparse_with_nas", "integer_rsparse_with_nas"), uns_types = c("scalar_character", "scalar_integer", "scalar_factor", "scalar_factor_ordered", "scalar_logical", "scalar_numeric", "scalar_character_with_nas", "scalar_integer_with_nas", "scalar_factor_with_nas", "scalar_factor_ordered_with_nas", "scalar_logical_with_nas", "scalar_numeric_with_nas", "vec_character", "vec_integer", "vec_factor", "vec_factor_ordered", "vec_logical", "vec_numeric", "vec_character_with_nas", "vec_integer_with_nas", "vec_factor_with_nas", "vec_factor_ordered_with_nas", "vec_logical_with_nas", "vec_numeric_with_nas", "df_character", "df_integer", "df_factor", "df_factor_ordered", "df_logical", "df_numeric", "df_character_with_nas", "df_integer_with_nas", "df_factor_with_nas", "df_factor_ordered_with_nas", "df_logical_with_nas", "df_numeric_with_nas", "mat_numeric_matrix", "mat_numeric_dense", "mat_numeric_csparse", "mat_numeric_rsparse", "mat_numeric_matrix_with_nas", "mat_numeric_dense_with_nas", "mat_numeric_csparse_with_nas", "mat_numeric_rsparse_with_nas", "mat_integer_matrix", "mat_integer_dense", "mat_integer_csparse", "mat_integer_rsparse", "mat_integer_matrix_with_nas", "mat_integer_dense_with_nas", "mat_integer_csparse_with_nas", "mat_integer_rsparse_with_nas", "list"), example = FALSE, format = c("list", "AnnData", "SingleCellExperiment", "Seurat") )
generate_dataset( n_obs = 10L, n_vars = 20L, x_type = "numeric_matrix", layer_types = c("numeric_matrix", "numeric_dense", "numeric_csparse", "numeric_rsparse", "numeric_matrix_with_nas", "numeric_dense_with_nas", "numeric_csparse_with_nas", "numeric_rsparse_with_nas", "integer_matrix", "integer_dense", "integer_csparse", "integer_rsparse", "integer_matrix_with_nas", "integer_dense_with_nas", "integer_csparse_with_nas", "integer_rsparse_with_nas"), obs_types = c("character", "integer", "factor", "factor_ordered", "logical", "numeric", "character_with_nas", "integer_with_nas", "factor_with_nas", "factor_ordered_with_nas", "logical_with_nas", "numeric_with_nas"), var_types = c("character", "integer", "factor", "factor_ordered", "logical", "numeric", "character_with_nas", "integer_with_nas", "factor_with_nas", "factor_ordered_with_nas", "logical_with_nas", "numeric_with_nas"), obsm_types = c("numeric_matrix", "numeric_dense", "numeric_csparse", "numeric_rsparse", "numeric_matrix_with_nas", "numeric_dense_with_nas", "numeric_csparse_with_nas", "numeric_rsparse_with_nas", "integer_matrix", "integer_dense", "integer_csparse", "integer_rsparse", "integer_matrix_with_nas", "integer_dense_with_nas", "integer_csparse_with_nas", "integer_rsparse_with_nas", "character", "integer", "factor", "factor_ordered", "logical", "numeric", "character_with_nas", "integer_with_nas", "factor_with_nas", "factor_ordered_with_nas", "logical_with_nas", "numeric_with_nas"), varm_types = c("numeric_matrix", "numeric_dense", "numeric_csparse", "numeric_rsparse", "numeric_matrix_with_nas", "numeric_dense_with_nas", "numeric_csparse_with_nas", "numeric_rsparse_with_nas", "integer_matrix", "integer_dense", "integer_csparse", "integer_rsparse", "integer_matrix_with_nas", "integer_dense_with_nas", "integer_csparse_with_nas", "integer_rsparse_with_nas", "character", "integer", "factor", "factor_ordered", "logical", "numeric", "character_with_nas", "integer_with_nas", "factor_with_nas", "factor_ordered_with_nas", "logical_with_nas", "numeric_with_nas"), obsp_types = c("numeric_matrix", "numeric_dense", "numeric_csparse", "numeric_rsparse", "numeric_matrix_with_nas", "numeric_dense_with_nas", "numeric_csparse_with_nas", "numeric_rsparse_with_nas", "integer_matrix", "integer_dense", "integer_csparse", "integer_rsparse", "integer_matrix_with_nas", "integer_dense_with_nas", "integer_csparse_with_nas", "integer_rsparse_with_nas"), varp_types = c("numeric_matrix", "numeric_dense", "numeric_csparse", "numeric_rsparse", "numeric_matrix_with_nas", "numeric_dense_with_nas", "numeric_csparse_with_nas", "numeric_rsparse_with_nas", "integer_matrix", "integer_dense", "integer_csparse", "integer_rsparse", "integer_matrix_with_nas", "integer_dense_with_nas", "integer_csparse_with_nas", "integer_rsparse_with_nas"), uns_types = c("scalar_character", "scalar_integer", "scalar_factor", "scalar_factor_ordered", "scalar_logical", "scalar_numeric", "scalar_character_with_nas", "scalar_integer_with_nas", "scalar_factor_with_nas", "scalar_factor_ordered_with_nas", "scalar_logical_with_nas", "scalar_numeric_with_nas", "vec_character", "vec_integer", "vec_factor", "vec_factor_ordered", "vec_logical", "vec_numeric", "vec_character_with_nas", "vec_integer_with_nas", "vec_factor_with_nas", "vec_factor_ordered_with_nas", "vec_logical_with_nas", "vec_numeric_with_nas", "df_character", "df_integer", "df_factor", "df_factor_ordered", "df_logical", "df_numeric", "df_character_with_nas", "df_integer_with_nas", "df_factor_with_nas", "df_factor_ordered_with_nas", "df_logical_with_nas", "df_numeric_with_nas", "mat_numeric_matrix", "mat_numeric_dense", "mat_numeric_csparse", "mat_numeric_rsparse", "mat_numeric_matrix_with_nas", "mat_numeric_dense_with_nas", "mat_numeric_csparse_with_nas", "mat_numeric_rsparse_with_nas", "mat_integer_matrix", "mat_integer_dense", "mat_integer_csparse", "mat_integer_rsparse", "mat_integer_matrix_with_nas", "mat_integer_dense_with_nas", "mat_integer_csparse_with_nas", "mat_integer_rsparse_with_nas", "list"), example = FALSE, format = c("list", "AnnData", "SingleCellExperiment", "Seurat") )
n_obs |
Number of observations to generate |
n_vars |
Number of variables to generate |
x_type |
Type of matrix to generate for X |
layer_types |
Types of matrices to generate for layers |
obs_types |
Types of vectors to generate for obs |
var_types |
Types of vectors to generate for var |
obsm_types |
Types of matrices to generate for obsm |
varm_types |
Types of matrices to generate for varm |
obsp_types |
Types of matrices to generate for obsp |
varp_types |
Types of matrices to generate for varp |
uns_types |
Types of objects to generate for uns |
example |
If |
format |
Object type to output, one of "list", "AnnData", "SingleCellExperiment", or "Seurat". |
Object containing the generated dataset as defined by output
dummy <- generate_dataset() ## Not run: dummy <- generate_dataset(format = "AnnData") dummy <- generate_dataset(format = "SingleCellExperiment") dummy <- generate_dataset(format = "Seurat") ## End(Not run)
dummy <- generate_dataset() ## Not run: dummy <- generate_dataset(format = "AnnData") dummy <- generate_dataset(format = "SingleCellExperiment") dummy <- generate_dataset(format = "Seurat") ## End(Not run)
Read data from a H5AD file
read_h5ad( path, to = c("InMemoryAnnData", "HDF5AnnData", "SingleCellExperiment", "Seurat"), mode = c("r", "r+", "a", "w", "w-", "x"), ... )
read_h5ad( path, to = c("InMemoryAnnData", "HDF5AnnData", "SingleCellExperiment", "Seurat"), mode = c("r", "r+", "a", "w", "w-", "x"), ... )
path |
Path to the H5AD file to read |
to |
The type of object to return. Must be one of: "InMemoryAnnData", "HDF5AnnData", "SingleCellExperiment", "Seurat" |
mode |
The mode to open the HDF5 file.
|
... |
Extra arguments provided to |
The object specified by to
h5ad_file <- system.file("extdata", "example.h5ad", package = "anndataR") # Read the H5AD as a SingleCellExperiment object if (requireNamespace("SingleCellExperiment", quietly = TRUE)) { sce <- read_h5ad(h5ad_file, to = "SingleCellExperiment") } # Read the H5AD as a Seurat object if (requireNamespace("SeuratObject", quietly = TRUE)) { seurat <- read_h5ad(h5ad_file, to = "Seurat") }
h5ad_file <- system.file("extdata", "example.h5ad", package = "anndataR") # Read the H5AD as a SingleCellExperiment object if (requireNamespace("SingleCellExperiment", quietly = TRUE)) { sce <- read_h5ad(h5ad_file, to = "SingleCellExperiment") } # Read the H5AD as a Seurat object if (requireNamespace("SeuratObject", quietly = TRUE)) { seurat <- read_h5ad(h5ad_file, to = "Seurat") }
to_Seurat()
converts an AnnData object to a Seurat object. Only one assay can be converted at a time.
Arguments are used to configure the conversion. If NULL
, the functions to_Seurat_guess_*
will be used to guess
the mapping.
to_Seurat( adata, assay_name = "RNA", layers_mapping = NULL, reduction_mapping = NULL, graph_mapping = NULL, misc_mapping = NULL )
to_Seurat( adata, assay_name = "RNA", layers_mapping = NULL, reduction_mapping = NULL, graph_mapping = NULL, misc_mapping = NULL )
adata |
An AnnData object to be converted |
assay_name |
Name of the assay to be created (default: "RNA"). |
layers_mapping |
A named list to map AnnData layers to Seurat layers. See section "Layer mapping" for more details. |
reduction_mapping |
A named list to map AnnData reductions to Seurat reductions. Each item in the list must be a named list with keys 'key', 'obsm', and 'varm'. See section "Reduction mapping" for more details. |
graph_mapping |
A named list to map AnnData graphs to Seurat graphs. Each item in the list must be a character vector of length 1. See section "Graph mapping" for more details. |
misc_mapping |
A named list to map miscellaneous data to the names of the data in the Seurat object. See section "Miscellaneous mapping" for more details. |
A Seurat object
A named list to map AnnData layers to Seurat layers. Each item in the list must be a character vector of length 1,
where the values correspond to the names of the layers in the AnnData object, and the names correspond
to the names of the layers in the resulting Seurat object. A value of NULL
corresponds to the AnnData X
slot.
Example: layers_mapping = list(counts = "counts", data = NULL, foo = "bar")
.
If NULL
, the internal function to_Seurat_guess_layers
will be used to guess the layer mapping as follows:
All AnnData layers are copied to Seurat layers by name.
A named list to map AnnData $obsm
and $varm
to Seurat reductions. Each item in the list must be a named list
with keys 'key'
, 'obsm'
, and 'varm'
.
Example: reduction_mapping = list(pca = list(key = "PC_", obsm = "X_pca", varm = "PCs"))
.
If NULL
, the internal function to_Seurat_guess_reductions
will be used to guess the reduction mapping as follows:
All $obsm
items starting with X_
are copied by name.
A named list mapping graph names to the names of the graphs in the AnnData object. Each item in the list must be a character vector of length 1. The values correspond to the names of the graphs in the resulting Seurat object, while the names correspond to the names of the graphs in the AnnData object.
Example: graph_mapping = list(nn = "connectivities")
.
If NULL
, the internal function to_Seurat_guess_graphs
will be used to guess the graph mapping as follows:
An obsp named connectivities
will be mapped to nn
.
Other graphs starting with connectivities_
are stripped of the prefix and copied by name.
A named list mapping miscellaneous data to the names of the data in the AnnData object. Each item in the list must be a vector with one or two elements. The first element must be one of: 'X', 'layers', 'obs', 'obsm', 'obsp', 'var', 'varm', 'varp', 'uns'. The second element is the name of the data in the corresponding slot. If the second element is not present, the whole slot as specified by the first element will be used.
Example: misc_mapping = list(uns = "uns", varp_neighbors = c("varp", "neighbors"))
.
If NULL
, the internal function to_Seurat_guess_misc
will be used to guess the miscellaneous mapping as follows:
If $uns
is defined, all values in $uns
are copied to the Seurat misc.
ad <- AnnData( X = matrix(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) ) to_Seurat(ad)
ad <- AnnData( X = matrix(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) ) to_Seurat(ad)
to_SingleCellExperiment()
converts an AnnData object
to a SingleCellExperiment object.
to_SingleCellExperiment( adata, assays_mapping = NULL, colData_mapping = NULL, rowData_mapping = NULL, reduction_mapping = NULL, colPairs_mapping = NULL, rowPairs_mapping = NULL, metadata_mapping = NULL )
to_SingleCellExperiment( adata, assays_mapping = NULL, colData_mapping = NULL, rowData_mapping = NULL, reduction_mapping = NULL, colPairs_mapping = NULL, rowPairs_mapping = NULL, metadata_mapping = NULL )
adata |
an AnnData object, e.g., InMemoryAnnData |
assays_mapping |
A named list mapping |
colData_mapping |
a named list mapping |
rowData_mapping |
a named list mapping |
reduction_mapping |
a named list mapping reduction names in
|
colPairs_mapping |
a named list mapping obsp names in |
rowPairs_mapping |
a named list mapping varp names in |
metadata_mapping |
a named list mapping uns names in |
to_SingleCellExperiment()
returns a SingleCellExperiment
representing the content of adata
.
if (interactive()) { ## useful when interacting with the SingleCellExperiment ! library(SingleCellExperiment) } ad <- 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) ) ## construct a SingleCellExperiment from an AnnData object sce <- to_SingleCellExperiment(ad) sce
if (interactive()) { ## useful when interacting with the SingleCellExperiment ! library(SingleCellExperiment) } ad <- 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) ) ## construct a SingleCellExperiment from an AnnData object sce <- to_SingleCellExperiment(ad) sce
Write an H5AD file
write_h5ad( object, path, compression = c("none", "gzip", "lzf"), mode = c("w-", "r", "r+", "a", "w", "x") )
write_h5ad( object, path, compression = c("none", "gzip", "lzf"), mode = c("w-", "r", "r+", "a", "w", "x") )
object |
The object to write, either a "SingleCellExperiment" or a "Seurat" object |
path |
Path of the file to write to |
compression |
The compression algorithm to use when writing the
HDF5 file. Can be one of |
mode |
The mode to open the HDF5 file.
|
path
invisibly
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) # Write a SingleCellExperiment as an H5AD if (requireNamespace("SingleCellExperiment", quietly = TRUE)) { ncells <- 100 counts <- matrix(rpois(20000, 5), ncol = ncells) logcounts <- log2(counts + 1) pca <- matrix(runif(ncells * 5), ncells) tsne <- matrix(rnorm(ncells * 2), ncells) sce <- SingleCellExperiment::SingleCellExperiment( assays = list(counts = counts, logcounts = logcounts), reducedDims = list(PCA = pca, tSNE = tsne) ) adata <- from_SingleCellExperiment(sce) h5ad_file <- tempfile(fileext = ".h5ad") adata$write_h5ad(h5ad_file) } # Write a Seurat as a H5AD if (requireNamespace("Seurat", quietly = TRUE)) { library(Seurat) counts <- matrix(1:15, 5L, 3L) dimnames(counts) <- list( LETTERS[1:5], letters[1:3] ) cell.metadata <- data.frame( row.names = letters[1:3], cell = 1:3 ) obj <- CreateSeuratObject(counts, meta.data = cell.metadata) gene.metadata <- data.frame( row.names = LETTERS[1:5], gene = 1:5 ) obj[["RNA"]] <- AddMetaData(GetAssay(obj), gene.metadata) adata <- from_Seurat(obj) h5ad_file <- tempfile(fileext = ".h5ad") adata$write_h5ad(h5ad_file) }
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) # Write a SingleCellExperiment as an H5AD if (requireNamespace("SingleCellExperiment", quietly = TRUE)) { ncells <- 100 counts <- matrix(rpois(20000, 5), ncol = ncells) logcounts <- log2(counts + 1) pca <- matrix(runif(ncells * 5), ncells) tsne <- matrix(rnorm(ncells * 2), ncells) sce <- SingleCellExperiment::SingleCellExperiment( assays = list(counts = counts, logcounts = logcounts), reducedDims = list(PCA = pca, tSNE = tsne) ) adata <- from_SingleCellExperiment(sce) h5ad_file <- tempfile(fileext = ".h5ad") adata$write_h5ad(h5ad_file) } # Write a Seurat as a H5AD if (requireNamespace("Seurat", quietly = TRUE)) { library(Seurat) counts <- matrix(1:15, 5L, 3L) dimnames(counts) <- list( LETTERS[1:5], letters[1:3] ) cell.metadata <- data.frame( row.names = letters[1:3], cell = 1:3 ) obj <- CreateSeuratObject(counts, meta.data = cell.metadata) gene.metadata <- data.frame( row.names = LETTERS[1:5], gene = 1:5 ) obj[["RNA"]] <- AddMetaData(GetAssay(obj), gene.metadata) adata <- from_Seurat(obj) h5ad_file <- tempfile(fileext = ".h5ad") adata$write_h5ad(h5ad_file) }