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 7 2 569 1 273 1 178 53
gene2 23 17 164 5 26 30 81 173 1
gene3 66 107 23 1 81 2 5 21 17
gene4 1 121 91 3 1 6 8 500 1
gene5 58 7 100 112 307 5 216 273 144
gene6 10 341 125 1 50 3 3 1 11
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 34 267 339 3 20 102 26 231
gene2 33 1 27 1 183 395 49 11
gene3 419 379 15 334 2 1 22 1
gene4 69 36 1 356 127 1 119 1
gene5 335 118 1 309 6 70 174 174
gene6 10 6 34 141 13 322 626 14
sample18 sample19 sample20
gene1 175 49 97
gene2 1 6 1
gene3 287 1 229
gene4 12 171 16
gene5 62 241 10
gene6 23 564 98
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 36.96606 -0.0409631 -1.9711520 -1.2225683 0
sample2 28.50380 0.2950463 0.1671848 1.2710449 0
sample3 77.98221 -1.0282046 -1.4236736 0.4726612 1
sample4 59.47101 -1.4239704 -0.1264069 -2.1995088 2
sample5 30.73332 -1.0728392 0.3776688 -1.2709749 0
sample6 35.89123 -1.4921051 -0.3662688 -0.3196610 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 131.4044 1.00006 1.8021251 0.179469 0.662994 235.286 242.256
gene2 67.3163 1.00007 0.4465512 0.504036 0.984680 203.986 210.956
gene3 99.2105 1.00005 0.1459134 0.702566 0.984680 215.985 222.955
gene4 57.2586 1.00005 0.0623288 0.802964 0.984680 198.152 205.122
gene5 109.6867 1.00007 0.2936885 0.587955 0.984680 248.947 255.917
gene6 121.9993 1.04764 0.1814121 0.670523 0.984680 228.523 235.541
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 131.4044 -0.725533 0.427993 -1.695196 0.0900382 0.282036 235.286
gene2 67.3163 1.050456 0.424515 2.474487 0.0133428 0.092982 203.986
gene3 99.2105 -1.058589 0.442473 -2.392440 0.0167368 0.092982 215.985
gene4 57.2586 -0.124854 0.400813 -0.311501 0.7554199 0.858432 198.152
gene5 109.6867 0.344821 0.370245 0.931332 0.3516821 0.666748 248.947
gene6 121.9993 1.146610 0.450536 2.544991 0.0109281 0.092982 228.523
BIC
<numeric>
gene1 242.256
gene2 210.956
gene3 222.955
gene4 205.122
gene5 255.917
gene6 235.541
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 131.4044 -1.485366 1.003650 -1.479965 0.1388827 0.547083 235.286
gene2 67.3163 -0.564355 0.983649 -0.573736 0.5661467 0.884604 203.986
gene3 99.2105 -2.466048 1.040076 -2.371026 0.0177388 0.177388 215.985
gene4 57.2586 1.339206 0.937572 1.428378 0.1531831 0.547083 198.152
gene5 109.6867 0.293603 0.868443 0.338080 0.7353029 0.972068 248.947
gene6 121.9993 -0.517977 1.061886 -0.487790 0.6256987 0.920145 228.523
BIC
<numeric>
gene1 242.256
gene2 210.956
gene3 222.955
gene4 205.122
gene5 255.917
gene6 235.541
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>
gene20 106.3079 1.00007 11.13793 0.000846273 0.0423137 203.841 210.811
gene30 124.0888 1.00016 9.51175 0.002045404 0.0511351 213.178 220.149
gene45 92.4817 1.16342 4.66226 0.039174751 0.6301610 211.393 218.526
gene23 40.8184 1.00021 3.48176 0.062113842 0.6301610 196.694 203.664
gene28 52.3966 1.00011 3.21656 0.072916317 0.6301610 205.819 212.790
gene49 87.2911 2.03428 6.52199 0.075619315 0.6301610 200.269 208.269
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: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.2 LTS
Matrix products: default
BLAS: /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so
LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.0
locale:
[1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
[3] LC_TIME=en_GB 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: America/New_York
tzcode source: system (glibc)
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
Di, Y, DW Schafer, JS Cumbie, and JH Chang. 2015. “NBPSeq: Negative Binomial Models for Rna-Sequencing Data.” R Package Version 0.3. 0, URL Http://CRAN. R-Project. Org/Package= NBPSeq.
Love, Michael I, Wolfgang Huber, and Simon Anders. 2014. “Moderated Estimation of Fold Change and Dispersion for Rna-Seq Data with Deseq2.” Genome Biology 15 (12): 550.
Robinson, Mark D, Davis J McCarthy, and Gordon K Smyth. 2010. “EdgeR: A Bioconductor Package for Differential Expression Analysis of Digital Gene Expression Data.” Bioinformatics 26 (1): 139–40.
Wood, Simon, and Maintainer Simon Wood. 2015. “Package ’Mgcv’.” R Package Version 1: 29.
Zhou, Yi-Hui, Kai Xia, and Fred A Wright. 2011. “A Powerful and Flexible Approach to the Analysis of Rna Sequence Count Data.” Bioinformatics 27 (19): 2672–8.