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 126 8 189 782 2 9 40 699 1
gene2 2 178 46 11 17 2 347 1 2
gene3 603 6 201 10 9 17 340 410 9
gene4 72 264 34 24 16 2 12 89 3
gene5 335 681 1 5 1 64 309 5 83
gene6 25 11 42 6 25 534 60 705 2
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 1 2 3 28 294 9 1 15
gene2 116 116 6 104 1 14 5 23
gene3 529 70 11 183 42 25 637 34
gene4 37 1 17 47 1002 331 36 1
gene5 1 30 13 389 404 1 20 404
gene6 32 3 1 15 9 12 77 1
sample18 sample19 sample20
gene1 38 130 165
gene2 21 82 4
gene3 15 81 1
gene4 5 3 18
gene5 67 370 1
gene6 32 41 64
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 74.27429 -0.32874314 -2.0933716 -1.2923141 0
sample2 23.68229 -0.52228202 0.2763394 -0.7240890 1
sample3 33.54247 -0.61337931 1.5867593 -0.4173168 2
sample4 38.08514 -0.04616288 0.4751899 -0.9184247 0
sample5 41.91163 -2.09773234 0.0964033 0.1111367 1
sample6 59.61730 0.78790988 -1.3104822 0.8607507 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 100.0667 1.00017 0.7607978 0.3831718 0.903027 220.890 227.860
gene2 55.0365 1.00006 0.0495338 0.8239566 0.936314 205.574 212.544
gene3 125.1557 1.00008 0.3545196 0.5516604 0.936314 246.831 253.801
gene4 94.3227 1.00005 3.2449846 0.0716577 0.534939 208.407 215.377
gene5 155.9123 1.00009 0.0171632 0.8961401 0.953341 240.188 247.158
gene6 51.0358 1.00005 0.0935062 0.7598387 0.936314 205.491 212.462
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 100.0667 0.372395 0.529284 0.703583 0.48169271 0.7769237 220.890
gene2 55.0365 0.771152 0.519383 1.484747 0.13761099 0.5973556 205.574
gene3 125.1557 0.384005 0.496199 0.773893 0.43899404 0.7769237 246.831
gene4 94.3227 0.471852 0.461832 1.021696 0.30692485 0.7769237 208.407
gene5 155.9123 1.879258 0.602169 3.120813 0.00180353 0.0450882 240.188
gene6 51.0358 -0.502264 0.452617 -1.109688 0.26713365 0.7769237 205.491
BIC
<numeric>
gene1 227.860
gene2 212.544
gene3 253.801
gene4 215.377
gene5 247.158
gene6 212.462
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 100.0667 -1.660849 1.115943 -1.488292 0.13667389 0.509385 220.890
gene2 55.0365 1.184576 1.103028 1.073931 0.28285362 0.543949 205.574
gene3 125.1557 -1.231892 1.047009 -1.176583 0.23936207 0.526503 246.831
gene4 94.3227 -2.677507 0.988069 -2.709839 0.00673159 0.112193 208.407
gene5 155.9123 -1.637134 1.262266 -1.296980 0.19463823 0.509385 240.188
gene6 51.0358 -0.853155 0.959268 -0.889381 0.37379833 0.622997 205.491
BIC
<numeric>
gene1 227.860
gene2 212.544
gene3 253.801
gene4 215.377
gene5 247.158
gene6 212.462
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
<numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
gene22 117.2649 1.00004 18.83915 1.44278e-05 0.000721388 232.803
gene44 117.9230 1.00006 12.98816 3.13511e-04 0.007837764 224.461
gene38 76.3417 1.00008 6.94292 8.41811e-03 0.140301774 219.402
gene25 56.4966 1.00008 4.00275 4.54423e-02 0.534939490 204.506
gene31 64.4542 1.00006 3.67541 5.52343e-02 0.534939490 214.789
gene4 94.3227 1.00005 3.24498 7.16577e-02 0.534939490 208.407
BIC
<numeric>
gene22 239.774
gene44 231.431
gene38 226.372
gene25 211.476
gene31 221.760
gene4 215.377
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.39.0
[3] NBAMSeq_1.21.0 SummarizedExperiment_1.35.0
[5] Biobase_2.65.0 GenomicRanges_1.57.0
[7] GenomeInfoDb_1.41.0 IRanges_2.39.0
[9] S4Vectors_0.43.0 BiocGenerics_0.51.0
[11] MatrixGenerics_1.17.0 matrixStats_1.3.0
loaded via a namespace (and not attached):
[1] KEGGREST_1.45.0 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.67.0 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.73.0 munsell_0.5.1 DESeq2_1.45.0
[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.31.0 cachem_1.0.8
[34] abind_1.4-5 nlme_3.1-164 genefilter_1.87.0
[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.5.0 magrittr_2.0.3
[49] S4Arrays_1.5.0 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_1.1.0 bit64_4.0.5 rmarkdown_2.26
[58] XVector_0.45.0 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.83.0 jsonlite_1.8.8
[73] R6_2.5.1 zlibbioc_1.51.0