Skip to contents

from_SingleCellExperiment() converts a SingleCellExperiment to an AnnData object.

Usage

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,
  ...
)

Arguments

sce

An object inheriting from SingleCellExperiment.

output_class

Name of the AnnData class. Must be one of "HDF5AnnData" or "InMemoryAnnData".

x_mapping

Name of the assay in sce to use as the X matrix in the AnnData object.

layers_mapping

A named list mapping assay names in sce to layers in the created AnnData object. The names of the list should be the names of the layers in the resulting AnnData object, and the values should be the names of the assays in the sce object.

obs_mapping

A named list mapping colData in sce to obs in the created AnnData object. The names of the list should be the names of the obs columns in the resulting AnnData object. The values of the list should be the names of the colData columns in sce.

var_mapping

A named list mapping rowData in sce to var in the created AnnData object. The names of the list should be the names of the var columns in the resulting AnnData object. The values of the list should be the names of the rowData columns in sce.

obsm_mapping

A named list mapping reducedDim in sce to obsm in the created AnnData object. The names of the list should be the names of the obsm in the resulting AnnData object. The values of the list should be a named list with as key the name of the obsm slot in the resulting AnnData object, and as value a list with the following elements

  • reducedDim

  • the name of the reducedDim in sce

varm_mapping

A named list mapping reducedDim in sce to varm in the created AnnData object. The names of the list should be the names of the varm in the resulting AnnData object. The values of the list should be a named list with as key the name of the varm slot in the resulting AnnData object, and as value a list with the following elements

  • reducedDim

  • the name of the reducedDim in sce, that is LinearEmbeddingMatrix of which you want the featureLoadings to end up in the varm slot

obsp_mapping

A named list mapping colPairs in sce to obsp in the created AnnData object. The names of the list should be the names of the obsp in the resulting AnnData object. The values of the list should be the names of the colPairs in sce.

varp_mapping

A named list mapping rowPairs in sce to varp in the created AnnData object. The names of the list should be the names of the varp in the resulting AnnData object. The values of the list should be the names of the rowPairs in sce.

uns_mapping

A named list mapping metadata in sce to uns in the created AnnData object. The names of the list should be the names of the uns in the resulting AnnData object. The values of the list should be the names of the metadata in sce.

...

Additional arguments to pass to the generator function.

Value

from_SingleCellExperiment() returns an AnnData object (e.g., InMemoryAnnData) representing the content of sce.

Examples

## construct an AnnData object from a SingleCellExperiment
library(SingleCellExperiment)
#> Loading required package: SummarizedExperiment
#> Loading required package: MatrixGenerics
#> Loading required package: matrixStats
#> 
#> Attaching package: ‘MatrixGenerics’
#> The following objects are masked from ‘package:matrixStats’:
#> 
#>     colAlls, colAnyNAs, colAnys, colAvgsPerRowSet, colCollapse,
#>     colCounts, colCummaxs, colCummins, colCumprods, colCumsums,
#>     colDiffs, colIQRDiffs, colIQRs, colLogSumExps, colMadDiffs,
#>     colMads, colMaxs, colMeans2, colMedians, colMins, colOrderStats,
#>     colProds, colQuantiles, colRanges, colRanks, colSdDiffs, colSds,
#>     colSums2, colTabulates, colVarDiffs, colVars, colWeightedMads,
#>     colWeightedMeans, colWeightedMedians, colWeightedSds,
#>     colWeightedVars, rowAlls, rowAnyNAs, rowAnys, rowAvgsPerColSet,
#>     rowCollapse, rowCounts, rowCummaxs, rowCummins, rowCumprods,
#>     rowCumsums, rowDiffs, rowIQRDiffs, rowIQRs, rowLogSumExps,
#>     rowMadDiffs, rowMads, rowMaxs, rowMeans2, rowMedians, rowMins,
#>     rowOrderStats, rowProds, rowQuantiles, rowRanges, rowRanks,
#>     rowSdDiffs, rowSds, rowSums2, rowTabulates, rowVarDiffs, rowVars,
#>     rowWeightedMads, rowWeightedMeans, rowWeightedMedians,
#>     rowWeightedSds, rowWeightedVars
#> Loading required package: GenomicRanges
#> Loading required package: stats4
#> Loading required package: BiocGenerics
#> 
#> Attaching package: ‘BiocGenerics’
#> The following object is masked from ‘package:SeuratObject’:
#> 
#>     intersect
#> The following objects are masked from ‘package:stats’:
#> 
#>     IQR, mad, sd, var, xtabs
#> The following objects are masked from ‘package:base’:
#> 
#>     Filter, Find, Map, Position, Reduce, anyDuplicated, aperm, append,
#>     as.data.frame, basename, cbind, colnames, dirname, do.call,
#>     duplicated, eval, evalq, get, grep, grepl, intersect, is.unsorted,
#>     lapply, mapply, match, mget, order, paste, pmax, pmax.int, pmin,
#>     pmin.int, rank, rbind, rownames, sapply, saveRDS, setdiff, table,
#>     tapply, union, unique, unsplit, which.max, which.min
#> Loading required package: S4Vectors
#> 
#> Attaching package: ‘S4Vectors’
#> The following object is masked from ‘package:utils’:
#> 
#>     findMatches
#> The following objects are masked from ‘package:base’:
#> 
#>     I, expand.grid, unname
#> Loading required package: IRanges
#> 
#> Attaching package: ‘IRanges’
#> The following object is masked from ‘package:sp’:
#> 
#>     %over%
#> Loading required package: GenomeInfoDb
#> Loading required package: Biobase
#> Welcome to Bioconductor
#> 
#>     Vignettes contain introductory material; view with
#>     'browseVignettes()'. To cite Bioconductor, see
#>     'citation("Biobase")', and for packages 'citation("pkgname")'.
#> 
#> Attaching package: ‘Biobase’
#> The following object is masked from ‘package:MatrixGenerics’:
#> 
#>     rowMedians
#> The following objects are masked from ‘package:matrixStats’:
#> 
#>     anyMissing, rowMedians
#> 
#> Attaching package: ‘SummarizedExperiment’
#> The following object is masked from ‘package:Seurat’:
#> 
#>     Assays
#> The following object is masked from ‘package:SeuratObject’:
#> 
#>     Assays
sce <- SingleCellExperiment(
  assays = list(counts = matrix(1:5, 5L, 3L)),
  colData = DataFrame(cell = 1:3),
  rowData = DataFrame(gene = 1:5)
)
from_SingleCellExperiment(sce, "InMemory")
#> AnnData object with n_obs × n_vars = 3 × 5
#>     obs: 'cell'
#>     var: 'gene'
#>     layers: 'counts'