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 4 6 155 26 1 2 1 1 1
gene2 22 1 48 48 1 108 19 2 31
gene3 835 43 28 14 3 4 308 23 58
gene4 11 205 7 100 3 1 1 423 142
gene5 39 1 541 2 169 188 262 46 1
gene6 4 26 612 1 7 1 298 5 17
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 1 1 34 17 93 17 210 1
gene2 496 46 18 4 7 356 25 13
gene3 102 389 34 3 77 1 3 93
gene4 62 47 8 20 70 15 118 4
gene5 29 345 1 628 53 25 114 78
gene6 56 23 10 1 69 1 105 162
sample18 sample19 sample20
gene1 47 191 5
gene2 15 76 6
gene3 1 139 10
gene4 11 511 2
gene5 1 16 30
gene6 26 10 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 77.16695 -2.9931128 -0.62400730 -1.0175335 2
sample2 33.19389 -0.5428027 -0.37128103 -0.9525296 2
sample3 22.91615 -0.5784937 0.44443170 -2.5807626 0
sample4 34.69313 -0.2349788 2.54787543 -0.4535905 0
sample5 53.12109 -1.3808090 -0.39424108 -0.6025536 0
sample6 48.00638 -0.3130268 -0.04284673 1.2385099 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 29.5124 1.00006 0.7447988 0.388189 0.808728 178.065 185.035
gene2 41.3324 1.00005 0.0337002 0.854472 0.968919 203.669 210.639
gene3 87.0233 1.00011 0.0278515 0.867759 0.968919 225.710 232.680
gene4 69.9674 1.00055 2.0484862 0.152685 0.477140 214.939 221.910
gene5 93.2368 1.00007 0.0515029 0.820580 0.968919 229.673 236.643
gene6 48.6789 1.00009 2.1155724 0.145816 0.477140 196.942 203.913
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 29.5124 0.6853972 0.364229 1.8817751 0.0598666 0.374166 178.065
gene2 41.3324 0.0232996 0.322011 0.0723567 0.9423181 0.969312 203.669
gene3 87.0233 -0.0149645 0.388979 -0.0384713 0.9693119 0.969312 225.710
gene4 69.9674 0.3540518 0.355460 0.9960371 0.3192321 0.577778 214.939
gene5 93.2368 -0.3081137 0.364491 -0.8453271 0.3979283 0.641820 229.673
gene6 48.6789 0.2165665 0.353156 0.6132325 0.5397226 0.769872 196.942
BIC
<numeric>
gene1 185.035
gene2 210.639
gene3 232.680
gene4 221.910
gene5 236.643
gene6 203.913
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 29.5124 -0.627151 1.08868 -0.576065 0.56457112 0.9220420 178.065
gene2 41.3324 0.446652 0.97463 0.458278 0.64675277 0.9421056 203.669
gene3 87.0233 1.751850 1.17499 1.490954 0.13597368 0.3999226 225.710
gene4 69.9674 1.730631 1.07096 1.615968 0.10610117 0.3789328 214.939
gene5 93.2368 -3.295256 1.11540 -2.954320 0.00313358 0.0313358 229.673
gene6 48.6789 -1.581761 1.05992 -1.492344 0.13560904 0.3999226 196.942
BIC
<numeric>
gene1 185.035
gene2 210.639
gene3 232.680
gene4 221.910
gene5 236.643
gene6 203.913
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 79.3169 1.00002 13.48661 0.000240205 0.0120102 205.216 212.186
gene33 53.3346 1.00022 4.78996 0.028669882 0.4710721 194.618 201.589
gene36 118.2441 1.00007 3.56999 0.058851441 0.4710721 220.179 227.150
gene35 103.4337 1.00017 3.52582 0.060482599 0.4710721 222.763 229.733
gene8 71.5521 1.00037 3.39813 0.065366652 0.4710721 223.740 230.711
gene17 67.4474 1.00005 3.03718 0.081387526 0.4710721 200.079 207.049
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.0 Patched (2024-04-24 r86482)
Platform: aarch64-apple-darwin20
Running under: macOS Ventura 13.6.6
Matrix products: default
BLAS: /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRblas.0.dylib
LAPACK: /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRlapack.dylib; LAPACK version 3.12.0
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.1 BiocParallel_1.37.1
[3] NBAMSeq_1.19.0 SummarizedExperiment_1.33.3
[5] Biobase_2.63.1 GenomicRanges_1.55.4
[7] GenomeInfoDb_1.39.14 IRanges_2.37.1
[9] S4Vectors_0.41.7 BiocGenerics_0.49.1
[11] MatrixGenerics_1.15.1 matrixStats_1.3.0
loaded via a namespace (and not attached):
[1] KEGGREST_1.43.1 gtable_0.3.5 xfun_0.43
[4] bslib_0.7.0 lattice_0.22-6 vctrs_0.6.5
[7] tools_4.4.0 generics_0.1.3 parallel_4.4.0
[10] RSQLite_2.3.6 tibble_3.2.1 fansi_1.0.6
[13] AnnotationDbi_1.65.2 highr_0.10 blob_1.2.4
[16] pkgconfig_2.0.3 Matrix_1.7-0 lifecycle_1.0.4
[19] GenomeInfoDbData_1.2.12 farver_2.1.1 compiler_4.4.0
[22] Biostrings_2.71.6 munsell_0.5.1 DESeq2_1.43.5
[25] codetools_0.2-20 htmltools_0.5.8.1 sass_0.4.9
[28] yaml_2.3.8 pillar_1.9.0 crayon_1.5.2
[31] jquerylib_0.1.4 DelayedArray_0.29.9 cachem_1.0.8
[34] abind_1.4-5 nlme_3.1-164 genefilter_1.85.1
[37] tidyselect_1.2.1 locfit_1.5-9.9 digest_0.6.35
[40] dplyr_1.1.4 labeling_0.4.3 splines_4.4.0
[43] fastmap_1.1.1 grid_4.4.0 colorspace_2.1-0
[46] cli_3.6.2 SparseArray_1.3.7 magrittr_2.0.3
[49] S4Arrays_1.3.7 survival_3.6-4 XML_3.99-0.16.1
[52] utf8_1.2.4 withr_3.0.0 scales_1.3.0
[55] UCSC.utils_0.99.7 bit64_4.0.5 rmarkdown_2.26
[58] XVector_0.43.1 httr_1.4.7 bit_4.0.5
[61] png_0.1-8 memoise_2.0.1 evaluate_0.23
[64] knitr_1.46 mgcv_1.9-1 rlang_1.1.3
[67] Rcpp_1.0.12 DBI_1.2.2 xtable_1.8-4
[70] glue_1.7.0 annotate_1.81.2 jsonlite_1.8.8
[73] R6_2.5.1 zlibbioc_1.49.3