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 3 3 5 1 3 5 360 7 51
gene2 1 472 298 3 2 345 3 322 56
gene3 2 4 7 3 55 192 1 3 49
gene4 2 157 9 2 3 598 38 1 15
gene5 58 1 1402 1 4 179 2 291 7
gene6 48 2 392 722 59 25 2 80 30
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 185 8 59 1 3 12 98 14
gene2 112 4 110 250 475 50 2 83
gene3 28 51 708 56 69 26 7 8
gene4 39 68 290 39 1 1 135 1
gene5 164 54 63 409 52 37 7 159
gene6 2 1 165 8 91 3 7 1
sample18 sample19 sample20
gene1 1 284 97
gene2 11 7 78
gene3 31 373 51
gene4 3 195 14
gene5 1 1 1
gene6 2 21 607
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 45.95569 0.8605322 -0.07537501 0.7998795 0
sample2 76.57352 0.3216642 1.22503599 -1.5682601 1
sample3 68.17386 0.7630422 0.19026692 -0.9469453 2
sample4 67.58064 -1.7951475 -0.87237910 0.1558409 0
sample5 59.68017 0.7823192 0.79753803 -0.6131603 2
sample6 64.12638 0.5452854 -1.01082960 -0.6096498 2
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 43.2998 1.00007 1.593277 0.2068561 0.737399 197.324 204.295
gene2 141.5531 1.00008 2.913397 0.0878597 0.488109 236.791 243.761
gene3 57.3542 1.00009 0.135897 0.7125162 0.876844 216.081 223.051
gene4 69.7815 1.00010 1.250297 0.2635504 0.742491 205.427 212.397
gene5 90.8988 1.00004 0.929780 0.3349417 0.742491 220.041 227.011
gene6 86.9726 1.00008 0.726844 0.3939199 0.742491 211.963 218.933
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 43.2998 0.4524963 0.540812 0.8366978 0.4027624 0.652739 197.324
gene2 141.5531 -0.3849312 0.551171 -0.6983875 0.4849349 0.713140 236.791
gene3 57.3542 0.0507267 0.523442 0.0969099 0.9227979 0.943190 216.081
gene4 69.7815 0.5691722 0.575478 0.9890426 0.3226423 0.597486 205.427
gene5 90.8988 1.3885037 0.592285 2.3443180 0.0190619 0.111337 220.041
gene6 86.9726 -1.1999353 0.548368 -2.1881927 0.0286556 0.130253 211.963
BIC
<numeric>
gene1 204.295
gene2 243.761
gene3 223.051
gene4 212.397
gene5 227.011
gene6 218.933
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 43.2998 0.83445758 0.895014 0.93234006 0.351161 0.605450 197.324
gene2 141.5531 -0.00120466 0.913607 -0.00131858 0.998948 0.998948 236.791
gene3 57.3542 0.86476573 0.865105 0.99960726 0.317501 0.587964 216.081
gene4 69.7815 1.54061473 0.948034 1.62506339 0.104149 0.325466 205.427
gene5 90.8988 -1.16152166 0.973004 -1.19374860 0.232576 0.531972 220.041
gene6 86.9726 0.17579468 0.910090 0.19316195 0.846832 0.901437 211.963
BIC
<numeric>
gene1 204.295
gene2 243.761
gene3 223.051
gene4 212.397
gene5 227.011
gene6 218.933
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 54.0188 1.00005 14.74381 0.000123048 0.00615238 206.944 213.915
gene9 89.4493 1.00005 10.88490 0.000969998 0.02424995 217.625 224.595
gene22 53.9364 1.00008 7.98217 0.004724271 0.07873785 202.806 209.776
gene42 71.6606 1.00011 7.10536 0.007693722 0.09617152 206.771 213.741
gene24 113.0432 1.00004 5.82089 0.015841448 0.15841448 232.232 239.202
gene8 37.8614 1.00005 5.39638 0.020183743 0.16409659 175.607 182.578
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 Patched (2026-04-24 r89963)
Platform: aarch64-apple-darwin23
Running under: macOS Tahoe 26.3.1
Matrix products: default
BLAS: /Library/Frameworks/R.framework/Versions/4.6/Resources/lib/libRblas.0.dylib
LAPACK: /Library/Frameworks/R.framework/Versions/4.6/Resources/lib/libRlapack.dylib; LAPACK version 3.12.1
locale:
[1] C/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
time zone: America/New_York
tzcode source: internal
attached base packages:
[1] stats4 stats graphics grDevices utils datasets methods
[8] base
other attached packages:
[1] ggplot2_4.0.3 BiocParallel_1.46.0
[3] NBAMSeq_1.28.0 SummarizedExperiment_1.42.0
[5] Biobase_2.72.0 GenomicRanges_1.64.0
[7] Seqinfo_1.2.0 IRanges_2.46.0
[9] S4Vectors_0.50.0 BiocGenerics_0.58.0
[11] generics_0.1.4 MatrixGenerics_1.24.0
[13] matrixStats_1.5.0
loaded via a namespace (and not attached):
[1] KEGGREST_1.52.0 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.74.0 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.2 lifecycle_1.0.5 compiler_4.6.0
[19] farver_2.1.2 Biostrings_2.80.0 DESeq2_1.52.0
[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.38.0 cachem_1.1.0
[31] abind_1.4-8 nlme_3.1-169 genefilter_1.94.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.12.0 magrittr_2.0.5 S4Arrays_1.12.0
[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.8.0
[52] rmarkdown_2.31 XVector_0.52.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.1
[64] xtable_1.8-8 glue_1.8.1 DBI_1.3.0
[67] annotate_1.90.0 jsonlite_2.0.0 R6_2.6.1