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      37       3     233      91      47      98       1      31     233
gene2     173       1       4       4     145       1       5      36      22
gene3      13      35     280       1       5     397       2      22       1
gene4      51      36       3       1       3     230     259       3      54
gene5       1      25       1      43     289       7      36      32     309
gene6       6      54      95      41       4      18      97      64      64
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1      113        9        9        1      169       29      301        9
gene2      126       62        3      186       31      237       47        1
gene3       91        5       48       30      490      120       72        2
gene4       55        5       62       15       21      379        6      618
gene5       34       17       19     1305      122        9        2      138
gene6        1       21        2        1        2       34       44      354
      sample18 sample19 sample20
gene1      145        3      623
gene2        1       10       13
gene3      146       12        8
gene4        1      402        2
gene5       11      224        1
gene6      558       10      835

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 59.78855 -0.1130351 -1.1258588 -2.0543204    1
sample2 40.58325  0.8857501 -1.0906763  0.4939448    0
sample3 30.10439  0.1833921  0.6427279 -0.1855790    0
sample4 41.14597  0.9050350  0.8845010  0.6910609    0
sample5 65.49492  0.4860724  1.4264302  1.0496904    2
sample6 61.12732  1.1176470  2.2566607  1.8188325    1

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   80.7274   1.00010 0.00226462 0.96291502 0.9825663   226.622   233.592
gene2   46.7694   1.00008 0.36676897 0.54486492 0.8681620   190.081   197.052
gene3   76.8399   1.00008 3.34815959 0.06728932 0.3364466   211.791   218.761
gene4  128.9087   1.13274 0.30151277 0.79628358 0.8684552   228.077   235.179
gene5  130.7393   1.00007 6.69544902 0.00966671 0.0690479   222.431   229.401
gene6   96.1041   1.00009 0.25337427 0.61473965 0.8681620   221.068   228.038

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   80.7274 -0.592137  0.634486 -0.933254  0.350689  0.892943   226.622
gene2   46.7694 -0.313193  0.596043 -0.525453  0.599268  0.907982   190.081
gene3   76.8399 -0.396409  0.645605 -0.614012  0.539208  0.892943   211.791
gene4  128.9087 -1.081177  0.809460 -1.335677  0.181655  0.684800   228.077
gene5  130.7393  0.220314  0.705447  0.312305  0.754809  0.916654   222.431
gene6   96.1041 -0.909008  0.687483 -1.322226  0.186093  0.684800   221.068
            BIC
      <numeric>
gene1   233.592
gene2   197.052
gene3   218.761
gene4   235.179
gene5   229.401
gene6   228.038

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   80.7274  0.289920   1.24261  0.233317 0.8155156  0.855558   226.622
gene2   46.7694  2.462774   1.16967  2.105526 0.0352455  0.440569   190.081
gene3   76.8399 -2.334418   1.27263 -1.834325 0.0666057  0.555048   211.791
gene4  128.9087 -0.913329   1.57375 -0.580352 0.5616772  0.855558   228.077
gene5  130.7393  0.886640   1.37853  0.643176 0.5201096  0.855558   222.431
gene6   96.1041 -0.348018   1.34582 -0.258591 0.7959509  0.855558   221.068
            BIC
      <numeric>
gene1   233.592
gene2   197.052
gene3   218.761
gene4   235.179
gene5   229.401
gene6   228.038

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>
gene11  119.6635   1.00006  12.71149 0.000363407 0.0181704   225.345   232.315
gene31   96.6180   1.00016   8.95796 0.002764802 0.0429987   223.666   230.637
gene8    85.2700   1.00010   8.65755 0.003258310 0.0429987   215.366   222.336
gene29  125.8419   1.00009   8.55863 0.003439892 0.0429987   225.853   232.823
gene10  114.2460   1.00007   7.84439 0.005098682 0.0467697   235.984   242.954
gene26   94.0304   1.00008   7.67128 0.005612363 0.0467697   225.340   232.311
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.1 Patched (2025-08-23 r88802)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.3 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.22-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.0               BiocParallel_1.43.4        
 [3] NBAMSeq_1.25.0              SummarizedExperiment_1.39.2
 [5] Biobase_2.69.1              GenomicRanges_1.61.5       
 [7] Seqinfo_0.99.2              IRanges_2.43.5             
 [9] S4Vectors_0.47.4            BiocGenerics_0.55.1        
[11] generics_0.1.4              MatrixGenerics_1.21.0      
[13] matrixStats_1.5.0          

loaded via a namespace (and not attached):
 [1] KEGGREST_1.49.1      gtable_0.3.6         xfun_0.53           
 [4] bslib_0.9.0          lattice_0.22-7       vctrs_0.6.5         
 [7] tools_4.5.1          parallel_4.5.1       tibble_3.3.0        
[10] AnnotationDbi_1.71.1 RSQLite_2.4.3        blob_1.2.4          
[13] pkgconfig_2.0.3      Matrix_1.7-4         RColorBrewer_1.1-3  
[16] S7_0.2.0             lifecycle_1.0.4      compiler_4.5.1      
[19] farver_2.1.2         Biostrings_2.77.2    DESeq2_1.49.4       
[22] codetools_0.2-20     htmltools_0.5.8.1    sass_0.4.10         
[25] yaml_2.3.10          pillar_1.11.1        crayon_1.5.3        
[28] jquerylib_0.1.4      DelayedArray_0.35.3  cachem_1.1.0        
[31] abind_1.4-8          nlme_3.1-168         genefilter_1.91.0   
[34] tidyselect_1.2.1     locfit_1.5-9.12      digest_0.6.37       
[37] dplyr_1.1.4          labeling_0.4.3       splines_4.5.1       
[40] fastmap_1.2.0        grid_4.5.1           cli_3.6.5           
[43] SparseArray_1.9.1    magrittr_2.0.4       S4Arrays_1.9.1      
[46] survival_3.8-3       dichromat_2.0-0.1    XML_3.99-0.19       
[49] withr_3.0.2          scales_1.4.0         bit64_4.6.0-1       
[52] rmarkdown_2.30       XVector_0.49.1       httr_1.4.7          
[55] bit_4.6.0            png_0.1-8            memoise_2.0.1       
[58] evaluate_1.0.5       knitr_1.50           mgcv_1.9-3          
[61] rlang_1.1.6          Rcpp_1.1.0           xtable_1.8-4        
[64] glue_1.8.0           DBI_1.2.3            annotate_1.87.0     
[67] 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.