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     314     220      78      26     191       1     119      13      48
gene2      54       2      91     375       1      15     119      15       6
gene3     504      82       3     114     504       4     228     148      23
gene4     560       4     398       1     158      51      74     169       1
gene5      53     184     400     244       2       1     121     221      56
gene6      24      45      62      15      20      80     149      89       2
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1       47       15        3      190      540      205       23        1
gene2        5       18       94        1      162       17        1        9
gene3        1       80        1       71       11        7      225        4
gene4       29      799        1        1      280       21        4        1
gene5       13      941       35      207        6        3      188        9
gene6      412      250        1      165       58       26        8      126
      sample18 sample19 sample20
gene1       18       44       22
gene2      141       19        1
gene3      138      143        5
gene4        1       31        1
gene5        1      632     1331
gene6      199        2       19

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.69074  1.7574283  0.35616797  0.6082323    1
sample2 22.38337 -0.1630608 -1.42379609  0.8317170    2
sample3 28.70908  0.3734551 -1.35523249 -0.2085717    1
sample4 41.75098  0.7894484  0.13706945 -1.2922998    1
sample5 34.96640  2.2162639  0.07763842  0.3669172    1
sample6 59.58079  0.1030320  0.58479007 -0.5326034    1

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   73.7650   1.00004 0.000393897 0.9850251  0.995067   232.351   239.321
gene2   49.9695   1.00005 0.106880717 0.7437918  0.995067   199.867   206.837
gene3   92.0971   1.00013 2.845960441 0.0916208  0.381753   230.677   237.647
gene4   89.7748   1.00007 0.006296591 0.9371847  0.995067   201.874   208.844
gene5  206.1077   1.00007 6.413708698 0.0113297  0.188829   255.317   262.287
gene6   59.5017   1.00005 2.020401125 0.1551896  0.566168   221.385   228.355

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   73.7650  0.4081456  0.372220  1.096517  0.272853  0.568443   232.351
gene2   49.9695  0.3869508  0.395083  0.979417  0.327374  0.607255   199.867
gene3   92.0971  0.3842441  0.402754  0.954041  0.340063  0.607255   230.677
gene4   89.7748 -0.0483954  0.405281 -0.119412  0.904949  0.961960   201.874
gene5  206.1077 -0.7318945  0.450414 -1.624937  0.104176  0.447010   255.317
gene6   59.5017 -0.5218413  0.333066 -1.566778  0.117167  0.447010   221.385
            BIC
      <numeric>
gene1   239.321
gene2   206.837
gene3   237.647
gene4   208.844
gene5   262.287
gene6   228.355

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
      <numeric>  <numeric> <numeric>  <numeric>   <numeric>  <numeric>
gene1   73.7650 -0.0925709   1.16000 -0.0798026 0.936394276 0.96250870
gene2   49.9695 -2.7321609   1.25342 -2.1797691 0.029274584 0.29274584
gene3   92.0971 -2.0140630   1.25840 -1.6005004 0.109487620 0.46685851
gene4   89.7748 -4.9805070   1.31662 -3.7827932 0.000155078 0.00775391
gene5  206.1077  1.1786030   1.41071  0.8354696 0.403453342 0.74713582
gene6   59.5017 -0.0640978   1.03328 -0.0620336 0.950536043 0.96250870
            AIC       BIC
      <numeric> <numeric>
gene1   232.351   239.321
gene2   199.867   206.837
gene3   230.677   237.647
gene4   201.874   208.844
gene5   255.317   262.287
gene6   221.385   228.355

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>
gene39   53.7132   1.00004   7.13295 0.00756862  0.188829   201.166   208.137
gene12   82.8388   1.00012   6.88767 0.00868050  0.188829   213.152   220.122
gene5   206.1077   1.00007   6.41371 0.01132973  0.188829   255.317   262.287
gene22   83.0278   1.00005   4.78418 0.02872720  0.313714   190.483   197.453
gene40  111.3719   1.00014   4.63275 0.03137138  0.313714   230.209   237.179
gene48   51.3919   1.00006   4.02046 0.04495386  0.329990   204.913   211.883
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-13 r84269)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 22.04.2 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       

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.34.0        
 [3] NBAMSeq_1.16.0              SummarizedExperiment_1.30.0
 [5] Biobase_2.60.0              GenomicRanges_1.52.0       
 [7] GenomeInfoDb_1.36.0         IRanges_2.34.0             
 [9] S4Vectors_0.38.0            BiocGenerics_0.46.0        
[11] MatrixGenerics_1.12.0       matrixStats_0.63.0         

loaded via a namespace (and not attached):
 [1] KEGGREST_1.40.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.62.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.68.0       munsell_0.5.0          
[25] DESeq2_1.40.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.26.0     cachem_1.0.7            nlme_3.1-162           
[37] genefilter_1.82.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.40.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.78.0        
[70] jsonlite_1.8.4          R6_2.5.1                zlibbioc_1.46.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.