SingleCellExperiment 1.10.1
The SingleCellExperiment
class is a lightweight Bioconductor container for storing and manipulating single-cell genomics data.
It extends the RangedSummarizedExperiment
class and follows similar conventions,
i.e., rows should represent features (genes, transcripts, genomic regions) and columns should represent cells.
It provides methods for storing dimensionality reduction results and data for alternative feature sets (e.g., synthetic spike-in transcripts, antibody-derived tags).
It is the central data structure for Bioconductor single-cell packages like scater and scran.
SingleCellExperiment
objects can be created via the constructor of the same name.
For example, if we have a count matrix in counts
, we can simply call:
library(SingleCellExperiment)
counts <- matrix(rpois(100, lambda = 10), ncol=10, nrow=10)
sce <- SingleCellExperiment(counts)
sce
## class: SingleCellExperiment
## dim: 10 10
## metadata(0):
## assays(1): ''
## rownames: NULL
## rowData names(0):
## colnames: NULL
## colData names(0):
## reducedDimNames(0):
## altExpNames(0):
In practice, it is often more useful to name the assay by passing in a named list:
sce <- SingleCellExperiment(list(counts=counts))
sce
## class: SingleCellExperiment
## dim: 10 10
## metadata(0):
## assays(1): counts
## rownames: NULL
## rowData names(0):
## colnames: NULL
## colData names(0):
## reducedDimNames(0):
## altExpNames(0):
It is similarly easy to set the column and row metadata by passing values to the appropriate arguments. We will not go into much detail here as most of this is covered by the SummarizedExperiment documentation, but to give an example:
pretend.cell.labels <- sample(letters, ncol(counts), replace=TRUE)
pretend.gene.lengths <- sample(10000, nrow(counts))
sce <- SingleCellExperiment(list(counts=counts),
colData=DataFrame(label=pretend.cell.labels),
rowData=DataFrame(length=pretend.gene.lengths),
metadata=list(study="GSE111111")
)
sce
## class: SingleCellExperiment
## dim: 10 10
## metadata(1): study
## assays(1): counts
## rownames: NULL
## rowData names(1): length
## colnames: NULL
## colData names(1): label
## reducedDimNames(0):
## altExpNames(0):
We can also construct a SingleCellExperiment
by coercing an existing (Ranged)SummarizedExperiment
object:
se <- SummarizedExperiment(list(counts=counts))
as(se, "SingleCellExperiment")
## class: SingleCellExperiment
## dim: 10 10
## metadata(0):
## assays(1): counts
## rownames: NULL
## rowData names(0):
## colnames: NULL
## colData names(0):
## reducedDimNames(0):
## altExpNames(0):
The set of operations that can be applied to a RangedSummarizedExperiment
are also applicable to any instance of a SingleCellExperiment
.
This includes access to assay data via assay()
, column metadata with colData()
, and so on.
Again, without going into too much detail:
dim(assay(sce))
## [1] 10 10
colnames(colData(sce))
## [1] "label"
colnames(rowData(sce))
## [1] "length"
To demonstrate the use of the class in the rest of the vignette, we will use the Allen data set from the scRNAseq package.
library(scRNAseq)
sce <- ReprocessedAllenData("tophat_counts")
sce
## class: SingleCellExperiment
## dim: 20816 379
## metadata(2): SuppInfo which_qc
## assays(1): tophat_counts
## rownames(20816): 0610007P14Rik 0610009B22Rik ... Zzef1 Zzz3
## rowData names(0):
## colnames(379): SRR2140028 SRR2140022 ... SRR2139341 SRR2139336
## colData names(22): NREADS NALIGNED ... Animal.ID passes_qc_checks_s
## reducedDimNames(0):
## altExpNames(1): ERCC
We compute log-transformed normalized expression values from the count matrix.
(We note that many of these steps can be performed as one-liners from the scater package,
but we will show them here in full to demonstrate the capabilities of the SingleCellExperiment
class.)
counts <- assay(sce, "tophat_counts")
libsizes <- colSums(counts)
size.factors <- libsizes/mean(libsizes)
logcounts(sce) <- log2(t(t(counts)/size.factors) + 1)
assayNames(sce)
## [1] "tophat_counts" "logcounts"
We obtain the PCA and t-SNE representations of the data and add them to the object with the reducedDims()<-
method.
Alternatively, we can representations one at a time with the reducedDim()<-
method (note the missing s
).
pca_data <- prcomp(t(logcounts(sce)), rank=50)
library(Rtsne)
set.seed(5252)
tsne_data <- Rtsne(pca_data$x[,1:50], pca = FALSE)
reducedDims(sce) <- list(PCA=pca_data$x, TSNE=tsne_data$Y)
sce
## class: SingleCellExperiment
## dim: 20816 379
## metadata(2): SuppInfo which_qc
## assays(2): tophat_counts logcounts
## rownames(20816): 0610007P14Rik 0610009B22Rik ... Zzef1 Zzz3
## rowData names(0):
## colnames(379): SRR2140028 SRR2140022 ... SRR2139341 SRR2139336
## colData names(22): NREADS NALIGNED ... Animal.ID passes_qc_checks_s
## reducedDimNames(2): PCA TSNE
## altExpNames(1): ERCC
The coordinates for all representations can be retrieved from a SingleCellExperiment
en masse with reducedDims()
or one at a time by name/index with reducedDim()
.
Each row of the coordinate matrix is assumed to correspond to a cell while each column represents a dimension.
reducedDims(sce)
## List of length 2
## names(2): PCA TSNE
reducedDimNames(sce)
## [1] "PCA" "TSNE"
head(reducedDim(sce, "PCA")[,1:2])
## PC1 PC2
## SRR2140028 -70.23670 33.78880
## SRR2140022 -41.00366 49.76633
## SRR2140055 -133.70763 7.68870
## SRR2140083 -36.26099 113.18806
## SRR2139991 -100.86439 15.94740
## SRR2140067 -97.71210 35.77881
head(reducedDim(sce, "TSNE")[,1:2])
## [,1] [,2]
## SRR2140028 -5.089887 -0.9326747
## SRR2140022 -5.021462 0.9825430
## SRR2140055 7.877810 -1.7014442
## SRR2140083 -3.446139 8.4331700
## SRR2139991 -4.332190 2.0750396
## SRR2140067 -5.229462 -3.1371659
Any subsetting by column of sce_sub
will also lead to subsetting of the dimensionality reduction results by cell.
This is convenient as it ensures our low-dimensional results are always synchronized with the gene expression data.
dim(reducedDim(sce, "PCA"))
## [1] 379 50
dim(reducedDim(sce[,1:10], "PCA"))
## [1] 10 50
In the SingleCellExperiment
, users can assign arbitrary names to entries of assays
.
To assist interoperability between packages, we provide some suggestions for what the names should be for particular types of data:
counts
: Raw count data, e.g., number of reads or transcripts for a particular gene.normcounts
: Normalized values on the same scale as the original counts.
For example, counts divided by cell-specific size factors that are centred at unity.logcounts
: Log-transformed counts or count-like values.
In most cases, this will be defined as log-transformed normcounts
, e.g., using log base 2 and a pseudo-count of 1.cpm
: Counts-per-million.
This is the read count for each gene in each cell, divided by the library size of each cell in millions.tpm
: Transcripts-per-million.
This is the number of transcripts for each gene in each cell, divided by the total number of transcripts in that cell (in millions).Each of these suggested names has an appropriate getter/setter method for convenient manipulation of the SingleCellExperiment
.
For example, we can take the (very specifically named) tophat_counts
name and assign it to counts
instead:
counts(sce) <- assay(sce, "tophat_counts")
sce
## class: SingleCellExperiment
## dim: 20816 379
## metadata(2): SuppInfo which_qc
## assays(3): tophat_counts logcounts counts
## rownames(20816): 0610007P14Rik 0610009B22Rik ... Zzef1 Zzz3
## rowData names(0):
## colnames(379): SRR2140028 SRR2140022 ... SRR2139341 SRR2139336
## colData names(22): NREADS NALIGNED ... Animal.ID passes_qc_checks_s
## reducedDimNames(2): PCA TSNE
## altExpNames(1): ERCC
dim(counts(sce))
## [1] 20816 379
This means that functions expecting count data can simply call counts()
without worrying about package-specific naming conventions.
Many scRNA-seq experiments contain sequencing data for multiple feature types beyond the endogenous genes:
Such features can be stored inside the SingleCellExperiment
via the concept of “alternative Experiments”.
These are nested SummarizedExperiment
instances that are guaranteed to have the same number and ordering of columns as the main SingleCellExperiment
itself.
Data for endogenous genes and other features can thus be kept separate - which is often desirable as they need to be processed differently - while still retaining the synchronization of operations on a single object.
To illustrate, consider the case of the spike-in transcripts in the Allen data.
The altExp()
method returns a self-contained SingleCellExperiment
instance containing only the spike-in transcripts.
altExp(sce)
## class: SingleCellExperiment
## dim: 92 379
## metadata(2): SuppInfo which_qc
## assays(1): tophat_counts
## rownames(92): ERCC-00002 ERCC-00003 ... ERCC-00170 ERCC-00171
## rowData names(0):
## colnames(379): SRR2140028 SRR2140022 ... SRR2139341 SRR2139336
## colData names(0):
## reducedDimNames(0):
## altExpNames(0):
Each alternative Experiment can have a different set of assays from the main SingleCellExperiment
.
This is useful in cases where the other feature types must be normalized or transformed differently.
Similarly, the alternative Experiments can have different rowData()
from the main object.
rowData(altExp(sce))$concentration <- runif(nrow(altExp(sce)))
rowData(altExp(sce))
## DataFrame with 92 rows and 1 column
## concentration
## <numeric>
## ERCC-00002 0.034970
## ERCC-00003 0.979286
## ERCC-00004 0.171510
## ERCC-00009 0.577606
## ERCC-00012 0.965111
## ... ...
## ERCC-00164 0.2945227
## ERCC-00165 0.4449455
## ERCC-00168 0.2821027
## ERCC-00170 0.0817114
## ERCC-00171 0.4178510
rowData(sce)
## DataFrame with 20816 rows and 0 columns
We provide the splitAltExps()
utility to easily split a SingleCellExperiment
into new alternative Experiments.
For example, if we wanted to split the RIKEN transcripts into a separate Experiment
- say, to ensure that they are not used in downstream analyses without explicitly throwing out the data -
we would do:
is.riken <- grepl("^[0-9]", rownames(sce))
sce <- splitAltExps(sce, ifelse(is.riken, "RIKEN", "gene"))
altExpNames(sce)
## [1] "ERCC" "RIKEN"
Alternatively, if we want to swap the main and alternative Experiments -
perhaps because the RIKEN transcripts were more interesting than expected, and we want to perform our analyses on them -
we can simply use swapAltExp()
to switch the RIKEN alternative Experiment with the gene expression data:
swapAltExp(sce, "RIKEN", saved="original")
## class: SingleCellExperiment
## dim: 611 379
## metadata(2): SuppInfo which_qc
## assays(3): tophat_counts logcounts counts
## rownames(611): 0610007P14Rik 0610009B22Rik ... 9930111J21Rik1
## 9930111J21Rik2
## rowData names(0):
## colnames(379): SRR2140028 SRR2140022 ... SRR2139341 SRR2139336
## colData names(22): NREADS NALIGNED ... Animal.ID passes_qc_checks_s
## reducedDimNames(2): PCA TSNE
## altExpNames(3): ERCC RIKEN original
The SingleCellExperiment
class provides the sizeFactors()
getter and setter methods,
to set and retrieve size factors from the colData
of the object.
Each size factor represents the scaling factor applied to a cell to normalize expression values prior to downstream comparisons,
e.g., to remove the effects of differences in library size and other cell-specific biases.
These methods are primarily intended for programmatic use in functions implementing normalization methods,
but users can also directly call this to inspect or define the size factors for their analysis.
# Making up some size factors and storing them:
sizeFactors(sce) <- 2^rnorm(ncol(sce))
summary(sizeFactors(sce))
## Min. 1st Qu. Median Mean 3rd Qu. Max.
## 0.09687 0.68484 1.05047 1.23057 1.49024 5.86064
# Deleting the size factors:
sizeFactors(sce) <- NULL
sizeFactors(sce)
## NULL
The colLabels()
getter and setters methods allow applications to set and retrieve cell labels from the colData
.
These labels can be derived from cluster annotations, assigned by classification algorithms, etc.
and are often used in downstream visualization and analyses.
While labels can be stored in any colData
field, the colLabels()
methods aim to provide some informal standardization
to a default location that downstream functions can search first when attempting to retrieve annotations.
# Making up some labels and storing them:
colLabels(sce) <- sample(letters, ncol(sce), replace=TRUE)
table(colLabels(sce))
##
## a b c d e f g h i j k l m n o p q r s t u v w x y z
## 18 11 9 10 23 16 19 16 19 16 8 18 14 15 16 12 17 18 8 18 14 17 14 14 14 5
# Deleting the labels:
colLabels(sce) <- NULL
colLabels(sce)
## NULL
sessionInfo()
## R version 4.0.0 (2020-04-24)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.11-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.11-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] parallel stats4 stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] Rtsne_0.15 scRNAseq_2.1.7
## [3] SingleCellExperiment_1.10.1 SummarizedExperiment_1.18.0
## [5] DelayedArray_0.14.0 matrixStats_0.56.0
## [7] Biobase_2.48.0 GenomicRanges_1.40.0
## [9] GenomeInfoDb_1.24.0 IRanges_2.22.1
## [11] S4Vectors_0.26.0 BiocGenerics_0.34.0
## [13] BiocStyle_2.16.0
##
## loaded via a namespace (and not attached):
## [1] Rcpp_1.0.4.6 lattice_0.20-41
## [3] assertthat_0.2.1 digest_0.6.25
## [5] mime_0.9 BiocFileCache_1.12.0
## [7] R6_2.4.1 RSQLite_2.2.0
## [9] evaluate_0.14 httr_1.4.1
## [11] pillar_1.4.3 zlibbioc_1.34.0
## [13] rlang_0.4.5 curl_4.3
## [15] blob_1.2.1 Matrix_1.2-18
## [17] rmarkdown_2.1 AnnotationHub_2.20.0
## [19] stringr_1.4.0 RCurl_1.98-1.2
## [21] bit_1.1-15.2 shiny_1.4.0.2
## [23] compiler_4.0.0 httpuv_1.5.2
## [25] xfun_0.13 pkgconfig_2.0.3
## [27] htmltools_0.4.0 tidyselect_1.0.0
## [29] tibble_3.0.1 GenomeInfoDbData_1.2.3
## [31] interactiveDisplayBase_1.26.0 bookdown_0.18
## [33] later_1.0.0 crayon_1.3.4
## [35] dplyr_0.8.5 dbplyr_1.4.3
## [37] bitops_1.0-6 rappdirs_0.3.1
## [39] grid_4.0.0 xtable_1.8-4
## [41] lifecycle_0.2.0 DBI_1.1.0
## [43] magrittr_1.5 stringi_1.4.6
## [45] XVector_0.28.0 promises_1.1.0
## [47] ellipsis_0.3.0 vctrs_0.2.4
## [49] tools_4.0.0 bit64_0.9-7
## [51] glue_1.4.0 BiocVersion_3.11.1
## [53] purrr_0.3.4 fastmap_1.0.1
## [55] yaml_2.2.1 AnnotationDbi_1.50.0
## [57] BiocManager_1.30.10 ExperimentHub_1.14.0
## [59] memoise_1.1.0 knitr_1.28