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      57     728     130     667     319      11      17       1      71
gene2       5       1     578      48     415      22     725     360       1
gene3       1       9      40      32      75     172      15     305       7
gene4       5       6      63      94       1     179      40     326      21
gene5      30       3       1     199       1       1       1      15      29
gene6      96      92     200       6      38       1      19       1     432
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1        3      468      295       54      248       16       42       10
gene2        4        1       45       24       25       41       19       19
gene3      127       33       19        2       12       59        1       16
gene4      306       13       73      137       48       99       20       24
gene5        1        6        1      282       81      309        1      131
gene6       35        8       13        1      235       49       21        5
      sample18 sample19 sample20
gene1       21       17      186
gene2       86        1        1
gene3       24      492        6
gene4       11      114       32
gene5      145        3      185
gene6       43      397       12

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.59802 -1.2992402  1.321249458 -0.3276219    0
sample2 26.98979  0.8984400  1.685156923 -1.3486892    1
sample3 22.58746 -0.3308876 -0.109758906 -1.8210088    0
sample4 75.00540 -0.2418883 -0.373472482  1.7404782    2
sample5 41.16414  0.1912184  0.009181675  0.5078990    1
sample6 78.99259 -0.2156855  0.820151155  0.6006739    2

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  109.7287   1.00006  1.9089112 0.167098279 0.4082453   246.968   253.938
gene2   96.4821   1.00009  0.0175710 0.894845864 0.9723407   218.202   225.172
gene3   65.2953   1.00007  0.0984228 0.753778294 0.9332676   213.156   220.127
gene4   62.5342   1.00008  0.6989357 0.403205869 0.6300092   226.316   233.286
gene5   57.9446   1.00018  3.0509096 0.080720651 0.3363360   198.373   205.343
gene6   78.9005   1.00047 11.3454542 0.000756459 0.0189115   208.806   215.776

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  109.7287  0.458575  0.358318  1.279799 0.2006159  0.496097   246.968
gene2   96.4821 -0.330121  0.420163 -0.785698 0.4320445  0.603215   218.202
gene3   65.2953  0.783219  0.364765  2.147189 0.0317782  0.372851   213.156
gene4   62.5342  0.342598  0.324920  1.054407 0.2916965  0.502245   226.316
gene5   57.9446 -0.804498  0.416495 -1.931590 0.0534101  0.381501   198.373
gene6   78.9005  0.392388  0.333488  1.176619 0.2393474  0.496097   208.806
            BIC
      <numeric>
gene1   253.938
gene2   225.172
gene3   220.127
gene4   233.286
gene5   205.343
gene6   215.776

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  109.7287  1.730610  0.958840  1.804900 0.0710904  0.444315   246.968
gene2   96.4821 -1.211532  1.126545 -1.075441 0.2821775  0.755018   218.202
gene3   65.2953 -0.486688  0.976726 -0.498285 0.6182833  0.792671   213.156
gene4   62.5342 -0.914234  0.870360 -1.050410 0.2935298  0.755018   226.316
gene5   57.9446  1.887251  1.097800  1.719121 0.0855924  0.466052   198.373
gene6   78.9005 -1.257213  0.900996 -1.395359 0.1629076  0.581813   208.806
            BIC
      <numeric>
gene1   253.938
gene2   225.172
gene3   220.127
gene4   233.286
gene5   205.343
gene6   215.776

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>
gene50   50.8881   1.00013  14.49045 0.000141681 0.00708404   200.590   207.561
gene6    78.9005   1.00047  11.34545 0.000756459 0.01891148   208.806   215.776
gene20   39.3753   1.00012   9.02666 0.002662254 0.04437089   196.025   202.995
gene31   81.9116   1.00011   6.55495 0.010465612 0.12776229   211.764   218.734
gene17   47.5932   1.00008   6.20092 0.012776229 0.12776229   203.450   210.420
gene40   61.0485   1.00008   3.96505 0.046467852 0.26892234   215.723   222.694
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.2 (2025-10-31)
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.1               BiocParallel_1.45.0        
 [3] NBAMSeq_1.27.0              SummarizedExperiment_1.41.0
 [5] Biobase_2.71.0              GenomicRanges_1.63.1       
 [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.51.1      gtable_0.3.6         xfun_0.56           
 [4] bslib_0.9.0          lattice_0.22-7       vctrs_0.7.1         
 [7] tools_4.5.2          parallel_4.5.2       AnnotationDbi_1.73.0
[10] RSQLite_2.4.5        blob_1.3.0           Matrix_1.7-4        
[13] RColorBrewer_1.1-3   S7_0.2.1             lifecycle_1.0.5     
[16] compiler_4.5.2       farver_2.1.2         Biostrings_2.79.4   
[19] DESeq2_1.51.6        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.37.0  cachem_1.1.0         abind_1.4-8         
[31] nlme_3.1-168         genefilter_1.93.0    locfit_1.5-9.12     
[34] digest_0.6.39        labeling_0.4.3       splines_4.5.2       
[37] maketools_1.3.2      fastmap_1.2.0        grid_4.5.2          
[40] cli_3.6.5            SparseArray_1.11.10  S4Arrays_1.11.1     
[43] survival_3.8-6       XML_3.99-0.20        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.51          
[55] mgcv_1.9-4           rlang_1.1.7          Rcpp_1.1.1          
[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.