Contents

1 Installation

The recommended way to install the zinbwave package is via Bioconductor.

source("https://bioconductor.org/biocLite.R")
biocLite("zinbwave")

Note that zinbwave requires R (>=3.4) and Bioconductor (>=3.6).

2 Introduction

This vignette provides an introductory example on how to work with the zinbwave package, which implements the ZINB-WaVE method proposed in (Risso et al. 2018).

First, let’s load the packages and set serial computations.

library(zinbwave)
library(scRNAseq)
library(matrixStats)
library(magrittr)
library(ggplot2)
library(biomaRt)

# Register BiocParallel Serial Execution
BiocParallel::register(BiocParallel::SerialParam())

2.1 The ZINB-WaVE model

ZINB-WaVE is a general and flexible model for the analysis of high-dimensional zero-inflated count data, such as those recorded in single-cell RNA-seq assays. Given \(n\) samples (typically, \(n\) single cells) and \(J\) features (typically, \(J\) genes) that can be counted for each sample, we denote with \(Y_{ij}\) the count of feature \(j\) (\(j=1,\ldots,J\)) for sample \(i\) (\(i=1,\ldots,n\)). To account for various technical and biological effects, typical of single-cell sequencing technologies, we model \(Y_{ij}\) as a random variable following a zero-inflated negative binomial (ZINB) distribution with parameters \(\mu_{ij}\), \(\theta_{ij}\), and \(\pi_{ij}\), and consider the following regression models for the parameters:

\[\begin{align} \label{eq:model1} \ln(\mu_{ij}) &= \left( X\beta_\mu + (V\gamma_\mu)^\top + W\alpha_\mu + O_\mu\right)_{ij}\,,\\ \label{eq:model2} \text{logit}(\pi_{ij}) &= \left(X\beta_\pi + (V\gamma_\pi)^\top + W\alpha_\pi + O_\pi\right)_{ij} \,, \\ \label{eq:model3} \ln(\theta_{ij}) &= \zeta_j \,, \end{align}\]

.

where the elements of the regression models are as follows.

  • \(X\) is a known \(n \times M\) matrix corresponding to \(M\) cell-level covariates and \({\bf \beta}=(\beta_\mu,\beta_\pi)\) its associated \(M \times J\) matrices of regression parameters. \(X\) can typically include covariates that induce variation of interest, such as cell types, or covariates that induce unwanted variation, such as batch or quality control (QC) measures. By default, it includes only a constant column of ones, \({\bf 1}_n\), to account for gene-specific intercepts.
  • \(V\) is a known \(J \times L\) matrix corresponding to \(J\) gene-level covariates, such as gene length or GC-content, and \({\bf \gamma} = (\gamma_\mu , \gamma_\pi)\) its associated \(L\times n\) matrices of regression parameters. By default, \(V\) only includes a constant column of ones, \({\bf 1}_J\), to account for cell-specific intercepts, such as size factors representing differences in library sizes.
  • \(W\) is an unobserved \(n \times K\) matrix corresponding to \(K\) unknown cell-level covariates, which could be of “unwanted variation” or of interest (such as cell type), and \({\bf \alpha} = (\alpha_\mu,\alpha_{\pi})\) its associated \(K \times J\) matrices of regression parameters.
  • \(O_\mu\) and \(O_\pi\) are known \(n \times J\) matrices of offsets.
  • \(\zeta\in\mathbb{R}^J\) is a vector of gene-specific dispersion parameters on the log scale.

2.2 Example dataset

To illustrate the methodology, we will make use of the Fluidigm C1 dataset of (Pollen et al. 2014). The data consist of 65 cells, each sequenced at high and low depth. The data are publicly available as part of the scRNAseq package, in the form of a SummarizedExperiment object.

data("fluidigm")
fluidigm
## class: SummarizedExperiment 
## dim: 26255 130 
## metadata(3): sample_info clusters which_qc
## assays(4): tophat_counts cufflinks_fpkm rsem_counts rsem_tpm
## rownames(26255): A1BG A1BG-AS1 ... ZZEF1 ZZZ3
## rowData names(0):
## colnames(130): SRR1275356 SRR1274090 ... SRR1275366 SRR1275261
## colData names(28): NREADS NALIGNED ... Cluster1 Cluster2
table(colData(fluidigm)$Coverage_Type)
## 
## High  Low 
##   65   65

3 Gene filtering

First, we filter out the lowly expressed genes, by removing those genes that do not have at least 5 reads in at least 5 samples.

filter <- rowSums(assay(fluidigm)>5)>5
table(filter)
## filter
## FALSE  TRUE 
## 16127 10128
fluidigm <- fluidigm[filter,]

This leaves us with 10128 genes.

We next identify the 100 most variable genes, which will be the input of our ZINB-WaVE procedure. Although we apply ZINB-WaVE to only these genes primarily for computational reasons, it is generally a good idea to focus on a subset of highly-variable genes, in order to remove transcriptional noise and focus on the more biologically meaningful signals. However, at least 1,000 genes are probably needed for real analyses.

assay(fluidigm) %>% log1p %>% rowVars -> vars
names(vars) <- rownames(fluidigm)
vars <- sort(vars, decreasing = TRUE)
head(vars)
##  IGFBPL1    STMN2     EGR1   ANP32E    CENPF     LDHA 
## 13.06109 12.24748 11.90608 11.67819 10.83797 10.72307
fluidigm <- fluidigm[names(vars)[1:100],]

Before proceeding, we rename the first assay of fluidigm “counts” to avoid needing to specify which assay we should use for the zinbwave workflow. This is an optional step.

assayNames(fluidigm)[1] <- "counts"

4 ZINB-WaVE

The easiest way to obtain the low-dimensional representation of the data with ZINB-WaVE is to use the zinbwave function. This function takes as input a SummarizedExperiment object and returns a SingleCellExperiment object.

fluidigm_zinb <- zinbwave(fluidigm, K = 2, epsilon=1000)

By default, the zinbwave function fits a ZINB model with \(X = {\bf 1}_n\) and \(V = {\bf 1}_J\). In this case, the model is a factor model akin to principal component analysis (PCA), where \(W\) is a factor matrix and \(\alpha_\mu\) and \(\alpha_\pi\) are loading matrices. By default, the epsilon parameter is set to the number of genes. We empirically found that a high epsilon is often required to obtained a good low-level representation. See ?zinbModel for details. Here we set epsilon=1000.

The parameter \(K\) controls how many latent variables we want to infer from the data. \(W\) is stored in the reducedDim slot of the object. (See the SingleCellExperiment vignette for details).

In this case, as we specified \(K=2\), we can visualize the resulting \(W\) matrix in a simple plot, color-coded by cell-type.

W <- reducedDim(fluidigm_zinb)

data.frame(W, bio=colData(fluidigm)$Biological_Condition,
           coverage=colData(fluidigm)$Coverage_Type) %>%
    ggplot(aes(W1, W2, colour=bio, shape=coverage)) + geom_point() + 
    scale_color_brewer(type = "qual", palette = "Set1") + theme_classic()

4.1 Adding covariates

The ZINB-WaVE model is more general than PCA, allowing the inclusion of additional sample and gene-level covariates that might help to infer the unknown factors.

4.1.1 Sample-level covariates

Typically, one could include batch information as sample-level covariate, to account for batch effects. Here, we illustrate this capability by including the coverage (high or low) as a sample-level covariate.

The column Coverage_Type in the colData of fluidigm contains the coverage information. We can specify a design matrix that includes an intercept and an indicator variable for the coverage, by using the formula interface of zinbFit.

fluidigm_cov <- zinbwave(fluidigm, K=2, X="~Coverage_Type", epsilon=1000)
W <- reducedDim(fluidigm_cov)

data.frame(W, bio=colData(fluidigm)$Biological_Condition,
           coverage=colData(fluidigm)$Coverage_Type) %>%
    ggplot(aes(W1, W2, colour=bio, shape=coverage)) + geom_point() + 
    scale_color_brewer(type = "qual", palette = "Set1") + theme_classic()

In this case, the inferred \(W\) matrix is essentially the same with or without covariates, indicating that the scaling factor included in the model (the \(\gamma\) parameters associated with the intercept of \(V\)) are enough to achieve a good low-dimensional representation of the data.

4.1.2 Gene-level covariates

Analogously, we can include gene-level covariates, as columns of \(V\). Here, we illustrate this capability by including gene length and GC-content.

We use the biomaRt package to compute gene length and GC-content.

mart <- useMart("ensembl")
mart <- useDataset("hsapiens_gene_ensembl", mart = mart)
bm <- getBM(attributes=c('hgnc_symbol', 'start_position',
                         'end_position', 'percentage_gene_gc_content'),
            filters = 'hgnc_symbol',
            values = rownames(fluidigm),
            mart = mart)

bm$length <- bm$end_position - bm$start_position
len <- tapply(bm$length, bm$hgnc_symbol, mean)
len <- len[rownames(fluidigm)]
gcc <- tapply(bm$percentage_gene_gc_content, bm$hgnc_symbol, mean)
gcc <- gcc[rownames(fluidigm)]

We then include the gene-level information as rowData in the fluidigm object.

rowData(fluidigm) <- data.frame(gccontent = gcc, length = len)
fluidigm_gcc <- zinbwave(fluidigm, K=2, V="~gccontent + log(length)", epsilon=1000)

5 t-SNE representation

A t-SNE representation of the data can be obtained by computing the cell distances in the reduced space and running the t-SNE algorithm on the distance.

set.seed(93024)

library(Rtsne)
W <- reducedDim(fluidigm_cov)
tsne_data <- Rtsne(W, pca = FALSE, perplexity=10, max_iter=5000)

data.frame(Dim1=tsne_data$Y[,1], Dim2=tsne_data$Y[,2], 
           bio=colData(fluidigm)$Biological_Condition,
           coverage=colData(fluidigm)$Coverage_Type) %>%
    ggplot(aes(Dim1, Dim2, colour=bio, shape=coverage)) + geom_point() + 
    scale_color_brewer(type = "qual", palette = "Set1") + theme_classic()

6 Normalized values and deviance residuals

Sometimes it is useful to have normalized values for visualization and residuals for model evaluation. Both quantities can be computed with the zinbwave() function.

fluidigm_norm <- zinbwave(fluidigm, K=2, epsilon=1000, normalizedValues=TRUE,
                    residuals = TRUE)

The fluidigm_norm object includes normalized values and residuals as additional assays.

fluidigm_norm
## class: SingleCellExperiment 
## dim: 100 130 
## metadata(3): sample_info clusters which_qc
## assays(7): counts cufflinks_fpkm ... residuals weights
## rownames(100): IGFBPL1 STMN2 ... SRSF7 FAM117B
## rowData names(0):
## colnames(130): SRR1275356 SRR1274090 ... SRR1275366 SRR1275261
## colData names(28): NREADS NALIGNED ... Cluster1 Cluster2
## reducedDimNames(1): zinbwave
## spikeNames(0):

7 The zinbFit function

The zinbwave function is a user-friendly function to obtain the low-dimensional representation of the data, and optionally the normalized values and residuals from the model.

However, it is sometimes useful to store all the parameter estimates and the value of the likelihood. The zinbFit function allows the user to create an object of class zinbModel that can be used to store all the parameter estimates and have greater control on the results.

zinb <- zinbFit(fluidigm, K=2, epsilon=1000)

As with zinbwave, by default, the zinbFit function fits a ZINB model with \(X = {\bf 1}_n\) and \(V = {\bf 1}_J\).

If a user has run zinbFit and wants to obtain normalized values or the low-dimensional representation of the data in a SingleCellExperiment format, they can pass the zinbModel object to zinbwave to avoid repeating all the computations.

fluidigm_zinb <- zinbwave(fluidigm, fitted_model = zinb, K = 2, epsilon=1000)

8 Differential Expression

The zinbwave package can be used to compute observational weights to “unlock” bulk RNA-seq tools for single-cell applications, as illustrated in (Van den Berge et al. 2018).

Since version 1.1.5, zinbwave computes the observational weights by default. See the man page of zinbwave. The weights are stored in an assay named weights and can be accessed with the following call.

weights <- assay(fluidigm_zinb, "weights")

Note that in this example, the value of the penalty parameter epsilon was set at 1000, although we do not recommend this for differential expression analysis in real applications. Our evaluations have shown that a value of epsilon=1e12 gives good performance across a range of datasets, although this number is still arbitrary. In general, values between 1e6 and 1e13 give best performances.

8.1 Differential expression with edgeR

Once we have the observational weights, we can use them in edgeR to perform differential expression. Specifically, we use a moderated F-test in which the denominator residual degrees of freedom are adjusted by the extent of zero inflation (see (Van den Berge et al. 2018) for details).

Here, we compare NPC to GW16. Note that we start from only 100 genes for computational reasons, but in real analyses we would use all the expressed genes.

library(edgeR)

dge <- DGEList(assay(fluidigm_zinb))
dge <- calcNormFactors(dge)

design <- model.matrix(~Biological_Condition, data = colData(fluidigm))
dge$weights <- weights
dge <- estimateDisp(dge, design)
fit <- glmFit(dge, design)

lrt <- glmWeightedF(fit, coef = 3)
topTags(lrt)
## Coefficient:  Biological_ConditionGW21+3 
##              logFC   logCPM       LR       PValue   padjFilter          FDR
## VIM      -4.768620 13.21770 47.43920 2.378498e-10 2.378498e-08 2.378498e-08
## FOS      -5.314301 14.50175 37.39214 8.940291e-09 4.470146e-07 4.470146e-07
## USP47    -3.900611 13.37159 29.91900 2.216395e-07 7.387983e-06 7.387983e-06
## PTN      -3.190171 13.22777 22.68395 4.691691e-06 1.172923e-04 1.172923e-04
## MIR100HG  2.388503 14.26684 18.25761 3.515603e-05 6.633045e-04 6.633045e-04
## NNAT     -2.062983 13.60868 17.90557 3.979827e-05 6.633045e-04 6.633045e-04
## SPARC    -3.202981 13.23880 16.08226 1.047911e-04 1.497016e-03 1.497016e-03
## SFRP1    -3.406073 13.01421 14.45986 2.437270e-04 3.046588e-03 3.046588e-03
## EGR1     -2.658631 14.93923 13.58175 3.193159e-04 3.547954e-03 3.547954e-03
## ST8SIA1  -3.338410 13.35889 12.64125 5.336245e-04 5.336245e-03 5.336245e-03

8.2 Differential expression with DESeq2

Analogously, we can use the weights in a DESeq2 analysis by using observation-level weights in the parameter estimation steps. In this case, there is no need to pass the weights to DESeq2 since they are already in the weights assay of the object.

library(DESeq2)

dds <- DESeqDataSet(fluidigm_zinb, design = ~ Biological_Condition)

dds <- DESeq(dds, sfType="poscounts", useT=TRUE, minmu=1e-6)
res <- lfcShrink(dds, contrast=c("Biological_Condition", "NPC", "GW16"),
                 type = "normal")
head(res)
## log2 fold change (MAP): Biological_Condition NPC vs GW16 
## Wald test p-value: Biological Condition NPC vs GW16 
## DataFrame with 6 rows and 6 columns
##                 baseMean    log2FoldChange             lfcSE
##                <numeric>         <numeric>         <numeric>
## IGFBPL1 2054.40279669602 -8.29482247680305 0.705612098941605
## STMN2   2220.07799968962 -10.0742847458669 0.769583053245998
## EGR1    1342.46496311301   -6.853900032812 0.662688086199466
## ANP32E   806.98309206076  1.99090585716592 0.502369874616181
## CENPF   255.632254758972    1.371091032346 0.553771652069953
## LDHA    311.760456182291  2.36714571054835 0.578056211747908
##                      stat               pvalue                 padj
##                 <numeric>            <numeric>            <numeric>
## IGFBPL1 -8.14358136637719 3.83754418626486e-16 2.01976009803414e-15
## STMN2   -9.37586868275918 6.86084508352232e-21 5.27757314117102e-20
## EGR1    -10.3437075893327 2.31663359990025e-18 1.44789599993765e-17
## ANP32E   3.96659788431347 0.000130773249669056 0.000297211931066037
## CENPF    3.04960073891026  0.00319116449522177  0.00580211726403959
## LDHA     4.91244091643819 4.25376629157773e-06 1.11203558492587e-05

Note that DESeq2’s default normalization procedure is based on geometric means of counts, which are zero for genes with at least one zero count. This greatly limits the number of genes that can be used for normalization in scRNA-seq applications. We therefore use the normalization method suggested in the phyloseq package, which calculates the geometric mean for a gene by only using its positive counts, so that genes with zero counts could still be used for normalization purposes. The phyloseq normalization procedure can be applied by setting the argument type equal to poscounts in DESeq.

For UMI data, for which the expected counts may be very low, the likelihood ratio test implemented in nbinomLRT should be used. For other protocols (i.e., non-UMI), the Wald test in nbinomWaldTest can be used, with null distribution a t-distribution with degrees of freedom corrected by the observational weights. In both cases, we recommend the minimum expected count to be set to a small value (e.g., minmu=1e-6).

9 Using zinbwave with Seurat

The factors inferred in the zinbwave model can be added as one of the low dimensional data representations in the Seurat object, for instance to find subpopulations using Seurat’s cluster analysis method.

Note that the following workflow has only been tested with Seurat’s version 2.3.0.

Here we create a simple Seurat object with the raw data. Please, refer to the Seurat’s vignettes for a typical analysis, which includes filtering, normalization, etc.

library(Seurat)

seu <- CreateSeuratObject(raw.data = counts(fluidigm_zinb))

We can then add our zinbwave factors in the Seurat object.

seu <- SetDimReduction(object = seu, reduction.type = "zinbwave", 
                       slot = "cell.embeddings",
                       new.data = reducedDim(fluidigm_zinb, "zinbwave"))
seu <- SetDimReduction(object = seu, reduction.type = "zinbwave", slot = "key",
                       new.data = "zinbwave")

Finally, we can use the zinbwave factors for cluster analysis.

seu <- FindClusters(object = seu, reduction.type = "zinbwave", 
                    dims.use = 1:2, #this should match K
                    resolution = 0.6, print.output = 0, save.SNN = TRUE)

10 A note on performance and parallel computing

The zinbwave package uses the BiocParallel package to allow for parallel computing. Here, we used the register command to ensure that the vignette runs with serial computations.

However, in real datasets, parallel computations can speed up the computations dramatically, in the presence of many genes and/or many cells.

There are two ways of allowing parallel computations in zinbwave. The first is to register() a parallel back-end (see ?BiocParallel::register for details). Alternatively, one can pass a BPPARAM object to zinbwave and zinbFit, e.g.

library(BiocParallel)
zinb_res <- zinbwave(fluidigm, K=2, BPPARAM=MulticoreParam(2))

We found that MulticoreParam() may have some performance issues on Mac; hence, we recommend DoparParam() when working on Mac.

11 Session Info

sessionInfo()
## R version 3.5.0 (2018-04-23)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.4 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] Seurat_2.3.0                Matrix_1.2-14              
##  [3] cowplot_0.9.2               DESeq2_1.20.0              
##  [5] edgeR_3.22.0                limma_3.36.0               
##  [7] Rtsne_0.13                  biomaRt_2.36.0             
##  [9] ggplot2_2.2.1               magrittr_1.5               
## [11] scRNAseq_1.5.0              zinbwave_1.2.0             
## [13] SingleCellExperiment_1.2.0  SummarizedExperiment_1.10.0
## [15] DelayedArray_0.6.0          BiocParallel_1.14.0        
## [17] matrixStats_0.53.1          Biobase_2.40.0             
## [19] GenomicRanges_1.32.0        GenomeInfoDb_1.16.0        
## [21] IRanges_2.14.0              S4Vectors_0.18.0           
## [23] BiocGenerics_0.26.0         BiocStyle_2.8.0            
## 
## loaded via a namespace (and not attached):
##   [1] R.utils_2.6.0          tidyselect_0.2.4       RSQLite_2.1.0         
##   [4] AnnotationDbi_1.42.0   htmlwidgets_1.2        ranger_0.9.0          
##   [7] grid_3.5.0             trimcluster_0.1-2      munsell_0.4.3         
##  [10] ica_1.0-1              codetools_0.2-15       withr_2.1.2           
##  [13] colorspace_1.3-2       knitr_1.20             rstudioapi_0.7        
##  [16] geometry_0.3-6         pspline_1.0-18         ROCR_1.0-7            
##  [19] dtw_1.18-1             robustbase_0.93-0      dimRed_0.1.0          
##  [22] labeling_0.3           lars_1.2               GenomeInfoDbData_1.1.0
##  [25] mnormt_1.5-5           bit64_0.9-7            rprojroot_1.3-2       
##  [28] ipred_0.9-6            xfun_0.1               diptest_0.75-7        
##  [31] R6_2.2.2               VGAM_1.0-5             locfit_1.5-9.1        
##  [34] flexmix_2.3-14         DRR_0.0.3              bitops_1.0-6          
##  [37] assertthat_0.2.0       SDMTools_1.1-221       scales_0.5.0          
##  [40] nnet_7.3-12            gtable_0.2.0           ddalpha_1.3.3         
##  [43] timeDate_3043.102      rlang_0.2.0            CVST_0.2-1            
##  [46] genefilter_1.62.0      scatterplot3d_0.3-41   RcppRoll_0.2.2        
##  [49] splines_3.5.0          lazyeval_0.2.1         ModelMetrics_1.1.0    
##  [52] acepack_1.4.1          broom_0.4.4            checkmate_1.8.5       
##  [55] yaml_2.1.18            reshape2_1.4.3         abind_1.4-5           
##  [58] backports_1.1.2        Hmisc_4.1-1            caret_6.0-79          
##  [61] tools_3.5.0            lava_1.6.1             bookdown_0.7          
##  [64] psych_1.8.3.3          gplots_3.0.1           RColorBrewer_1.1-2    
##  [67] proxy_0.4-22           stabledist_0.7-1       ggridges_0.5.0        
##  [70] Rcpp_0.12.16           plyr_1.8.4             base64enc_0.1-3       
##  [73] progress_1.1.2         zlibbioc_1.26.0        purrr_0.2.4           
##  [76] RCurl_1.95-4.10        prettyunits_1.0.2      rpart_4.1-13          
##  [79] pbapply_1.3-4          zoo_1.8-1              sfsmisc_1.1-2         
##  [82] cluster_2.0.7-1        data.table_1.10.4-3    lmtest_0.9-36         
##  [85] RANN_2.5.1             mvtnorm_1.0-7          fitdistrplus_1.0-9    
##  [88] gsl_1.9-10.3           evaluate_0.10.1        xtable_1.8-2          
##  [91] XML_3.98-1.11          mclust_5.4             gridExtra_2.3         
##  [94] compiler_3.5.0         tibble_1.4.2           KernSmooth_2.23-15    
##  [97] R.oo_1.22.0            htmltools_0.3.6        segmented_0.5-3.0     
## [100] pcaPP_1.9-73           Formula_1.2-2          snow_0.4-2            
## [103] tidyr_0.8.0            geneplotter_1.58.0     tclust_1.3-1          
## [106] lubridate_1.7.4        DBI_0.8                diffusionMap_1.1-0    
## [109] magic_1.5-8            MASS_7.3-50            fpc_2.1-11            
## [112] R.methodsS3_1.7.1      gdata_2.18.0           metap_0.9             
## [115] bindr_0.1.1            gower_0.1.2            igraph_1.2.1          
## [118] pkgconfig_2.0.1        sn_1.5-2               numDeriv_2016.8-1     
## [121] foreign_0.8-70         recipes_0.1.2          foreach_1.4.4         
## [124] annotate_1.58.0        XVector_0.20.0         prodlim_2018.04.18    
## [127] stringr_1.3.0          digest_0.6.15          tsne_0.1-3            
## [130] copula_0.999-18        ADGofTest_0.3          softImpute_1.4        
## [133] rmarkdown_1.9          htmlTable_1.11.2       kernlab_0.9-26        
## [136] gtools_3.5.0           modeltools_0.2-21      nlme_3.1-137          
## [139] bindrcpp_0.2.2         pillar_1.2.2           lattice_0.20-35       
## [142] httr_1.3.1             DEoptimR_1.0-8         survival_2.42-3       
## [145] glue_1.2.0             FNN_1.1                png_0.1-7             
## [148] prabclus_2.2-6         iterators_1.0.9        glmnet_2.0-16         
## [151] bit_1.1-12             mixtools_1.1.0         class_7.3-14          
## [154] stringi_1.1.7          blob_1.1.1             doSNOW_1.0.16         
## [157] latticeExtra_0.6-28    caTools_1.17.1         memoise_1.1.0         
## [160] dplyr_0.7.4            irlba_2.3.2            ape_5.1

References

Pollen, Alex A, Tomasz J Nowakowski, Joe Shuga, Xiaohui Wang, Anne A Leyrat, Jan H Lui, Nianzhen Li, et al. 2014. “Low-coverage single-cell mRNA sequencing reveals cellular heterogeneity and activated signaling pathways in developing cerebral cortex.” Nature Biotechnology 32 (10):1053–8.

Risso, D, F Perraudeau, S Gribkova, S Dudoit, and Vert JP. 2018. “A General and Flexible Method for Signal Extraction from Single-Cell RNA-Seq Data.” Nature Communications 9:284.

Van den Berge, Koen, Fanny Perraudeau, Charlotte Soneson, Michael I Love, Davide Risso, Jean-Philippe Vert, Mark D Robinson, Sandrine Dudoit, and Lieven Clement. 2018. “Observation Weights to Unlock Bulk Rna-Seq Tools for Zero Inflation and Single-Cell Applications.” bioRxiv. Cold Spring Harbor Laboratory, 250126.