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 184 66 139 1 34 432 1 158 18
gene2 1 620 192 29 2 47 67 362 836
gene3 3 26 149 333 361 25 235 68 253
gene4 26 5 1 1 1 29 71 27 20
gene5 1 147 73 237 505 131 932 73 214
gene6 3 278 117 228 59 7 201 5 369
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 40 9 4 161 57 40 113 9
gene2 50 21 161 25 191 41 18 1
gene3 161 138 12 34 3 190 1 24
gene4 3 277 287 269 1 7 14 7
gene5 7 61 6 213 31 4 785 102
gene6 3 13 1 69 290 11 35 1
sample18 sample19 sample20
gene1 407 44 5
gene2 85 8 239
gene3 13 72 1
gene4 436 7 9
gene5 1 380 12
gene6 50 72 52
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 20.31646 -1.6606170 1.5752713 -0.01890820 2
sample2 62.10704 1.2470952 0.2358893 0.16123225 2
sample3 32.93064 0.7546723 1.2934512 -0.16475713 2
sample4 23.39285 1.6755144 1.3896307 -1.96937892 1
sample5 40.35618 0.3460260 -1.2403966 1.31629931 2
sample6 53.00657 0.1456312 0.8628724 0.05175416 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 76.6234 1.00004 1.815781 0.17782907 0.4507364 225.281 232.251
gene2 137.3865 1.00006 7.674797 0.00560168 0.0466807 233.924 240.894
gene3 79.5262 1.00012 0.204024 0.65158080 0.9049733 232.733 239.703
gene4 44.2820 1.00011 0.122278 0.72675935 0.9449164 200.379 207.349
gene5 137.6483 1.00030 0.344015 0.55793570 0.8399440 253.944 260.914
gene6 89.6196 1.00006 1.676411 0.19541792 0.4507364 218.649 225.620
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 76.6234 -0.576590 0.356652 -1.616677 0.1059481 0.481582 225.281
gene2 137.3865 0.529586 0.361943 1.463176 0.1434193 0.551613 233.924
gene3 79.5262 0.979159 0.384580 2.546049 0.0108950 0.181583 232.733
gene4 44.2820 -0.933481 0.402859 -2.317144 0.0204959 0.204959 200.379
gene5 137.6483 0.248385 0.410002 0.605815 0.5446376 0.839424 253.944
gene6 89.6196 -0.271223 0.358038 -0.757527 0.4487340 0.775635 218.649
BIC
<numeric>
gene1 232.251
gene2 240.894
gene3 239.703
gene4 207.349
gene5 260.914
gene6 225.620
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 76.6234 -0.742060 1.02861 -0.721421 0.470650 0.739713 225.281
gene2 137.3865 1.335780 1.02611 1.301793 0.192987 0.577180 233.924
gene3 79.5262 -0.344811 1.09317 -0.315425 0.752439 0.917609 232.733
gene4 44.2820 -1.862268 1.16448 -1.599225 0.109771 0.552737 200.379
gene5 137.6483 -0.886031 1.18322 -0.748828 0.453961 0.739713 253.944
gene6 89.6196 1.593699 1.03102 1.545753 0.122164 0.555292 218.649
BIC
<numeric>
gene1 232.251
gene2 240.894
gene3 239.703
gene4 207.349
gene5 260.914
gene6 225.620
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>
gene43 61.8037 1.00026 16.21644 0.000057179 0.00285895 211.399 218.369
gene30 79.0847 1.00006 10.64270 0.001105634 0.02074594 195.776 202.747
gene28 85.8542 1.00006 10.31284 0.001321829 0.02074594 183.449 190.419
gene46 71.4307 1.00008 9.53628 0.002015424 0.02074594 202.130 209.100
gene50 38.5462 1.00004 9.48287 0.002074594 0.02074594 181.931 188.902
gene2 137.3865 1.00006 7.67480 0.005601678 0.04668065 233.924 240.894
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.5.0 beta (2025-04-02 r88102)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.2 LTS
Matrix products: default
BLAS: /home/biocbuild/bbs-3.22-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_3.5.2 BiocParallel_1.43.0
[3] NBAMSeq_1.25.0 SummarizedExperiment_1.39.0
[5] Biobase_2.69.0 GenomicRanges_1.61.0
[7] GenomeInfoDb_1.45.0 IRanges_2.43.0
[9] S4Vectors_0.47.0 BiocGenerics_0.55.0
[11] generics_0.1.3 MatrixGenerics_1.21.0
[13] matrixStats_1.5.0
loaded via a namespace (and not attached):
[1] KEGGREST_1.49.0 gtable_0.3.6 xfun_0.52
[4] bslib_0.9.0 lattice_0.22-7 vctrs_0.6.5
[7] tools_4.5.0 parallel_4.5.0 tibble_3.2.1
[10] AnnotationDbi_1.71.0 RSQLite_2.3.9 blob_1.2.4
[13] pkgconfig_2.0.3 Matrix_1.7-3 lifecycle_1.0.4
[16] GenomeInfoDbData_1.2.14 farver_2.1.2 compiler_4.5.0
[19] Biostrings_2.77.0 munsell_0.5.1 DESeq2_1.49.0
[22] codetools_0.2-20 htmltools_0.5.8.1 sass_0.4.10
[25] yaml_2.3.10 pillar_1.10.2 crayon_1.5.3
[28] jquerylib_0.1.4 DelayedArray_0.35.0 cachem_1.1.0
[31] abind_1.4-8 nlme_3.1-168 genefilter_1.91.0
[34] tidyselect_1.2.1 locfit_1.5-9.12 digest_0.6.37
[37] dplyr_1.1.4 labeling_0.4.3 splines_4.5.0
[40] fastmap_1.2.0 grid_4.5.0 colorspace_2.1-1
[43] cli_3.6.4 SparseArray_1.9.0 magrittr_2.0.3
[46] S4Arrays_1.9.0 survival_3.8-3 XML_3.99-0.18
[49] withr_3.0.2 scales_1.3.0 UCSC.utils_1.5.0
[52] bit64_4.6.0-1 rmarkdown_2.29 XVector_0.49.0
[55] httr_1.4.7 bit_4.6.0 png_0.1-8
[58] memoise_2.0.1 evaluate_1.0.3 knitr_1.50
[61] mgcv_1.9-3 rlang_1.1.6 Rcpp_1.0.14
[64] xtable_1.8-4 glue_1.8.0 DBI_1.2.3
[67] annotate_1.87.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.