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 [aut] (<https://orcid.org/0000-0001-5949-2191>, bschilder), Chananchida Sang-aram [ctb] (<https://orcid.org/0000-0002-0922-0822>, csangara) |
Maintainer: | Robrecht Cannoodt <[email protected]> |
License: | MIT + file LICENSE |
Version: | 0.99.0 |
Built: | 2024-09-14 05:54:46 UTC |
Source: | https://github.com/scverse/anndataR |
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 = shape )
AnnData( X = NULL, obs = NULL, var = NULL, layers = NULL, obsm = NULL, varm = NULL, obsp = NULL, varp = NULL, uns = NULL, shape = shape )
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.
from_Seurat( seurat_obj, output_class = c("InMemoryAnnData", "HDF5AnnData"), assay = NULL, X = "counts", ... )
from_Seurat( seurat_obj, output_class = c("InMemoryAnnData", "HDF5AnnData"), assay = NULL, X = "counts", ... )
seurat_obj |
An object inheriting from Seurat. |
output_class |
Name of the AnnData class. Must be one of |
assay |
Assay to be converted. If NULL, |
X |
Which of 'counts', 'data', or 'scale.data' will be used for X. By default, 'counts' will be used (if it is not empty), followed by 'data', then 'scale.data'. The remaining non-empty slots will be stored in different layers. |
... |
Additional arguments passed to the generator function. |
For more information on the functionality of an AnnData object, see anndataR-package.
from_SingleCellExperiment()
converts a
SingleCellExperiment to an AnnData object.
from_SingleCellExperiment( sce, output_class = c("InMemory", "HDF5AnnData"), ... )
from_SingleCellExperiment( sce, output_class = c("InMemory", "HDF5AnnData"), ... )
sce |
An object inheriting from SingleCellExperiment. |
output_class |
Name of the AnnData class. Must be one of |
... |
Additional arguments passed to the generator function. See the "Details" section for more information on which parameters |
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") }
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") write_h5ad(adata, 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) ) h5ad_file <- tempfile(fileext = ".h5ad") write_h5ad(sce, h5ad_file) } # Write a Seurat as a H5AD if (requireNamespace("SeuratObject", quietly = TRUE)) { # TODO: uncomment this code when the seurat converter is fixed # counts <- matrix(1:15, 3L, 5L) # dimnames(counts) <- list( # letters[1:3], # LETTERS[1:5] # ) # gene.metadata <- data.frame( # row.names = LETTERS[1:5], # gene = 1:5 # ) # obj <- SeuratObject::CreateSeuratObject(counts, meta.data = gene.metadata) # cell.metadata <- data.frame( # row.names = letters[1:3], # cell = 1:3 # ) # obj <- SeuratObject::AddMetaData(obj, cell.metadata) # # h5ad_file <- tempfile(fileext = ".h5ad") # write_h5ad(obj, 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") write_h5ad(adata, 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) ) h5ad_file <- tempfile(fileext = ".h5ad") write_h5ad(sce, h5ad_file) } # Write a Seurat as a H5AD if (requireNamespace("SeuratObject", quietly = TRUE)) { # TODO: uncomment this code when the seurat converter is fixed # counts <- matrix(1:15, 3L, 5L) # dimnames(counts) <- list( # letters[1:3], # LETTERS[1:5] # ) # gene.metadata <- data.frame( # row.names = LETTERS[1:5], # gene = 1:5 # ) # obj <- SeuratObject::CreateSeuratObject(counts, meta.data = gene.metadata) # cell.metadata <- data.frame( # row.names = letters[1:3], # cell = 1:3 # ) # obj <- SeuratObject::AddMetaData(obj, cell.metadata) # # h5ad_file <- tempfile(fileext = ".h5ad") # write_h5ad(obj, h5ad_file) }