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 63 51 1 1 1 59 28 6
gene2 24 298 116 1 63 141 41 175 87
gene3 1 100 6 1 4 193 7 164 262
gene4 89 34 74 253 93 18 1 75 8
gene5 1 1 58 19 24 339 19 62 2
gene6 27 524 87 3 193 13 538 41 266
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 1 3 1 110 2 14 3 1
gene2 36 149 1 331 2 93 9 129
gene3 11 4 1 44 213 66 100 38
gene4 49 62 10 14 3 30 122 64
gene5 1 1 79 536 23 5 268 80
gene6 9 96 287 210 2 61 44 69
sample18 sample19 sample20
gene1 1 252 2
gene2 292 6 213
gene3 95 36 3
gene4 2 16 9
gene5 1 98 20
gene6 3 147 132
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 43.10920 0.2726211 -0.2549687 -1.3950122 0
sample2 42.58649 -0.7243301 -0.3191373 0.3616755 1
sample3 22.25855 -1.3578757 -0.1776224 0.6299800 0
sample4 56.19190 -0.7083012 -0.5424075 -1.5691176 0
sample5 40.19584 0.1119836 -1.7301169 0.1534927 0
sample6 47.54824 -1.6438241 -2.0753446 1.6084766 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 27.4310 1.00005 4.162602 0.0413411 0.172254 163.304 170.275
gene2 100.3415 1.00007 0.693298 0.4050619 0.823670 235.902 242.873
gene3 47.4098 1.00011 0.195486 0.6584270 0.840060 212.721 219.691
gene4 51.0118 1.00008 2.091807 0.1480910 0.411364 208.066 215.036
gene5 58.7841 1.00005 0.191631 0.6616413 0.840060 204.578 211.548
gene6 129.5763 1.00008 1.176630 0.2780996 0.662142 243.446 250.417
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 27.4310 -1.881073 0.550155 -3.419170 0.000628124 0.0222734 163.304
gene2 100.3415 -0.282401 0.473662 -0.596208 0.551036277 0.6719955 235.902
gene3 47.4098 0.334259 0.520719 0.641919 0.520925879 0.6719955 212.721
gene4 51.0118 -0.267819 0.447428 -0.598573 0.549457546 0.6719955 208.066
gene5 58.7841 -1.146769 0.512308 -2.238436 0.025192629 0.1259631 204.578
gene6 129.5763 -0.281757 0.469958 -0.599536 0.548815785 0.6719955 243.446
BIC
<numeric>
gene1 170.275
gene2 242.873
gene3 219.691
gene4 215.036
gene5 211.548
gene6 250.417
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 27.4310 -0.755928 1.050079 -0.719877 0.4716008 0.758044 163.304
gene2 100.3415 1.245788 0.895712 1.390835 0.1642754 0.410689 235.902
gene3 47.4098 0.394419 0.984574 0.400599 0.6887157 0.845329 212.721
gene4 51.0118 0.755892 0.846193 0.893286 0.3717041 0.665280 208.066
gene5 58.7841 -0.960163 0.970318 -0.989535 0.3224016 0.665280 204.578
gene6 129.5763 -1.769738 0.890084 -1.988283 0.0467804 0.272979 243.446
BIC
<numeric>
gene1 170.275
gene2 242.873
gene3 219.691
gene4 215.036
gene5 211.548
gene6 250.417
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
<numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
gene16 111.2252 1.00010 21.24697 4.31629e-06 0.000215815 213.559
gene8 53.0768 1.00007 6.84529 8.89135e-03 0.116864266 186.283
gene31 55.4110 1.00003 6.40008 1.14133e-02 0.116864266 202.947
gene22 88.5019 1.00018 6.37614 1.15813e-02 0.116864266 218.341
gene7 51.8110 1.00006 6.15125 1.31353e-02 0.116864266 197.548
gene41 37.1682 1.00008 6.03522 1.40237e-02 0.116864266 187.539
BIC
<numeric>
gene16 220.529
gene8 193.254
gene31 209.917
gene22 225.312
gene7 204.518
gene41 194.509
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 Patched (2025-04-21 r88169)
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.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 RColorBrewer_1.1-3
[16] lifecycle_1.0.4 GenomeInfoDbData_1.2.14 compiler_4.5.0
[19] farver_2.1.2 Biostrings_2.77.0 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.1 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 cli_3.6.5
[43] SparseArray_1.9.0 magrittr_2.0.3 S4Arrays_1.9.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.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