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     269     133      66       1      28     314     130       3       1
gene2       2      49     650       2     106       1       7      26     110
gene3      41     159       1    1156      51      10      61       4       1
gene4       7       1     170       1       2     361      56       7     520
gene5     192      36      98       1      80       6       1     557       1
gene6      42      53       2      19     329      93       8       2       1
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1        1       32       60      300       66        3       27      125
gene2       14       74      335        5       16       27      113        1
gene3      116        5       31        1        1      521      541       62
gene4        1      103        4       23      111      192        2       16
gene5        3      235      487        5        1      327       55       97
gene6        6      385       69      215      727      604        1        1
      sample18 sample19 sample20
gene1        6       67       16
gene2       53       18      666
gene3       13      115       19
gene4        2       79      116
gene5      128        1        1
gene6        4      211        3

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 27.21477  0.7638499  1.8793190 -0.5151936    1
sample2 69.83544  0.1505509  2.1620222 -0.8835014    1
sample3 66.05655  0.1320033 -1.2333532  0.5249466    1
sample4 48.30346 -0.2549355  0.6483301  2.1115087    2
sample5 26.82218  1.2718286  0.5079491 -0.3745474    0
sample6 41.92866  0.9454574 -0.3607100 -0.6655164    2

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   59.8501   1.00006  2.490994 0.1145201  0.396779   215.451   222.422
gene2  103.8130   1.00010  2.589710 0.1075949  0.396779   223.353   230.323
gene3  147.8854   1.00005  5.513094 0.0188796  0.134854   220.537   227.507
gene4   79.1549   1.00003  2.269795 0.1319223  0.396779   219.182   226.153
gene5  105.8756   1.00009  0.616135 0.4325418  0.786433   220.280   227.250
gene6  102.2639   1.00010  0.530771 0.4663447  0.786433   228.658   235.628

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   59.8501  0.487324  0.332255  1.466717 0.1424530  0.339174   215.451
gene2  103.8130 -0.647702  0.374029 -1.731690 0.0833289  0.339174   223.353
gene3  147.8854  0.329537  0.363807  0.905803 0.3650402  0.550750   220.537
gene4   79.1549 -0.361155  0.410629 -0.879517 0.3791208  0.550750   219.182
gene5  105.8756 -0.558978  0.406086 -1.376501 0.1686667  0.375700   220.280
gene6  102.2639 -0.889558  0.418959 -2.123259 0.0337322  0.296966   228.658
            BIC
      <numeric>
gene1   222.422
gene2   230.323
gene3   227.507
gene4   226.153
gene5   227.250
gene6   235.628

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   59.8501  1.110497  0.984359  1.128141 0.259260219 0.6822637   215.451
gene2  103.8130 -3.779279  1.102106 -3.429142 0.000605493 0.0192984   223.353
gene3  147.8854 -1.880281  1.056893 -1.779064 0.075229227 0.4179402   220.537
gene4   79.1549  0.385887  1.205686  0.320056 0.748925772 0.9017678   219.182
gene5  105.8756 -1.135945  1.191846 -0.953097 0.340541074 0.7328485   220.280
gene6  102.2639 -2.036401  1.230133 -1.655431 0.097837119 0.4891856   228.658
            BIC
      <numeric>
gene1   222.422
gene2   230.323
gene3   227.507
gene4   226.153
gene5   227.250
gene6   235.628

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>
gene10  140.3054   1.00011  10.25348 0.00136508 0.0682539   225.194   232.165
gene37   46.6320   1.00005   7.03048 0.00801560 0.1348542   199.554   206.524
gene8    32.4476   1.00004   5.90227 0.01512609 0.1348542   185.531   192.501
gene48   36.0030   1.79379   8.24554 0.01525087 0.1348542   169.344   177.104
gene16   68.4084   1.00004   5.81031 0.01593599 0.1348542   202.318   209.288
gene39   61.6792   1.00035   5.72256 0.01678315 0.1348542   203.959   210.929
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.46.0        
 [3] NBAMSeq_1.28.0              SummarizedExperiment_1.42.0
 [5] Biobase_2.72.0              GenomicRanges_1.64.0       
 [7] Seqinfo_1.2.0               IRanges_2.46.0             
 [9] S4Vectors_0.50.1            BiocGenerics_0.58.1        
[11] generics_0.1.4              MatrixGenerics_1.24.0      
[13] matrixStats_1.5.0           rmarkdown_2.31             

loaded via a namespace (and not attached):
 [1] KEGGREST_1.52.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.74.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.80.1   
[19] DESeq2_1.52.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.38.2  cachem_1.1.0         abind_1.4-8         
[31] nlme_3.1-169         genefilter_1.94.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.12.2   S4Arrays_1.12.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.52.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.90.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.