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 86 2 15 2 2 253 203 34 1
gene2 58 31 11 25 88 149 2 41 20
gene3 3 1 110 19 60 39 51 38 551
gene4 8 3 587 242 2 161 206 21 191
gene5 163 39 91 25 5 1 26 92 2
gene6 6 3 33 312 47 71 17 213 1
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 96 3 233 1 328 1 391 413
gene2 1 338 2 2 19 65 135 91
gene3 242 1004 25 8 293 11 1 43
gene4 1 1 20 7 42 185 35 4
gene5 1 275 20 13 654 10 1 21
gene6 37 20 1 2 60 265 1 88
sample18 sample19 sample20
gene1 463 309 53
gene2 88 112 576
gene3 2 50 1
gene4 360 69 2
gene5 457 77 206
gene6 2 17 117
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 72.28458 0.6519901 -0.7124010 -1.27666549 2
sample2 38.18440 -1.5560521 -0.7883353 -0.50988764 0
sample3 78.03110 -0.2719248 -0.3702485 1.82829721 0
sample4 45.16483 0.3847899 1.2381355 0.24274274 1
sample5 54.31841 -2.7341289 -0.9211927 0.10598287 0
sample6 39.70366 -1.0945870 -0.1296539 0.01002405 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 120.2233 1.00008 2.404772 0.120988 0.465340 229.438 236.408
gene2 66.2243 1.00006 0.028360 0.866434 0.904858 222.305 229.275
gene3 94.8374 1.00004 0.484840 0.486265 0.813723 222.607 229.578
gene4 96.2173 1.00041 0.244189 0.621266 0.813723 226.021 232.991
gene5 75.7973 1.00007 1.980188 0.159384 0.550438 219.694 226.664
gene6 63.8863 1.00003 1.058605 0.303553 0.661357 204.916 211.886
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 120.2233 1.0436069 0.344035 3.033434 0.00241788 0.0604469 229.438
gene2 66.2243 0.2012330 0.297098 0.677329 0.49819710 0.8303285 222.305
gene3 94.8374 -0.0450286 0.342453 -0.131488 0.89538903 0.9705086 222.607
gene4 96.2173 0.1801169 0.362934 0.496279 0.61969723 0.8852818 226.021
gene5 75.7973 0.0491319 0.321220 0.152954 0.87843433 0.9705086 219.694
gene6 63.8863 -0.2879550 0.318746 -0.903398 0.36631452 0.7894912 204.916
BIC
<numeric>
gene1 236.408
gene2 229.275
gene3 229.578
gene4 232.991
gene5 226.664
gene6 211.886
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 120.2233 -0.247235 0.959850 -0.257577 0.79673348 0.885259 229.438
gene2 66.2243 0.102425 0.843959 0.121363 0.90340388 0.921841 222.305
gene3 94.8374 -2.841032 0.977485 -2.906471 0.00365531 0.158402 222.607
gene4 96.2173 -0.139244 1.031161 -0.135036 0.89258316 0.921841 226.021
gene5 75.7973 2.267965 0.909787 2.492852 0.01267216 0.158402 219.694
gene6 63.8863 0.489355 0.901678 0.542716 0.58732549 0.793683 204.916
BIC
<numeric>
gene1 236.408
gene2 229.275
gene3 229.578
gene4 232.991
gene5 226.664
gene6 211.886
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>
gene8 35.9545 1.00006 6.31924 0.0119458 0.366456 181.774 188.744
gene24 104.3405 1.00011 5.60849 0.0178887 0.366456 238.253 245.223
gene42 133.0468 1.00016 4.96582 0.0258792 0.366456 227.839 234.809
gene15 86.0546 1.00006 4.74903 0.0293164 0.366456 218.454 225.424
gene14 49.5655 1.00008 4.10970 0.0426561 0.406607 210.345 217.316
gene36 25.8818 1.00005 3.68180 0.0550215 0.406607 180.332 187.302
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.1.0 beta (2021-05-03 r80259)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 20.04.2 LTS
Matrix products: default
BLAS: /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
LAPACK: /home/biocbuild/bbs-3.14-bioc/R/lib/libRlapack.so
locale:
[1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
[3] LC_TIME=en_GB 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
attached base packages:
[1] parallel stats4 stats graphics grDevices utils datasets
[8] methods base
other attached packages:
[1] ggplot2_3.3.3 BiocParallel_1.27.0
[3] NBAMSeq_1.9.0 SummarizedExperiment_1.23.0
[5] Biobase_2.53.0 GenomicRanges_1.45.0
[7] GenomeInfoDb_1.29.0 IRanges_2.27.0
[9] S4Vectors_0.31.0 BiocGenerics_0.39.0
[11] MatrixGenerics_1.5.0 matrixStats_0.58.0
loaded via a namespace (and not attached):
[1] httr_1.4.2 sass_0.4.0 bit64_4.0.5
[4] jsonlite_1.7.2 splines_4.1.0 bslib_0.2.5.1
[7] assertthat_0.2.1 highr_0.9 blob_1.2.1
[10] GenomeInfoDbData_1.2.6 yaml_2.2.1 pillar_1.6.1
[13] RSQLite_2.2.7 lattice_0.20-44 glue_1.4.2
[16] digest_0.6.27 RColorBrewer_1.1-2 XVector_0.33.0
[19] colorspace_2.0-1 htmltools_0.5.1.1 Matrix_1.3-3
[22] DESeq2_1.33.0 XML_3.99-0.6 pkgconfig_2.0.3
[25] genefilter_1.75.0 zlibbioc_1.39.0 purrr_0.3.4
[28] xtable_1.8-4 scales_1.1.1 tibble_3.1.2
[31] annotate_1.71.0 mgcv_1.8-35 KEGGREST_1.33.0
[34] farver_2.1.0 generics_0.1.0 ellipsis_0.3.2
[37] withr_2.4.2 cachem_1.0.5 survival_3.2-11
[40] magrittr_2.0.1 crayon_1.4.1 memoise_2.0.0
[43] evaluate_0.14 fansi_0.4.2 nlme_3.1-152
[46] tools_4.1.0 lifecycle_1.0.0 stringr_1.4.0
[49] locfit_1.5-9.4 munsell_0.5.0 DelayedArray_0.19.0
[52] AnnotationDbi_1.55.0 Biostrings_2.61.0 compiler_4.1.0
[55] jquerylib_0.1.4 rlang_0.4.11 grid_4.1.0
[58] RCurl_1.98-1.3 labeling_0.4.2 bitops_1.0-7
[61] rmarkdown_2.8 gtable_0.3.0 DBI_1.1.1
[64] R6_2.5.0 knitr_1.33 dplyr_1.0.6
[67] fastmap_1.1.0 bit_4.0.4 utf8_1.2.1
[70] stringi_1.6.2 Rcpp_1.0.6 vctrs_0.3.8
[73] geneplotter_1.71.0 png_0.1-7 tidyselect_1.1.1
[76] xfun_0.23
Di, Y, DW Schafer, JS Cumbie, and JH Chang. 2015. “NBPSeq: Negative Binomial Models for Rna-Sequencing Data.” R Package Version 0.3. 0, URL Http://CRAN. R-Project. Org/Package= NBPSeq.
Love, Michael I, Wolfgang Huber, and Simon Anders. 2014. “Moderated Estimation of Fold Change and Dispersion for Rna-Seq Data with Deseq2.” Genome Biology 15 (12): 550.
Robinson, Mark D, Davis J McCarthy, and Gordon K Smyth. 2010. “EdgeR: A Bioconductor Package for Differential Expression Analysis of Digital Gene Expression Data.” Bioinformatics 26 (1): 139–40.
Wood, Simon, and Maintainer Simon Wood. 2015. “Package ’Mgcv’.” R Package Version 1: 29.
Zhou, Yi-Hui, Kai Xia, and Fred A Wright. 2011. “A Powerful and Flexible Approach to the Analysis of Rna Sequence Count Data.” Bioinformatics 27 (19): 2672–8.