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     447       6     189     173     660       6      57      14      59
gene2       1       1      13      19       2      57     171      39      16
gene3      83       2       3      26       4      70      94     203     166
gene4     188      63     108       3     317     174       2      10     190
gene5      46       1      10      82       2       2      11       1       2
gene6      64      66     852     151      78      88       4       1      30
      sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1      511      123      181        1        4      308       67      121
gene2        8       50      101       25       51      284       13      121
gene3      114      553      261       41        9       82       49      398
gene4      106      129      279       56      250      466       39        3
gene5      119       34        2       48        1      547        8        5
gene6       33        1      376        2        1      204        8      509
      sample18 sample19 sample20
gene1        5        9      174
gene2       31       27        1
gene3       45      106      103
gene4      399      204       17
gene5       26      141      132
gene6       16        3     1103

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 30.66316  0.9197954 -0.67261583 -1.223635570    2
sample2 69.57554  1.5787574 -1.66481889  0.221199611    1
sample3 34.17805 -0.3004315  1.19648809  0.003652822    0
sample4 67.37383 -1.2610021 -0.84742589  0.871489633    1
sample5 34.64948 -2.1675244 -1.30741001 -0.288487532    2
sample6 47.88003 -0.5407355  0.07943524 -0.185664564    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  128.4559   1.00035 1.7873865 0.18124213 0.4491028   246.360   253.330
gene2   39.1156   1.00007 0.3961528 0.52910008 0.7585839   203.051   210.021
gene3   89.9361   1.00005 0.0826494 0.77384184 0.8657102   239.651   246.621
gene4  129.6968   1.00008 3.5926601 0.05804429 0.3040951   253.006   259.976
gene5   38.5948   1.00004 1.1605420 0.28138026 0.5210746   182.995   189.965
gene6  111.6350   1.00014 7.9796805 0.00473321 0.0788868   230.549   237.519

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  128.4559 -0.2811095  0.287376 -0.978195  0.327978  0.729785   246.360
gene2   39.1156 -0.1891584  0.266038 -0.711020  0.477072  0.822538   203.051
gene3   89.9361  0.4073337  0.253912  1.604229  0.108664  0.417937   239.651
gene4  129.6968 -0.2213903  0.271676 -0.814907  0.415126  0.741296   253.006
gene5   38.5948  0.3153874  0.263116  1.198664  0.230658  0.576646   182.995
gene6  111.6350  0.0736847  0.307060  0.239968  0.810355  0.920858   230.549
            BIC
      <numeric>
gene1   253.330
gene2   210.021
gene3   246.621
gene4   259.976
gene5   189.965
gene6   237.519

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  128.4559  0.287457  0.846404  0.339621 0.73414176  0.872092   246.360
gene2   39.1156 -0.779832  0.784121 -0.994530 0.31996506  0.734333   203.051
gene3   89.9361 -0.781760  0.746888 -1.046689 0.29524299  0.734333   239.651
gene4  129.6968  0.776996  0.801765  0.969106 0.33249217  0.734333   253.006
gene5   38.5948  2.462611  0.806275  3.054306 0.00225582  0.112791   182.995
gene6  111.6350 -2.040184  0.903882 -2.257134 0.02399970  0.194605   230.549
            BIC
      <numeric>
gene1   253.330
gene2   210.021
gene3   246.621
gene4   259.976
gene5   189.965
gene6   237.519

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>
gene46   85.9535   1.00006  11.80775 0.000590149 0.0295075   212.849   219.819
gene19  100.8409   1.00005   8.42182 0.003708768 0.0788868   225.344   232.315
gene6   111.6350   1.00014   7.97968 0.004733210 0.0788868   230.549   237.519
gene10   99.5487   1.00006   7.05229 0.007919015 0.0796559   220.348   227.318
gene43  104.3044   1.00016   7.04252 0.007965585 0.0796559   216.586   223.557
gene22   62.2711   1.00003   6.23009 0.012563784 0.1046982   205.659   212.629
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 beta (2024-04-15 r86425)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 22.04.4 LTS

Matrix products: default
BLAS:   /home/biocbuild/bbs-3.19-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.0               BiocParallel_1.37.1        
 [3] NBAMSeq_1.19.0              SummarizedExperiment_1.33.3
 [5] Biobase_2.63.1              GenomicRanges_1.55.4       
 [7] GenomeInfoDb_1.39.14        IRanges_2.37.1             
 [9] S4Vectors_0.41.6            BiocGenerics_0.49.1        
[11] MatrixGenerics_1.15.1       matrixStats_1.3.0          

loaded via a namespace (and not attached):
 [1] KEGGREST_1.43.0         gtable_0.3.4            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.65.2    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.71.5       munsell_0.5.1           DESeq2_1.43.5          
[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.29.9     cachem_1.0.8           
[34] abind_1.4-5             nlme_3.1-164            genefilter_1.85.1      
[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.3.5       magrittr_2.0.3         
[49] S4Arrays_1.3.7          survival_3.5-8          XML_3.99-0.16.1        
[52] utf8_1.2.4              withr_3.0.0             scales_1.3.0           
[55] UCSC.utils_0.99.7       bit64_4.0.5             rmarkdown_2.26         
[58] XVector_0.43.1          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.81.2         jsonlite_1.8.8         
[73] R6_2.5.1                zlibbioc_1.49.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.