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 2 1 220 7 2 680 31 102 3
gene2 27 33 163 38 2 127 38 17 11
gene3 239 2 114 169 58 1 191 16 104
gene4 111 12 14 11 1 428 140 61 2
gene5 2 17 46 2 10 62 696 12 7
gene6 33 5 602 22 378 1 4 55 9
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 163 1 39 107 16 236 2 33
gene2 44 447 44 1 3 1 2 377
gene3 1 1 113 133 942 50 1 206
gene4 47 63 1 109 9 132 16 197
gene5 3 1 9 1 6 33 190 141
gene6 160 54 252 635 4 34 7 39
sample18 sample19 sample20
gene1 1 6 8
gene2 58 283 173
gene3 14 323 14
gene4 26 6 144
gene5 53 24 9
gene6 18 3 240
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 55.98011 1.63751205 -1.0997371 -1.4039810 1
sample2 34.13251 -0.30089995 -1.0545540 0.8582900 2
sample3 35.45191 1.06459769 -0.8992885 1.5800864 1
sample4 23.35149 0.36402180 -0.2432098 0.8107233 0
sample5 63.07076 1.18663274 0.1002464 -1.4819245 1
sample6 30.20543 0.08313703 -1.1940850 1.0330631 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 59.7293 1.00013 0.877717 0.3488883 0.625390 203.922 210.892
gene2 89.4233 1.00007 0.152211 0.6965176 0.735514 221.021 227.991
gene3 130.2202 1.00005 0.187691 0.6648821 0.722698 234.136 241.106
gene4 66.8300 1.00009 0.371514 0.5422801 0.646065 221.517 228.488
gene5 93.6086 1.00012 3.020610 0.0822465 0.316333 201.164 208.134
gene6 96.5862 1.00005 0.795308 0.3725290 0.625390 226.890 233.860
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 59.7293 0.393588 0.561928 0.700425 0.48366208 0.7205455 203.922
gene2 89.4233 0.119313 0.512014 0.233027 0.81574030 0.9083350 221.021
gene3 130.2202 -0.284248 0.559786 -0.507781 0.61160723 0.8047464 234.136
gene4 66.8300 -0.194596 0.520177 -0.374096 0.70833308 0.8854163 221.517
gene5 93.6086 -1.505587 0.555575 -2.709962 0.00672908 0.0560757 201.164
gene6 96.5862 0.500417 0.512661 0.976118 0.32900615 0.6580123 226.890
BIC
<numeric>
gene1 210.892
gene2 227.991
gene3 241.106
gene4 228.488
gene5 208.134
gene6 233.860
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
<numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
gene1 59.7293 1.02057436 1.45589 0.70099645 0.48330523 0.8417310
gene2 89.4233 -1.99492602 1.32965 -1.50033795 0.13352688 0.3709080
gene3 130.2202 -4.86829784 1.55995 -3.12080897 0.00180355 0.0300592
gene4 66.8300 1.16090818 1.34767 0.86141666 0.38900860 0.7463774
gene5 93.6086 0.00418844 1.42983 0.00292932 0.99766274 0.9976627
gene6 96.5862 -2.83852425 1.38238 -2.05336560 0.04003713 0.1662133
AIC BIC
<numeric> <numeric>
gene1 203.922 210.892
gene2 221.021 227.991
gene3 234.136 241.106
gene4 221.517 228.488
gene5 201.164 208.134
gene6 226.890 233.860
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>
gene15 93.1031 1.00004 17.30777 3.27278e-05 0.00163639 204.871 211.841
gene45 119.5809 1.00010 10.01116 1.55687e-03 0.02693888 214.326 221.297
gene13 63.7973 1.00004 9.94162 1.61633e-03 0.02693888 209.394 216.364
gene31 94.9188 1.00006 6.08822 1.36145e-02 0.15559727 226.186 233.156
gene46 41.8465 1.00004 5.85245 1.55597e-02 0.15559727 184.418 191.388
gene19 140.9082 1.00007 5.26594 2.17545e-02 0.16285409 236.081 243.052
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.1 (2024-06-14)
Platform: aarch64-apple-darwin20
Running under: macOS Ventura 13.6.7
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.40.0
[3] NBAMSeq_1.22.0 SummarizedExperiment_1.36.0
[5] Biobase_2.66.0 GenomicRanges_1.58.0
[7] GenomeInfoDb_1.42.0 IRanges_2.40.0
[9] S4Vectors_0.44.0 BiocGenerics_0.52.0
[11] MatrixGenerics_1.18.0 matrixStats_1.4.1
loaded via a namespace (and not attached):
[1] KEGGREST_1.46.0 gtable_0.3.6 xfun_0.49
[4] bslib_0.8.0 lattice_0.22-6 vctrs_0.6.5
[7] tools_4.4.1 generics_0.1.3 parallel_4.4.1
[10] RSQLite_2.3.7 tibble_3.2.1 fansi_1.0.6
[13] AnnotationDbi_1.68.0 blob_1.2.4 pkgconfig_2.0.3
[16] Matrix_1.7-1 lifecycle_1.0.4 GenomeInfoDbData_1.2.13
[19] farver_2.1.2 compiler_4.4.1 Biostrings_2.74.0
[22] munsell_0.5.1 DESeq2_1.46.0 codetools_0.2-20
[25] htmltools_0.5.8.1 sass_0.4.9 yaml_2.3.10
[28] pillar_1.9.0 crayon_1.5.3 jquerylib_0.1.4
[31] DelayedArray_0.32.0 cachem_1.1.0 abind_1.4-8
[34] nlme_3.1-166 genefilter_1.88.0 tidyselect_1.2.1
[37] locfit_1.5-9.10 digest_0.6.37 dplyr_1.1.4
[40] labeling_0.4.3 splines_4.4.1 fastmap_1.2.0
[43] grid_4.4.1 colorspace_2.1-1 cli_3.6.3
[46] SparseArray_1.6.0 magrittr_2.0.3 S4Arrays_1.6.0
[49] survival_3.7-0 XML_3.99-0.17 utf8_1.2.4
[52] withr_3.0.2 scales_1.3.0 UCSC.utils_1.2.0
[55] bit64_4.5.2 rmarkdown_2.29 XVector_0.46.0
[58] httr_1.4.7 bit_4.5.0 png_0.1-8
[61] memoise_2.0.1 evaluate_1.0.1 knitr_1.49
[64] mgcv_1.9-1 rlang_1.1.4 Rcpp_1.0.13-1
[67] DBI_1.2.3 xtable_1.8-4 glue_1.8.0
[70] annotate_1.84.0 jsonlite_1.8.9 R6_2.5.1
[73] zlibbioc_1.52.0