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 3 5 2 2 61 13 212 6 107
gene2 1 10 106 1 22 1 335 2 36
gene3 7 98 18 1 30 64 41 35 17
gene4 41 6 1 60 16 143 1 11 4
gene5 10 3 319 1 1 312 43 3 22
gene6 1 81 449 4 362 4 19 74 532
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 273 2 71 30 43 14 1 1
gene2 5 1 120 6 16 2 7 805
gene3 55 120 6 59 1 11 14 1
gene4 572 31 1 344 58 12 20 330
gene5 26 109 5 8 26 76 200 166
gene6 267 35 69 59 243 188 21 632
sample18 sample19 sample20
gene1 155 289 799
gene2 64 95 1
gene3 308 313 10
gene4 27 7 4
gene5 139 340 12
gene6 11 95 4
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 36.45066 -0.02413811 1.2940361 0.57054527 0
sample2 58.61363 0.05496400 0.8194974 0.81123820 1
sample3 63.32145 -0.00449474 0.3444390 0.69527733 1
sample4 60.62626 0.26733796 0.1804601 0.04029040 2
sample5 44.44742 0.77958638 0.5657089 0.09589169 1
sample6 75.39852 0.54212324 0.8143508 0.97449953 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 90.0788 1.00143 0.000371608 0.9952494 0.995249 217.150 224.121
gene2 80.7395 1.00025 1.114631449 0.2910806 0.748239 192.597 199.568
gene3 42.4583 1.00008 0.016321577 0.8985639 0.916902 206.056 213.026
gene4 81.8683 1.00031 0.201885112 0.6537538 0.806844 216.570 223.540
gene5 80.2695 1.00006 5.572949511 0.0182452 0.210294 218.095 225.066
gene6 145.3781 1.00005 0.154988158 0.6938859 0.806844 239.810 246.781
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 90.0788 -0.124537 0.399229 -0.311943 0.75508363 0.8850563 217.150
gene2 80.7395 -1.069347 0.388345 -2.753602 0.00589435 0.0982391 192.597
gene3 42.4583 -0.772035 0.340774 -2.265536 0.02347981 0.1230501 206.056
gene4 81.8683 -0.440174 0.413360 -1.064868 0.28693568 0.5313624 216.570
gene5 80.2695 -0.705700 0.351417 -2.008156 0.04462673 0.1716413 218.095
gene6 145.3781 -0.740655 0.334364 -2.215113 0.02675232 0.1230501 239.810
BIC
<numeric>
gene1 224.121
gene2 199.568
gene3 213.026
gene4 223.540
gene5 225.066
gene6 246.781
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 90.0788 1.4089402 1.058595 1.3309531 0.1832044 0.436201 217.150
gene2 80.7395 0.3354875 1.035418 0.3240117 0.7459292 0.925890 192.597
gene3 42.4583 -0.2760074 0.900537 -0.3064920 0.7592301 0.925890 206.056
gene4 81.8683 0.0691645 1.093739 0.0632368 0.9495779 0.961397 216.570
gene5 80.2695 -2.0069057 0.932441 -2.1523139 0.0313726 0.235583 218.095
gene6 145.3781 0.8913475 0.887066 1.0048265 0.3149804 0.627228 239.810
BIC
<numeric>
gene1 224.121
gene2 199.568
gene3 213.026
gene4 223.540
gene5 225.066
gene6 246.781
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>
gene26 110.2313 1.00010 13.22154 0.00027713 0.0138565 205.076 212.046
gene34 55.0388 1.00023 9.13066 0.00252013 0.0630032 205.670 212.640
gene40 106.3222 1.00015 6.35095 0.01174566 0.1957609 228.787 235.757
gene5 80.2695 1.00006 5.57295 0.01824521 0.2102943 218.095 225.066
gene33 118.1699 1.00015 5.17247 0.02297736 0.2102943 240.353 247.323
gene35 54.7496 1.00019 5.00944 0.02523531 0.2102943 207.447 214.418
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.1 (2022-06-23)
Platform: x86_64-pc-linux-gnu (64-bit)
Running under: Ubuntu 20.04.5 LTS
Matrix products: default
BLAS: /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
LAPACK: /home/biocbuild/bbs-3.16-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.6 BiocParallel_1.32.0
[3] NBAMSeq_1.14.0 SummarizedExperiment_1.28.0
[5] Biobase_2.58.0 GenomicRanges_1.50.0
[7] GenomeInfoDb_1.34.0 IRanges_2.32.0
[9] S4Vectors_0.36.0 BiocGenerics_0.44.0
[11] MatrixGenerics_1.10.0 matrixStats_0.62.0
loaded via a namespace (and not attached):
[1] httr_1.4.4 sass_0.4.2 bit64_4.0.5
[4] jsonlite_1.8.3 splines_4.2.1 bslib_0.4.0
[7] assertthat_0.2.1 highr_0.9 blob_1.2.3
[10] GenomeInfoDbData_1.2.9 yaml_2.3.6 pillar_1.8.1
[13] RSQLite_2.2.18 lattice_0.20-45 glue_1.6.2
[16] digest_0.6.30 RColorBrewer_1.1-3 XVector_0.38.0
[19] colorspace_2.0-3 htmltools_0.5.3 Matrix_1.5-1
[22] DESeq2_1.38.0 XML_3.99-0.12 pkgconfig_2.0.3
[25] genefilter_1.80.0 zlibbioc_1.44.0 xtable_1.8-4
[28] scales_1.2.1 tibble_3.1.8 annotate_1.76.0
[31] mgcv_1.8-41 KEGGREST_1.38.0 farver_2.1.1
[34] generics_0.1.3 withr_2.5.0 cachem_1.0.6
[37] cli_3.4.1 survival_3.4-0 magrittr_2.0.3
[40] crayon_1.5.2 memoise_2.0.1 evaluate_0.17
[43] fansi_1.0.3 nlme_3.1-160 tools_4.2.1
[46] lifecycle_1.0.3 stringr_1.4.1 locfit_1.5-9.6
[49] munsell_0.5.0 DelayedArray_0.24.0 AnnotationDbi_1.60.0
[52] Biostrings_2.66.0 compiler_4.2.1 jquerylib_0.1.4
[55] rlang_1.0.6 grid_4.2.1 RCurl_1.98-1.9
[58] labeling_0.4.2 bitops_1.0-7 rmarkdown_2.17
[61] gtable_0.3.1 codetools_0.2-18 DBI_1.1.3
[64] R6_2.5.1 knitr_1.40 dplyr_1.0.10
[67] fastmap_1.1.0 bit_4.0.4 utf8_1.2.2
[70] stringi_1.7.8 parallel_4.2.1 Rcpp_1.0.9
[73] vctrs_0.5.0 geneplotter_1.76.0 png_0.1-7
[76] tidyselect_1.2.0 xfun_0.34
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.