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 301 76 7 1 523 22 268 114 59
gene2 130 262 23 13 97 325 6 68 204
gene3 454 168 1 22 2 148 13 352 3
gene4 34 234 51 2 208 54 8 193 28
gene5 11 468 27 1 1 75 15 1 899
gene6 4 56 156 186 24 1 119 93 1
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 16 149 158 115 19 5 9 28
gene2 1 1 5 7 55 7 222 68
gene3 7 13 64 11 1 2 106 45
gene4 206 19 46 1 178 319 229 775
gene5 1 37 14 46 35 54 255 26
gene6 121 66 90 272 3 13 128 6
sample18 sample19 sample20
gene1 53 29 1
gene2 1 10 1
gene3 84 50 14
gene4 202 273 9
gene5 152 8 269
gene6 1 165 1
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.67279 1.33351808 -1.815513148 0.26370699 0
sample2 35.36804 -0.83583927 -0.790266252 -2.43361280 1
sample3 37.25504 -1.08571507 0.007593744 0.05564406 1
sample4 64.83742 -0.17374276 -0.240979010 0.55838657 2
sample5 36.95051 -0.00803284 -1.752551858 -3.50006262 2
sample6 53.36819 -1.97145100 1.339464343 0.48793539 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 65.1262 1.00012 0.19779512 0.6565476 0.820908 226.604 233.574
gene2 64.7553 1.00005 0.42405558 0.5149499 0.820908 211.279 218.250
gene3 54.9511 1.00018 0.00381677 0.9518291 0.954606 205.547 212.517
gene4 134.6692 1.00004 0.27619095 0.5992381 0.820908 252.458 259.428
gene5 121.4099 1.00017 2.87008131 0.0902861 0.410392 226.122 233.092
gene6 65.7235 1.00009 1.10861883 0.2924203 0.595862 217.974 224.945
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 65.1262 0.2338915 0.350751 0.666830 0.50488066 0.8948568 226.604
gene2 64.7553 -0.7554963 0.391295 -1.930757 0.05351306 0.3869572 211.279
gene3 54.9511 -1.1164631 0.357613 -3.121984 0.00179637 0.0875513 205.547
gene4 134.6692 -0.0459519 0.373952 -0.122882 0.90220079 0.9206131 252.458
gene5 121.4099 0.2506278 0.413575 0.606003 0.54451302 0.8948568 226.122
gene6 65.7235 -0.4201158 0.388017 -1.082727 0.27892973 0.8203816 217.974
BIC
<numeric>
gene1 233.574
gene2 218.250
gene3 212.517
gene4 259.428
gene5 233.092
gene6 224.945
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 65.1262 0.663021 0.824177 0.804465 0.4211287 0.638651 226.604
gene2 64.7553 -1.648574 0.915091 -1.801541 0.0716177 0.275453 211.279
gene3 54.9511 -1.334259 0.814570 -1.637991 0.1014235 0.338078 205.547
gene4 134.6692 -1.622770 0.880143 -1.843757 0.0652185 0.271744 252.458
gene5 121.4099 -1.261539 0.976117 -1.292406 0.1962166 0.545046 226.122
gene6 65.7235 0.954959 0.910726 1.048569 0.2943766 0.574153 217.974
BIC
<numeric>
gene1 233.574
gene2 218.250
gene3 212.517
gene4 259.428
gene5 233.092
gene6 224.945
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>
gene39 156.1090 1.00011 15.78134 7.07333e-05 0.00353667 207.235 214.205
gene12 105.8409 1.00005 12.34760 4.41854e-04 0.01104636 213.983 220.953
gene38 37.8851 1.00005 10.05135 1.52300e-03 0.02179965 189.188 196.158
gene44 141.3677 1.00018 9.80395 1.74397e-03 0.02179965 238.575 245.545
gene22 106.2786 1.00011 8.91235 2.83456e-03 0.02834562 226.302 233.272
gene23 22.6642 1.00004 7.53431 6.05443e-03 0.05045359 168.987 175.957
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.3.1 (2023-06-16)
Platform: aarch64-apple-darwin20 (64-bit)
Running under: macOS Ventura 13.6.1
Matrix products: default
BLAS: /Library/Frameworks/R.framework/Versions/4.3-arm64/Resources/lib/libRblas.0.dylib
LAPACK: /Library/Frameworks/R.framework/Versions/4.3-arm64/Resources/lib/libRlapack.dylib; LAPACK version 3.11.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.4.2 BiocParallel_1.36.0
[3] NBAMSeq_1.18.0 SummarizedExperiment_1.32.0
[5] Biobase_2.62.0 GenomicRanges_1.54.1
[7] GenomeInfoDb_1.38.0 IRanges_2.36.0
[9] S4Vectors_0.40.1 BiocGenerics_0.48.1
[11] MatrixGenerics_1.14.0 matrixStats_1.0.0
loaded via a namespace (and not attached):
[1] KEGGREST_1.42.0 gtable_0.3.3 xfun_0.39
[4] bslib_0.5.0 lattice_0.21-8 vctrs_0.6.3
[7] tools_4.3.1 bitops_1.0-7 generics_0.1.3
[10] parallel_4.3.1 RSQLite_2.3.1 AnnotationDbi_1.64.1
[13] tibble_3.2.1 fansi_1.0.4 highr_0.10
[16] blob_1.2.4 pkgconfig_2.0.3 Matrix_1.6-0
[19] lifecycle_1.0.3 GenomeInfoDbData_1.2.10 farver_2.1.1
[22] compiler_4.3.1 Biostrings_2.70.1 munsell_0.5.0
[25] DESeq2_1.42.0 codetools_0.2-19 htmltools_0.5.5
[28] sass_0.4.6 RCurl_1.98-1.12 yaml_2.3.7
[31] pillar_1.9.0 crayon_1.5.2 jquerylib_0.1.4
[34] DelayedArray_0.28.0 cachem_1.0.8 abind_1.4-5
[37] nlme_3.1-162 genefilter_1.84.0 tidyselect_1.2.0
[40] locfit_1.5-9.8 digest_0.6.33 dplyr_1.1.2
[43] labeling_0.4.2 splines_4.3.1 fastmap_1.1.1
[46] grid_4.3.1 colorspace_2.1-0 cli_3.6.1
[49] SparseArray_1.2.0 magrittr_2.0.3 S4Arrays_1.2.0
[52] survival_3.5-5 XML_3.99-0.14 utf8_1.2.3
[55] withr_2.5.0 scales_1.2.1 bit64_4.0.5
[58] rmarkdown_2.23 XVector_0.42.0 httr_1.4.6
[61] bit_4.0.5 png_0.1-8 memoise_2.0.1
[64] evaluate_0.21 knitr_1.43 mgcv_1.9-0
[67] rlang_1.1.1 Rcpp_1.0.11 DBI_1.1.3
[70] xtable_1.8-4 glue_1.6.2 annotate_1.80.0
[73] jsonlite_1.8.7 R6_2.5.1 zlibbioc_1.48.0