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 19 15 32 5 3 256 206 73 152
gene2 1 29 29 5 1 710 244 5 19
gene3 106 126 1 72 56 48 127 1 3
gene4 27 13 1124 226 2 5 310 1 296
gene5 373 34 1 1 57 4 12 1 1
gene6 5 58 1 1 183 107 81 1268 228
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 2 1 27 62 309 83 21 85
gene2 113 1 52 2 9 1 12 13
gene3 41 50 560 56 44 1 15 5
gene4 303 8 210 3 3 238 1 8
gene5 8 2 18 13 34 11 473 14
gene6 67 273 128 110 79 135 346 245
sample18 sample19 sample20
gene1 634 47 6
gene2 192 11 1
gene3 2 10 42
gene4 791 162 1
gene5 316 65 457
gene6 33 16 8
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 27.72766 1.5486634 0.1701477 -1.1598660 0
sample2 55.32085 0.1250517 2.3123090 0.1174632 2
sample3 57.62663 0.3538677 0.4820402 0.5747464 2
sample4 28.95509 -0.2606925 1.2875045 0.4958720 1
sample5 58.89602 -1.9020619 -0.1418307 -0.7944310 2
sample6 47.17494 0.9479394 0.2035777 -0.2471525 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:
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 92.6186 1.00016 1.1866979 0.276048 0.761134 232.136 239.106
gene2 89.9338 1.00010 0.0863429 0.769049 0.967222 199.202 206.172
gene3 62.5782 1.00014 0.2371985 0.626346 0.967222 216.614 223.584
gene4 172.5156 1.00038 0.0528435 0.818361 0.967222 239.867 246.837
gene5 92.2048 1.00011 0.0480500 0.826715 0.967222 214.935 221.905
gene6 114.1850 1.00023 0.7144809 0.398002 0.904550 250.043 257.013
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 92.6186 0.5756334 0.535281 1.0753863 0.28220182 0.6413678 232.136
gene2 89.9338 1.5076918 0.600984 2.5087051 0.01211746 0.0865533 199.202
gene3 62.5782 -0.2331701 0.546836 -0.4263985 0.66981748 0.9089865 216.614
gene4 172.5156 1.8076984 0.666389 2.7126760 0.00667423 0.0777521 239.867
gene5 92.2048 0.5816079 0.605222 0.9609821 0.33656118 0.6731224 214.935
gene6 114.1850 0.0291906 0.505099 0.0577918 0.95391447 0.9671968 250.043
BIC
<numeric>
gene1 239.106
gene2 206.172
gene3 223.584
gene4 246.837
gene5 221.905
gene6 257.013
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 92.6186 0.773449 1.06342 0.727320 0.467030 0.686809 232.136
gene2 89.9338 3.090815 1.20817 2.558252 0.010520 0.131500 199.202
gene3 62.5782 -1.503919 1.08568 -1.385238 0.165980 0.593131 216.614
gene4 172.5156 1.450402 1.31918 1.099469 0.271563 0.686236 239.867
gene5 92.2048 0.361470 1.20275 0.300536 0.763768 0.844027 214.935
gene6 114.1850 -0.394863 1.00374 -0.393392 0.694030 0.835123 250.043
BIC
<numeric>
gene1 239.106
gene2 206.172
gene3 223.584
gene4 246.837
gene5 221.905
gene6 257.013
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>
gene17 61.0425 1.00003 14.56663 0.000135698 0.00678491 172.904 179.874
gene8 100.5595 1.00014 10.06430 0.001513774 0.03784435 216.902 223.872
gene9 116.1673 1.00005 6.18501 0.012887067 0.21478446 227.740 234.711
gene23 136.1298 1.00007 5.61030 0.017858099 0.22322623 205.764 212.734
gene40 111.8965 1.00013 5.12068 0.023663514 0.23663514 241.243 248.213
gene42 48.7367 1.00010 3.91510 0.047873535 0.31673854 198.574 205.544
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.2 (2025-10-31)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.3 LTS
Matrix products: default
BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
locale:
[1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
[3] LC_TIME=en_US.UTF-8 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: Etc/UTC
tzcode source: system (glibc)
attached base packages:
[1] stats4 stats graphics grDevices utils datasets methods
[8] base
other attached packages:
[1] ggplot2_4.0.1 BiocParallel_1.44.0
[3] NBAMSeq_1.26.0 SummarizedExperiment_1.40.0
[5] Biobase_2.70.0 GenomicRanges_1.62.1
[7] Seqinfo_1.0.0 IRanges_2.44.0
[9] S4Vectors_0.48.0 BiocGenerics_0.56.0
[11] generics_0.1.4 MatrixGenerics_1.22.0
[13] matrixStats_1.5.0 rmarkdown_2.30
loaded via a namespace (and not attached):
[1] KEGGREST_1.50.0 gtable_0.3.6 xfun_0.56
[4] bslib_0.10.0 lattice_0.22-7 vctrs_0.7.1
[7] tools_4.5.2 parallel_4.5.2 AnnotationDbi_1.72.0
[10] RSQLite_2.4.5 blob_1.3.0 Matrix_1.7-4
[13] RColorBrewer_1.1-3 S7_0.2.1 lifecycle_1.0.5
[16] compiler_4.5.2 farver_2.1.2 Biostrings_2.78.0
[19] DESeq2_1.50.2 codetools_0.2-20 htmltools_0.5.9
[22] sys_3.4.3 buildtools_1.0.0 sass_0.4.10
[25] yaml_2.3.12 crayon_1.5.3 jquerylib_0.1.4
[28] DelayedArray_0.36.0 cachem_1.1.0 abind_1.4-8
[31] nlme_3.1-168 genefilter_1.92.0 locfit_1.5-9.12
[34] digest_0.6.39 labeling_0.4.3 splines_4.5.2
[37] maketools_1.3.2 fastmap_1.2.0 grid_4.5.2
[40] cli_3.6.5 SparseArray_1.10.8 S4Arrays_1.10.1
[43] survival_3.8-6 XML_3.99-0.20 withr_3.0.2
[46] scales_1.4.0 bit64_4.6.0-1 XVector_0.50.0
[49] httr_1.4.7 bit_4.6.0 png_0.1-8
[52] memoise_2.0.1 evaluate_1.0.5 knitr_1.51
[55] mgcv_1.9-4 rlang_1.1.7 Rcpp_1.1.1
[58] xtable_1.8-4 glue_1.8.0 DBI_1.2.3
[61] annotate_1.88.0 jsonlite_2.0.0 R6_2.6.1