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     303     102      10      23     271      20      11     116    1390
gene2      31     128       1      16     104      38     146     313    1190
gene3     132     253       2     185      75      25       2       2       7
gene4       5      47      88      12      18      25      27       9     113
gene5      16       9     254      99      14     306     258      64      39
gene6      12      83      56     136       1      16       1       5     348
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1       10        2        2      161        1      184     1188       19
gene2        3        1      277       10      406        7       92      105
gene3       57        2        1       25      120        1      343       22
gene4       63      196       29      177        6       61        1       29
gene5       10       54      323      148        5        3        1       50
gene6        2      204       81      121        3       12      254      103
      sample18 sample19 sample20
gene1       50       94     1996
gene2        2       43        1
gene3       50      426       31
gene4        1       63       25
gene5       72       61       31
gene6       17       56      411

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 54.99922  0.15374288  0.3285383 -0.5506622    1
sample2 65.15015 -0.84206363  0.1730560 -0.7765419    2
sample3 78.23617 -2.39732742  0.4769664 -0.2584027    2
sample4 53.16097 -0.95691514  1.5124306 -1.0397088    2
sample5 69.60082 -1.27269699 -0.1918755 -0.7732594    1
sample6 45.42110 -0.09838307 -0.2462131  0.2219302    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:

  • 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  235.6258   1.00010 3.9673686 0.0463983  0.289990   252.370   259.340
gene2  100.7962   1.00027 0.0161532 0.8998364  0.970936   233.735   240.705
gene3   77.4360   1.00010 2.2637231 0.1324445  0.509980   212.488   219.458
gene4   38.7618   1.00006 1.2670842 0.2603314  0.723143   207.681   214.651
gene5   70.8317   1.00006 0.2446006 0.6209807  0.887115   228.250   235.220
gene6   71.4947   1.00009 6.0517736 0.0138985  0.138985   212.798   219.769

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  235.6258 -0.290251  0.543785 -0.533761  0.593507  0.909959   252.370
gene2  100.7962 -0.317272  0.563806 -0.562731  0.573618  0.909959   233.735
gene3   77.4360  0.321846  0.491651  0.654622  0.512711  0.909959   212.488
gene4   38.7618  0.426054  0.440947  0.966225  0.333932  0.909959   207.681
gene5   70.8317 -0.499894  0.462914 -1.079885  0.280193  0.909959   228.250
gene6   71.4947 -0.433687  0.422043 -1.027589  0.304143  0.909959   212.798
            BIC
      <numeric>
gene1   259.340
gene2   240.705
gene3   219.458
gene4   214.651
gene5   235.220
gene6   219.769

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  235.6258  0.917232   1.45150  0.631921 5.27438e-01 0.71275472   252.370
gene2  100.7962 -1.401967   1.49829 -0.935710 3.49423e-01 0.62396876   233.735
gene3   77.4360  5.522962   1.40094  3.942327 8.06947e-05 0.00403474   212.488
gene4   38.7618  0.307070   1.17263  0.261865 7.93425e-01 0.82964162   207.681
gene5   70.8317 -1.829585   1.22914 -1.488507 1.36617e-01 0.49683558   228.250
gene6   71.4947  0.995074   1.13658  0.875501 3.81301e-01 0.65741590   212.798
            BIC
      <numeric>
gene1   259.340
gene2   240.705
gene3   219.458
gene4   214.651
gene5   235.220
gene6   219.769

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>
gene16   90.7758   1.00005  10.29983 0.00133090 0.0537591   218.703   225.673
gene11   93.5924   1.00049   9.42237 0.00215037 0.0537591   224.683   231.654
gene20   33.0865   1.00003   8.07279 0.00449446 0.0749077   181.216   188.186
gene18   75.5262   1.00007   6.72130 0.00952906 0.1191133   223.518   230.488
gene6    71.4947   1.00009   6.05177 0.01389853 0.1389853   212.798   219.769
gene41  110.2157   1.00008   5.28241 0.02155082 0.1795901   219.535   226.505
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.5.1 (2025-06-13)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.3 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_4.0.0               BiocParallel_1.45.0        
 [3] NBAMSeq_1.27.0              SummarizedExperiment_1.41.0
 [5] Biobase_2.71.0              GenomicRanges_1.63.0       
 [7] Seqinfo_1.1.0               IRanges_2.45.0             
 [9] S4Vectors_0.49.0            BiocGenerics_0.57.0        
[11] generics_0.1.4              MatrixGenerics_1.23.0      
[13] matrixStats_1.5.0           rmarkdown_2.30             

loaded via a namespace (and not attached):
 [1] KEGGREST_1.49.2      gtable_0.3.6         xfun_0.54           
 [4] bslib_0.9.0          lattice_0.22-7       vctrs_0.6.5         
 [7] tools_4.5.1          parallel_4.5.1       AnnotationDbi_1.73.0
[10] RSQLite_2.4.3        blob_1.2.4           Matrix_1.7-4        
[13] RColorBrewer_1.1-3   S7_0.2.0             lifecycle_1.0.4     
[16] compiler_4.5.1       farver_2.1.2         Biostrings_2.79.1   
[19] DESeq2_1.51.0        codetools_0.2-20     htmltools_0.5.8.1   
[22] sys_3.4.3            buildtools_1.0.0     sass_0.4.10         
[25] yaml_2.3.10          crayon_1.5.3         jquerylib_0.1.4     
[28] DelayedArray_0.37.0  cachem_1.1.0         abind_1.4-8         
[31] nlme_3.1-168         genefilter_1.91.0    locfit_1.5-9.12     
[34] digest_0.6.37        labeling_0.4.3       splines_4.5.1       
[37] maketools_1.3.2      fastmap_1.2.0        grid_4.5.1          
[40] cli_3.6.5            SparseArray_1.11.1   S4Arrays_1.11.0     
[43] survival_3.8-3       XML_3.99-0.19        withr_3.0.2         
[46] scales_1.4.0         bit64_4.6.0-1        XVector_0.51.0      
[49] httr_1.4.7           bit_4.6.0            png_0.1-8           
[52] memoise_2.0.1        evaluate_1.0.5       knitr_1.50          
[55] mgcv_1.9-3           rlang_1.1.6          Rcpp_1.1.0          
[58] xtable_1.8-4         glue_1.8.0           DBI_1.2.3           
[61] annotate_1.89.0      jsonlite_2.0.0       R6_2.6.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.