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       5     195      16      38       1     258       1      81      52
gene2     574     119      26      52       1      24       8     175     204
gene3       2      78     199       2      43       1      48       1       1
gene4       1     320       3       2      70       5      96       1     163
gene5     352     209     277      13     274      28     118      54      23
gene6     113     211      96     116      42       2      24      36       2
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1      180      119       45       40      100      170        5       10
gene2       23       74       17      133       20        8       78        1
gene3      227        8        7        7        6      109       48      419
gene4       21        1       60       49      209      101        4        3
gene5        9        8       11        1        7      130      269        9
gene6      485       34        8        4        1        8       32        8
      sample18 sample19 sample20
gene1        1      533       12
gene2       36       35        9
gene3        1        8      111
gene4      132       11      191
gene5      228      159        1
gene6      277       25        4

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 27.58401  0.6501462 -1.4634808  0.2495924    2
sample2 52.96869 -1.1945369 -0.4601252  0.1985181    0
sample3 77.50772 -1.1041742 -1.5669763 -1.9828937    0
sample4 25.84204  0.5408499  0.9292463 -0.4723164    1
sample5 76.86154 -0.6434235 -0.7285971  1.4055751    0
sample6 22.98632 -0.6297166  0.7319863 -0.6933401    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   79.3411   1.00006  3.481707 0.06206167 0.2216488   215.490   222.461
gene2   65.4688   1.00005  4.086654 0.04323461 0.1712771   217.997   224.968
gene3   53.3500   1.00006  8.953142 0.00277099 0.0435694   198.590   205.561
gene4   53.0708   1.00003  2.013322 0.15592877 0.4331355   204.645   211.615
gene5   77.4277   1.00003  0.120025 0.72906045 0.8477447   231.593   238.563
gene6   62.1240   1.00016  0.327752 0.56725582 0.7509692   217.823   224.794

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   79.3411  0.2558010  0.385711  0.6631927  0.507207  0.845345   215.490
gene2   65.4688 -0.2975635  0.384450 -0.7739980  0.438932  0.783807   217.997
gene3   53.3500  0.0433835  0.444812  0.0975321  0.922304  0.944133   198.590
gene4   53.0708 -0.3643095  0.410775 -0.8868825  0.375142  0.755291   204.645
gene5   77.4277 -0.3457732  0.417767 -0.8276708  0.407857  0.755291   231.593
gene6   62.1240 -0.6820386  0.447316 -1.5247347  0.127325  0.432541   217.823
            BIC
      <numeric>
gene1   222.461
gene2   224.968
gene3   205.561
gene4   211.615
gene5   238.563
gene6   224.794

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   79.3411 -3.6295911   1.10358 -3.2889143 0.00100575 0.0251437   215.490
gene2   65.4688  0.7374638   1.08279  0.6810789 0.49582155 0.7753749   217.997
gene3   53.3500  2.3222829   1.23593  1.8789762 0.06024775 0.3347097   198.590
gene4   53.0708 -1.8110701   1.19411 -1.5166721 0.12934950 0.4974981   204.645
gene5   77.4277  0.6393265   1.17726  0.5430625 0.58708674 0.8260584   231.593
gene6   62.1240 -0.0287999   1.26000 -0.0228571 0.98176424 0.9817642   217.823
            BIC
      <numeric>
gene1   222.461
gene2   224.968
gene3   205.561
gene4   211.615
gene5   238.563
gene6   224.794

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>
gene43   62.5764   1.00011  12.21163 0.000475519 0.0215939   196.188   203.159
gene34   38.0694   1.00002  11.09945 0.000863755 0.0215939   185.753   192.723
gene3    53.3500   1.00006   8.95314 0.002770992 0.0435694   198.590   205.561
gene23   38.6098   1.00003   8.53465 0.003485553 0.0435694   183.161   190.131
gene11   88.4396   1.00003   7.30672 0.006870682 0.0613455   213.484   220.454
gene30   96.0069   1.00007   7.18312 0.007361460 0.0613455   210.314   217.284
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.0 beta (2024-04-15 r86425 ucrt)
Platform: x86_64-w64-mingw32/x64
Running under: Windows Server 2022 x64 (build 20348)

Matrix products: default


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.1               BiocParallel_1.38.0        
 [3] NBAMSeq_1.20.0              SummarizedExperiment_1.34.0
 [5] Biobase_2.64.0              GenomicRanges_1.56.0       
 [7] GenomeInfoDb_1.40.0         IRanges_2.38.0             
 [9] S4Vectors_0.42.0            BiocGenerics_0.50.0        
[11] MatrixGenerics_1.16.0       matrixStats_1.3.0          

loaded via a namespace (and not attached):
 [1] KEGGREST_1.44.0         gtable_0.3.5            xfun_0.43              
 [4] bslib_0.7.0             lattice_0.22-6          vctrs_0.6.5            
 [7] tools_4.4.0             generics_0.1.3          parallel_4.4.0         
[10] RSQLite_2.3.6           tibble_3.2.1            fansi_1.0.6            
[13] AnnotationDbi_1.66.0    highr_0.10              blob_1.2.4             
[16] pkgconfig_2.0.3         Matrix_1.7-0            lifecycle_1.0.4        
[19] GenomeInfoDbData_1.2.12 farver_2.1.1            compiler_4.4.0         
[22] Biostrings_2.72.0       munsell_0.5.1           DESeq2_1.44.0          
[25] codetools_0.2-20        snow_0.4-4              htmltools_0.5.8.1      
[28] sass_0.4.9              yaml_2.3.8              pillar_1.9.0           
[31] crayon_1.5.2            jquerylib_0.1.4         DelayedArray_0.30.0    
[34] cachem_1.0.8            abind_1.4-5             nlme_3.1-164           
[37] genefilter_1.86.0       tidyselect_1.2.1        locfit_1.5-9.9         
[40] digest_0.6.35           dplyr_1.1.4             labeling_0.4.3         
[43] splines_4.4.0           fastmap_1.1.1           grid_4.4.0             
[46] colorspace_2.1-0        cli_3.6.2               SparseArray_1.4.0      
[49] magrittr_2.0.3          S4Arrays_1.4.0          survival_3.6-4         
[52] XML_3.99-0.16.1         utf8_1.2.4              withr_3.0.0            
[55] scales_1.3.0            UCSC.utils_1.0.0        bit64_4.0.5            
[58] rmarkdown_2.26          XVector_0.44.0          httr_1.4.7             
[61] bit_4.0.5               png_0.1-8               memoise_2.0.1          
[64] evaluate_0.23           knitr_1.46              mgcv_1.9-1             
[67] rlang_1.1.3             Rcpp_1.0.12             DBI_1.2.2              
[70] xtable_1.8-4            glue_1.7.0              annotate_1.82.0        
[73] jsonlite_1.8.8          R6_2.5.1                zlibbioc_1.50.0        

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.