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      11     311     213      27       5     374      81       6      17
gene2       1       8     124      64       1       1      95     280     100
gene3      66      16     125      50       2      96      68     171      13
gene4       6       1       1      76      28     291       1      30     166
gene5       6       7       1       1       9      41       1     238      17
gene6      41     148       8      29       8     160      82       2     302
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1       98      101       22      126        4        1      108       77
gene2      158        1      206      184       16        3        5        2
gene3       10        1      183        4        1       35        9       61
gene4      106      676      247        8        1        1        6       93
gene5       91       22       70        5        1       35      288      237
gene6       17      159      254      170        6       85       46      208
      sample18 sample19 sample20
gene1       18        3       65
gene2      189       50       18
gene3      142       83      103
gene4       37        3      155
gene5        4       30      465
gene6       10       69        2

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 22.20976 -0.6861990 -0.46831152 -1.13695678    1
sample2 51.48763 -0.2329253  1.38927429 -0.09283967    1
sample3 41.91769 -0.4907470  1.68299231 -0.45066115    1
sample4 57.49310  0.7780957 -0.41440236  0.64666067    0
sample5 34.94452 -0.4807003 -1.86523946  0.11369116    0
sample6 36.39463 -0.4595581 -0.05138813  0.12378722    2

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   60.4409   1.00008  1.355475 0.244330151 0.4072169   217.769   224.739
gene2   68.0942   1.00004 12.133842 0.000495205 0.0125431   209.872   216.842
gene3   50.1512   1.00007  0.060392 0.806021162 0.8976087   216.196   223.167
gene4   70.7558   1.78759  3.843925 0.201767279 0.3602987   216.246   224.001
gene5   59.1123   1.54425  1.674631 0.325896316 0.4792593   205.701   213.213
gene6   72.2962   1.00005  1.022798 0.311882309 0.4725490   231.047   238.017

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   60.4409 -0.606706  0.472863 -1.283049 0.1994750 0.4533523   217.769
gene2   68.0942  1.380261  0.549856  2.510223 0.0120655 0.0861821   209.872
gene3   50.1512  0.172402  0.503254  0.342575 0.7319181 0.8827780   216.196
gene4   70.7558 -0.554656  0.618544 -0.896712 0.3698727 0.6776470   216.246
gene5   59.1123  0.251822  0.556411  0.452583 0.6508489 0.8452655   205.701
gene6   72.2962 -0.727265  0.511348 -1.422252 0.1549531 0.3873828   231.047
            BIC
      <numeric>
gene1   224.739
gene2   216.842
gene3   223.167
gene4   224.001
gene5   213.213
gene6   238.017

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   60.4409  0.652802  0.819469  0.7966165  0.425674  0.625991   217.769
gene2   68.0942  0.857123  0.950514  0.9017469  0.367191  0.583869   209.872
gene3   50.1512  1.194465  0.871285  1.3709239  0.170399  0.405711   216.196
gene4   70.7558  0.061954  1.050754  0.0589615  0.952983  0.952983   216.246
gene5   59.1123  0.817297  0.949819  0.8604766  0.389526  0.590191   205.701
gene6   72.2962 -0.432981  0.874472 -0.4951342  0.620505  0.738697   231.047
            BIC
      <numeric>
gene1   224.739
gene2   216.842
gene3   223.167
gene4   224.001
gene5   213.213
gene6   238.017

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>
gene2    68.0942   1.00004  12.13384 0.000495205 0.0125431   209.872   216.842
gene43   66.7524   1.00006  12.10937 0.000501724 0.0125431   199.390   206.360
gene17   87.0458   1.00005   9.37105 0.002204833 0.0344369   219.077   226.047
gene27  129.5501   1.00005   8.96357 0.002754954 0.0344369   242.456   249.426
gene25   53.9172   1.00004   8.45000 0.003651689 0.0365169   193.746   200.716
gene13   29.8594   1.00012   7.45604 0.006324872 0.0527073   179.124   186.094
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 Under development (unstable) (2022-10-25 r83175)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 22.04.1 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.17-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       

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

other attached packages:
 [1] ggplot2_3.3.6               BiocParallel_1.33.0        
 [3] NBAMSeq_1.15.0              SummarizedExperiment_1.29.0
 [5] Biobase_2.59.0              GenomicRanges_1.51.0       
 [7] GenomeInfoDb_1.35.0         IRanges_2.33.0             
 [9] S4Vectors_0.37.0            BiocGenerics_0.45.0        
[11] MatrixGenerics_1.11.0       matrixStats_0.62.0         

loaded via a namespace (and not attached):
 [1] tidyselect_1.2.0       farver_2.1.1           dplyr_1.0.10          
 [4] blob_1.2.3             Biostrings_2.67.0      bitops_1.0-7          
 [7] fastmap_1.1.0          RCurl_1.98-1.9         XML_3.99-0.12         
[10] digest_0.6.30          lifecycle_1.0.3        survival_3.4-0        
[13] KEGGREST_1.39.0        RSQLite_2.2.18         magrittr_2.0.3        
[16] genefilter_1.81.0      compiler_4.3.0         rlang_1.0.6           
[19] sass_0.4.2             tools_4.3.0            utf8_1.2.2            
[22] yaml_2.3.6             knitr_1.40             labeling_0.4.2        
[25] bit_4.0.4              DelayedArray_0.25.0    RColorBrewer_1.1-3    
[28] withr_2.5.0            grid_4.3.0             fansi_1.0.3           
[31] xtable_1.8-4           colorspace_2.0-3       scales_1.2.1          
[34] cli_3.4.1              rmarkdown_2.17         crayon_1.5.2          
[37] generics_0.1.3         httr_1.4.4             DBI_1.1.3             
[40] cachem_1.0.6           stringr_1.4.1          zlibbioc_1.45.0       
[43] splines_4.3.0          assertthat_0.2.1       parallel_4.3.0        
[46] AnnotationDbi_1.61.0   XVector_0.39.0         vctrs_0.5.0           
[49] Matrix_1.5-1           jsonlite_1.8.3         geneplotter_1.77.0    
[52] bit64_4.0.5            locfit_1.5-9.6         jquerylib_0.1.4       
[55] annotate_1.77.0        glue_1.6.2             codetools_0.2-18      
[58] stringi_1.7.8          gtable_0.3.1           munsell_0.5.0         
[61] tibble_3.1.8           pillar_1.8.1           htmltools_0.5.3       
[64] GenomeInfoDbData_1.2.9 R6_2.5.1               evaluate_0.17         
[67] lattice_0.20-45        highr_0.9              png_0.1-7             
[70] memoise_2.0.1          bslib_0.4.0            Rcpp_1.0.9            
[73] nlme_3.1-160           mgcv_1.8-41            DESeq2_1.39.0         
[76] xfun_0.34              pkgconfig_2.0.3       

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.