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 29 1040 3 351 1 1 84 171 34
gene2 106 472 390 281 81 996 10 121 147
gene3 413 161 6 113 27 59 202 357 2
gene4 5 26 53 19 564 23 1 155 1
gene5 10 39 1 1 1 48 23 612 13
gene6 186 1 1 3 182 1 141 94 101
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 35 274 88 163 269 1 19 5
gene2 143 123 1 2 6 2 409 61
gene3 50 15 1 1 1 340 5 61
gene4 1 37 20 193 57 1 493 61
gene5 15 2 2 2 115 11 7 379
gene6 2 87 1 236 452 3 1 261
sample18 sample19 sample20
gene1 116 288 11
gene2 325 171 5
gene3 53 8 51
gene4 268 16 7
gene5 239 65 1
gene6 5 7 3
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 57.01996 -0.13897305 0.9820239 -2.3117869 0
sample2 22.25832 -0.04656733 0.3061022 0.7687949 2
sample3 44.80903 -0.31575079 -0.9391511 -0.6626072 0
sample4 69.07252 2.07598523 2.3136288 0.3837749 2
sample5 62.86006 -0.45802181 -1.3645738 0.8396898 0
sample6 51.06175 -0.47314006 0.6905879 0.7114551 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 85.8243 1.00011 7.070031 0.00784358 0.1092815 229.843 236.813
gene2 171.3144 1.00003 0.958728 0.32754460 0.6157582 260.104 267.075
gene3 72.8575 1.00006 0.101358 0.75025056 0.8931554 222.219 229.189
gene4 87.8976 1.00006 0.381442 0.53684841 0.7456228 216.888 223.858
gene5 51.9033 1.00005 8.572105 0.00341532 0.0853829 191.192 198.163
gene6 70.5429 1.00010 2.591803 0.10744707 0.3837395 201.743 208.714
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 85.8243 0.2453994 0.385249 0.6369894 0.524132 0.714657 229.843
gene2 171.3144 0.6536084 0.457489 1.4286859 0.153095 0.505451 260.104
gene3 72.8575 0.4763072 0.421482 1.1300770 0.258444 0.605409 222.219
gene4 87.8976 -0.0179741 0.417014 -0.0431019 0.965620 0.966765 216.888
gene5 51.9033 -0.3237725 0.402786 -0.8038326 0.421494 0.679828 191.192
gene6 70.5429 0.5035128 0.464615 1.0837215 0.278488 0.605409 201.743
BIC
<numeric>
gene1 236.813
gene2 267.075
gene3 229.189
gene4 223.858
gene5 198.163
gene6 208.714
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 85.8243 0.920687 0.847657 1.086156 2.77410e-01 0.556556980 229.843
gene2 171.3144 -0.500449 1.000086 -0.500406 6.16789e-01 0.790755768 260.104
gene3 72.8575 0.571042 0.925849 0.616777 5.37382e-01 0.746363947 222.219
gene4 87.8976 -2.590706 0.910343 -2.845855 4.42923e-03 0.036910273 216.888
gene5 51.9033 -0.301397 0.896451 -0.336211 7.36712e-01 0.833400527 191.192
gene6 70.5429 -4.241619 0.985429 -4.304339 1.67485e-05 0.000837426 201.743
BIC
<numeric>
gene1 236.813
gene2 267.075
gene3 229.189
gene4 223.858
gene5 198.163
gene6 208.714
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>
gene30 75.9464 1.00015 11.40022 0.000735435 0.0367718 211.263 218.233
gene5 51.9033 1.00005 8.57211 0.003415317 0.0853829 191.192 198.163
gene1 85.8243 1.00011 7.07003 0.007843580 0.1092815 229.843 236.813
gene21 190.4171 1.00008 6.77913 0.009225911 0.1092815 230.879 237.849
gene7 80.6558 1.00011 6.47799 0.010928151 0.1092815 207.146 214.116
gene45 59.6353 1.00010 5.60821 0.017886074 0.1490506 200.062 207.032
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: x86_64-apple-darwin20
Running under: macOS Monterey 12.7.4
Matrix products: default
BLAS: /Library/Frameworks/R.framework/Versions/4.4-x86_64/Resources/lib/libRblas.0.dylib
LAPACK: /Library/Frameworks/R.framework/Versions/4.4-x86_64/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
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.