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      23      30      92       6       3     126      12     129       1
gene2      40      65       1       1      19     180      18      70     145
gene3      11       5      63     194     224      20     282     147     799
gene4     512       7       2       3     389       5     134      16       1
gene5     143      20       6      47      20       1       7     272       2
gene6       9     292       1      72     324       2     266     495       9
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1        1      550        3       21       91        4        2      211
gene2        2        6      165        8       23        5        1      117
gene3        8     1491        1      120       12        1       31      132
gene4        1       24       30       82       37        2        2       11
gene5       10      288      223        8     1140       56      512      105
gene6      235      241      168       26      143       59        1      134
      sample18 sample19 sample20
gene1        1      104       10
gene2       94       18        9
gene3        9        2        9
gene4       82        3       21
gene5        1        6       37
gene6        1        2      102

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 47.80916 -1.2025830 -1.7838450 -0.64103356    1
sample2 60.63673 -1.7494413  0.8675466 -1.85818205    1
sample3 43.73349  1.2382517 -0.2184462 -0.08237318    0
sample4 69.49117 -1.5701106 -1.7502891 -0.32889378    2
sample5 66.04341 -0.2695398  0.1446839 -0.39066974    2
sample6 66.93442  1.3510984 -0.3389946 -1.36065887    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   49.4800   1.00006 3.0550426 0.0805026  0.351480   190.030   197.001
gene2   42.5847   1.00005 0.0999655 0.7519882  0.883283   200.960   207.930
gene3  140.1690   1.00007 0.3445285 0.5572545  0.883283   227.461   234.431
gene4   46.6957   1.00010 0.1057019 0.7451861  0.883283   203.810   210.781
gene5  102.2155   1.00040 5.9903314 0.0143946  0.179933   231.640   238.611
gene6  104.9993   1.00021 1.4346728 0.2310764  0.525174   237.464   244.434

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   49.4800  0.763569  0.313281  2.437327 0.01479629 0.0882349   190.030
gene2   42.5847  1.121028  0.346958  3.231023 0.00123348 0.0282617   200.960
gene3  140.1690  1.093542  0.359822  3.039122 0.00237269 0.0296586   227.461
gene4   46.6957 -0.578297  0.403091 -1.434658 0.15138467 0.3983807   203.810
gene5  102.2155 -0.769893  0.395222 -1.948000 0.05141499 0.1606719   231.640
gene6  104.9993 -0.179747  0.386673 -0.464855 0.64203543 0.8229957   237.464
            BIC
      <numeric>
gene1   197.001
gene2   207.930
gene3   234.431
gene4   210.781
gene5   238.611
gene6   244.434

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   49.4800  -2.88178  0.983146 -2.931186 0.0033767  0.168835   190.030
gene2   42.5847   1.07252  1.073860  0.998756 0.3179131  0.662319   200.960
gene3  140.1690   1.53995  1.119216  1.375919 0.1688467  0.582952   227.461
gene4   46.6957   1.12347  1.266944  0.886759 0.3752089  0.721556   203.810
gene5  102.2155  -2.37423  1.238384 -1.917199 0.0552127  0.552127   231.640
gene6  104.9993   1.42148  1.212236  1.172612 0.2409514  0.582952   237.464
            BIC
      <numeric>
gene1   197.001
gene2   207.930
gene3   234.431
gene4   210.781
gene5   238.611
gene6   244.434

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>
gene21   78.5311   1.00006   9.78161 0.00176305 0.0881526   221.408   228.379
gene23   64.4750   1.00005   7.03745 0.00798530 0.1799331   209.306   216.276
gene10   32.2968   1.00005   6.19349 0.01282573 0.1799331   184.630   191.600
gene5   102.2155   1.00040   5.99033 0.01439465 0.1799331   231.640   238.611
gene39  121.2724   1.00013   5.45928 0.01947697 0.1947697   232.276   239.246
gene27   59.0689   1.00025   4.73751 0.02954236 0.2103366   217.156   224.127
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 RC (2024-04-16 r86468)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 22.04.4 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.20-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.5.1               BiocParallel_1.39.0        
 [3] NBAMSeq_1.21.0              SummarizedExperiment_1.35.0
 [5] Biobase_2.65.0              GenomicRanges_1.57.0       
 [7] GenomeInfoDb_1.41.0         IRanges_2.39.0             
 [9] S4Vectors_0.43.0            BiocGenerics_0.51.0        
[11] MatrixGenerics_1.17.0       matrixStats_1.3.0          

loaded via a namespace (and not attached):
 [1] KEGGREST_1.45.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.67.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.73.0       munsell_0.5.1           DESeq2_1.45.0          
[25] codetools_0.2-20        htmltools_0.5.8.1       sass_0.4.9             
[28] yaml_2.3.8              pillar_1.9.0            crayon_1.5.2           
[31] jquerylib_0.1.4         DelayedArray_0.31.0     cachem_1.0.8           
[34] abind_1.4-5             nlme_3.1-164            genefilter_1.87.0      
[37] tidyselect_1.2.1        locfit_1.5-9.9          digest_0.6.35          
[40] dplyr_1.1.4             labeling_0.4.3          splines_4.4.0          
[43] fastmap_1.1.1           grid_4.4.0              colorspace_2.1-0       
[46] cli_3.6.2               SparseArray_1.5.0       magrittr_2.0.3         
[49] S4Arrays_1.5.0          survival_3.6-4          XML_3.99-0.16.1        
[52] utf8_1.2.4              withr_3.0.0             scales_1.3.0           
[55] UCSC.utils_1.1.0        bit64_4.0.5             rmarkdown_2.26         
[58] XVector_0.45.0          httr_1.4.7              bit_4.0.5              
[61] png_0.1-8               memoise_2.0.1           evaluate_0.23          
[64] knitr_1.46              mgcv_1.9-1              rlang_1.1.3            
[67] Rcpp_1.0.12             DBI_1.2.2               xtable_1.8-4           
[70] glue_1.7.0              annotate_1.83.0         jsonlite_1.8.8         
[73] R6_2.5.1                zlibbioc_1.51.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.