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 66 1 3 142 1 1 3 41 7
gene2 3 162 44 95 5 2 116 18 2
gene3 48 11 6 15 2 2 278 1 27
gene4 163 48 178 354 144 1172 93 9 20
gene5 8 8 1 3 257 11 407 8 331
gene6 56 105 48 12 75 2 38 88 5
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 20 43 68 465 216 27 97 440
gene2 164 1 116 2 26 22 161 38
gene3 500 4 75 260 137 38 205 6
gene4 36 572 37 38 6 3 78 31
gene5 1 15 270 8 103 22 295 225
gene6 6 18 234 22 54 1 1 2
sample18 sample19 sample20
gene1 11 710 317
gene2 1 1 96
gene3 1 1 9
gene4 2 137 169
gene5 155 13 160
gene6 7 6 22
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 39.32583 -1.466817457 0.253565457 -0.4836906 0
sample2 64.96879 0.009883497 -0.005491285 0.3062376 2
sample3 48.43803 -1.379046346 0.245591511 -1.1989460 0
sample4 72.23607 0.091576885 -2.368038614 1.5361869 0
sample5 58.18745 -0.581581836 -0.183409295 1.2285452 2
sample6 68.01419 -0.710004188 0.553171134 -0.7510875 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 83.4950 1.00004 2.85330461 0.09119168 0.3039723 230.078 237.048
gene2 40.7426 1.00013 0.00143641 0.97128357 0.9896601 202.512 209.482
gene3 50.4769 1.00014 1.01512994 0.31365235 0.7025948 204.952 211.922
gene4 244.9670 1.00011 7.85478951 0.00507215 0.0634018 250.354 257.324
gene5 78.4993 1.00017 0.30093433 0.58329219 0.7495266 229.711 236.681
gene6 30.8898 1.00008 3.29121405 0.06966618 0.2631689 192.840 199.810
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 83.4950 -0.5531634 0.359970 -1.536692 0.12436883 0.518203 230.078
gene2 40.7426 0.1824976 0.346935 0.526028 0.59886894 0.893616 202.512
gene3 50.4769 0.0480085 0.351486 0.136587 0.89135712 0.938854 204.952
gene4 244.9670 -0.9428239 0.343721 -2.742996 0.00608815 0.101469 250.354
gene5 78.4993 -0.1652298 0.350587 -0.471295 0.63743019 0.893616 229.711
gene6 30.8898 0.1267646 0.310087 0.408803 0.68268420 0.893616 192.840
BIC
<numeric>
gene1 237.048
gene2 209.482
gene3 211.922
gene4 257.324
gene5 236.681
gene6 199.810
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 83.4950 -0.281249 0.902434 -0.311656 0.7553023 0.939681 230.078
gene2 40.7426 1.359620 0.868752 1.565026 0.1175769 0.576695 202.512
gene3 50.4769 1.165995 0.878958 1.326564 0.1846528 0.576695 204.952
gene4 244.9670 -1.211238 0.860779 -1.407141 0.1593857 0.576695 250.354
gene5 78.4993 -0.189416 0.881924 -0.214776 0.8299417 0.953234 229.711
gene6 30.8898 1.670768 0.779445 2.143536 0.0320701 0.446951 192.840
BIC
<numeric>
gene1 237.048
gene2 209.482
gene3 211.922
gene4 257.324
gene5 236.681
gene6 199.810
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>
gene16 70.9235 1.00010 14.17082 0.000166878 0.00494364 203.322 210.292
gene7 60.5208 1.00007 13.85170 0.000197746 0.00494364 206.712 213.682
gene32 115.5726 1.00007 9.83554 0.001712764 0.02854607 240.093 247.063
gene4 244.9670 1.00011 7.85479 0.005072146 0.06340182 250.354 257.324
gene10 110.1209 1.00006 6.93856 0.008438374 0.08438374 217.894 224.864
gene50 80.7912 1.00005 6.53371 0.010587286 0.08754528 216.427 223.397
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.1 (2024-06-14)
Platform: x86_64-apple-darwin20
Running under: macOS Monterey 12.7.6
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.40.0
[3] NBAMSeq_1.22.0 SummarizedExperiment_1.36.0
[5] Biobase_2.66.0 GenomicRanges_1.58.0
[7] GenomeInfoDb_1.42.0 IRanges_2.40.0
[9] S4Vectors_0.44.0 BiocGenerics_0.52.0
[11] MatrixGenerics_1.18.0 matrixStats_1.4.1
loaded via a namespace (and not attached):
[1] KEGGREST_1.46.0 gtable_0.3.6 xfun_0.48
[4] bslib_0.8.0 lattice_0.22-6 vctrs_0.6.5
[7] tools_4.4.1 generics_0.1.3 parallel_4.4.1
[10] RSQLite_2.3.7 tibble_3.2.1 fansi_1.0.6
[13] AnnotationDbi_1.68.0 highr_0.11 blob_1.2.4
[16] pkgconfig_2.0.3 Matrix_1.7-1 lifecycle_1.0.4
[19] GenomeInfoDbData_1.2.13 farver_2.1.2 compiler_4.4.1
[22] Biostrings_2.74.0 munsell_0.5.1 DESeq2_1.46.0
[25] codetools_0.2-20 htmltools_0.5.8.1 sass_0.4.9
[28] yaml_2.3.10 pillar_1.9.0 crayon_1.5.3
[31] jquerylib_0.1.4 DelayedArray_0.32.0 cachem_1.1.0
[34] abind_1.4-8 nlme_3.1-166 genefilter_1.88.0
[37] tidyselect_1.2.1 locfit_1.5-9.10 digest_0.6.37
[40] dplyr_1.1.4 labeling_0.4.3 splines_4.4.1
[43] fastmap_1.2.0 grid_4.4.1 colorspace_2.1-1
[46] cli_3.6.3 SparseArray_1.6.0 magrittr_2.0.3
[49] S4Arrays_1.6.0 survival_3.7-0 XML_3.99-0.17
[52] utf8_1.2.4 withr_3.0.2 scales_1.3.0
[55] UCSC.utils_1.2.0 bit64_4.5.2 rmarkdown_2.28
[58] XVector_0.46.0 httr_1.4.7 bit_4.5.0
[61] png_0.1-8 memoise_2.0.1 evaluate_1.0.1
[64] knitr_1.48 mgcv_1.9-1 rlang_1.1.4
[67] Rcpp_1.0.13 DBI_1.2.3 xtable_1.8-4
[70] glue_1.8.0 annotate_1.84.0 jsonlite_1.8.9
[73] R6_2.5.1 zlibbioc_1.52.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.