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 5 73 211 1 1 21 4 29 2
gene2 94 2 406 321 1 1 927 188 56
gene3 122 19 4 4 13 92 147 9 3
gene4 30 133 1 40 4 1 252 1 109
gene5 60 1 5 179 1 58 14 14 94
gene6 83 1 490 3 37 69 225 2 11
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 8 397 130 8 114 2 1 8
gene2 147 1 66 32 2 380 1 281
gene3 112 330 3 6 14 4 2 13
gene4 1 1 2 324 504 58 307 8
gene5 29 3 7 18 39 70 156 7
gene6 3 1 51 57 49 2 36 9
sample18 sample19 sample20
gene1 3 245 76
gene2 1 62 120
gene3 14 1 81
gene4 112 112 4
gene5 1 265 38
gene6 12 372 3
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 26.61320 -0.87991740 1.30022139 0.7591436 0
sample2 50.46816 -0.58531533 0.38352668 -1.4390917 2
sample3 27.22425 -0.34183604 1.05833006 0.2887099 0
sample4 23.86635 -0.78617637 -2.74945847 0.7555255 2
sample5 22.60160 0.02654447 0.04689304 1.2074754 1
sample6 34.10561 0.18027857 -0.97000841 -0.5046642 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 60.2456 1.00008 0.0338331 0.854154 0.978178 205.681 212.651
gene2 117.2252 1.00076 0.6142221 0.434159 0.657817 232.413 239.384
gene3 41.7175 1.00010 0.0199017 0.888032 0.978178 195.048 202.019
gene4 82.9041 1.00006 2.0529531 0.151933 0.544674 210.312 217.282
gene5 38.0813 1.00007 1.2487526 0.263817 0.544674 201.331 208.301
gene6 60.6378 1.00010 1.9720928 0.160274 0.544674 202.123 209.094
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 60.2456 -0.244169 0.462634 -0.527780 0.597652067 0.7470651 205.681
gene2 117.2252 0.613154 0.458428 1.337513 0.181055403 0.3848923 232.413
gene3 41.7175 0.839487 0.392062 2.141211 0.032257029 0.2016064 195.048
gene4 82.9041 -1.520444 0.444982 -3.416863 0.000633472 0.0223685 210.312
gene5 38.0813 -0.308562 0.377336 -0.817738 0.413506860 0.5907241 201.331
gene6 60.6378 -0.539059 0.413594 -1.303354 0.192454024 0.3848923 202.123
BIC
<numeric>
gene1 212.651
gene2 239.384
gene3 202.019
gene4 217.282
gene5 208.301
gene6 209.094
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 60.2456 -0.8003810 1.010147 -0.7923408 0.4281620 0.898202 205.681
gene2 117.2252 -1.4262200 0.996559 -1.4311443 0.1523889 0.634954 232.413
gene3 41.7175 -0.0106729 0.858454 -0.0124327 0.9900804 0.993753 195.048
gene4 82.9041 -0.4731953 0.926044 -0.5109859 0.6093610 0.952126 210.312
gene5 38.0813 -0.7132697 0.818456 -0.8714820 0.3834911 0.898202 201.331
gene6 60.6378 -2.1074616 0.878809 -2.3980877 0.0164809 0.206011 202.123
BIC
<numeric>
gene1 212.651
gene2 239.384
gene3 202.019
gene4 217.282
gene5 208.301
gene6 209.094
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>
gene11 118.1638 1.00010 6.73250 0.00947325 0.209200 223.214 230.184
gene12 49.9469 1.00010 6.40016 0.01142031 0.209200 204.402 211.373
gene38 34.1588 1.00005 6.23157 0.01255200 0.209200 178.839 185.809
gene27 53.0171 1.00007 5.13647 0.02343362 0.244152 201.422 208.392
gene42 104.1530 1.00004 5.06539 0.02441516 0.244152 229.687 236.657
gene26 73.6163 1.00003 4.51399 0.03362624 0.280219 206.450 213.420
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.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.4 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
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.