Typical RNA-seq call from DESeq2

Note: the typical RNA-seq workflow for users would be to call apeglm estimation from within the lfcShrink function from the DESeq2 package. The unevaluated code chunk shows how to obtain apeglm shrinkage estimates after running DESeq. See the DESeq2 vignette for more details. The lfcShrink wrapper function takes care of many details below, and unifies the interface for multiple shrinkage estimators. The coefficient to shrink can be specified either by name or by number (following the order in resultsNames(dds)). Be aware that DESeq2’s lfcShrink interface provides LFCs on the log2 scale, while apeglm provides coefficients on the natural log scale.

Example RNA-seq analysis

Here we show example code which mimics what will happen inside the lfcShrink function when using the apeglm method (Zhu, Ibrahim, and Love 2018).

Load a prepared SummarizedExperiment:

##                 SRR1039508 SRR1039509 SRR1039512 SRR1039513 SRR1039516
## ENSG00000000003        679        448        873        408       1138
## ENSG00000000005          0          0          0          0          0
## ENSG00000000419        467        515        621        365        587
## ENSG00000000457        260        211        263        164        245
## ENSG00000000460         60         55         40         35         78
## ENSG00000000938          0          0          2          0          1
##                 SRR1039517 SRR1039520 SRR1039521
## ENSG00000000003       1047        770        572
## ENSG00000000005          0          0          0
## ENSG00000000419        799        417        508
## ENSG00000000457        331        233        229
## ENSG00000000460         63         76         60
## ENSG00000000938          0          0          0

For demonstration, we will use 5000 genes of the airway dataset, the first from those genes with at least 10 counts across all samples.

First run a DESeq2 differential expression analysis (Love, Huber, and Anders 2014) (size factors, and dispersion estimates could similarly be estimated using edgeR):

Defining data and parameter objects necessary for apeglm. We must multiply the coefficients from DESeq2 by a factor, because apeglm provides natural log coefficients. Again, this would be handled inside of lfcShrink in DESeq2 for a typical RNA-seq analysis.

Running apeglm

Here apeglm on 5000 genes takes less than a minute on a laptop. It scales with number of genes, the number of samples and the number of variables in the design formula, where here we have 5 coefficients (one for the four cell cultures and one for the difference due to dexamethasone treatment).

We provide apeglm with the SummarizedExperiment although the function can also run on a matrix of counts or other observed data. We specify a coef as well as a threshold which we discuss below. Note that we multiple the threshold by log(2) to convert from log2 scale to natural log scale.

##    user  system elapsed 
##  25.296   0.000  25.296
## List of 7
##  $ no.shrink            : int [1:4] 1 2 3 4
##  $ prior.mean           : num 0
##  $ prior.scale          : num 0.299
##  $ prior.df             : num 1
##  $ prior.no.shrink.mean : num 0
##  $ prior.no.shrink.scale: num 15
##  $ prior.var            : num 0.0894

There are faster implementations of apeglm specifically for negative binomial likelihoods. The version nbinomR is ~5 times faster than the default general.

##    user  system elapsed 
##   5.308   0.000   5.308

The version nbinomCR is ~10 times faster than the default general.

##    user  system elapsed 
##   2.596   0.000   2.597

The version nbinomC returns only the MAP coefficients and can be ~50-100 times faster than the default general. The MAP coefficients are the same as returned by nbinomCR above, we just skip the calculation of posterior SD. A variant of nbinomC is nbinomC* which includes random starts.

##    user  system elapsed 
##   0.163   0.000   0.162

Among other output, we have the estimated coefficients attached to the ranges of the SummarizedExperiment used as input:

## [1] "CompressedGRangesList"
## attr(,"package")
## [1] "GenomicRanges"
## DataFrame with 5000 rows and 5 columns
##                       X.Intercept.         cellN061011         cellN080611
##                          <numeric>           <numeric>           <numeric>
## ENSG00000000003     6.645656858832  0.0607881068431617   0.229672330993859
## ENSG00000000419   6.23452224381182  -0.084999062909147 -0.0151414731605632
## ENSG00000000457   5.45295046291592   0.060167029688361 -0.0556310344682003
## ENSG00000000460   3.77047988584934   0.545203896906746   0.261356100290922
## ENSG00000000938 -0.314689549938794   -4.17748981892207  -0.539678670081039
## ...                            ...                 ...                 ...
## ENSG00000122756  -3.62985171699038    4.03910910191739   -2.00591343576517
## ENSG00000122778   5.23905962019645  0.0800614479140051    0.78571824141725
## ENSG00000122779   5.66860926529946 -0.0569747440556599  0.0501292974549099
## ENSG00000122783   6.01534721390051  0.0368918236863277  -0.157828298844423
## ENSG00000122786   9.78251386175167   -0.18234928440018  0.0942307382694629
##                          cellN61311                dextrt
##                           <numeric>             <numeric>
## ENSG00000000003  -0.152039407731068    -0.264701484867069
## ENSG00000000419 -0.0476622243521883     0.115132140591595
## ENSG00000000457  0.0484622689335007   0.00907000351291652
## ENSG00000000460   0.349429848485217   -0.0425238103250921
## ENSG00000000938   -4.20084641372255  -0.00836574242756413
## ...                             ...                   ...
## ENSG00000122756     2.9579939014709 -0.000701743005250805
## ENSG00000122778   0.240098604814436    -0.285457153652284
## ENSG00000122779 -0.0441315171481841    -0.109945276886632
## ENSG00000122783   -0.29739986153585    0.0690120376649482
## ENSG00000122786   -0.12028002880088     0.523514860323912

We can compare the coefficients from apeglm with those from DESeq2. apeglm provides coefficients on the natural log scale, so we must convert to log2 scale by multiplying by log2(exp(1)). Note that DESeq2’s lfcShrink function converts apeglm coefficients to the log2 scale internally.

##    user  system elapsed 
##   1.972   0.000   1.973

Here we plot apeglm estimators against DESeq2:

Here we plot MLE, DESeq2 and apeglm estimators against the mean of normalized counts:

Specific coefficients

Note that p-values and FSR define different events, and are not on the same scale. An FSR of 0.5 means that the estimated sign is as bad as random guess.

The s-value was proposed by Stephens (2016), as a statistic giving the aggregate false sign rate for tests with equal or lower s-value than the one considered. We recommend using a lower threshold on s-values than typically used for adjusted p-values, for example one might be interested in sets with 0.01 or 0.005 aggregate FSR.

More scrutiny can be applied by using an LFC threshold greater than zero, and asking for the probability of a “false-sign-or-small” (FSOS) event: that the effect size is not further from zero in distance than the threshold amount. We can run the svalue function on these per-gene probabilities to produce s-values that bound the FSOS rate for sets of genes. By specifying threshold=log(2) * 1 above, apeglm will then output a vector thresh in the results list that gives the per-gene probabilities of false-sign-or-small events.

Modeling zero-inflated counts

We have created a separate GitHub repository giving an example of how the apeglm estimator can be used for Zero-Inflated Negative Binomial data. This approach uses the zinbwave method and Bioconductor package to estimate the probability of each zero belonging to the excess zero component. We compare using a Negative Binomial likelihood with the excess zeros down-weighted and using a Zero-Inflated Negative Binomial likelihood. These two approaches with apeglm perform similarly but we note that the first approach involves less additional code and is faster to compute.

https://github.com/mikelove/zinbwave-apeglm

Modeling ratios of counts

We also show an short example using an alternative likelihood to the negative binomial. Suppose we have allele-specific counts for n=20 vs 20 samples across 5000 genes. We can define a binomial model and test for the allelic balance across groups of samples.

Here we will simulate allele counts from our existing dataset for demonstration. We spike in 10 genes with strong allelic imbalance (instead of an allelic ratio close to 0.5, these will have a ratio of 0.75).

We define a beta-binomial likelihood function which uses the total counts as a parameter, and the logit function as a link:

We first need to estimate MLE coefficients and standard errors.

##    user  system elapsed 
##  18.735   0.000  18.734

##    user  system elapsed 
##  32.142   0.000  32.148

In the apeglm plot, we color in red the genes with a low aggregate probability of false-sign-or-small (FSOS) events (s-value < .01), where we’ve again defined “small” on the log odds scale using the threshold argument above.

## [1] 1.098612
## 
## TRUE 
##   10

Acknowledgments

We have benefited in the development of apeglm from feedback or contributions from the following individuals:

Wolfgang Huber, Cecile Le Sueur, Charlotte Soneson

Session Info

## R version 3.6.0 (2019-04-26)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.9-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.9-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] emdbook_1.3.11              apeglm_1.6.0               
##  [3] DESeq2_1.24.0               airway_1.3.0               
##  [5] SummarizedExperiment_1.14.0 DelayedArray_0.10.0        
##  [7] BiocParallel_1.18.0         matrixStats_0.54.0         
##  [9] Biobase_2.44.0              GenomicRanges_1.36.0       
## [11] GenomeInfoDb_1.20.0         IRanges_2.18.0             
## [13] S4Vectors_0.22.0            BiocGenerics_0.30.0        
## 
## loaded via a namespace (and not attached):
##  [1] bit64_0.9-7            splines_3.6.0          Formula_1.2-3         
##  [4] assertthat_0.2.1       latticeExtra_0.6-28    blob_1.1.1            
##  [7] GenomeInfoDbData_1.2.1 yaml_2.2.0             numDeriv_2016.8-1     
## [10] pillar_1.3.1           RSQLite_2.1.1          backports_1.1.4       
## [13] lattice_0.20-38        glue_1.3.1             bbmle_1.0.20          
## [16] digest_0.6.18          RColorBrewer_1.1-2     XVector_0.24.0        
## [19] checkmate_1.9.1        colorspace_1.4-1       htmltools_0.3.6       
## [22] Matrix_1.2-17          plyr_1.8.4             XML_3.98-1.19         
## [25] pkgconfig_2.0.2        genefilter_1.66.0      zlibbioc_1.30.0       
## [28] purrr_0.3.2            xtable_1.8-4           scales_1.0.0          
## [31] htmlTable_1.13.1       tibble_2.1.1           annotate_1.62.0       
## [34] ggplot2_3.1.1          nnet_7.3-12            lazyeval_0.2.2        
## [37] survival_2.44-1.1      magrittr_1.5           crayon_1.3.4          
## [40] memoise_1.1.0          evaluate_0.13          MASS_7.3-51.4         
## [43] foreign_0.8-71         tools_3.6.0            data.table_1.12.2     
## [46] stringr_1.4.0          locfit_1.5-9.1         munsell_0.5.0         
## [49] cluster_2.0.9          AnnotationDbi_1.46.0   compiler_3.6.0        
## [52] rlang_0.3.4            grid_3.6.0             RCurl_1.95-4.12       
## [55] rstudioapi_0.10        htmlwidgets_1.3        bitops_1.0-6          
## [58] base64enc_0.1-3        rmarkdown_1.12         gtable_0.3.0          
## [61] codetools_0.2-16       DBI_1.0.0              R6_2.4.0              
## [64] gridExtra_2.3          knitr_1.22             dplyr_0.8.0.1         
## [67] bit_1.1-14             Hmisc_4.2-0            stringi_1.4.3         
## [70] Rcpp_1.0.1             geneplotter_1.62.0     rpart_4.1-15          
## [73] acepack_1.4.1          coda_0.19-2            tidyselect_0.2.5      
## [76] xfun_0.6

References

Love, Michael I., Wolfgang Huber, and Simon Anders. 2014. “Moderated estimation of fold change and dispersion for RNA-seq data with DESeq2.” Genome Biology 15 (12):550. https://doi.org/10.1186/s13059-014-0550-8.

Stephens, Matthew. 2016. “False Discovery Rates: A New Deal.” Biostatistics 18 (2). https://doi.org/10.1093/biostatistics/kxw041.

Zhu, Anqi, Joseph G. Ibrahim, and Michael I. Love. 2018. “Heavy-Tailed Prior Distributions for Sequence Count Data: Removing the Noise and Preserving Large Differences.” Bioinformatics. https://doi.org/10.1093/bioinformatics/bty895.