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      20      14       3     107      35       5     587      12       3
gene2       8       3     154       5     210      11       2       1      18
gene3       2       2      95      30      33      66      69       7       1
gene4      64      10     109      14     237     196       2      90       2
gene5      80       3       1       5      14      14       1      55      64
gene6       1      73      32     137      45      21       1       2       5
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1       18      159       23        2        9        1       29      116
gene2      250        1       23        1      354        6        1      194
gene3      137       28      127      614       24       95        2       18
gene4        8      115      146       94       85      198      355      189
gene5      139       24      197        1       12       60       68        3
gene6        2        2       48        8        8      499      125       53
      sample18 sample19 sample20
gene1      433        3       26
gene2      794        4        1
gene3       54       39      214
gene4      501        1       72
gene5        6        3        1
gene6       95      120      173

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 32.47681 -0.3393877 -0.86450450 -0.4914680    2
sample2 34.42932 -0.6072906  0.85222425  0.4293722    0
sample3 79.68896  1.3679244 -0.07571716  0.5567646    2
sample4 45.95901  1.9026196 -1.68535593  1.5158835    0
sample5 62.60046 -0.1886570  0.01124564 -0.4382059    2
sample6 43.03864  0.0831485  1.01891319 -1.4635382    0

design is a formula which specifies how to model the samples. Compared with other packages performing DE analysis including DESeq2 (Love et al. 2014), edgeR (Robinson et al. 2010), NBPSeq (Di et al. 2015) and BBSeq (Zhou et al. 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 et al. 2014), edgeR (Robinson et al. 2010), NBPSeq (Di et al. 2015) or BBSeq (Zhou et al. 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   77.1933   1.00008  0.4573032 0.498885227 0.7976577   213.515   220.485
gene2   96.9756   1.00005 15.2087127 0.000095752 0.0047876   200.830   207.800
gene3   59.7984   1.00007  0.0275812 0.868228668 0.9391381   213.324   220.294
gene4  107.9808   1.00006  4.6040361 0.031901604 0.3987700   238.006   244.976
gene5   34.4374   1.00007  0.0646142 0.799503143 0.9296548   192.983   199.953
gene6   60.2583   1.00005  3.7148929 0.053933141 0.5046917   202.575   209.545

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
      <numeric>  <numeric> <numeric> <numeric>   <numeric>   <numeric>
gene1   77.1933 -0.0396853  0.388279 -0.102208 9.18592e-01 0.937338307
gene2   96.9756 -1.7756695  0.415250 -4.276145 1.90158e-05 0.000950789
gene3   59.7984  0.7941877  0.312254  2.543406 1.09778e-02 0.091481288
gene4  107.9808 -0.3742565  0.324668 -1.152737 2.49018e-01 0.655311208
gene5   34.4374 -0.4695188  0.391844 -1.198228 2.30828e-01 0.655311208
gene6   60.2583  0.2024720  0.332061  0.609743 5.42032e-01 0.875250279
            AIC       BIC
      <numeric> <numeric>
gene1   213.515   220.485
gene2   200.830   207.800
gene3   213.324   220.294
gene4   238.006   244.976
gene5   192.983   199.953
gene6   202.575   209.545

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   77.1933 -0.536345  0.825479 -0.649738  0.515861  0.679161   213.515
gene2   96.9756  0.350364  0.871408  0.402067  0.687635  0.781403   200.830
gene3   59.7984  0.520310  0.665959  0.781294  0.434629  0.679161   213.324
gene4  107.9808  0.527289  0.687685  0.766759  0.443225  0.679161   238.006
gene5   34.4374 -0.294672  0.830506 -0.354810  0.722732  0.803035   192.983
gene6   60.2583  0.680824  0.673514  1.010854  0.312086  0.679161   202.575
            BIC
      <numeric>
gene1   220.485
gene2   207.800
gene3   220.294
gene4   244.976
gene5   199.953
gene6   209.545

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>
gene2    96.9756   1.00005  15.20871 0.000095752 0.00478760   200.830   207.800
gene29   59.2691   1.00006  12.91414 0.000326255 0.00815638   210.081   217.051
gene35   86.6472   1.00006   6.53827 0.010559752 0.17599587   224.740   231.710
gene4   107.9808   1.00006   4.60404 0.031901604 0.39877005   238.006   244.976
gene6    60.2583   1.00005   3.71489 0.053933141 0.50469171   202.575   209.545
gene50   76.0783   1.00009   3.33838 0.067706578 0.50469171   228.171   235.141
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.6.0 (2026-04-24)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.4 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=en_US.UTF-8    
 [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.3               BiocParallel_1.47.0        
 [3] NBAMSeq_1.29.0              SummarizedExperiment_1.43.0
 [5] Biobase_2.73.1              GenomicRanges_1.65.0       
 [7] Seqinfo_1.3.0               IRanges_2.47.2             
 [9] S4Vectors_0.51.3            BiocGenerics_0.59.6        
[11] generics_0.1.4              MatrixGenerics_1.25.0      
[13] matrixStats_1.5.0           rmarkdown_2.31             

loaded via a namespace (and not attached):
 [1] KEGGREST_1.53.0      gtable_0.3.6         xfun_0.57           
 [4] bslib_0.11.0         lattice_0.22-9       vctrs_0.7.3         
 [7] tools_4.6.0          parallel_4.6.0       AnnotationDbi_1.75.0
[10] RSQLite_3.53.1       blob_1.3.0           Matrix_1.7-5        
[13] RColorBrewer_1.1-3   S7_0.2.2             lifecycle_1.0.5     
[16] compiler_4.6.0       farver_2.1.2         Biostrings_2.81.2   
[19] DESeq2_1.53.0        codetools_0.2-20     htmltools_0.5.9     
[22] sys_3.4.3            buildtools_1.0.0     sass_0.4.10         
[25] yaml_2.3.12          crayon_1.5.3         jquerylib_0.1.4     
[28] DelayedArray_0.39.3  cachem_1.1.0         abind_1.4-8         
[31] nlme_3.1-169         genefilter_1.95.0    locfit_1.5-9.12     
[34] digest_0.6.39        labeling_0.4.3       splines_4.6.0       
[37] maketools_1.3.2      fastmap_1.2.0        grid_4.6.0          
[40] cli_3.6.6            SparseArray_1.13.2   S4Arrays_1.13.0     
[43] survival_3.8-6       XML_3.99-0.23        withr_3.0.2         
[46] scales_1.4.0         bit64_4.8.2          XVector_0.53.0      
[49] httr_1.4.8           bit_4.6.0            png_0.1-9           
[52] memoise_2.0.1        evaluate_1.0.5       knitr_1.51          
[55] mgcv_1.9-4           rlang_1.2.0          Rcpp_1.1.1-1.1      
[58] xtable_1.8-8         glue_1.8.1           DBI_1.3.0           
[61] annotate_1.91.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.