To install and load NBAMSeq
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.
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     530      15      94     248       1      37     126      67      77
gene2      11      12       3     201      52       4       4     110      17
gene3       9       7       1       1      32     220     337      91     422
gene4       3       3       5     151     206       2       1     262       3
gene5       3      53     180      64     314      12       2       1      91
gene6      83       2      19       1       1      97    1298       1      26
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1        1       20        2      363       29        1       38      174
gene2       54      586       17        9      199       35        1      127
gene3        1        1        5      126        3        9      106        2
gene4       16      404       57       10      212       38        1       79
gene5      149      294      622        2       11        1      280      239
gene6      135      127        1        5      135       93       51      138
      sample18 sample19 sample20
gene1      242       60      701
gene2      227       46       89
gene3      435        2        1
gene4        1        6        1
gene5       66       45       12
gene6       77       30      266colData 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 68.02791 -0.6161603 -0.1301888  1.3166746    2
sample2 75.51706 -0.7063917  0.8411168 -0.1783757    0
sample3 61.63397 -0.3797644  0.9230199  1.4793823    2
sample4 25.85360  0.3205278  0.1807490 -0.5231973    1
sample5 32.68332 -0.3750781  0.8562689  0.9641655    2
sample6 45.64452  1.7378885  0.3242683  0.7408498    2design 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:
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:
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 var4Differential expression analysis can be performed by NBAMSeq function:
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:
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.
DataFrame with 6 rows and 7 columns
       baseMean       edf        stat     pvalue      padj       AIC       BIC
      <numeric> <numeric>   <numeric>  <numeric> <numeric> <numeric> <numeric>
gene1  150.1351   1.00010 0.294831937 0.58722738  0.863570   241.187   248.157
gene2   73.9266   1.00010 1.509395927 0.21928369  0.580796   221.000   227.970
gene3   60.6195   1.00010 2.073674804 0.14990659  0.499689   205.165   212.135
gene4   53.5216   1.00011 0.829544740 0.36251367  0.697570   205.482   212.452
gene5  109.4553   1.00038 0.000360007 0.99164934  0.991649   234.712   241.682
gene6  103.9510   1.00006 6.687226291 0.00971908  0.161985   227.629   234.599For linear continuous covariates, base mean, estimated coefficient, standard error, test statistics, p-value, and adjusted p-value will be returned.
DataFrame with 6 rows and 8 columns
       baseMean         coef        SE         stat    pvalue      padj
      <numeric>    <numeric> <numeric>    <numeric> <numeric> <numeric>
gene1  150.1351  0.491522029  0.457555  1.074237102  0.282716  0.579238
gene2   73.9266  0.645825550  0.403676  1.599863035  0.109629  0.416217
gene3   60.6195  0.000119174  0.514245  0.000231746  0.999815  0.999815
gene4   53.5216 -0.197023397  0.506856 -0.388717073  0.697485  0.830340
gene5  109.4553  0.022854550  0.469582  0.048669998  0.961182  0.999815
gene6  103.9510  0.201451690  0.468924  0.429604134  0.667484  0.814004
            AIC       BIC
      <numeric> <numeric>
gene1   241.187   248.157
gene2   221.000   227.970
gene3   205.165   212.135
gene4   205.482   212.452
gene5   234.712   241.682
gene6   227.629   234.599For discrete covariates, the contrast argument should be specified. e.g.  contrast = c("var4", "2", "0") means comparing level 2 vs. level 0 in var4.
DataFrame with 6 rows and 8 columns
       baseMean      coef        SE       stat    pvalue      padj       AIC
      <numeric> <numeric> <numeric>  <numeric> <numeric> <numeric> <numeric>
gene1  150.1351  1.316918  1.089958  1.2082284 0.2269594  0.667528   241.187
gene2   73.9266  0.526804  0.961311  0.5480056 0.5836881  0.775198   221.000
gene3   60.6195 -0.905719  1.217645 -0.7438287 0.4569801  0.771474   205.165
gene4   53.5216  0.080167  1.209217  0.0662967 0.9471416  0.986606   205.482
gene5  109.4553 -2.501405  1.119020 -2.2353536 0.0253941  0.211618   234.712
gene6  103.9510  1.223695  1.116206  1.0962988 0.2729480  0.682370   227.629
            BIC
      <numeric>
gene1   248.157
gene2   227.970
gene3   212.135
gene4   212.452
gene5   241.682
gene6   234.599We 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>
gene50   21.2811   1.00003   8.66888 0.00323827  0.161914   168.367   175.338
gene47  163.6726   1.00018   6.97203 0.00829817  0.161985   250.624   257.595
gene6   103.9510   1.00006   6.68723 0.00971908  0.161985   227.629   234.599
gene15  131.0714   1.00006   4.82839 0.02801147  0.238817   206.091   213.061
gene34  107.1092   1.00011   4.74397 0.02942056  0.238817   214.321   221.291
gene16  101.9969   1.00006   4.53552 0.03321386  0.238817   220.562   227.532library(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))R version 4.1.0 (2021-05-18)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 20.04.2 LTS
Matrix products: default
BLAS:   /home/biocbuild/bbs-3.13-bioc/R/lib/libRblas.so
LAPACK: /home/biocbuild/bbs-3.13-bioc/R/lib/libRlapack.so
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       
attached base packages:
[1] parallel  stats4    stats     graphics  grDevices utils     datasets 
[8] methods   base     
other attached packages:
 [1] ggplot2_3.3.3               BiocParallel_1.26.0        
 [3] NBAMSeq_1.8.0               SummarizedExperiment_1.22.0
 [5] Biobase_2.52.0              GenomicRanges_1.44.0       
 [7] GenomeInfoDb_1.28.0         IRanges_2.26.0             
 [9] S4Vectors_0.30.0            BiocGenerics_0.38.0        
[11] MatrixGenerics_1.4.0        matrixStats_0.58.0         
loaded via a namespace (and not attached):
 [1] httr_1.4.2             sass_0.4.0             bit64_4.0.5           
 [4] jsonlite_1.7.2         splines_4.1.0          bslib_0.2.5.1         
 [7] assertthat_0.2.1       highr_0.9              blob_1.2.1            
[10] GenomeInfoDbData_1.2.6 yaml_2.2.1             pillar_1.6.1          
[13] RSQLite_2.2.7          lattice_0.20-44        glue_1.4.2            
[16] digest_0.6.27          RColorBrewer_1.1-2     XVector_0.32.0        
[19] colorspace_2.0-1       htmltools_0.5.1.1      Matrix_1.3-3          
[22] DESeq2_1.32.0          XML_3.99-0.6           pkgconfig_2.0.3       
[25] genefilter_1.74.0      zlibbioc_1.38.0        purrr_0.3.4           
[28] xtable_1.8-4           scales_1.1.1           tibble_3.1.2          
[31] annotate_1.70.0        mgcv_1.8-35            KEGGREST_1.32.0       
[34] farver_2.1.0           generics_0.1.0         ellipsis_0.3.2        
[37] withr_2.4.2            cachem_1.0.5           survival_3.2-11       
[40] magrittr_2.0.1         crayon_1.4.1           memoise_2.0.0         
[43] evaluate_0.14          fansi_0.4.2            nlme_3.1-152          
[46] tools_4.1.0            lifecycle_1.0.0        stringr_1.4.0         
[49] locfit_1.5-9.4         munsell_0.5.0          DelayedArray_0.18.0   
[52] AnnotationDbi_1.54.0   Biostrings_2.60.0      compiler_4.1.0        
[55] jquerylib_0.1.4        rlang_0.4.11           grid_4.1.0            
[58] RCurl_1.98-1.3         labeling_0.4.2         bitops_1.0-7          
[61] rmarkdown_2.8          gtable_0.3.0           DBI_1.1.1             
[64] R6_2.5.0               knitr_1.33             dplyr_1.0.6           
[67] fastmap_1.1.0          bit_4.0.4              utf8_1.2.1            
[70] stringi_1.6.2          Rcpp_1.0.6             vctrs_0.3.8           
[73] geneplotter_1.70.0     png_0.1-7              tidyselect_1.1.1      
[76] xfun_0.23             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.