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      53      26      10       1       3      24      36       9     717
gene2      15       1      24       1       3      87       4     165     114
gene3       1     144      86      61       8       1       1     102      23
gene4      63       1      89       1       4       1      15      26       7
gene5       1       2       5      73      14      79       1      11     159
gene6     493      16     227      53      13       2       1      73      29
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1      144      148        1      230       14       29        6        5
gene2        3       23       50       13      271        6      880       13
gene3      275        4       14      408        1        1       11      673
gene4        4      132       56        9        8       25       12      153
gene5      611        6        7        5       53        1       63       73
gene6      153      308      144       19      109       54      125      256
      sample18 sample19 sample20
gene1        1        1      142
gene2        3      118       75
gene3       16        2      209
gene4      126       21       63
gene5      320       70       35
gene6        1        2       40

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 76.51054 -0.6496998  0.2927702 -1.5273647    0
sample2 26.37295 -0.9217091 -0.9537890  1.5869146    0
sample3 60.35212 -2.3911142 -0.2198217 -0.2261561    1
sample4 48.53997 -0.3861458  0.2509984 -0.4601455    0
sample5 54.50890  0.4179040 -0.3389644 -1.1644001    2
sample6 24.44823  1.6967647  0.7766880 -0.1418359    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
      <numeric> <numeric>   <numeric>   <numeric>  <numeric> <numeric>
gene1   64.1040   1.34270 1.17354e+00 5.50273e-01 0.70547827   202.708
gene2   82.5379   1.00020 4.21507e-01 5.16361e-01 0.68572125   213.991
gene3   92.6232   1.00006 1.02758e+00 3.10793e-01 0.51798758   211.195
gene4   40.2564   1.00007 1.85418e+00 1.73322e-01 0.33331126   200.916
gene5   78.0204   1.00009 1.30103e-05 9.99348e-01 0.99934824   212.934
gene6  116.8766   1.00006 1.75100e+01 2.92882e-05 0.00146441   219.870
            BIC
      <numeric>
gene1   210.019
gene2   220.962
gene3   218.165
gene4   207.886
gene5   219.905
gene6   226.840

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   64.1040 -0.0910321  0.478085 -0.1904100 0.8489879  0.943320   202.708
gene2   82.5379  0.8700025  0.498938  1.7437102 0.0812096  0.290034   213.991
gene3   92.6232  0.2917115  0.520528  0.5604150 0.5751964  0.835619   211.195
gene4   40.2564  0.2587779  0.473784  0.5461934 0.5849330  0.835619   200.916
gene5   78.0204  0.9774854  0.520471  1.8780800 0.0603702  0.251543   212.934
gene6  116.8766 -0.0167402  0.398504 -0.0420075 0.9664927  0.982245   219.870
            BIC
      <numeric>
gene1   210.019
gene2   220.962
gene3   218.165
gene4   207.886
gene5   219.905
gene6   226.840

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   64.1040 -2.90309267  1.185422 -2.44899533 0.0143255  0.160638   202.708
gene2   82.5379 -0.00798848  1.207086 -0.00661798 0.9947197  0.994720   213.991
gene3   92.6232 -0.71786831  1.281608 -0.56013088 0.5753902  0.799153   211.195
gene4   40.2564  0.29066022  1.149747  0.25280360 0.8004200  0.889356   200.916
gene5   78.0204  1.03388791  1.242819  0.83188931 0.4054714  0.652694   212.934
gene6  116.8766 -2.34348371  0.973425 -2.40746276 0.0160638  0.160638   219.870
            BIC
      <numeric>
gene1   210.019
gene2   220.962
gene3   218.165
gene4   207.886
gene5   219.905
gene6   226.840

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>
gene6   116.8766   1.00006  17.51001 2.92882e-05 0.00146441   219.870   226.840
gene36   72.3637   2.09155  18.84096 4.98531e-04 0.01246327   186.942   194.999
gene45   58.9811   1.00005   8.55268 3.45180e-03 0.05752997   200.370   207.340
gene37   86.2045   1.00006   7.76642 5.32574e-03 0.06389474   218.340   225.310
gene30  117.5830   1.00010   7.22448 7.19852e-03 0.06389474   226.923   233.893
gene47  117.4313   1.00005   7.11031 7.66737e-03 0.06389474   236.993   243.963
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.1 (2023-06-16)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 22.04.3 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.4               BiocParallel_1.36.0        
 [3] NBAMSeq_1.18.0              SummarizedExperiment_1.32.0
 [5] Biobase_2.62.0              GenomicRanges_1.54.1       
 [7] GenomeInfoDb_1.38.0         IRanges_2.36.0             
 [9] S4Vectors_0.40.1            BiocGenerics_0.48.1        
[11] MatrixGenerics_1.14.0       matrixStats_1.0.0          

loaded via a namespace (and not attached):
 [1] KEGGREST_1.42.0         gtable_0.3.4            xfun_0.41              
 [4] bslib_0.5.1             lattice_0.22-5          vctrs_0.6.4            
 [7] tools_4.3.1             bitops_1.0-7            generics_0.1.3         
[10] parallel_4.3.1          RSQLite_2.3.2           AnnotationDbi_1.64.0   
[13] tibble_3.2.1            fansi_1.0.5             highr_0.10             
[16] blob_1.2.4              pkgconfig_2.0.3         Matrix_1.6-1.1         
[19] lifecycle_1.0.3         GenomeInfoDbData_1.2.11 farver_2.1.1           
[22] compiler_4.3.1          Biostrings_2.70.1       munsell_0.5.0          
[25] DESeq2_1.42.0           codetools_0.2-19        htmltools_0.5.6.1      
[28] sass_0.4.7              RCurl_1.98-1.12         yaml_2.3.7             
[31] pillar_1.9.0            crayon_1.5.2            jquerylib_0.1.4        
[34] DelayedArray_0.28.0     cachem_1.0.8            abind_1.4-5            
[37] nlme_3.1-163            genefilter_1.84.0       tidyselect_1.2.0       
[40] locfit_1.5-9.8          digest_0.6.33           dplyr_1.1.3            
[43] labeling_0.4.3          splines_4.3.1           fastmap_1.1.1          
[46] grid_4.3.1              colorspace_2.1-0        cli_3.6.1              
[49] SparseArray_1.2.0       magrittr_2.0.3          S4Arrays_1.2.0         
[52] survival_3.5-7          XML_3.99-0.14           utf8_1.2.4             
[55] withr_2.5.2             scales_1.2.1            bit64_4.0.5            
[58] rmarkdown_2.25          XVector_0.42.0          httr_1.4.7             
[61] bit_4.0.5               png_0.1-8               memoise_2.0.1          
[64] evaluate_0.23           knitr_1.45              mgcv_1.9-0             
[67] rlang_1.1.1             Rcpp_1.0.11             DBI_1.1.3              
[70] xtable_1.8-4            glue_1.6.2              annotate_1.80.0        
[73] jsonlite_1.8.7          R6_2.5.1                zlibbioc_1.48.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.