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 257 310 517 19 1 1 1 84 92
gene2 9 1 19 1 368 2 56 1 482
gene3 142 77 93 200 74 39 22 99 38
gene4 19 4 4 234 1 20 42 98 136
gene5 65 2 165 126 381 9 2 5 58
gene6 61 347 10 5 1 125 28 24 156
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 3 171 33 6 2 6 4 2
gene2 15 1 133 8 258 9 634 73
gene3 434 1 4 43 149 2 78 29
gene4 3 23 223 30 49 240 43 97
gene5 1 1 20 194 2 8 2 203
gene6 5 52 6 17 118 642 12 2
sample18 sample19 sample20
gene1 4 166 1
gene2 160 37 41
gene3 30 3 16
gene4 3 2 4
gene5 17 33 72
gene6 15 1 6
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.39437 0.3720971 0.03496827 1.0638518 1
sample2 69.34883 0.3363310 0.13003263 1.3202320 1
sample3 57.40232 -1.7681632 -1.15926239 0.4378728 1
sample4 46.36819 -0.4022975 -0.53442071 1.9860639 1
sample5 35.46529 -0.6202689 -0.22168022 0.2726405 2
sample6 75.79146 1.4774305 -0.59126508 -0.3102660 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 95.7599 1.00047 0.7874128 0.37501290 0.6122610 207.638 214.609
gene2 101.7836 1.00020 8.3111210 0.00394832 0.0394832 218.130 225.100
gene3 76.1389 1.00006 0.0613751 0.80445718 0.9420134 228.761 235.732
gene4 59.9637 1.43022 1.3346453 0.36153995 0.6122610 215.617 223.015
gene5 75.7577 1.00013 9.0843970 0.00258105 0.0322631 207.086 214.056
gene6 106.8134 1.00016 0.2284826 0.63285565 0.8327048 211.926 218.896
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 95.7599 -0.485851 0.497646 -0.976299 0.3289162 0.7827395 207.638
gene2 101.7836 -1.073688 0.439253 -2.444351 0.0145113 0.0761843 218.130
gene3 76.1389 -0.244427 0.402559 -0.607184 0.5437291 0.8458008 228.761
gene4 59.9637 -0.348149 0.466233 -0.746727 0.4552284 0.8211566 215.617
gene5 75.7577 -0.401881 0.433682 -0.926673 0.3540965 0.7827395 207.086
gene6 106.8134 1.178036 0.438972 2.683628 0.0072828 0.0438384 211.926
BIC
<numeric>
gene1 214.609
gene2 225.100
gene3 235.732
gene4 223.015
gene5 214.056
gene6 218.896
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 95.7599 -1.805477 1.69771 -1.063475 0.2875664 0.586177 207.638
gene2 101.7836 -3.463739 1.49139 -2.322491 0.0202065 0.112258 218.130
gene3 76.1389 0.716611 1.37287 0.521982 0.6016828 0.810264 228.761
gene4 59.9637 -1.589462 1.46270 -1.086662 0.2771861 0.586177 215.617
gene5 75.7577 -3.481129 1.47764 -2.355864 0.0184797 0.112258 207.086
gene6 106.8134 -1.757074 1.49100 -1.178455 0.2386153 0.586177 211.926
BIC
<numeric>
gene1 214.609
gene2 225.100
gene3 235.732
gene4 223.015
gene5 214.056
gene6 218.896
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>
gene24 68.0425 1.00005 15.40773 8.73068e-05 0.00436534 209.640 216.611
gene16 188.2261 1.00018 12.69317 3.67656e-04 0.00688247 225.245 232.215
gene34 108.4120 1.00013 12.47469 4.12948e-04 0.00688247 235.345 242.315
gene5 75.7577 1.00013 9.08440 2.58105e-03 0.03226308 207.086 214.056
gene2 101.7836 1.00020 8.31112 3.94832e-03 0.03948318 218.130 225.100
gene14 89.8530 1.00011 7.91500 4.90673e-03 0.04088944 206.014 212.984
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 beta (2024-04-15 r86425)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 22.04.4 LTS
Matrix products: default
BLAS: /home/biocbuild/bbs-3.19-bioc/R/lib/libRblas.so
LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
locale:
[1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
[3] LC_TIME=en_US.UTF-8 LC_COLLATE=en_US.UTF-8
[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
time zone: America/New_York
tzcode source: system (glibc)
attached base packages:
[1] stats4 stats graphics grDevices utils datasets methods
[8] base
other attached packages:
[1] ggplot2_3.5.1 BiocParallel_1.38.0
[3] NBAMSeq_1.20.0 SummarizedExperiment_1.34.0
[5] Biobase_2.64.0 GenomicRanges_1.56.0
[7] GenomeInfoDb_1.40.0 IRanges_2.38.0
[9] S4Vectors_0.42.0 BiocGenerics_0.50.0
[11] MatrixGenerics_1.16.0 matrixStats_1.3.0
loaded via a namespace (and not attached):
[1] KEGGREST_1.44.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.66.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.72.0 munsell_0.5.1 DESeq2_1.44.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.30.0 cachem_1.0.8
[34] abind_1.4-5 nlme_3.1-164 genefilter_1.86.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.4.0 magrittr_2.0.3
[49] S4Arrays_1.4.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.0.0 bit64_4.0.5 rmarkdown_2.26
[58] XVector_0.44.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.82.0 jsonlite_1.8.8
[73] R6_2.5.1 zlibbioc_1.50.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.