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 3 53 12 1 56 31 63 63
gene2 11 547 81 57 290 105 1 1 3
gene3 73 4 1 29 1 44 42 1 267
gene4 9 1 22 334 52 487 187 17 214
gene5 5 57 36 2 5 33 17 18 106
gene6 186 33 171 45 34 17 483 9 35
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 194 5 54 136 1 194 131 240
gene2 1 74 80 22 79 36 1 59
gene3 1 193 15 11 2 280 598 183
gene4 1 46 66 138 16 467 3 84
gene5 21 126 56 24 65 1 275 69
gene6 94 4 35 333 79 26 12 7
sample18 sample19 sample20
gene1 13 1 179
gene2 3 153 38
gene3 12 28 1
gene4 34 1 1
gene5 32 15 88
gene6 60 1 2
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.66254 0.01265036 0.9934029 -1.7142257 1
sample2 34.66720 -0.10435286 0.3137609 0.5980733 0
sample3 56.25247 -0.65644423 1.7424033 1.0073511 1
sample4 78.24353 1.43103207 0.8477051 1.3220032 2
sample5 51.59431 -0.64136458 0.8038664 1.1385558 1
sample6 53.37846 0.08941207 0.6622098 0.3194802 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 62.4659 1.00006 0.0493945 0.824246 0.945626 217.952 224.922
gene2 70.1680 1.00017 0.2069966 0.649419 0.942503 209.409 216.380
gene3 86.0899 1.00008 1.5093814 0.219295 0.577093 206.729 213.699
gene4 97.3626 1.00006 1.7793996 0.182247 0.577093 225.370 232.340
gene5 47.8611 1.00011 0.1829719 0.668920 0.942503 211.169 218.139
gene6 63.2627 1.00019 1.6538605 0.198580 0.577093 218.465 225.436
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 62.4659 0.5455755 0.531811 1.025882 3.04947e-01 0.58643716 217.952
gene2 70.1680 -0.6500294 0.506925 -1.282299 1.99738e-01 0.48434398 209.409
gene3 86.0899 2.2221952 0.562946 3.947437 7.89922e-05 0.00394961 206.729
gene4 97.3626 1.4552846 0.566246 2.570059 1.01681e-02 0.08473427 225.370
gene5 47.8611 -0.0738979 0.474452 -0.155754 8.76227e-01 0.94676709 211.169
gene6 63.2627 0.4354026 0.482650 0.902109 3.66999e-01 0.65535566 218.465
BIC
<numeric>
gene1 224.922
gene2 216.380
gene3 213.699
gene4 232.340
gene5 218.139
gene6 225.436
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 62.4659 -0.473506 1.055792 -0.448484 0.6538040 0.797322 217.952
gene2 70.1680 -2.186699 0.994064 -2.199756 0.0278242 0.231868 209.409
gene3 86.0899 2.145727 1.141301 1.880072 0.0600983 0.241070 206.729
gene4 97.3626 1.783394 1.143254 1.559929 0.1187767 0.371177 225.370
gene5 47.8611 1.396437 0.943920 1.479401 0.1390331 0.408921 211.169
gene6 63.2627 -0.868513 0.962588 -0.902268 0.3669145 0.611524 218.465
BIC
<numeric>
gene1 224.922
gene2 216.380
gene3 213.699
gene4 232.340
gene5 218.139
gene6 225.436
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>
gene19 96.5271 1.00021 9.22165 0.00239504 0.119752 229.806 236.776
gene43 50.9227 1.00007 5.16122 0.02310830 0.345641 206.526 213.496
gene37 49.3237 1.00003 4.90794 0.02673575 0.345641 190.119 197.089
gene35 85.1063 1.00009 4.70989 0.03000369 0.345641 217.393 224.364
gene32 66.0225 1.00241 4.32738 0.03777153 0.345641 197.776 204.749
gene22 52.7951 1.00008 4.11411 0.04253480 0.345641 199.708 206.679
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.4.2 (2024-10-31)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.1 LTS
Matrix products: default
BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
locale:
[1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
[3] LC_TIME=en_US.UTF-8 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: Etc/UTC
tzcode source: system (glibc)
attached base packages:
[1] stats4 stats graphics grDevices utils datasets methods
[8] base
other attached packages:
[1] ggplot2_3.5.1 BiocParallel_1.41.0
[3] NBAMSeq_1.23.0 SummarizedExperiment_1.37.0
[5] Biobase_2.67.0 GenomicRanges_1.59.1
[7] GenomeInfoDb_1.43.4 IRanges_2.41.2
[9] S4Vectors_0.45.2 BiocGenerics_0.53.5
[11] generics_0.1.3 MatrixGenerics_1.19.1
[13] matrixStats_1.5.0 rmarkdown_2.29
loaded via a namespace (and not attached):
[1] KEGGREST_1.47.0 gtable_0.3.6 xfun_0.50
[4] bslib_0.8.0 lattice_0.22-6 vctrs_0.6.5
[7] tools_4.4.2 parallel_4.4.2 tibble_3.2.1
[10] AnnotationDbi_1.69.0 RSQLite_2.3.9 blob_1.2.4
[13] pkgconfig_2.0.3 Matrix_1.7-2 lifecycle_1.0.4
[16] GenomeInfoDbData_1.2.13 farver_2.1.2 compiler_4.4.2
[19] Biostrings_2.75.3 munsell_0.5.1 DESeq2_1.47.1
[22] codetools_0.2-20 htmltools_0.5.8.1 sys_3.4.3
[25] buildtools_1.0.0 sass_0.4.9 yaml_2.3.10
[28] pillar_1.10.1 crayon_1.5.3 jquerylib_0.1.4
[31] DelayedArray_0.33.4 cachem_1.1.0 abind_1.4-8
[34] nlme_3.1-167 genefilter_1.89.0 locfit_1.5-9.10
[37] digest_0.6.37 labeling_0.4.3 splines_4.4.2
[40] maketools_1.3.1 fastmap_1.2.0 grid_4.4.2
[43] colorspace_2.1-1 cli_3.6.3 SparseArray_1.7.4
[46] magrittr_2.0.3 S4Arrays_1.7.1 survival_3.8-3
[49] XML_3.99-0.18 withr_3.0.2 scales_1.3.0
[52] UCSC.utils_1.3.1 bit64_4.6.0-1 XVector_0.47.2
[55] httr_1.4.7 bit_4.5.0.1 png_0.1-8
[58] memoise_2.0.1 evaluate_1.0.3 knitr_1.49
[61] mgcv_1.9-1 rlang_1.1.5 Rcpp_1.0.14
[64] xtable_1.8-4 glue_1.8.0 DBI_1.2.3
[67] annotate_1.85.0 jsonlite_1.8.9 R6_2.5.1