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      67    1008       1      74     291      11     101      20      39
gene2       1     286      71      60     213     261       3      20      92
gene3       1       1       1       4     153      17      19      44       3
gene4      15       8    1246      12      77     166       3       4       4
gene5     116      16      41       7      34     279     107     772       2
gene6       1     142       1     299      10      20     170       1       8
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1        1      153        1       13        1        1      111      165
gene2        2       77        1        1        1       48      183      209
gene3        5       22        1       44      128        1       20      206
gene4      111      354        2      493       65      169        2        6
gene5      122        1        1        1        1       22        1       37
gene6      208       44        1       87       36       22       17      137
      sample18 sample19 sample20
gene1       65      124      189
gene2        1       77       30
gene3        6       71      315
gene4        5      398      161
gene5       31        2       38
gene6       38        1       42

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 20.16715  0.2184479  0.1793180  0.10342530    2
sample2 28.67469 -0.5348193 -0.2253440  1.96910886    0
sample3 31.76493 -1.3643465 -0.5393728  0.01762547    2
sample4 67.48672 -1.4457682  1.9760324  1.14831552    2
sample5 29.61630  1.2858705 -1.1045828 -1.13049729    0
sample6 29.88034  0.5723151 -0.5904115  1.67435753    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   83.4015   1.00058 6.20576e-05 0.99882756  0.998828   228.837   235.808
gene2   81.4488   1.00008 2.08382e+00 0.14890072  0.553475   217.382   224.353
gene3   40.9799   1.00005 1.00272e+00 0.31671103  0.719798   195.484   202.454
gene4   96.6560   1.00006 8.93141e+00 0.00280399  0.140199   230.531   237.501
gene5   90.4174   1.00004 4.67954e-01 0.49395286  0.873316   204.182   211.152
gene6   51.4623   1.00012 7.37499e+00 0.00661694  0.165424   208.368   215.339

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   83.4015  0.1077510  0.379196  0.2841563 0.7762906  0.953695   228.837
gene2   81.4488  0.7928548  0.381580  2.0778191 0.0377260  0.188630   217.382
gene3   40.9799  0.3057559  0.363647  0.8408039 0.4004578  0.767597   195.484
gene4   96.6560  0.2340405  0.359351  0.6512867 0.5148614  0.858102   230.531
gene5   90.4174  0.9291528  0.384709  2.4152071 0.0157263  0.112330   204.182
gene6   51.4623 -0.0167848  0.349853 -0.0479769 0.9617346  0.968665   208.368
            BIC
      <numeric>
gene1   235.808
gene2   224.353
gene3   202.454
gene4   237.501
gene5   211.152
gene6   215.339

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   83.4015 -0.863572  0.883408 -0.977546 0.3282987  0.820747   228.837
gene2   81.4488 -0.647128  0.885823 -0.730539 0.4650611  0.868908   217.382
gene3   40.9799 -1.483732  0.846367 -1.753061 0.0795916  0.720496   195.484
gene4   96.6560  1.334247  0.837774  1.592609 0.1112479  0.720496   230.531
gene5   90.4174 -0.538872  0.899820 -0.598866 0.5492620  0.868908   204.182
gene6   51.4623 -1.004732  0.818355 -1.227746 0.2195423  0.784080   208.368
            BIC
      <numeric>
gene1   235.808
gene2   224.353
gene3   202.454
gene4   237.501
gene5   211.152
gene6   215.339

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>
gene4    96.6560   1.00006   8.93141 0.00280399  0.140199   230.531   237.501
gene6    51.4623   1.00012   7.37499 0.00661694  0.165424   208.368   215.339
gene18   67.2896   1.00012   5.47153 0.01933368  0.322228   210.267   217.237
gene26   90.6064   1.00005   4.81565 0.02820648  0.352581   226.111   233.081
gene20   62.5938   1.00011   3.82949 0.05037651  0.479528   203.800   210.770
gene45  185.1473   1.00007   3.60711 0.05754332  0.479528   251.780   258.751
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.3.0 RC (2023-04-18 r84287)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 22.04.2 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.18-bioc/R/lib/libRblas.so 
LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.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_3.4.2               BiocParallel_1.35.0        
 [3] NBAMSeq_1.17.0              SummarizedExperiment_1.31.0
 [5] Biobase_2.61.0              GenomicRanges_1.53.0       
 [7] GenomeInfoDb_1.37.0         IRanges_2.35.0             
 [9] S4Vectors_0.39.0            BiocGenerics_0.47.0        
[11] MatrixGenerics_1.13.0       matrixStats_0.63.0         

loaded via a namespace (and not attached):
 [1] KEGGREST_1.41.0         gtable_0.3.3            xfun_0.39              
 [4] bslib_0.4.2             lattice_0.21-8          vctrs_0.6.2            
 [7] tools_4.3.0             bitops_1.0-7            generics_0.1.3         
[10] parallel_4.3.0          tibble_3.2.1            fansi_1.0.4            
[13] AnnotationDbi_1.63.0    RSQLite_2.3.1           highr_0.10             
[16] blob_1.2.4              pkgconfig_2.0.3         Matrix_1.5-4           
[19] lifecycle_1.0.3         GenomeInfoDbData_1.2.10 farver_2.1.1           
[22] compiler_4.3.0          Biostrings_2.69.0       munsell_0.5.0          
[25] DESeq2_1.41.0           codetools_0.2-19        htmltools_0.5.5        
[28] sass_0.4.5              RCurl_1.98-1.12         yaml_2.3.7             
[31] crayon_1.5.2            pillar_1.9.0            jquerylib_0.1.4        
[34] DelayedArray_0.27.0     cachem_1.0.7            nlme_3.1-162           
[37] genefilter_1.83.0       tidyselect_1.2.0        locfit_1.5-9.7         
[40] digest_0.6.31           dplyr_1.1.2             labeling_0.4.2         
[43] splines_4.3.0           fastmap_1.1.1           grid_4.3.0             
[46] colorspace_2.1-0        cli_3.6.1               magrittr_2.0.3         
[49] survival_3.5-5          XML_3.99-0.14           utf8_1.2.3             
[52] withr_2.5.0             scales_1.2.1            bit64_4.0.5            
[55] rmarkdown_2.21          XVector_0.41.0          httr_1.4.5             
[58] bit_4.0.5               png_0.1-8               memoise_2.0.1          
[61] evaluate_0.20           knitr_1.42              mgcv_1.8-42            
[64] rlang_1.1.0             Rcpp_1.0.10             xtable_1.8-4           
[67] glue_1.6.2              DBI_1.1.3               annotate_1.79.0        
[70] jsonlite_1.8.4          R6_2.5.1                zlibbioc_1.47.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.