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 321 160 310 14 167 118 5 264 77
gene2 736 2 5 32 134 3 115 1013 116
gene3 1 112 1 27 53 289 26 201 33
gene4 90 48 6 52 53 192 258 356 42
gene5 81 12 17 45 86 3 43 2 98
gene6 242 122 18 38 7 5 149 614 27
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 237 45 127 75 5 2 1 309
gene2 24 1 61 596 166 14 6 2
gene3 2 1 1 6 32 20 1 1
gene4 23 115 140 1 781 41 7 2
gene5 191 61 334 7 18 7 54 1
gene6 3 17 7 114 149 11 9 9
sample18 sample19 sample20
gene1 316 4 271
gene2 1 7 25
gene3 153 177 1067
gene4 47 48 11
gene5 92 1 568
gene6 1 273 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 65.53931 -1.41484246 -0.4286336 0.8210214 0
sample2 67.40928 -0.22382841 1.1054109 0.7785891 2
sample3 72.15319 0.02647715 -0.9673179 -0.5929070 1
sample4 36.94218 -1.15264117 -0.3431567 -0.2390098 2
sample5 75.84780 0.69098219 -1.2337943 0.3538889 2
sample6 70.66884 -1.45000600 -0.1383750 -0.9596376 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 128.5210 1.00005 0.0567017 0.811899471 0.88249942 251.524 258.494
gene2 159.8815 1.00003 14.3845115 0.000148782 0.00743909 217.286 224.256
gene3 102.0467 1.85582 2.5322205 0.284216863 0.47369477 214.139 221.961
gene4 96.3421 1.00006 1.3272192 0.249323484 0.44522051 226.962 233.932
gene5 69.7584 1.00004 2.9155560 0.087737044 0.33533234 219.650 226.620
gene6 92.4897 1.00005 1.7019983 0.192032137 0.40773716 214.266 221.236
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 128.5210 -0.441345 0.458029 -0.963572 0.3352604 0.718358 251.524
gene2 159.8815 0.855387 0.448642 1.906612 0.0565708 0.404077 217.286
gene3 102.0467 -0.490912 0.525848 -0.933564 0.3505287 0.718358 214.139
gene4 96.3421 -0.183150 0.379806 -0.482219 0.6296503 0.852578 226.962
gene5 69.7584 -0.254720 0.439487 -0.579585 0.5621945 0.826757 219.650
gene6 92.4897 0.518804 0.432105 1.200642 0.2298899 0.638583 214.266
BIC
<numeric>
gene1 258.494
gene2 224.256
gene3 221.961
gene4 233.932
gene5 226.620
gene6 221.236
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 128.5210 0.0540399 1.17315 0.0460639 0.9632593 0.963259 251.524
gene2 159.8815 -2.7387434 1.14412 -2.3937473 0.0166772 0.166732 217.286
gene3 102.0467 2.0432624 1.32708 1.5396638 0.1236423 0.364844 214.139
gene4 96.3421 0.4522329 0.97891 0.4619759 0.6440986 0.805123 226.962
gene5 69.7584 -0.4903549 1.12470 -0.4359857 0.6628471 0.808350 219.650
gene6 92.4897 -1.3231486 1.10595 -1.1963939 0.2315429 0.445275 214.266
BIC
<numeric>
gene1 258.494
gene2 224.256
gene3 221.961
gene4 233.932
gene5 226.620
gene6 221.236
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>
gene2 159.8815 1.00003 14.38451 0.000148782 0.00743909 217.286 224.256
gene23 78.3548 1.00004 9.35953 0.002219253 0.05548133 207.935 214.905
gene27 33.7972 1.00006 7.19056 0.007330558 0.10880910 182.010 188.980
gene15 75.3277 1.00004 6.88315 0.008704728 0.10880910 217.003 223.973
gene26 91.5546 1.00004 6.13869 0.013229514 0.11858628 215.383 222.353
gene8 48.7895 1.00006 6.00996 0.014230354 0.11858628 209.388 216.359
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.2.0 RC (2022-04-19 r82224)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 20.04.4 LTS
Matrix products: default
BLAS: /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
LAPACK: /home/biocbuild/bbs-3.15-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] stats4 stats graphics grDevices utils datasets methods
[8] base
other attached packages:
[1] ggplot2_3.3.5 BiocParallel_1.30.0
[3] NBAMSeq_1.12.0 SummarizedExperiment_1.26.0
[5] Biobase_2.56.0 GenomicRanges_1.48.0
[7] GenomeInfoDb_1.32.0 IRanges_2.30.0
[9] S4Vectors_0.34.0 BiocGenerics_0.42.0
[11] MatrixGenerics_1.8.0 matrixStats_0.62.0
loaded via a namespace (and not attached):
[1] httr_1.4.2 sass_0.4.1 bit64_4.0.5
[4] jsonlite_1.8.0 splines_4.2.0 bslib_0.3.1
[7] assertthat_0.2.1 highr_0.9 blob_1.2.3
[10] GenomeInfoDbData_1.2.8 yaml_2.3.5 pillar_1.7.0
[13] RSQLite_2.2.12 lattice_0.20-45 glue_1.6.2
[16] digest_0.6.29 RColorBrewer_1.1-3 XVector_0.36.0
[19] colorspace_2.0-3 htmltools_0.5.2 Matrix_1.4-1
[22] DESeq2_1.36.0 XML_3.99-0.9 pkgconfig_2.0.3
[25] genefilter_1.78.0 zlibbioc_1.42.0 purrr_0.3.4
[28] xtable_1.8-4 scales_1.2.0 tibble_3.1.6
[31] annotate_1.74.0 mgcv_1.8-40 KEGGREST_1.36.0
[34] farver_2.1.0 generics_0.1.2 ellipsis_0.3.2
[37] withr_2.5.0 cachem_1.0.6 cli_3.3.0
[40] survival_3.3-1 magrittr_2.0.3 crayon_1.5.1
[43] memoise_2.0.1 evaluate_0.15 fansi_1.0.3
[46] nlme_3.1-157 tools_4.2.0 lifecycle_1.0.1
[49] stringr_1.4.0 locfit_1.5-9.5 munsell_0.5.0
[52] DelayedArray_0.22.0 AnnotationDbi_1.58.0 Biostrings_2.64.0
[55] compiler_4.2.0 jquerylib_0.1.4 rlang_1.0.2
[58] grid_4.2.0 RCurl_1.98-1.6 labeling_0.4.2
[61] bitops_1.0-7 rmarkdown_2.14 gtable_0.3.0
[64] DBI_1.1.2 R6_2.5.1 knitr_1.38
[67] dplyr_1.0.8 fastmap_1.1.0 bit_4.0.4
[70] utf8_1.2.2 stringi_1.7.6 parallel_4.2.0
[73] Rcpp_1.0.8.3 vctrs_0.4.1 geneplotter_1.74.0
[76] png_0.1-7 tidyselect_1.1.2 xfun_0.30
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.