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 90 406 38 271 119 88 47 1295 1
gene2 134 141 1 15 468 62 2 7 22
gene3 165 53 141 176 1 23 45 121 2
gene4 74 73 46 3 1 98 5 128 744
gene5 663 125 11 40 1 37 519 2 201
gene6 5 5 1 41 1 6 174 32 2
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 1 450 40 157 49 10 1 58
gene2 196 55 121 65 1 927 3 222
gene3 3 1 61 192 117 1281 1 1
gene4 398 54 93 81 1 7 6 115
gene5 18 1 204 36 1 170 3 274
gene6 477 301 95 2 1 124 7 3
sample18 sample19 sample20
gene1 2 22 64
gene2 156 9 71
gene3 2 100 8
gene4 65 4 56
gene5 1 49 13
gene6 58 2 48
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 57.32107 -0.6853969 -1.5421821 0.3352787 2
sample2 62.26834 0.5617706 -1.1061032 -1.5320630 0
sample3 55.28553 -0.2017167 0.8937113 2.4003258 1
sample4 67.48336 0.7136253 -1.2438200 -1.2035347 1
sample5 22.90693 -1.7536426 -1.1405534 0.3556751 0
sample6 32.03116 1.5485640 -0.3297687 -1.6364384 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 156.0128 1.00019 4.29228 0.03830614 0.191531 245.428 252.399
gene2 103.5116 1.00009 6.19014 0.01285210 0.102446 226.567 233.537
gene3 94.4109 1.00007 9.18645 0.00243978 0.040663 220.150 227.121
gene4 70.0585 1.00005 6.72255 0.00952395 0.102446 221.148 228.118
gene5 80.6402 1.00016 10.11317 0.00147503 0.040663 222.714 229.684
gene6 58.8544 1.00008 1.54503 0.21391566 0.538125 205.526 212.497
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 156.0128 -0.0072981 0.444766 -0.0164089 0.98690822 0.993996 245.428
gene2 103.5116 -0.2369971 0.356386 -0.6650016 0.50604944 0.766742 226.567
gene3 94.4109 1.2305044 0.445315 2.7632214 0.00572339 0.133241 220.150
gene4 70.0585 0.6585564 0.361292 1.8227806 0.06833663 0.310621 221.148
gene5 80.6402 0.2166535 0.422408 0.5129014 0.60802031 0.844473 222.714
gene6 58.8544 0.1470811 0.453130 0.3245894 0.74549188 0.924676 205.526
BIC
<numeric>
gene1 252.399
gene2 233.537
gene3 227.121
gene4 228.118
gene5 229.684
gene6 212.497
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 156.0128 -0.4971640 1.045751 -0.4754131 0.63449252 0.9492159 245.428
gene2 103.5116 -2.7442848 0.845270 -3.2466375 0.00116777 0.0194628 226.567
gene3 94.4109 -0.0498924 1.032432 -0.0483252 0.96145709 0.9674059 220.150
gene4 70.0585 -2.1795518 0.853567 -2.5534642 0.01066572 0.1066572 221.148
gene5 80.6402 0.5150649 0.991415 0.5195251 0.60339461 0.9492159 222.714
gene6 58.8544 0.0434989 1.064533 0.0408620 0.96740593 0.9674059 205.526
BIC
<numeric>
gene1 252.399
gene2 233.537
gene3 227.121
gene4 228.118
gene5 229.684
gene6 212.497
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>
gene5 80.6402 1.00016 10.11317 0.00147503 0.040663 222.714 229.684
gene9 45.4798 1.00005 9.60305 0.00194345 0.040663 183.569 190.539
gene3 94.4109 1.00007 9.18645 0.00243978 0.040663 220.150 227.121
gene4 70.0585 1.00005 6.72255 0.00952395 0.102446 221.148 228.118
gene22 67.3830 1.00011 6.53992 0.01055810 0.102446 214.666 221.636
gene2 103.5116 1.00009 6.19014 0.01285210 0.102446 226.567 233.537
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 RC (2025-04-04 r88129)
Platform: aarch64-apple-darwin20
Running under: macOS Ventura 13.7.5
Matrix products: default
BLAS: /Library/Frameworks/R.framework/Versions/4.5-arm64/Resources/lib/libRblas.0.dylib
LAPACK: /Library/Frameworks/R.framework/Versions/4.5-arm64/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_3.5.2 BiocParallel_1.42.0
[3] NBAMSeq_1.24.1 SummarizedExperiment_1.38.1
[5] Biobase_2.68.0 GenomicRanges_1.60.0
[7] GenomeInfoDb_1.44.0 IRanges_2.42.0
[9] S4Vectors_0.46.0 BiocGenerics_0.54.0
[11] generics_0.1.3 MatrixGenerics_1.20.0
[13] matrixStats_1.5.0
loaded via a namespace (and not attached):
[1] KEGGREST_1.48.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.70.0 RSQLite_2.3.11 blob_1.2.4
[13] pkgconfig_2.0.3 Matrix_1.7-3 RColorBrewer_1.1-3
[16] lifecycle_1.0.4 GenomeInfoDbData_1.2.14 compiler_4.5.0
[19] farver_2.1.2 Biostrings_2.76.0 DESeq2_1.48.1
[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.34.1 cachem_1.1.0
[31] abind_1.4-8 nlme_3.1-168 genefilter_1.90.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 cli_3.6.5
[43] SparseArray_1.8.0 magrittr_2.0.3 S4Arrays_1.8.0
[46] survival_3.8-3 dichromat_2.0-0.1 XML_3.99-0.18
[49] withr_3.0.2 scales_1.4.0 UCSC.utils_1.4.0
[52] bit64_4.6.0-1 rmarkdown_2.29 XVector_0.48.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.86.0 jsonlite_2.0.0 R6_2.6.1