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:

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       5      99      29      29      54      53      28       4
gene2     146     193       7       9       1       1      30     188     124
gene3       3     449      57     270       1       9      21      37       1
gene4     188       4     119       3     103       6       5     214     372
gene5     153     117       4     162     197      37      91      48      14
gene6      98      11     143      49       9      28      41      10     288
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1      223       28       11      121        7        1      102       71
gene2      176       66      171      510        1      283      332        2
gene3        1        1       44      166       16       21      216       88
gene4       64        1        1        1       96       73       10        9
gene5       37       61        1      259       32      199      115       45
gene6      161      219        9        1       20       73      231      172
      sample18 sample19 sample20
gene1        2        8       17
gene2        6      121       16
gene3       53       36        3
gene4        1        1       45
gene5      151        3      133
gene6        1      102       52

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 71.40647  0.01096858  0.02200244 -0.6255840    0
sample2 26.24122 -1.28048917 -1.09283057  2.2348423    1
sample3 59.93755  0.93849016 -0.19313324  1.0771340    1
sample4 56.12270 -0.09659377 -0.48703777  0.2230210    0
sample5 27.83534 -0.49305519 -0.14944696  2.3650277    1
sample6 47.37790  0.09499500 -0.57060670 -0.5317726    0

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:

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.

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   31.1704   1.00007 0.061275363 0.8046351  0.894039   205.150   212.120
gene2   75.8792   1.00007 1.100418272 0.2942220  0.565812   231.732   238.702
gene3   56.4444   1.00004 0.000161856 0.9913795  0.991380   208.627   215.597
gene4   51.1687   1.00006 0.377355790 0.5390663  0.753449   202.698   209.668
gene5   66.7720   1.00006 1.548428381 0.2134010  0.533503   231.035   238.005
gene6   61.3641   1.00007 5.175695662 0.0229159  0.143224   222.447   229.417

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   31.1704  0.174203  0.393255  0.442978 0.6577817  0.925486   205.150
gene2   75.8792 -0.160408  0.451654 -0.355157 0.7224720  0.925486   231.732
gene3   56.4444 -0.560189  0.443100 -1.264250 0.2061401  0.675714   208.627
gene4   51.1687  0.835532  0.483239  1.729022 0.0838052  0.465584   202.698
gene5   66.7720 -0.260274  0.360275 -0.722431 0.4700298  0.783383   231.035
gene6   61.3641 -0.050458  0.375624 -0.134331 0.8931406  0.992378   222.447
            BIC
      <numeric>
gene1   212.120
gene2   238.702
gene3   215.597
gene4   209.668
gene5   238.005
gene6   229.417

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   31.1704  0.0548536  0.836950  0.0655399 0.94774412 0.9934832   205.150
gene2   75.8792 -0.1630824  0.962504 -0.1694355 0.86545413 0.9626678   231.732
gene3   56.4444 -0.8788901  0.945953 -0.9291057 0.35283430 0.8525097   208.627
gene4   51.1687  1.8729886  1.030708  1.8171860 0.06918864 0.3797730   202.698
gene5   66.7720 -1.5863456  0.771717 -2.0556058 0.03982053 0.2844324   231.035
gene6   61.3641  2.2053204  0.806768  2.7335259 0.00626602 0.0626602   222.447
            BIC
      <numeric>
gene1   212.120
gene2   238.702
gene3   215.597
gene4   209.668
gene5   238.005
gene6   229.417

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>
gene45   93.0939   1.00006  19.59653 1.00781e-05 0.000503905   217.146
gene36   74.2207   1.00004   8.01519 4.63939e-03 0.079630738   222.650
gene40   58.6414   1.00014   7.96324 4.77784e-03 0.079630738   200.458
gene25  107.3228   1.00004   7.17727 7.38519e-03 0.092314858   243.709
gene26   33.1684   1.00009   6.39151 1.14698e-02 0.112030861   195.240
gene20   91.5436   1.00008   6.11086 1.34437e-02 0.112030861   231.091
             BIC
       <numeric>
gene45   224.116
gene36   229.621
gene40   207.429
gene25   250.679
gene26   202.211
gene20   238.061
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 RC (2026-04-17 r89917)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.4 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.23-bioc/R/lib/libRblas.so 
LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.0

locale:
 [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
 [3] LC_TIME=en_GB              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: America/New_York
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          

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       tibble_3.3.1        
[10] AnnotationDbi_1.74.0 RSQLite_2.4.6        blob_1.3.0          
[13] pkgconfig_2.0.3      Matrix_1.7-5         RColorBrewer_1.1-3  
[16] S7_0.2.2             lifecycle_1.0.5      compiler_4.6.0      
[19] farver_2.1.2         Biostrings_2.80.0    DESeq2_1.52.0       
[22] codetools_0.2-20     htmltools_0.5.9      sass_0.4.10         
[25] yaml_2.3.12          crayon_1.5.3         pillar_1.11.1       
[28] jquerylib_0.1.4      DelayedArray_0.38.0  cachem_1.1.0        
[31] abind_1.4-8          nlme_3.1-169         genefilter_1.94.0   
[34] tidyselect_1.2.1     locfit_1.5-9.12      digest_0.6.39       
[37] dplyr_1.2.1          labeling_0.4.3       splines_4.6.0       
[40] fastmap_1.2.0        grid_4.6.0           cli_3.6.6           
[43] SparseArray_1.12.0   magrittr_2.0.5       S4Arrays_1.12.0     
[46] survival_3.8-6       dichromat_2.0-0.1    XML_3.99-0.23       
[49] withr_3.0.2          scales_1.4.0         bit64_4.8.0         
[52] rmarkdown_2.31       XVector_0.52.0       httr_1.4.8          
[55] bit_4.6.0            otel_0.2.0           png_0.1-9           
[58] memoise_2.0.1        evaluate_1.0.5       knitr_1.51          
[61] mgcv_1.9-4           rlang_1.2.0          Rcpp_1.1.1-1.1      
[64] xtable_1.8-8         glue_1.8.1           DBI_1.3.0           
[67] 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–8.