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      23       3      42     340       9      22      11      17      51
gene2     937     546    1355       1      81      23      19     691       2
gene3       1       3      61     100      29       4      17      27     105
gene4      10     348     111      33     312     110      60      56       1
gene5      24       1      38       3       2     129      56       5      78
gene6       2      47      38       1     290       5       1      22      90
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1       70       84       13      241      110       34        3        1
gene2       61      180       92      131       10       27        1        1
gene3       43      273        1      273        2        1       16       23
gene4       54        2      892       32      303      496      140       94
gene5      198        2      191        1       11      164       95       62
gene6      123       27      131      226       36       23        1      515
      sample18 sample19 sample20
gene1        4      143       99
gene2      425       13       56
gene3        1        1        2
gene4      280      179      289
gene5      406       29        1
gene6       46      156       18

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 69.04681  1.1134977  1.2191987 -0.39846583    2
sample2 33.78697  0.3604733 -0.8851305  0.39951712    1
sample3 53.21970  1.2106192  0.9559055 -0.14925302    1
sample4 32.95667  1.3124109 -1.9681201 -0.07401552    1
sample5 21.25676 -1.7294309 -1.3545982  0.52285935    0
sample6 52.93654  0.5127478 -1.3759760  1.51451944    1

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   50.4485   1.00005  0.284171 0.5940713  0.750214   217.252   224.222
gene2  172.4241   1.00014  0.274869 0.6001714  0.750214   246.633   253.604
gene3   39.1822   1.00006  0.925402 0.3360977  0.579479   192.496   199.466
gene4  139.0595   1.00007  3.274000 0.0704111  0.220035   260.723   267.693
gene5   56.1488   1.00004  2.787454 0.0950168  0.250044   213.227   220.197
gene6  102.0912   1.00006  0.654775 0.4184617  0.653846   219.275   226.245

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   50.4485  0.8718479  0.515219  1.6921889  0.090610  0.411863   217.252
gene2  172.4241 -0.3139281  0.594388 -0.5281533  0.597393  0.846463   246.633
gene3   39.1822 -0.0243157  0.588711 -0.0413032  0.967054  0.979182   192.496
gene4  139.0595 -0.1602507  0.509193 -0.3147148  0.752978  0.855657   260.723
gene5   56.1488 -0.8516221  0.588689 -1.4466410  0.147997  0.437105   213.227
gene6  102.0912 -0.3991684  0.511003 -0.7811468  0.434716  0.724527   219.275
            BIC
      <numeric>
gene1   224.222
gene2   253.604
gene3   199.466
gene4   267.693
gene5   220.197
gene6   226.245

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   50.4485  0.0834528   1.07843  0.0773839 0.93831815 0.9774147   217.252
gene2  172.4241  1.5560327   1.24020  1.2546581 0.20960292 0.5254011   246.633
gene3   39.1822 -3.2309866   1.23216 -2.6222090 0.00873618 0.0581765   192.496
gene4  139.0595 -0.4283711   1.06255 -0.4031550 0.68683419 0.9242040   260.723
gene5   56.1488  1.5052023   1.16383  1.2933196 0.19590055 0.5254011   213.227
gene6  102.0912 -2.4778750   1.06038 -2.3367845 0.01945039 0.1080577   219.275
            BIC
      <numeric>
gene1   224.222
gene2   253.604
gene3   199.466
gene4   267.693
gene5   220.197
gene6   226.245

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>
gene13   50.4506   1.00005  19.63187 9.94847e-06 0.000497424   201.880
gene37  107.5952   1.00007  12.16388 4.87476e-04 0.012186909   225.613
gene23   48.3546   1.00004   9.17235 2.45803e-03 0.040967109   183.316
gene7   120.8753   1.00007   7.29218 6.92927e-03 0.082501068   240.491
gene17   85.6520   1.00031   6.68189 9.76520e-03 0.082501068   226.828
gene35  106.8198   1.00013   6.52658 1.06348e-02 0.082501068   213.337
             BIC
       <numeric>
gene13   208.850
gene37   232.584
gene23   190.286
gene7    247.461
gene17   233.799
gene35   220.307
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.0            BiocGenerics_0.58.0        
[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.10.0         lattice_0.22-9       vctrs_0.7.3         
 [7] tools_4.6.0          parallel_4.6.0       AnnotationDbi_1.74.0
[10] RSQLite_2.4.6        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.0   
[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.1  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.0   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.0          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.