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 25 229 1 48 8 615 188 98 167
gene2 98 9 3 5 29 67 123 10 1
gene3 103 199 146 10 9 9 113 1 45
gene4 164 1 30 163 192 421 82 2 5
gene5 799 1 95 50 1 26 50 39 86
gene6 785 285 1 50 76 22 76 16 455
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 459 85 1 1 129 1 81 69
gene2 55 18 1 200 7 150 8 48
gene3 14 48 1 38 139 231 71 174
gene4 7 73 38 279 3 19 5 43
gene5 6 20 86 1 325 7 3 16
gene6 148 165 20 111 8 1 3 210
sample18 sample19 sample20
gene1 8 177 26
gene2 50 1183 148
gene3 60 50 51
gene4 28 35 31
gene5 429 19 70
gene6 4 12 124
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 38.51204 -0.6678354 1.6709069 -0.3413777 1
sample2 71.21410 -0.4988147 -0.9298130 -1.0556452 2
sample3 65.22731 0.3528193 -0.3233422 -0.5270079 0
sample4 45.75614 -0.5392666 -1.0856483 0.6210935 0
sample5 75.73885 -1.1519942 0.3344813 1.6515987 1
sample6 76.96546 0.1413459 -0.8102947 1.5693547 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 116.0184 1.00009 2.7391315 0.0979262 0.325777 231.878 238.848
gene2 93.1707 1.00011 0.3637875 0.5464682 0.758984 224.546 231.517
gene3 76.3340 1.00006 4.0484705 0.0442157 0.280692 224.961 231.931
gene4 82.6128 1.00011 0.1763165 0.6747481 0.843435 210.056 217.027
gene5 78.6076 1.00011 0.0177276 0.8943970 0.956187 221.948 228.919
gene6 81.5067 1.00009 2.6275680 0.1050478 0.325777 229.439 236.409
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 116.0184 1.15576440 0.583077 1.98218106 0.0474590 0.379112 231.878
gene2 93.1707 0.00482234 0.605040 0.00797029 0.9936407 0.993641 224.546
gene3 76.3340 0.16161598 0.506391 0.31915239 0.7496110 0.903428 224.961
gene4 82.6128 -0.89696252 0.478130 -1.87598044 0.0606580 0.379112 210.056
gene5 78.6076 0.99099287 0.593078 1.67093165 0.0947352 0.394730 221.948
gene6 81.5067 -0.93341441 0.539728 -1.72941660 0.0837346 0.380612 229.439
BIC
<numeric>
gene1 238.848
gene2 231.517
gene3 231.931
gene4 217.027
gene5 228.919
gene6 236.409
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 116.0184 -1.705392 0.951935 -1.791500 0.0732132 0.291052 231.878
gene2 93.1707 -0.874424 0.990694 -0.882638 0.3774316 0.767800 224.546
gene3 76.3340 -0.648339 0.828196 -0.782833 0.4337250 0.767800 224.961
gene4 82.6128 -1.755842 0.783702 -2.240445 0.0250620 0.208850 210.056
gene5 78.6076 0.447886 0.970178 0.461653 0.6443299 0.894903 221.948
gene6 81.5067 1.052208 0.881423 1.193760 0.2325717 0.646033 229.439
BIC
<numeric>
gene1 238.848
gene2 231.517
gene3 231.931
gene4 217.027
gene5 228.919
gene6 236.409
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>
gene45 55.9752 1.00005 11.20844 0.000814655 0.0407327 187.873 194.843
gene25 133.1563 1.00009 7.96924 0.004759913 0.1031898 240.473 247.443
gene24 91.1278 1.00018 7.49489 0.006191387 0.1031898 208.166 215.136
gene26 52.5734 1.00004 6.51574 0.010694199 0.1336775 196.061 203.031
gene50 102.7870 1.00004 5.66189 0.017342079 0.1734208 220.883 227.853
gene31 53.3687 1.00004 4.32287 0.037610153 0.2806924 190.491 197.461
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-pc-linux-gnu
Running under: Ubuntu 24.04.1 LTS
Matrix products: default
BLAS: /home/biocbuild/bbs-3.20-bioc/R/lib/libRblas.so
LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0
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
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.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.