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       1       1      24      54     113      63      27       2      12
gene2      30     138      99      85      20       3     210      26       3
gene3     157     136      58       4      16      65      16       3       2
gene4       2      12       8       8       1       1      96     208       2
gene5       1     118      20      10       1       6      14       4      15
gene6      61      70      58     132      41     279      48       4      11
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1       97       11       59        3        7       60       93        6
gene2       37       79       93       17       55       44        1       10
gene3       19      134       19        1      399       35       84       17
gene4      584        6        9        2        5      126       13       48
gene5        1       50       77       20        1      406       16        1
gene6       10       88      161       82      119      128        5       21
      sample18 sample19 sample20
gene1      147       56       27
gene2        1       77       82
gene3      138      313      597
gene4        6        2       30
gene5       11       33        1
gene6       13       84        7

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 73.80921 -0.5296293 -0.38644692  0.3314055    0
sample2 62.64547 -1.4221333  0.72979185 -0.7813876    2
sample3 42.21892 -0.1342837 -0.02289044 -0.1576072    1
sample4 41.15586  0.5441802  0.46856961 -0.1227523    1
sample5 32.68910 -0.8832228  1.15048656 -0.5047426    1
sample6 37.19406 -1.1186365 -0.53091742  0.5848132    1

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   37.9490   1.00007 0.7814047 0.3767691  0.819063   204.296   211.266
gene2   50.5655   1.00006 2.0634587 0.1508802  0.538858   211.916   218.886
gene3  129.5485   1.00007 3.1144645 0.0776163  0.388082   229.087   236.057
gene4   49.0171   1.00003 0.0255126 0.8731673  0.914205   182.852   189.822
gene5   33.1165   1.00011 0.5601260 0.4542135  0.825876   180.221   187.191
gene6   67.1276   1.00009 0.9154465 0.3387066  0.771265   226.631   233.601

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   37.9490 -0.375223  0.341397  -1.09908 2.71734e-01 0.55241945   204.296
gene2   50.5655 -0.404034  0.318444  -1.26877 2.04522e-01 0.52845201   211.916
gene3  129.5485 -1.399246  0.344755  -4.05867 4.93539e-05 0.00246769   229.087
gene4   49.0171 -0.571546  0.347666  -1.64395 1.00186e-01 0.38533010   182.852
gene5   33.1165 -0.461115  0.368958  -1.24978 2.11381e-01 0.52845201   180.221
gene6   67.1276 -0.264273  0.312682  -0.84518 3.98011e-01 0.62189141   226.631
            BIC
      <numeric>
gene1   211.266
gene2   218.886
gene3   236.057
gene4   189.822
gene5   187.191
gene6   233.601

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   37.9490  0.991322   1.15017  0.861891 0.3887474  0.716305   204.296
gene2   50.5655  0.544486   1.07866  0.504782 0.6137122  0.870925   211.916
gene3  129.5485  1.984795   1.15507  1.718328 0.0857368  0.389713   229.087
gene4   49.0171  0.504952   1.19958  0.420939 0.6737996  0.870925   182.852
gene5   33.1165  3.198216   1.25235  2.553770 0.0106564  0.133205   180.221
gene6   67.1276  0.459218   1.05367  0.435827 0.6629623  0.870925   226.631
            BIC
      <numeric>
gene1   211.266
gene2   218.886
gene3   236.057
gene4   189.822
gene5   187.191
gene6   233.601

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>
gene22  119.0806   1.00005  17.76521 2.57896e-05 0.00128948   204.776   211.747
gene44   94.8223   1.00009  15.66249 7.53275e-05 0.00188319   215.551   222.521
gene32   45.8879   1.00005  10.47680 1.20926e-03 0.02015431   189.474   196.445
gene25   81.1212   1.00006   6.40270 1.13983e-02 0.14247837   202.247   209.218
gene46  107.5987   1.00006   5.49553 1.90712e-02 0.19071162   201.358   208.329
gene26  158.5995   1.00004   4.90149 2.68364e-02 0.22363702   220.513   227.483
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 Under development (unstable) (2023-10-22 r85388)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 22.04.3 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.19-bioc/R/lib/libRblas.so 
LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0

locale:
 [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
 [3] LC_TIME=en_GB              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       

time zone: America/New_York
tzcode source: system (glibc)

attached base packages:
[1] stats4    stats     graphics  grDevices utils     datasets  methods  
[8] base     

other attached packages:
 [1] ggplot2_3.4.4               BiocParallel_1.37.0        
 [3] NBAMSeq_1.19.0              SummarizedExperiment_1.33.0
 [5] Biobase_2.63.0              GenomicRanges_1.55.1       
 [7] GenomeInfoDb_1.39.0         IRanges_2.37.0             
 [9] S4Vectors_0.41.1            BiocGenerics_0.49.1        
[11] MatrixGenerics_1.15.0       matrixStats_1.0.0          

loaded via a namespace (and not attached):
 [1] KEGGREST_1.43.0         gtable_0.3.4            xfun_0.41              
 [4] bslib_0.5.1             lattice_0.22-5          vctrs_0.6.4            
 [7] tools_4.4.0             bitops_1.0-7            generics_0.1.3         
[10] parallel_4.4.0          RSQLite_2.3.2           AnnotationDbi_1.65.0   
[13] tibble_3.2.1            fansi_1.0.5             highr_0.10             
[16] blob_1.2.4              pkgconfig_2.0.3         Matrix_1.6-1.1         
[19] lifecycle_1.0.3         GenomeInfoDbData_1.2.11 farver_2.1.1           
[22] compiler_4.4.0          Biostrings_2.71.1       munsell_0.5.0          
[25] DESeq2_1.43.0           codetools_0.2-19        htmltools_0.5.6.1      
[28] sass_0.4.7              RCurl_1.98-1.12         yaml_2.3.7             
[31] pillar_1.9.0            crayon_1.5.2            jquerylib_0.1.4        
[34] DelayedArray_0.29.0     cachem_1.0.8            abind_1.4-5            
[37] nlme_3.1-163            genefilter_1.85.0       tidyselect_1.2.0       
[40] locfit_1.5-9.8          digest_0.6.33           dplyr_1.1.3            
[43] labeling_0.4.3          splines_4.4.0           fastmap_1.1.1          
[46] grid_4.4.0              colorspace_2.1-0        cli_3.6.1              
[49] SparseArray_1.3.0       magrittr_2.0.3          S4Arrays_1.3.0         
[52] survival_3.5-7          XML_3.99-0.14           utf8_1.2.4             
[55] withr_2.5.2             scales_1.2.1            bit64_4.0.5            
[58] rmarkdown_2.25          XVector_0.43.0          httr_1.4.7             
[61] bit_4.0.5               png_0.1-8               memoise_2.0.1          
[64] evaluate_0.23           knitr_1.45              mgcv_1.9-0             
[67] rlang_1.1.1             Rcpp_1.0.11             DBI_1.1.3              
[70] xtable_1.8-4            glue_1.6.2              annotate_1.81.0        
[73] jsonlite_1.8.7          R6_2.5.1                zlibbioc_1.49.0        

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.