Installation

To install and load NBAMSeq

if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
BiocManager::install("NBAMSeq")
library(NBAMSeq)

Introduction

High-throughput sequencing experiments followed by differential expression analysis is a widely used approach to detect genomic biomarkers. A fundamental step in differential expression analysis is to model the association between gene counts and covariates of interest. NBAMSeq is a flexible statistical model based on the generalized additive model and allows for information sharing across genes in variance estimation. Specifically, we model the logarithm of mean gene counts as sums of smooth functions with the smoothing parameters and coefficients estimated simultaneously by a nested iteration. The variance is estimated by the Bayesian shrinkage approach to fully exploit the information across all genes.

The workflow of NBAMSeq contains three main steps:

Here we illustrate each of these steps respectively.

Data input

Users are expected to provide three parts of input, i.e. countData, colData, and design.

countData is a matrix of gene counts generated by RNASeq experiments.

## An example of countData
n = 50  ## n stands for number of genes
m = 20   ## m stands for sample size
countData = matrix(rnbinom(n*m, mu=100, size=1/3), ncol = m) + 1
mode(countData) = "integer"
colnames(countData) = paste0("sample", 1:m)
rownames(countData) = paste0("gene", 1:n)
head(countData)
      sample1 sample2 sample3 sample4 sample5 sample6 sample7 sample8 sample9
gene1      42      12      64     276       6     299       1       9      86
gene2     610       1      99     806     138      54       6       6       3
gene3     214      16       9     235      19       1      12       6      32
gene4      85     289       9     422       2     133     221       2       2
gene5     505     140      11       1      23      21       6      14      65
gene6       3      15     254      58      12     347      85       4      13
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1       72       83        1       13        1        1      142       20
gene2      554        8      420       13        6        1        1        1
gene3        1      116       17      546      105        3       12     1063
gene4      106       29        1      471       11        1       49       60
gene5      938       41       88        1      218        8      301       16
gene6      108        7       35      281     1316       40        9       29
      sample18 sample19 sample20
gene1        2      122        2
gene2       90        7        1
gene3       36       47       21
gene4       68      657        1
gene5      235        5        1
gene6        2        8       19

colData is a data frame which contains the covariates of samples. The sample order in colData should match the sample order in countData.

## An example of colData
pheno = runif(m, 20, 80)
var1 = rnorm(m)
var2 = rnorm(m)
var3 = rnorm(m)
var4 = as.factor(sample(c(0,1,2), m, replace = TRUE))
colData = data.frame(pheno = pheno, var1 = var1, var2 = var2,
    var3 = var3, var4 = var4)
rownames(colData) = paste0("sample", 1:m)
head(colData)
           pheno        var1       var2       var3 var4
sample1 79.79865 -0.06153226  0.9076525 -0.5013926    2
sample2 65.33880  0.98638962 -0.7226681  2.1574485    0
sample3 46.79358  0.47751731 -0.5039563 -1.3706989    2
sample4 52.01470 -0.13822694  0.6972845  1.5473712    2
sample5 35.14663  0.82205264  1.2476779  1.2422261    1
sample6 66.79806 -0.65425491 -0.5837493 -0.2428002    2

design is a formula which specifies how to model the samples. Compared with other packages performing DE analysis including DESeq2 (Love, Huber, and Anders 2014), edgeR (Robinson, McCarthy, and Smyth 2010), NBPSeq (Di et al. 2015) and BBSeq (Zhou, Xia, and Wright 2011), NBAMSeq supports the nonlinear model of covariates via mgcv (Wood and Wood 2015). To indicate the nonlinear covariate in the model, users are expected to use s(variable_name) in the design formula. In our example, if we would like to model pheno as a nonlinear covariate, the design formula should be:

design = ~ s(pheno) + var1 + var2 + var3 + var4

Several notes should be made regarding the design formula:

We then construct the NBAMSeqDataSet using countData, colData, and design:

gsd = NBAMSeqDataSet(countData = countData, colData = colData, design = design)
gsd
class: NBAMSeqDataSet 
dim: 50 20 
metadata(1): fitted
assays(1): counts
rownames(50): gene1 gene2 ... gene49 gene50
rowData names(0):
colnames(20): sample1 sample2 ... sample19 sample20
colData names(5): pheno var1 var2 var3 var4

Differential expression analysis

Differential expression analysis can be performed by NBAMSeq function:

gsd = NBAMSeq(gsd)

Several other arguments in NBAMSeq function are available for users to customize the analysis.

library(BiocParallel)
gsd = NBAMSeq(gsd, parallel = TRUE)

Pulling out DE results

Results of DE analysis can be pulled out by results function. For continuous covariates, the name argument should be specified indicating the covariate of interest. For nonlinear continuous covariates, base mean, effective degrees of freedom (edf), test statistics, p-value, and adjusted p-value will be returned.

res1 = results(gsd, name = "pheno")
head(res1)
DataFrame with 6 rows and 7 columns
       baseMean       edf      stat     pvalue      padj       AIC       BIC
      <numeric> <numeric> <numeric>  <numeric> <numeric> <numeric> <numeric>
gene1   50.6332   1.00003 8.0173123 0.00463400 0.0579250   197.202   204.172
gene2  137.0070   1.00005 0.0176062 0.89462202 0.9502937   215.172   222.143
gene3   88.0940   1.00006 0.1540755 0.69473999 0.8892969   218.179   225.149
gene4   99.0527   1.00004 2.6706575 0.10222217 0.3006534   226.258   233.228
gene5  120.0785   1.00004 6.8599270 0.00881689 0.0734741   229.064   236.034
gene6   96.9594   1.00002 3.5402373 0.05990633 0.1953475   228.219   235.189

For linear continuous covariates, base mean, estimated coefficient, standard error, test statistics, p-value, and adjusted p-value will be returned.

res2 = results(gsd, name = "var1")
head(res2)
DataFrame with 6 rows and 8 columns
       baseMean      coef        SE      stat     pvalue      padj       AIC
      <numeric> <numeric> <numeric> <numeric>  <numeric> <numeric> <numeric>
gene1   50.6332 -0.642995  0.382280 -1.681999 0.09256908 0.2892784   197.202
gene2  137.0070 -0.740851  0.464692 -1.594285 0.11087232 0.3260951   215.172
gene3   88.0940 -0.394955  0.376514 -1.048979 0.29418783 0.5883757   218.179
gene4   99.0527 -1.384980  0.431403 -3.210408 0.00132547 0.0331367   226.258
gene5  120.0785 -0.129879  0.429034 -0.302725 0.76209913 0.9406000   229.064
gene6   96.9594 -0.283982  0.406718 -0.698227 0.48503506 0.6978130   228.219
            BIC
      <numeric>
gene1   204.172
gene2   222.143
gene3   225.149
gene4   233.228
gene5   236.034
gene6   235.189

For discrete covariates, the contrast argument should be specified. e.g.  contrast = c("var4", "2", "0") means comparing level 2 vs. level 0 in var4.

res3 = results(gsd, contrast = c("var4", "2", "0"))
head(res3)
DataFrame with 6 rows and 8 columns
       baseMean      coef        SE      stat      pvalue      padj       AIC
      <numeric> <numeric> <numeric> <numeric>   <numeric> <numeric> <numeric>
gene1   50.6332  2.924086   1.12270  2.604501 0.009200810  0.115010   197.202
gene2  137.0070  4.606767   1.37670  3.346248 0.000819132  0.019993   215.172
gene3   88.0940 -1.633997   1.09340 -1.494412 0.135067947  0.350687   218.179
gene4   99.0527 -0.583349   1.21398 -0.480526 0.630853432  0.789194   226.258
gene5  120.0785  0.489929   1.24737  0.392769 0.694490329  0.789194   229.064
gene6   96.9594 -0.642245   1.18273 -0.543020 0.587116334  0.789194   228.219
            BIC
      <numeric>
gene1   204.172
gene2   222.143
gene3   225.149
gene4   233.228
gene5   236.034
gene6   235.189

Visualization

We suggest two approaches to visualize the nonlinear associations. The first approach is to plot the smooth components of a fitted negative binomial additive model by plot.gam function in mgcv (Wood and Wood 2015). This can be done by calling makeplot function and passing in NBAMSeqDataSet object. Users are expected to provide the phenotype of interest in phenoname argument and gene of interest in genename argument.

## assuming we are interested in the nonlinear relationship between gene10's 
## expression and "pheno"
makeplot(gsd, phenoname = "pheno", genename = "gene10", main = "gene10")

In addition, to explore the nonlinear association of covariates, it is also instructive to look at log normalized counts vs. variable scatter plot. Below we show how to produce such plot.

## here we explore the most significant nonlinear association
res1 = res1[order(res1$pvalue),]
topgene = rownames(res1)[1]  
sf = getsf(gsd)  ## get the estimated size factors
## divide raw count by size factors to obtain normalized counts
countnorm = t(t(countData)/sf) 
head(res1)
DataFrame with 6 rows and 7 columns
        baseMean       edf      stat      pvalue      padj       AIC       BIC
       <numeric> <numeric> <numeric>   <numeric> <numeric> <numeric> <numeric>
gene31   84.9508   1.00008  17.89424 2.37181e-05 0.0011859   209.366   216.336
gene37   94.5447   1.00008   8.07125 4.49956e-03 0.0579250   212.226   219.197
gene13   84.4954   1.00011   8.06034 4.52787e-03 0.0579250   212.581   219.551
gene1    50.6332   1.00003   8.01731 4.63400e-03 0.0579250   197.202   204.172
gene15   98.4627   1.00003   7.46910 6.27837e-03 0.0627837   220.124   227.094
gene5   120.0785   1.00004   6.85993 8.81689e-03 0.0734741   229.064   236.034
library(ggplot2)
setTitle = topgene
df = data.frame(pheno = pheno, logcount = log2(countnorm[topgene,]+1))
ggplot(df, aes(x=pheno, y=logcount))+geom_point(shape=19,size=1)+
    geom_smooth(method='loess')+xlab("pheno")+ylab("log(normcount + 1)")+
    annotate("text", x = max(df$pheno)-5, y = max(df$logcount)-1, 
    label = paste0("edf: ", signif(res1[topgene,"edf"],digits = 4)))+
    ggtitle(setTitle)+
    theme(text = element_text(size=10), plot.title = element_text(hjust = 0.5))

Session info

sessionInfo()
R version 4.0.3 (2020-10-10)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 18.04.5 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.12-bioc/R/lib/libRblas.so
LAPACK: /home/biocbuild/bbs-3.12-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] ggplot2_3.3.2               BiocParallel_1.24.0        
 [3] NBAMSeq_1.6.1               SummarizedExperiment_1.20.0
 [5] Biobase_2.50.0              GenomicRanges_1.42.0       
 [7] GenomeInfoDb_1.26.0         IRanges_2.24.0             
 [9] S4Vectors_0.28.0            BiocGenerics_0.36.0        
[11] MatrixGenerics_1.2.0        matrixStats_0.57.0         

loaded via a namespace (and not attached):
 [1] Rcpp_1.0.5             locfit_1.5-9.4         lattice_0.20-41       
 [4] digest_0.6.27          R6_2.5.0               RSQLite_2.2.1         
 [7] evaluate_0.14          httr_1.4.2             pillar_1.4.6          
[10] zlibbioc_1.36.0        rlang_0.4.8            annotate_1.68.0       
[13] blob_1.2.1             Matrix_1.2-18          rmarkdown_2.5         
[16] labeling_0.4.2         splines_4.0.3          geneplotter_1.68.0    
[19] stringr_1.4.0          RCurl_1.98-1.2         bit_4.0.4             
[22] munsell_0.5.0          DelayedArray_0.16.0    compiler_4.0.3        
[25] xfun_0.18              pkgconfig_2.0.3        mgcv_1.8-33           
[28] htmltools_0.5.0        tidyselect_1.1.0       tibble_3.0.4          
[31] GenomeInfoDbData_1.2.4 XML_3.99-0.5           withr_2.3.0           
[34] crayon_1.3.4           dplyr_1.0.2            bitops_1.0-6          
[37] grid_4.0.3             nlme_3.1-150           xtable_1.8-4          
[40] gtable_0.3.0           lifecycle_0.2.0        DBI_1.1.0             
[43] magrittr_1.5           scales_1.1.1           stringi_1.5.3         
[46] farver_2.0.3           XVector_0.30.0         genefilter_1.72.0     
[49] ellipsis_0.3.1         vctrs_0.3.4            generics_0.0.2        
[52] RColorBrewer_1.1-2     tools_4.0.3            bit64_4.0.5           
[55] glue_1.4.2             DESeq2_1.30.0          purrr_0.3.4           
[58] survival_3.2-7         yaml_2.2.1             AnnotationDbi_1.52.0  
[61] colorspace_1.4-1       memoise_1.1.0          knitr_1.30            

References

Di, Y, DW Schafer, JS Cumbie, and JH Chang. 2015. “NBPSeq: Negative Binomial Models for Rna-Sequencing Data.” R Package Version 0.3. 0, URL Http://CRAN. R-Project. Org/Package= NBPSeq.

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.

Robinson, Mark D, Davis J McCarthy, and Gordon K Smyth. 2010. “EdgeR: A Bioconductor Package for Differential Expression Analysis of Digital Gene Expression Data.” Bioinformatics 26 (1): 139–40.

Wood, Simon, and Maintainer Simon Wood. 2015. “Package ’Mgcv’.” R Package Version 1: 29.

Zhou, Yi-Hui, Kai Xia, and Fred A Wright. 2011. “A Powerful and Flexible Approach to the Analysis of Rna Sequence Count Data.” Bioinformatics 27 (19): 2672–8.