Using scran to perform basic analyses of single-cell RNA-seq data

Package: scran
Author: Aaron Lun (alun@wehi.edu.au)
Compilation date: 2018-08-07

Introduction

Single-cell RNA sequencing (scRNA-seq) is a widely used technique for profiling gene expression in individual cells. This allows molecular biology to be studied at a resolution that cannot be matched by bulk sequencing of cell populations. Common analyses include detection of highly variable and correlated genes across cells, or assignment of cells to cell cycle phases. Cell-specific biases also need to be normalized in a manner that is robust to low counts and technical noise. The scran package implements methods to perform these analyses. This vignette provides brief descriptions of these methods and some toy examples to demonstrate their use.

Setting up the data

We start off with a count matrix where each row is a gene and each column is a cell. These can be obtained by mapping read sequences to a reference genome, and then counting the number of reads mapped to the exons of each gene. (See, for example, the Rsubread package to do both of these tasks.) Alternatively, pseudo-alignment methods can be used to quantify the abundance of each transcript in each cell. For simplicity, though, we'll just simulate some counts here from a negative binomial distribution.

ngenes <- 10000
ncells <- 200
mu <- 2^runif(ngenes, -1, 5)
gene.counts <- matrix(rnbinom(ngenes*ncells, mu=mu, size=10), nrow=ngenes)

We add some arbitrary Ensembl gene IDs to give the impression that this is real (mouse) data.

library(org.Mm.eg.db)
all.ensembl <- unique(toTable(org.Mm.egENSEMBL)$ensembl_id)
rownames(gene.counts) <- sample(all.ensembl, ngenes)

We also have a set of counts for spike-in transcripts. These are appended to the counts for the endogenous genes. In practice, the reads should have been mapped to the spike-in transcipts by including the spike-in sequences in the genome index.

nspikes <- 100
ncells <- 200
mu <- 2^runif(nspikes, -1, 5)
spike.counts <- matrix(rnbinom(nspikes*ncells, mu=mu, size=10), nrow=nspikes)
rownames(spike.counts) <- paste0("ERCC-", seq_len(nspikes))
all.counts <- rbind(gene.counts, spike.counts)

Finally, we construct a SingleCellExperiment object to store all of the data. We also indicate which rows correspond to spike-in transcripts. This is done through the calculateQCMetrics method from scater, which takes a named list of sets of control genes. We indicate which sets of controls are spike-ins using the setSpike setter function. (In this case, there is only one control set, so the process may seem more complicated than necessary. The usefulness of this setup becomes more obvious when multiple control sets are present.) This information can be easily extracted later on using the isSpike, spikes and whichSpike methods.

library(scran)
sce <- SingleCellExperiment(list(counts=all.counts))
isSpike(sce, "MySpike") <- grep("^ERCC", rownames(sce))

This is simulated data, so we assume that quality control has already been applied to remove low-quality cells or low-abundance genes. Check out the scater and cellity packages for more details. Also see https://www.bioconductor.org/help/workflows/simpleSingleCell/ for a workflow where all these steps are used in real data analyses.

Cell cycle phase assignment

We use a pre-defined classifier to assign cells into their cell cycle phases. This classifier was constructed from a training data set by identifying pairs of genes where the difference in expression within each pair changed sign across phases. Thus, by examining the sign of the difference in test data, the phase to which the cell belongs can be identified. Classifiers for human and mouse data are provided with the package – for other systems, classifiers can be constructed from a training set using the sandbag function.

mm.pairs <- readRDS(system.file("exdata", "mouse_cycle_markers.rds", package="scran"))

The classification itself is done using the cyclone function, given the count data and the trained classifier. This yields a number of scores representing the consistency of the signs with each phase.

assigned <- cyclone(sce, pairs=mm.pairs)
head(assigned$scores)
##      G1     S   G2M
## 1 0.755 0.184 0.400
## 2 0.495 0.348 0.496
## 3 0.441 0.399 0.524
## 4 0.975 0.393 0.030
## 5 0.681 0.342 0.282
## 6 0.476 0.360 0.265

Cells are considered to be in G1 phase, if the G1 score is above 0.5 and the G2/M score is below 0.5; to be in G2/M phase, if the G2/M score is above 0.5 and the G1 score is below 0.5; to be in S phase, if both scores are below 0.5; and to be unknown, if both scores are above 0.5. Despite the availability of a S score, it tends to be more accurate to assign cells based on the G1 and G2/M scores only.

table(assigned$phases)
## 
##  G1 G2M   S 
## 151  20  29

Note that it is generally best practice to perform cell cycle phase assignment before filtering out low-abundance genes.

Normalizing cell-specific biases

Based on the gene counts

Cell-specific biases are normalized using the computeSumFactors method, which implements the deconvolution strategy for scRNA-seq normalization. This computes size factors that are used to scale the counts in each cell. The assumption is that most genes are not differentially expressed (DE) between cells, such that any differences in expression across the majority of genes represents some technical bias that should be removed.

sce <- computeSumFactors(sce)
summary(sizeFactors(sce))
##    Min. 1st Qu.  Median    Mean 3rd Qu.    Max. 
##  0.9808  0.9951  1.0000  1.0000  1.0045  1.0217

For larger data sets, clustering should be performed with the quickCluster function before normalization. Briefly, cells are grouped into clusters of similar expression; normalization is applied within each cluster to compute size factors for each cell; and the factors are rescaled by normalization between clusters. This reduces the risk of violating the above assumption when many genes are DE between clusters in a heterogeneous population.

larger.sce <- SingleCellExperiment(list(counts=cbind(all.counts, all.counts, all.counts)))
clusters <- quickCluster(larger.sce)
larger.sce <- computeSumFactors(larger.sce, cluster=clusters)

Note that computeSumFactors will automatically remove low-abundance genes, which provides some protection against zero or negative size factor estimates. We also assume that quality control on the cells has already been performed, as low-quality cells with few expressed genes can often have negative size factor estimates.

Based on the spike-in counts

An alternative approach is to normalize based on the spike-in counts. The idea is that the same quantity of spike-in RNA was added to each cell prior to library preparation. Size factors are computed to scale the counts such that the total coverage of the spike-in transcripts is equal across cells. The main practical difference is that spike-in normalization preserves differences in total RNA content between cells, whereas computeSumFactors and other non-DE methods do not.

sce2 <- computeSpikeFactors(sce)
summary(sizeFactors(sce2))
##    Min. 1st Qu.  Median    Mean 3rd Qu.    Max. 
##  0.8240  0.9605  1.0003  1.0000  1.0367  1.1914

Even if we decide to use the deconvolution size factors, it is strongly recommended to compute a separate set of size factors for the spike-ins. This is because the spike-ins are not affected by total mRNA content. Using the deconvolution size factors will over-normalize the spike-in counts, whereas the spike-in size factors are more appropriate. To obtain the latter without overwriting the former, we set general.use=FALSE in our call to computeSpikeFactors. This means that the spike-in-based size factors will be computed and stored in the SingleCellExperiment object, but will only be used by the spike-in transcripts. (Obviously, if the spike-in size factors were already being used for normalization, e.g., in sce2, then this extra step is unnecessary.)

sce <- computeSpikeFactors(sce, general.use=FALSE)

Computing normalized expression values

Normalized expression values are calculated using the normalize method from scater. This will use the deconvolution size factors for the endogenous genes, and the spike-in-based size factors for the spike-in transcripts. Each expression value can be interpreted as a log-transformed “normalized count”, and can be used in downstream applications like clustering or dimensionality reduction.

sce <- normalize(sce)

Detecting highly variable genes

Highly variable genes (HVGs) are detected by decomposing the total variance of each gene into its biological and technical components. This avoids prioritizing low-abundance genes that have large variances due to technical noise. First, we fit a mean-variance trend to the normalized log-expression values with trendVar. By default, this done using only the spike-in transcripts, as these should only exhibit technical noise.

fit <- trendVar(sce, parametric=TRUE)

The fitted value of the trend is used as an estimate of the technical component. We subtract the fitted value from the total variance to obtain the biological component for each gene. HVGs are defined as the top set of genes with the largest biological components.

decomp <- decomposeVar(sce, fit)
top.hvgs <- order(decomp$bio, decreasing=TRUE)
head(decomp[top.hvgs,])
## DataFrame with 6 rows and 6 columns
##                                mean             total               bio
##                           <numeric>         <numeric>         <numeric>
## ENSMUSG00000098831  2.6201162964783 0.759087170538644 0.253133889889331
## ENSMUSG00000097172  2.5268038098624 0.734292658610173 0.215273273077961
## ENSMUSG00000025532 2.03838187826058 0.807680468225129 0.209202594963209
## ENSMUSG00000041959 1.86646293629607  0.81369037279531 0.198495550432575
## ENSMUSG00000034820 2.58461437246598  0.70290243979077 0.192211172264215
## ENSMUSG00000032356 2.64685619493326 0.691532346737567 0.188935841423244
##                                 tech              p.value
##                            <numeric>            <numeric>
## ENSMUSG00000098831 0.505953280649313 6.15797009727912e-06
## ENSMUSG00000097172 0.519019385532213 0.000105426336938057
## ENSMUSG00000025532  0.59847787326192 0.000736414387447314
## ENSMUSG00000041959 0.615194822362734  0.00154663982133305
## ENSMUSG00000034820 0.510691267526556 0.000339272947467935
## ENSMUSG00000032356 0.502596505314323  0.00034386234724774
##                                   FDR
##                             <numeric>
## ENSMUSG00000098831 0.0307898504863956
## ENSMUSG00000097172 0.0969784950318998
## ENSMUSG00000025532  0.222265246185973
## ENSMUSG00000041959  0.291818834213784
## ENSMUSG00000034820  0.163743974879876
## ENSMUSG00000032356  0.163743974879876

We examined this in more detail by constructing a mean-variance plot. Here, the black points represent the endogenous genes; the red points represent spike-in transcripts; and the red line represents the mean-variance trend fitted to the spike-ins.

plot(decomp$mean, decomp$total, xlab="Mean log-expression", ylab="Variance")
o <- order(decomp$mean)
lines(decomp$mean[o], decomp$tech[o], col="red", lwd=2)
points(fit$mean, fit$var, col="red", pch=16)

If spike-ins are absent or of poor quality, an alternative is to fit the trend to the gene variances directly with use.spikes=FALSE. This assumes that technical noise is the major contributor to the variance of most genes in the data set, such that the trend still represents the technical component. The resulting fit can then be used in decomposeVar as described above.

alt.fit <- trendVar(sce, use.spikes=FALSE) 
alt.decomp <- decomposeVar(sce, alt.fit)

If the data set already contains some uninteresting substructure (e.g., batch effects), we can block on this by setting the design argument in trendVar. This ensures that the substructure does not inflate the variance estimates. For example, if the cells were prepared in two separate batches, we could construct a design matrix incorporating this information with model.matrix and pass it to trendVar. The same design will also be used in decomposeVar.

batch <- rep(c("1", "2"), each=100)
design <- model.matrix(~batch)
alt.fit2 <- trendVar(sce, design=design)
alt.decomp2 <- decomposeVar(sce, alt.fit)

Alternative approaches to identifying HVGs are implemented in the DM and technicalCV2 functions. These are based on the coefficient of variation for count data, which provides more power for rare subpopulations but is also more sensitive to outliers. Also see the improvedCV2 function, which is intended as a more stable counterpart of technicalCV2.

Detecting correlated genes

Another useful procedure is to identify significant pairwise correlations between pairs of HVGs. The idea is to distinguish between HVGs caused by random stochasticity, and those that are driving systematic heterogeneity, e.g., between subpopulations. Correlations are computed in the correlatePairs method using a slightly modified version of Spearman's rho. Testing is performed against the null hypothesis of independent genes, using a permutation method in correlateNull to construct a null distribution.

null.dist <- correlateNull(ncol(sce))
# Only using the first 200 genes as a demonstration.
cor.pairs <- correlatePairs(sce, subset.row=top.hvgs[1:200], null.dist=null.dist)
head(cor.pairs)
## DataFrame with 6 rows and 6 columns
##                gene1              gene2                rho
##          <character>        <character>          <numeric>
## 1 ENSMUSG00000116408 ENSMUSG00000047230 -0.263427585689642
## 2 ENSMUSG00000041360 ENSMUSG00000042817  0.261203030075752
## 3 ENSMUSG00000066279 ENSMUSG00000086552 -0.258202955073877
## 4 ENSMUSG00000030793 ENSMUSG00000074665  0.255447386184655
## 5 ENSMUSG00000062061 ENSMUSG00000028551 -0.254217355433886
## 6 ENSMUSG00000028015 ENSMUSG00000019718 -0.246612165304133
##                p.value               FDR   limited
##              <numeric>         <numeric> <logical>
## 1 0.000181999818000182 0.940701499637991     FALSE
## 2 0.000181999818000182 0.940701499637991     FALSE
## 3 0.000253999746000254 0.940701499637991     FALSE
## 4  0.00027999972000028 0.940701499637991     FALSE
## 5 0.000311999688000312 0.940701499637991     FALSE
## 6 0.000435999564000436 0.940701499637991     FALSE

As with variance estimation, if uninteresting substructure is present, this should be blocked on using the design argument in both correlateNull and correlatePairs.

null.dist2 <- correlateNull(design=design, iter=1e5) # fewer iterations, to speed it up.
cor.pairs2 <- correlatePairs(sce, subset.row=top.hvgs[1:200], 
    null.dist=null.dist2, design=design)

The pairs can be used for choosing marker genes in experimental validation, and to construct gene-gene association networks. In other situations, the pairs may not be of direct interest - rather, we just want to know whether a gene is correlated with any other gene. This is often the case if we are to select a set of correlated HVGs for use in downstream steps like clustering or dimensionality reduction. To do so, we set per.gene=TRUE to compute a single set of statistics for each gene, rather than for each pair.

cor.genes <- correlatePairs(sce, subset.row=top.hvgs[1:200], 
    null.dist=null.dist, per.gene=TRUE)

Significant correlations are defined at a false discovery rate (FDR) threshold of, e.g., 5%. Note that the p-values are calculated by permutation and will have a lower bound. If there were insufficient permutation iterations, a warning will be issued suggesting that more iterations be performed.

Converting to other formats

The SingleCellExperiment object can be easily converted into other formats using the convertTo method. This allows analyses to be performed using other pipelines and packages. For example, if DE analyses were to be performed using edgeR, the count data in sce could be used to construct a DGEList.

y <- convertTo(sce, type="edgeR")

By default, rows corresponding to spike-in transcripts are dropped when get.spikes=FALSE. As such, the rows of y may not correspond directly to the rows of sce – users should match by row name to ensure correct cross-referencing between objects. Normalization factors are also automatically computed from the size factors.

The same conversion strategy roughly applies to the other supported formats. DE analyses can be performed using DESeq2 by converting the object to a DESeqDataSet. Cells can be ordered on pseudotime with monocle by converting the object to a CellDataSet (in this case, normalized unlogged expression values are stored).

Graphical interaction with the data

scran provides several functions for exploring scRNA-seq data:

These functions will generate and load Shiny apps for interactive data exploration.

Summary

This vignette describes the main functions in the scran package for basic analysis of single-cell RNA-seq data. We cover normalization, cell cycle phase assignment, HVG detection and correlation testing. Conversion to other formats can also be performed in preparation for analyses with other packages in the Bioconductor project. Further information can be obtained by examining the documentation for each function (e.g., ?convertTo); reading the workflow at https://www.bioconductor.org/help/workflows/simpleSingleCell; or asking for help on the Bioconductor support site (please read the posting guide beforehand).

Session information

sessionInfo()
## R version 3.5.1 Patched (2018-07-12 r74967)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.5 LTS
## 
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.7-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.7-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] org.Mm.eg.db_3.6.0          AnnotationDbi_1.42.1       
##  [3] scran_1.8.4                 SingleCellExperiment_1.2.0 
##  [5] SummarizedExperiment_1.10.1 DelayedArray_0.6.4         
##  [7] matrixStats_0.54.0          Biobase_2.40.0             
##  [9] GenomicRanges_1.32.6        GenomeInfoDb_1.16.0        
## [11] IRanges_2.14.10             S4Vectors_0.18.3           
## [13] BiocGenerics_0.26.0         BiocParallel_1.14.2        
## [15] BiocStyle_2.8.2             knitr_1.20                 
## 
## loaded via a namespace (and not attached):
##  [1] bitops_1.0-6             bit64_0.9-7             
##  [3] rprojroot_1.3-2          dynamicTreeCut_1.63-1   
##  [5] tools_3.5.1              backports_1.1.2         
##  [7] R6_2.2.2                 DT_0.4                  
##  [9] vipor_0.4.5              DBI_1.0.0               
## [11] lazyeval_0.2.1           colorspace_1.3-2        
## [13] tidyselect_0.2.4         gridExtra_2.3           
## [15] bit_1.1-14               compiler_3.5.1          
## [17] scales_0.5.0             stringr_1.3.1           
## [19] digest_0.6.15            rmarkdown_1.10          
## [21] XVector_0.20.0           scater_1.8.3            
## [23] pkgconfig_2.0.1          htmltools_0.3.6         
## [25] limma_3.36.2             htmlwidgets_1.2         
## [27] rlang_0.2.1              RSQLite_2.1.1           
## [29] FNN_1.1.2                shiny_1.1.0             
## [31] DelayedMatrixStats_1.2.0 bindr_0.1.1             
## [33] dplyr_0.7.6              RCurl_1.95-4.11         
## [35] magrittr_1.5             GenomeInfoDbData_1.1.0  
## [37] Matrix_1.2-14            Rcpp_0.12.18            
## [39] ggbeeswarm_0.6.0         munsell_0.5.0           
## [41] Rhdf5lib_1.2.1           viridis_0.5.1           
## [43] stringi_1.2.4            yaml_2.2.0              
## [45] edgeR_3.22.3             zlibbioc_1.26.0         
## [47] rhdf5_2.24.0             plyr_1.8.4              
## [49] grid_3.5.1               blob_1.1.1              
## [51] promises_1.0.1           shinydashboard_0.7.0    
## [53] crayon_1.3.4             lattice_0.20-35         
## [55] locfit_1.5-9.1           pillar_1.3.0            
## [57] igraph_1.2.2             rjson_0.2.20            
## [59] reshape2_1.4.3           glue_1.3.0              
## [61] evaluate_0.11            data.table_1.11.4       
## [63] httpuv_1.4.5             gtable_0.2.0            
## [65] purrr_0.2.5              assertthat_0.2.0        
## [67] ggplot2_3.0.0            mime_0.5                
## [69] xtable_1.8-2             later_0.7.3             
## [71] viridisLite_0.3.0        tibble_1.4.2            
## [73] memoise_1.1.0            beeswarm_0.2.3          
## [75] tximport_1.8.0           bindrcpp_0.2.2          
## [77] statmod_1.4.30