To install and load NBAMSeq
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:
Step 1: Data input using NBAMSeqDataSet;
Step 2: Differential expression (DE) analysis using NBAMSeq function;
Step 3: Pulling out DE results using results function.
Here we illustrate each of these steps respectively.
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 5 6 125 17 130 44 11 8 80
gene2 8 1 92 3 3 53 126 54 6
gene3 438 172 315 107 152 157 8 111 1
gene4 528 6 67 81 475 24 107 9 8
gene5 189 108 11 11 1 3 142 119 171
gene6 1 1 1 3 206 47 1 76 398
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 242 15 6 1 1 8 55 1
gene2 508 1 149 78 22 136 600 1
gene3 92 36 1773 94 21 28 21 16
gene4 2 9 36 2 9 11 303 706
gene5 280 61 349 143 19 7 96 72
gene6 749 26 1 331 139 103 397 1
sample18 sample19 sample20
gene1 53 47 1
gene2 27 1 2
gene3 622 24 1
gene4 29 698 6
gene5 25 3 1
gene6 923 19 241
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 61.46436 0.48488642 0.44497161 -0.4563828 1
sample2 73.87543 -0.51295750 -0.36834933 -0.1866260 1
sample3 73.90741 0.05099374 0.18196949 -1.4098181 0
sample4 62.69165 -0.09175661 -0.00233061 -0.7547526 0
sample5 37.00102 -0.16130247 -0.81531746 -0.3064635 0
sample6 69.70738 -0.74887019 1.61315887 -0.4502291 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:
Several notes should be made regarding the design formula:
multiple nonlinear covariates are supported, e.g. design = ~ s(pheno) + s(var1) + var2 + var3 + var4;
the nonlinear covariate cannot be a discrete variable, e.g. design = ~ s(pheno) + var1 + var2 + var3 + s(var4) as var4 is a factor, and it makes no sense to model a factor as nonlinear;
at least one nonlinear covariate should be provided in design. If all covariates are assumed to have linear effect on gene count, use DESeq2 (Love, Huber, and Anders 2014), edgeR (Robinson, McCarthy, and Smyth 2010), NBPSeq (Di et al. 2015) or BBSeq (Zhou, Xia, and Wright 2011) instead. e.g. design = ~ pheno + var1 + var2 + var3 + var4 is not supported in NBAMSeq;
design matrix is not supported.
We then construct the NBAMSeqDataSet using countData, colData, and design:
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 can be performed by NBAMSeq function:
Several other arguments in NBAMSeq function are available for users to customize the analysis.
gamma argument can be used to control the smoothness of the nonlinear function. Higher gamma means the nonlinear function will be more smooth. See the gamma argument of gam function in mgcv (Wood and Wood 2015) for details. Default gamma is 2.5;
fitlin is either TRUE or FALSE indicating whether linear model should be fitted after fitting the nonlinear model;
parallel is either TRUE or FALSE indicating whether parallel should be used. e.g. Run NBAMSeq with parallel = TRUE:
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.
DataFrame with 6 rows and 7 columns
baseMean edf stat pvalue padj AIC BIC
<numeric> <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
gene1 32.8470 1.00011 0.186055 0.66633142 0.8235049 190.323 197.293
gene2 73.5573 1.00042 0.555048 0.45692563 0.7648026 214.238 221.209
gene3 174.6556 1.00011 3.519373 0.06066466 0.2333256 253.414 260.384
gene4 107.8852 1.00010 0.347691 0.55545188 0.7935027 222.088 229.058
gene5 75.7354 1.00004 0.608087 0.43551530 0.7648026 227.350 234.320
gene6 118.2652 1.00005 8.230969 0.00411899 0.0514873 233.593 240.563
For linear continuous covariates, base mean, estimated coefficient, standard error, test statistics, p-value, and adjusted p-value will be returned.
DataFrame with 6 rows and 8 columns
baseMean coef SE stat pvalue padj AIC
<numeric> <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
gene1 32.8470 -0.281894 0.548373 -0.514054 0.6072141 0.759018 190.323
gene2 73.5573 -1.280929 0.672754 -1.904009 0.0569090 0.304808 214.238
gene3 174.6556 0.792067 0.612348 1.293490 0.1958416 0.425743 253.414
gene4 107.8852 0.712152 0.531107 1.340882 0.1799587 0.411976 222.088
gene5 75.7354 0.338665 0.588235 0.575731 0.5647969 0.754633 227.350
gene6 118.2652 -1.535893 0.674420 -2.277354 0.0227651 0.284563 233.593
BIC
<numeric>
gene1 197.293
gene2 221.209
gene3 260.384
gene4 229.058
gene5 234.320
gene6 240.563
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.
DataFrame with 6 rows and 8 columns
baseMean coef SE stat pvalue padj AIC
<numeric> <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
gene1 32.8470 -2.168969 0.926189 -2.341821 0.0191899 0.0959495 190.323
gene2 73.5573 -0.919149 1.119646 -0.820928 0.4116871 0.6432610 214.238
gene3 174.6556 -1.375356 1.015699 -1.354099 0.1757049 0.3514098 253.414
gene4 107.8852 -2.682472 0.888597 -3.018773 0.0025380 0.0423001 222.088
gene5 75.7354 -0.623436 0.976338 -0.638545 0.5231188 0.6985940 227.350
gene6 118.2652 2.747998 1.118874 2.456038 0.0140478 0.0959495 233.593
BIC
<numeric>
gene1 197.293
gene2 221.209
gene3 260.384
gene4 229.058
gene5 234.320
gene6 240.563
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>
gene31 20.7009 1.00004 10.37745 0.00127592 0.0514873 155.444 162.414
gene17 162.7021 1.00015 9.18713 0.00243817 0.0514873 216.002 222.973
gene49 129.4705 1.00009 8.71790 0.00315152 0.0514873 213.789 220.759
gene6 118.2652 1.00005 8.23097 0.00411899 0.0514873 233.593 240.563
gene38 83.5410 1.00015 5.11568 0.02373053 0.1748130 219.798 226.768
gene34 67.1388 1.00011 4.70071 0.03016560 0.1748130 215.743 222.713
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))R version 4.6.0 RC (2026-04-17 r89917)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.4 LTS
Matrix products: default
BLAS: /home/biocbuild/bbs-3.24-bioc/R/lib/libRblas.so
LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.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_4.0.2 BiocParallel_1.45.0
[3] NBAMSeq_1.27.0 SummarizedExperiment_1.41.1
[5] Biobase_2.71.0 GenomicRanges_1.63.2
[7] Seqinfo_1.1.0 IRanges_2.45.0
[9] S4Vectors_0.49.2 BiocGenerics_0.57.1
[11] generics_0.1.4 MatrixGenerics_1.23.0
[13] matrixStats_1.5.0
loaded via a namespace (and not attached):
[1] KEGGREST_1.51.1 gtable_0.3.6 xfun_0.57
[4] bslib_0.10.0 lattice_0.22-9 vctrs_0.7.3
[7] tools_4.6.0 parallel_4.6.0 tibble_3.3.1
[10] AnnotationDbi_1.73.1 RSQLite_2.4.6 blob_1.3.0
[13] pkgconfig_2.0.3 Matrix_1.7-5 RColorBrewer_1.1-3
[16] S7_0.2.1-1 lifecycle_1.0.5 compiler_4.6.0
[19] farver_2.1.2 Biostrings_2.79.5 DESeq2_1.51.7
[22] codetools_0.2-20 htmltools_0.5.9 sass_0.4.10
[25] yaml_2.3.12 crayon_1.5.3 pillar_1.11.1
[28] jquerylib_0.1.4 DelayedArray_0.37.1 cachem_1.1.0
[31] abind_1.4-8 nlme_3.1-169 genefilter_1.93.0
[34] tidyselect_1.2.1 locfit_1.5-9.12 digest_0.6.39
[37] dplyr_1.2.1 labeling_0.4.3 splines_4.6.0
[40] fastmap_1.2.0 grid_4.6.0 cli_3.6.6
[43] SparseArray_1.11.13 magrittr_2.0.5 S4Arrays_1.11.1
[46] survival_3.8-6 dichromat_2.0-0.1 XML_3.99-0.23
[49] withr_3.0.2 scales_1.4.0 bit64_4.6.0-1
[52] rmarkdown_2.31 XVector_0.51.0 httr_1.4.8
[55] bit_4.6.0 otel_0.2.0 png_0.1-9
[58] memoise_2.0.1 evaluate_1.0.5 knitr_1.51
[61] mgcv_1.9-4 rlang_1.2.0 Rcpp_1.1.1-1
[64] xtable_1.8-8 glue_1.8.1 DBI_1.3.0
[67] annotate_1.89.0 jsonlite_2.0.0 R6_2.6.1
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.