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 198 38 200 10 2 547 1 203 44
gene2 221 5 35 38 57 11 2 503 4
gene3 27 17 562 48 1 37 374 2 16
gene4 11 5 2 8 2 153 110 1 469
gene5 248 214 25 3 35 21 24 2 73
gene6 382 1 24 47 65 189 5 22 115
sample10 sample11 sample12 sample13 sample14 sample15 sample16 sample17
gene1 87 25 1 2 16 2 1177 305
gene2 12 28 1 1 2 1 4 33
gene3 1 70 15 1 27 66 8 38
gene4 7 5 12 214 1 737 1 262
gene5 157 2 503 5 151 1 5 27
gene6 317 109 1 1 233 48 159 315
sample18 sample19 sample20
gene1 1 129 1
gene2 86 460 5
gene3 1 140 18
gene4 22 35 55
gene5 1 33 7
gene6 323 83 43
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.54951 0.04488663 -0.71294090 0.0640629 1
sample2 49.92429 0.83912488 0.13588385 0.1130017 2
sample3 67.91117 0.96257435 -0.74942510 -0.7990847 1
sample4 62.99399 -0.87713156 -0.92663976 0.4771915 1
sample5 41.74307 0.90722197 -0.01735581 -1.5746186 2
sample6 73.50637 -1.17876044 -1.76796116 -1.1692036 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 124.8252 1.00041 0.409282 0.5224953 0.761648 222.550 229.521
gene2 78.9734 1.00008 0.907578 0.3407855 0.669120 199.872 206.842
gene3 94.3556 1.00006 1.157500 0.2820112 0.669120 212.541 219.512
gene4 79.2734 1.00005 6.030721 0.0140618 0.100441 206.715 213.685
gene5 57.6009 1.00005 4.540161 0.0331124 0.206952 198.279 205.249
gene6 89.5401 1.00015 0.842993 0.3586933 0.669120 241.729 248.700
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 124.8252 1.3165099 0.554546 2.374031 0.0175951 0.0799776 222.550
gene2 78.9734 1.4065711 0.512671 2.743614 0.0060767 0.0412761 199.872
gene3 94.3556 1.0349204 0.534193 1.937353 0.0527022 0.1756741 212.541
gene4 79.2734 -0.7677040 0.490859 -1.564000 0.1178176 0.3100464 206.715
gene5 57.6009 0.4761483 0.425916 1.117940 0.2635927 0.4881346 198.279
gene6 89.5401 -0.0587845 0.473015 -0.124276 0.9010966 0.9586134 241.729
BIC
<numeric>
gene1 229.521
gene2 206.842
gene3 219.512
gene4 213.685
gene5 205.249
gene6 248.700
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 124.8252 -2.023208 1.103220 -1.833911 0.0666672 0.3703733 222.550
gene2 78.9734 1.841593 1.063241 1.732055 0.0832637 0.3784713 199.872
gene3 94.3556 0.965255 1.096824 0.880045 0.3788349 0.7892393 212.541
gene4 79.2734 -3.049029 1.004993 -3.033881 0.0024143 0.0402383 206.715
gene5 57.6009 1.214025 0.877904 1.382867 0.1667055 0.4473678 198.279
gene6 89.5401 0.140884 0.973044 0.144787 0.8848792 0.9351798 241.729
BIC
<numeric>
gene1 229.521
gene2 206.842
gene3 219.512
gene4 213.685
gene5 205.249
gene6 248.700
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 81.7933 1.00012 12.04693 0.000519048 0.0163565 216.732 223.703
gene46 75.9742 1.00008 11.22271 0.000808667 0.0163565 217.922 224.893
gene50 41.4218 1.00005 10.79412 0.001018701 0.0163565 188.058 195.028
gene7 154.7134 1.00011 10.33237 0.001308517 0.0163565 237.272 244.243
gene31 61.4991 1.00009 6.79512 0.009145027 0.0914503 206.649 213.619
gene13 51.1736 1.00005 6.25324 0.012400312 0.1004411 201.091 208.061
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 Under development (unstable) (2024-10-21 r87258)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.1 LTS
Matrix products: default
BLAS: /home/biocbuild/bbs-3.21-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.41.0
[3] NBAMSeq_1.23.0 SummarizedExperiment_1.37.0
[5] Biobase_2.67.0 GenomicRanges_1.59.0
[7] GenomeInfoDb_1.43.0 IRanges_2.41.0
[9] S4Vectors_0.45.0 BiocGenerics_0.53.0
[11] MatrixGenerics_1.19.0 matrixStats_1.4.1
loaded via a namespace (and not attached):
[1] KEGGREST_1.47.0 gtable_0.3.6 xfun_0.48
[4] bslib_0.8.0 lattice_0.22-6 vctrs_0.6.5
[7] tools_4.5.0 generics_0.1.3 parallel_4.5.0
[10] RSQLite_2.3.7 tibble_3.2.1 fansi_1.0.6
[13] AnnotationDbi_1.69.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.5.0
[22] Biostrings_2.75.0 munsell_0.5.1 DESeq2_1.47.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.33.0 cachem_1.1.0
[34] abind_1.4-8 nlme_3.1-166 genefilter_1.89.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.5.0
[43] fastmap_1.2.0 grid_4.5.0 colorspace_2.1-1
[46] cli_3.6.3 SparseArray_1.7.0 magrittr_2.0.3
[49] S4Arrays_1.7.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.3.0 bit64_4.5.2 rmarkdown_2.28
[58] XVector_0.47.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.85.0 jsonlite_1.8.9
[73] R6_2.5.1 zlibbioc_1.53.0