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      36       1      24      42     176     128     254      61      76
gene2       7     154     251       1      10     329       6       1      43
gene3      41      46       1     400     232      12     221      32      94
gene4    1302     138       1       7      79     129      15      23       1
gene5       1      28       1     136       1       2       1      67       1
gene6     145       1       1       6     104     177      53      98     225
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1       11       27      191        1       97       39        1       70
gene2        1       13        2        1       11        1        1       48
gene3      140        1       37       17        3        3       17      205
gene4      101        1        8      167       13        3       41       21
gene5        1        3        2       30       16        1       14        5
gene6       22       27       13      198      444      683      164      415
      sample18 sample19 sample20
gene1      258        9       44
gene2       76        7      255
gene3        8       53        3
gene4      134       31        7
gene5        8        6      114
gene6        1        5        3

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 35.29267 -0.231470167 -1.1528658  0.05801744    0
sample2 30.05642 -0.121775976  1.8556718 -0.30875485    2
sample3 46.71960  1.131782772 -0.5696251  1.70792103    1
sample4 42.19015 -0.334950747 -1.1720077 -0.22341174    0
sample5 56.83200 -0.405029061 -1.4332656  0.10465325    0
sample6 48.75202  0.001871494 -1.7051015 -0.26506865    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:

  • 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   67.8088   1.00005  0.319907 0.57169808 0.9162836   227.508   234.478
gene2   65.4074   1.00038  0.795235 0.37263536 0.8986182   200.381   207.352
gene3   67.1238   1.00007  0.878490 0.34868188 0.8986182   212.158   219.129
gene4   83.2885   1.00006  8.838854 0.00294972 0.0721591   216.026   222.996
gene5   20.1673   1.00005  1.122337 0.28943942 0.8986182   164.883   171.853
gene6  100.5017   1.00008  0.217798 0.64078607 0.9162836   229.774   236.745

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   67.8088 -0.491513  0.689601 -0.712750  0.476001  0.610257   227.508
gene2   65.4074  0.953351  0.854170  1.116114  0.264373  0.583781   200.381
gene3   67.1238 -0.532217  0.692255 -0.768817  0.442002  0.610257   212.158
gene4   83.2885  0.132785  0.744304  0.178402  0.858407  0.894174   216.026
gene5   20.1673 -0.746732  0.815579 -0.915585  0.359884  0.599206   164.883
gene6  100.5017 -1.145177  0.723532 -1.582758  0.113477  0.378255   229.774
            BIC
      <numeric>
gene1   234.478
gene2   207.352
gene3   219.129
gene4   222.996
gene5   171.853
gene6   236.745

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   67.8088  0.599749  0.935989  0.640765 0.5216752  0.724231   227.508
gene2   65.4074  1.985626  1.159109  1.713062 0.0867011  0.312381   200.381
gene3   67.1238 -0.789314  0.901335 -0.875716 0.3811846  0.614814   212.158
gene4   83.2885  0.511225  1.005608  0.508374 0.6111912  0.724231   216.026
gene5   20.1673  0.970726  1.071794  0.905702 0.3650936  0.614814   164.883
gene6  100.5017 -1.443829  0.934282 -1.545388 0.1222523  0.382039   229.774
            BIC
      <numeric>
gene1   234.478
gene2   207.352
gene3   219.129
gene4   222.996
gene5   171.853
gene6   236.745

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>
gene45  131.3899   1.00018  14.84635 0.000116552 0.0058276   228.734   235.705
gene4    83.2885   1.00006   8.83885 0.002949721 0.0721591   216.026   222.996
gene20   59.6507   1.00008   8.14052 0.004329546 0.0721591   203.183   210.153
gene13  102.1395   1.00006   6.01169 0.014215635 0.1776954   222.264   229.234
gene22  118.8957   1.00008   4.10980 0.042647941 0.4264794   245.888   252.859
gene21   91.8227   1.00015   3.27993 0.070169785 0.5251234   211.427   218.397
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.4.1 (2024-06-14)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.1 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_3.5.1               BiocParallel_1.41.0        
 [3] NBAMSeq_1.23.0              SummarizedExperiment_1.35.5
 [5] Biobase_2.67.0              GenomicRanges_1.57.2       
 [7] GenomeInfoDb_1.41.2         IRanges_2.39.2             
 [9] S4Vectors_0.43.2            BiocGenerics_0.53.0        
[11] MatrixGenerics_1.17.1       matrixStats_1.4.1          
[13] rmarkdown_2.28             

loaded via a namespace (and not attached):
 [1] KEGGREST_1.45.1         gtable_0.3.6            xfun_0.48              
 [4] bslib_0.8.0             lattice_0.22-6          vctrs_0.6.5            
 [7] tools_4.4.1             parallel_4.4.1          RSQLite_2.3.7          
[10] tibble_3.2.1            fansi_1.0.6             AnnotationDbi_1.69.0   
[13] highr_0.11              blob_1.2.4              pkgconfig_2.0.3        
[16] Matrix_1.7-1            lifecycle_1.0.4         GenomeInfoDbData_1.2.13
[19] farver_2.1.2            compiler_4.4.1          Biostrings_2.75.0      
[22] munsell_0.5.1           DESeq2_1.47.0           codetools_0.2-20       
[25] htmltools_0.5.8.1       sys_3.4.3               buildtools_1.0.0       
[28] sass_0.4.9              yaml_2.3.10             pillar_1.9.0           
[31] crayon_1.5.3            jquerylib_0.1.4         DelayedArray_0.33.1    
[34] cachem_1.1.0            abind_1.4-8             nlme_3.1-166           
[37] genefilter_1.87.0       locfit_1.5-9.10         digest_0.6.37          
[40] labeling_0.4.3          splines_4.4.1           maketools_1.3.1        
[43] fastmap_1.2.0           grid_4.4.1              colorspace_2.1-1       
[46] cli_3.6.3               SparseArray_1.5.45      magrittr_2.0.3         
[49] S4Arrays_1.5.11         survival_3.7-0          XML_3.99-0.17          
[52] utf8_1.2.4              withr_3.0.2             scales_1.3.0           
[55] UCSC.utils_1.1.0        bit64_4.5.2             XVector_0.45.0         
[58] httr_1.4.7              bit_4.5.0               png_0.1-8              
[61] memoise_2.0.1           evaluate_1.0.1          knitr_1.48             
[64] mgcv_1.9-1              rlang_1.1.4             Rcpp_1.0.13            
[67] DBI_1.2.3               xtable_1.8-4            glue_1.8.0             
[70] annotate_1.85.0         jsonlite_1.8.9          R6_2.5.1               
[73] zlibbioc_1.51.2        

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.