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 95 75 25 15 1 13 9 122 78
gene2 5 4 1 17 14 32 179 1 1
gene3 2 112 49 29 737 41 189 52 13
gene4 12 4 60 426 502 3 315 88 379
gene5 40 2 10 21 18 24 1 78 63
gene6 455 55 155 482 317 1 64 31 1
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 1 505 1 53 18 1 297 180
gene2 109 2 1 142 1 45 1 3
gene3 3 66 433 491 1 44 34 1
gene4 229 93 2 10 33 7 269 60
gene5 137 59 154 45 38 3 309 426
gene6 1 176 7 1 2 1 21 1
sample18 sample19 sample20
gene1 117 42 70
gene2 42 425 46
gene3 21 41 7
gene4 124 256 4
gene5 24 6 773
gene6 282 2 584
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.34509 1.3950252 -0.86961843 1.559718431 2
sample2 38.82249 0.6909785 -1.14926454 -0.005511607 2
sample3 22.28316 2.3473260 0.88241085 -0.624859041 2
sample4 52.22037 0.5876980 -0.77637298 -0.376341252 0
sample5 33.31987 -1.6537328 0.44375437 -0.788551373 1
sample6 50.81314 1.5055999 -0.05318027 0.281177932 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 68.2967 1.00006 1.3829747 0.239630 0.499367 225.925 232.896
gene2 42.5276 1.00013 2.0997861 0.147363 0.449838 194.148 201.118
gene3 87.9044 1.00012 0.3816055 0.536811 0.813349 231.892 238.862
gene4 91.8586 1.00012 0.6792965 0.409926 0.732010 240.280 247.250
gene5 87.6959 1.00053 0.0211743 0.885951 0.914422 232.614 239.585
gene6 88.3517 1.00013 1.9435339 0.163294 0.449838 218.065 225.035
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 68.2967 0.476100 0.454574 1.047353 0.2949370 0.641442 225.925
gene2 42.5276 -0.409998 0.501045 -0.818286 0.4131938 0.645615 194.148
gene3 87.9044 -0.553788 0.458570 -1.207640 0.2271859 0.641442 231.892
gene4 91.8586 -0.388925 0.403980 -0.962733 0.3356814 0.641442 240.280
gene5 87.6959 -0.621959 0.440800 -1.410978 0.1582510 0.565182 232.614
gene6 88.3517 0.925041 0.494261 1.871563 0.0612671 0.340373 218.065
BIC
<numeric>
gene1 232.896
gene2 201.118
gene3 238.862
gene4 247.250
gene5 239.585
gene6 225.035
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 68.2967 -1.200968 1.12625 -1.0663387 0.286271 0.624520 225.925
gene2 42.5276 -0.840691 1.24075 -0.6775676 0.498046 0.690239 194.148
gene3 87.9044 0.319122 1.13285 0.2816989 0.778174 0.849418 231.892
gene4 91.8586 -0.907549 0.99979 -0.9077401 0.364016 0.624520 240.280
gene5 87.6959 -0.838869 1.09044 -0.7692931 0.441719 0.675674 232.614
gene6 88.3517 -0.113209 1.22939 -0.0920857 0.926630 0.943969 218.065
BIC
<numeric>
gene1 232.896
gene2 201.118
gene3 238.862
gene4 247.250
gene5 239.585
gene6 225.035
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>
gene46 62.2706 1.00009 13.20936 0.000278748 0.0139374 205.376 212.346
gene38 131.1233 1.00006 8.53648 0.003481852 0.0654754 225.018 231.988
gene33 74.7865 1.00007 8.31729 0.003928524 0.0654754 199.203 206.173
gene31 61.9187 1.00138 7.31453 0.006850542 0.0856318 209.838 216.810
gene7 68.4228 1.00021 5.97739 0.014491945 0.1294200 214.149 221.119
gene24 130.4832 1.00005 5.85593 0.015530405 0.1294200 207.197 214.167
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 r88126)
Platform: aarch64-apple-darwin20
Running under: macOS Ventura 13.7.1
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.0 SummarizedExperiment_1.38.0
[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.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.76.0 munsell_0.5.1 DESeq2_1.48.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.34.0 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 colorspace_2.1-1
[43] cli_3.6.4 SparseArray_1.8.0 magrittr_2.0.3
[46] S4Arrays_1.8.0 survival_3.8-3 XML_3.99-0.18
[49] withr_3.0.2 scales_1.3.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