NBAMSeq: Negative Binomial Additive Model for RNA-Seq Data

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:

  • Step 1: Data input using NBAMSeqDataSet;

  • Step 2: Differential expression (DE) analysis using NBAMSeq function;

  • Step 3: Pulling out DE results using results function.

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      75       1     247     185       3      75      12      49      14
gene2       1       1       1       2       1       2     236       5      28
gene3      93       1     353       2      10     140      61       1      15
gene4      40     406      12       2      63       6     179       2    1006
gene5      50      21       3     197     499      67     375     354      11
gene6     598      35      36       5     277       1     320      23     163
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1      121        1        4       35       71      148      293       11
gene2        7       27      377      317       38       24        3      240
gene3        3        1      160        4       84       40      193        1
gene4       80        4        1      176      196       13        6       26
gene5      716       38      233       30        1       17       44       27
gene6       99        3        1      120      841       57      132      104
      sample18 sample19 sample20
gene1       16      105      359
gene2        9      143        5
gene3       25        1        3
gene4        2       21       89
gene5      211      158        1
gene6      187      317       13

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 20.85488  0.1085909  0.2051840  0.40686804    1
sample2 47.58432  0.3254707 -0.8259428 -0.74382905    2
sample3 39.79753 -2.2383103  0.9718586 -0.04329676    0
sample4 34.40999  0.7158546  1.9093966 -0.33332164    2
sample5 59.29838 -0.1153544  0.7358736  0.30800894    1
sample6 43.62370  0.7276623 -0.9591835  0.29886542    0

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:

  • multiple nonlinear covariates are supported, e.g. design = ~ s(pheno) + s(var1) + var2 + var3 + var4;

  • the nonlinear covariate cannot be a discrete variable, e.g.  design = ~ s(pheno) + var1 + var2 + var3 + s(var4) as var4 is a factor, and it makes no sense to model a factor as nonlinear;

  • at least one nonlinear covariate should be provided in design. If all covariates are assumed to have linear effect on gene count, use DESeq2 (Love, Huber, and Anders 2014), edgeR (Robinson, McCarthy, and Smyth 2010), NBPSeq (Di et al. 2015) or BBSeq (Zhou, Xia, and Wright 2011) instead. e.g.  design = ~ pheno + var1 + var2 + var3 + var4 is not supported in NBAMSeq;

  • design matrix is not supported.

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.

  • gamma argument can be used to control the smoothness of the nonlinear function. Higher gamma means the nonlinear function will be more smooth. See the gamma argument of gam function in mgcv (Wood and Wood 2015) for details. Default gamma is 2.5;

  • fitlin is either TRUE or FALSE indicating whether linear model should be fitted after fitting the nonlinear model;

  • parallel is either TRUE or FALSE indicating whether parallel should be used. e.g. Run NBAMSeq with parallel = TRUE:

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   78.3164   1.00005 4.541112321 0.0330959  0.247227   227.123   234.093
gene2   52.8920   1.00008 2.606172309 0.1064819  0.266205   191.950   198.920
gene3   44.2394   1.00003 1.399510119 0.2368137  0.408299   193.577   200.548
gene4   85.5592   1.00019 0.000438537 0.9876592  0.987659   227.186   234.156
gene5  119.0688   1.00011 2.034182654 0.1538067  0.334362   247.133   254.103
gene6  112.2625   1.00007 1.708511543 0.1912353  0.380210   246.140   253.110

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   78.3164 -0.148164  0.389395 -0.380497 0.7035763  0.818112   227.123
gene2   52.8920 -0.867486  0.404573 -2.144201 0.0320168  0.371130   191.950
gene3   44.2394 -0.209042  0.412816 -0.506379 0.6125906  0.802403   193.577
gene4   85.5592  0.504634  0.460337  1.096229 0.2729786  0.718365   227.186
gene5  119.0688  0.351033  0.414236  0.847422 0.3967597  0.781411   247.133
gene6  112.2625 -0.453094  0.390400 -1.160592 0.2458080  0.682800   246.140
            BIC
      <numeric>
gene1   234.093
gene2   198.920
gene3   200.548
gene4   234.156
gene5   254.103
gene6   253.110

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   78.3164 -0.0168959  0.892656 -0.0189277 0.9848988  0.984899   227.123
gene2   52.8920 -0.4667111  0.917745 -0.5085409 0.6110741  0.787825   191.950
gene3   44.2394 -1.1700651  0.943216 -1.2405060 0.2147883  0.488155   193.577
gene4   85.5592  1.1174465  1.053434  1.0607657 0.2887964  0.501321   227.186
gene5  119.0688  1.1198112  0.947839  1.1814355 0.2374298  0.494645   247.133
gene6  112.2625  1.7345973  0.896471  1.9349177 0.0530004  0.304585   246.140
            BIC
      <numeric>
gene1   234.093
gene2   198.920
gene3   200.548
gene4   234.156
gene5   254.103
gene6   253.110

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
       <numeric> <numeric> <numeric>   <numeric>   <numeric> <numeric>
gene37   68.7427   1.00008  19.55376 1.03566e-05 0.000517832   200.482
gene48   48.9653   1.00008  13.41799 2.49517e-04 0.005270027   202.319
gene14  178.2826   1.00008  12.97303 3.16202e-04 0.005270027   237.479
gene16  162.6701   1.00007   8.01430 4.64286e-03 0.058035693   251.915
gene50   72.0162   1.00004   6.93278 8.46660e-03 0.084666001   214.962
gene1    78.3164   1.00005   4.54111 3.30959e-02 0.247227469   227.123
             BIC
       <numeric>
gene37   207.453
gene48   209.289
gene14   244.449
gene16   258.885
gene50   221.932
gene1    234.093
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.4.1 (2024-06-14)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04 LTS

Matrix products: default
BLAS:   /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3 
LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so;  LAPACK version 3.12.0

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       

time zone: Etc/UTC
tzcode source: system (glibc)

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

other attached packages:
 [1] ggplot2_3.5.1               BiocParallel_1.39.0        
 [3] NBAMSeq_1.21.0              SummarizedExperiment_1.35.1
 [5] Biobase_2.65.0              GenomicRanges_1.57.1       
 [7] GenomeInfoDb_1.41.1         IRanges_2.39.2             
 [9] S4Vectors_0.43.2            BiocGenerics_0.51.0        
[11] MatrixGenerics_1.17.0       matrixStats_1.3.0          
[13] rmarkdown_2.27             

loaded via a namespace (and not attached):
 [1] KEGGREST_1.45.1         gtable_0.3.5            xfun_0.46              
 [4] bslib_0.8.0             lattice_0.22-6          vctrs_0.6.5            
 [7] tools_4.4.1             parallel_4.4.1          RSQLite_2.3.7          
[10] tibble_3.2.1            fansi_1.0.6             AnnotationDbi_1.67.0   
[13] highr_0.11              blob_1.2.4              pkgconfig_2.0.3        
[16] Matrix_1.7-0            lifecycle_1.0.4         GenomeInfoDbData_1.2.12
[19] farver_2.1.2            compiler_4.4.1          Biostrings_2.73.1      
[22] munsell_0.5.1           DESeq2_1.45.3           codetools_0.2-20       
[25] htmltools_0.5.8.1       sys_3.4.2               buildtools_1.0.0       
[28] sass_0.4.9              yaml_2.3.10             pillar_1.9.0           
[31] crayon_1.5.3            jquerylib_0.1.4         DelayedArray_0.31.11   
[34] cachem_1.1.0            abind_1.4-5             nlme_3.1-166           
[37] genefilter_1.87.0       locfit_1.5-9.10         digest_0.6.36          
[40] labeling_0.4.3          splines_4.4.1           maketools_1.3.0        
[43] fastmap_1.2.0           grid_4.4.1              colorspace_2.1-1       
[46] cli_3.6.3               SparseArray_1.5.31      magrittr_2.0.3         
[49] S4Arrays_1.5.7          survival_3.7-0          XML_3.99-0.17          
[52] utf8_1.2.4              withr_3.0.1             scales_1.3.0           
[55] UCSC.utils_1.1.0        bit64_4.0.5             XVector_0.45.0         
[58] httr_1.4.7              bit_4.0.5               png_0.1-8              
[61] memoise_2.0.1           evaluate_0.24.0         knitr_1.48             
[64] mgcv_1.9-1              rlang_1.1.4             Rcpp_1.0.13            
[67] DBI_1.2.3               xtable_1.8-4            glue_1.7.0             
[70] annotate_1.83.0         jsonlite_1.8.8          R6_2.5.1               
[73] zlibbioc_1.51.1        

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–78.