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     185     224      21       1     208     455       8      57       9
gene2      15     255       1       2       6     124      49     481       1
gene3      21       9     359      10       3       2       9      43       3
gene4     154     124       1     116     152     357     146     284      56
gene5     179       2       1     294      44      12     201      12      11
gene6       2     226     160      19     329      46     142      26      21
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1      141      186       68        1        2        1      132      375
gene2       11        1      175      200       32       18       52        2
gene3       19      267       96        7        1      723       58      224
gene4       65        5        1        1        4       29        2        1
gene5       91      459       44      741        8       11        1        3
gene6      344        1       19       51        1        3      747       80
      sample18 sample19 sample20
gene1       67       69      116
gene2      227        1        4
gene3       19      323      156
gene4      942        4       90
gene5        4      307        2
gene6      228      136        6

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 51.05974  0.1796586 -0.4944812 -0.2978065    1
sample2 60.28843  0.4860911 -0.9966152  1.2176327    1
sample3 41.76345  1.1783782 -0.7291128  1.4032860    0
sample4 63.63253 -0.2544556  1.1673852  2.4063905    0
sample5 21.49585  1.0871780  0.0645908 -1.2198269    1
sample6 44.28833 -1.6371411 -0.4868558 -0.2289371    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   94.1884   1.00007  0.069150 0.79262623  0.910036   229.534   236.504
gene2   58.8858   1.00004  1.046289 0.30639690  0.589225   197.056   204.027
gene3  121.3131   1.00004  0.299729 0.58407809  0.847352   226.950   233.920
gene4  111.1213   1.00008  6.872240 0.00876037  0.104032   224.468   231.438
gene5  124.0881   1.00006  0.190059 0.66295615  0.847352   226.581   233.551
gene6  109.5878   1.00007  0.192337 0.66105018  0.847352   237.443   244.413

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   94.1884  0.570262  0.329588  1.730227 0.08358976  0.563886   229.534
gene2   58.8858 -0.386718  0.328944 -1.175635 0.23974095  0.673275   197.056
gene3  121.3131  1.130666  0.369804  3.057473 0.00223212  0.111606   226.950
gene4  111.1213 -0.222158  0.387631 -0.573118 0.56656502  0.885258   224.468
gene5  124.0881  0.162123  0.414033  0.391569 0.69537653  0.922975   226.581
gene6  109.5878  0.275666  0.372303  0.740436 0.45903567  0.791441   237.443
            BIC
      <numeric>
gene1   236.504
gene2   204.027
gene3   233.920
gene4   231.438
gene5   233.551
gene6   244.413

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   94.1884  0.765390   1.15265  0.6640292 0.50667163 0.6682748   229.534
gene2   58.8858  3.042428   1.13934  2.6703472 0.00757729 0.0757729   197.056
gene3  121.3131 -1.981081   1.28739 -1.5388386 0.12384370 0.4128123   226.950
gene4  111.1213 -2.460492   1.39317 -1.7661140 0.07737673 0.2976028   224.468
gene5  124.0881  0.131364   1.44005  0.0912214 0.92731666 0.9462415   226.581
gene6  109.5878  1.325091   1.29545  1.0228820 0.30636365 0.5470779   237.443
            BIC
      <numeric>
gene1   236.504
gene2   204.027
gene3   233.920
gene4   231.438
gene5   233.551
gene6   244.413

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>
gene34   78.1512   1.00003  11.91768 0.000556155 0.0155125   189.025   195.995
gene46  138.2990   1.00006  11.71412 0.000620498 0.0155125   243.769   250.740
gene4   111.1213   1.00008   6.87224 0.008760369 0.1040324   224.468   231.438
gene11   75.1442   1.00007   6.79070 0.009166309 0.1040324   222.928   229.898
gene30  201.4262   1.00006   6.36904 0.011619546 0.1040324   234.809   241.780
gene36   72.6573   1.00005   6.19944 0.012782615 0.1040324   205.738   212.708
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.0 (2025-04-11 ucrt)
Platform: x86_64-w64-mingw32/x64
Running under: Windows Server 2022 x64 (build 20348)

Matrix products: default
  LAPACK version 3.12.1

locale:
[1] LC_COLLATE=C                          
[2] LC_CTYPE=English_United States.utf8   
[3] LC_MONETARY=English_United States.utf8
[4] LC_NUMERIC=C                          
[5] LC_TIME=English_United States.utf8    

time zone: America/New_York
tzcode source: internal

attached base packages:
[1] stats4    stats     graphics  grDevices utils     datasets  methods  
[8] base     

other attached packages:
 [1] ggplot2_3.5.2               BiocParallel_1.43.3        
 [3] NBAMSeq_1.25.0              SummarizedExperiment_1.39.0
 [5] Biobase_2.69.0              GenomicRanges_1.61.0       
 [7] GenomeInfoDb_1.45.4         IRanges_2.43.0             
 [9] S4Vectors_0.47.0            BiocGenerics_0.55.0        
[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.0      gtable_0.3.6         xfun_0.52           
 [4] bslib_0.9.0          lattice_0.22-7       vctrs_0.6.5         
 [7] tools_4.5.0          parallel_4.5.0       tibble_3.2.1        
[10] AnnotationDbi_1.71.0 RSQLite_2.4.0        blob_1.2.4          
[13] pkgconfig_2.0.3      Matrix_1.7-3         RColorBrewer_1.1-3  
[16] lifecycle_1.0.4      compiler_4.5.0       farver_2.1.2        
[19] Biostrings_2.77.1    DESeq2_1.49.1        codetools_0.2-20    
[22] snow_0.4-4           htmltools_0.5.8.1    sass_0.4.10         
[25] yaml_2.3.10          pillar_1.10.2        crayon_1.5.3        
[28] jquerylib_0.1.4      DelayedArray_0.35.1  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.0       
[40] fastmap_1.2.0        grid_4.5.0           cli_3.6.5           
[43] SparseArray_1.9.0    magrittr_2.0.3       S4Arrays_1.9.1      
[46] survival_3.8-3       dichromat_2.0-0.1    XML_3.99-0.18       
[49] withr_3.0.2          scales_1.4.0         UCSC.utils_1.5.0    
[52] bit64_4.6.0-1        rmarkdown_2.29       XVector_0.49.0      
[55] httr_1.4.7           bit_4.6.0            png_0.1-8           
[58] memoise_2.0.1        evaluate_1.0.3       knitr_1.50          
[61] mgcv_1.9-3           rlang_1.1.6          Rcpp_1.0.14         
[64] xtable_1.8-4         glue_1.8.0           DBI_1.2.3           
[67] annotate_1.87.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.