PhenoPath models gene expression expression \(y\) in terms of a latent pathway score (pseudotime) \(z\). Uniquely, the evolution of genes along the trajectory isn’t common to each gene but can be perturbed by an additional sample-specific covariate \(\beta\). For example, this could be the mutational status of each sample or a drug that each sample was exposed to.
This software infers both the latent pathway scores \(z_n\) and interaction coefficients \(\beta_{pg}\) for samples \(n = 1, \ldots, N\), covariates \(p = 1, \ldots, P\) and genes \(G = 1, \ldots, G\).
Inference is performed using co-ordinate ascent mean field variational inference. This attempts to find a set of approximating distributions \(q(\mathbf{\theta}) = \prod_i q_i(\theta_i)\) for each variable \(i\) by minimising the KL divergence between the KL divergence between the variational distribution and the true posterior. For a good overview of variational inference see Blei, Kucukelbir, and McAuliffe (2016).
We can simulate data according to the PhenoPath model via a call to simulate_phenopath()
:
set.seed(123L)
sim <- simulate_phenopath()
This returns a list with four entries:
print(str(sim))
## List of 4
## $ parameters:Classes 'tbl_df', 'tbl' and 'data.frame': 40 obs. of 4 variables:
## ..$ alpha : num [1:40] -1 1 1 1 -1 1 -1 -1 -1 -1 ...
## ..$ lambda: num [1:40] 0 0 0 0 0 0 0 0 0 0 ...
## ..$ beta : num [1:40] 0 0 0 0 0 0 0 0 0 0 ...
## ..$ regime: chr [1:40] "de" "de" "de" "de" ...
## $ y : num [1:100, 1:40] -0.999 -1.001 -0.999 -0.998 -1.001 ...
## $ x : num [1:100] 1 1 1 1 1 1 1 1 1 1 ...
## $ z : num [1:100] -0.713 -0.3512 1.6085 -0.0218 0.0426 ...
## NULL
parameters
is a data frame with the simulated parameters, with a column for each of the parameters \(\alpha\), \(\beta\) and \(\lambda\), and a row for each gene. There is an additional column specifying from which regime that gene has been simulated (see ?simulate_phenopath
for details).y
is the \(N \times G\) matrix of gene expressionx
is the \(N\)-length vector of covariatesz
is the true latent pseudotimesBy default this simulates the model for \(N= 100\) cells and \(G=40\) genes.
For 8 representative genes we can visualise what the expression looks like over pseudotime:
df_gex <- tbl_df(sim$y[,c(1,3,11,13,21,23,31,33)]) %>%
mutate(x = factor(sim[['x']]), z = sim[['z']]) %>%
gather(gene, expression, -x, -z)
ggplot(df_gex, aes(x = z, y = expression, color = x)) +
geom_point() +
facet_wrap(~ gene, nrow = 2) +
scale_color_brewer(palette = "Set1")
We see for the first two genes there is differential expression only, genes 3 and 4 show a pseudotime dependence, genes 5 and 6 show pseudotime-covariate interactions (but marginally no differential expression), while genes 7 and 8 show differential expression, pseudotime dependence and pseudotime-covariate interactions.
We can further plot this in PCA space, coloured by both covariate and pseudotime to get an idea of how these look in the reduced dimension:
pca_df <- tbl_df(prcomp(sim$y)$x[,1:2]) %>%
mutate(x = factor(sim[['x']]), z = sim[['z']])
ggplot(pca_df, aes(x = PC1, y = PC2, color = x)) +
geom_point() + scale_colour_brewer(palette = "Set1")
ggplot(pca_df, aes(x = PC1, y = PC2, color = z)) +
geom_point()
PhenoPath fits models with a call to the phenopath
function, which requires at least an expression matrix y
and a covariate vector x
. The expression data should represent something comparable to logged counts, such as \(log_2(\text{TPM}+1)\). Note that buy default PhenoPath centre-scales the input expression.
For use with ExpressionSets
see the section on input formats. For this example we choose an ELBO tolerance of 1e-6
and ELBO calculations thinned by 40
. A discussion of how to control variational inference can be found below.
fit <- phenopath(sim$y, sim$x, elbo_tol = 1e-6, thin = 40)
## Iteration ELBO Change (%)
## [ 40 ] -371.296314122487 Inf
## [ 80 ] -370.787133339349 0.00343310714797981
## [ 120 ] -370.468673891154 0.00214903088060264
## [ 160 ] -370.26850032826 0.00135154329031753
## [ 200 ] -370.142208271187 0.000852996863437216
## [ 240 ] -370.062304913902 0.000539796651971746
## [ 280 ] -370.011642753417 0.00034230112401092
## [ 320 ] -369.979467866318 0.000217409950371418
## [ 360 ] -369.959007967484 0.000138257877181887
## [ 400 ] -369.945984685818 8.80079944417858e-05
## [ 440 ] -369.937688586131 5.60641693425983e-05
## [ 480 ] -369.932400584315 3.57362710521528e-05
## [ 520 ] -369.929028355195 2.27897033087284e-05
## [ 560 ] -369.926877027796 1.45388692502608e-05
## [ 600 ] -369.925504170121 9.27793339182394e-06
## [ 640 ] -369.924627881691 5.92207414394565e-06
## [ 680 ] -369.924068446346 3.78074443590521e-06
## [ 720 ] -369.923711241635 2.41404308536804e-06
## [ 760 ] -369.923483136216 1.54157163292751e-06
## [ 800 ] -369.92333745799 9.84516327119203e-07
print(fit)
## PhenoPath fit with 100 cells and 40 genes
The phenopath
function will print progress on iterations, ELBO, and % change in ELBO that may be turned off by setting verbose = FALSE
in the call.
Once the model has been fit it is important to check convergence with a call to plot_elbo(fit)
to ensure the ELBO is approximately flat:
plot_elbo(fit)
The posterior mean estimates of the pseudotimes \(z\) sit in fit$m_z
that can be extracted using the trajectory
function. We can visualise these compared to both the true pseudotimes and the first principal component of the data:
qplot(sim$z, trajectory(fit)) +
xlab("True z") + ylab("Phenopath z")
qplot(sim$z, pca_df$PC1) +
xlab("True z") + ylab("PC1")
We see that this has high correlation with the true pseudotimes:
cor(sim$z, trajectory(fit))
## [1] 0.9975191
Next, we’re interested in the interactions between the latent space and the covariates. There are three functions to help us here:
interaction_effects
retrieves the posterior interaction effect sizesinteraction_sds
retrieves the posterior interaction standard deviationssignificant_interactions
applies a Bayesian significant test to the interaction parametersNote that if \(P=1\) (ie there is only one covariate) each of these will return a vector, while if \(P>1\) then a matrix is returned.
Alternatively, one can call the interactions
function that returns a data frame with the following entries:
feature
The feature (usually gene)covariate
The covariate, specified from the order originally supplied to the call to phenopath
interaction_effect_size
The effect size of the interaction (\(\beta\) from the statistical model)significant
Boolean for whether the interaction effect is significantly different from 0chi
The precision of the ARD prior on \(\beta\)pathway_loading
The pathway loading \(\lambda\), showing the overall effect for each gene marginalised over the covariateIn our simulated data above, the first 20 genes were simulated with no interaction effect while the latter 20 were simulated with interaction effects. We can pull out the interaction effect sizes, standard deviations, and significance test results into a data frame and plot:
gene_names <- paste0("gene", seq_len(ncol(fit$m_beta)))
df_beta <- data_frame(beta = interaction_effects(fit),
beta_sd = interaction_sds(fit),
is_sig = significant_interactions(fit),
gene = gene_names)
df_beta$gene <- fct_relevel(df_beta$gene, gene_names)
ggplot(df_beta, aes(x = gene, y = beta, color = is_sig)) +
geom_point() +
geom_errorbar(aes(ymin = beta - 2 * beta_sd, ymax = beta + 2 * beta_sd)) +
theme(axis.text.x = element_text(angle = 90, hjust = 1),
axis.title.x = element_blank()) +
ylab(expression(beta)) +
scale_color_brewer(palette = "Set2", name = "Significant")
A typical analysis might follow by graphing the largest effect size; we can do this as follows:
which_largest <- which.max(df_beta$beta)
df_large <- data_frame(
y = sim[['y']][, which_largest],
x = factor(sim[['x']]),
z = sim[['z']]
)
ggplot(df_large, aes(x = z, y = y, color = x)) +
geom_point() +
scale_color_brewer(palette = "Set1") +
stat_smooth()
## `geom_smooth()` using method = 'loess'
ExpressionSet
as inputAlternatively you might have expression values in an ExpressionSet
. For single-cell data it is highly recommended to use the scater package (McCarthy et al. (2017)) in which case data is stored in a class derived from ExpressionSets
called SCESet
s.
We’ll first construct an example SCESet
using our previous simulation data:
suppressPackageStartupMessages(library(scater))
exprs_mat <- t(sim$y)
pdata <- Biobase::AnnotatedDataFrame(data.frame(x = sim$x))
sce <- newSCESet(exprs_mat, phenoData = pdata)
sce
## SCESet (storageMode: lockedEnvironment)
## assayData: 40 features, 100 samples
## element names: exprs
## protocolData: none
## phenoData
## sampleNames: 1 2 ... 100 (100 total)
## varLabels: x
## varMetadata: labelDescription
## featureData: none
## experimentData: use 'experimentData(object)'
## Annotation:
Note that PhenoPath will use whatever is in the exprs
slot of an ExpressionSet
as gene expression values.
We can then pass the \(x\) covariates to PhenoPath in three ways:
pData(sce)
to usepData(sce)
For our example, these three look like
fit <- phenopath(sce, sim$x) # 1
fit <- phenopath(sce, "x") # 2
fit <- phenopath(sce, ~ x) # 3
Note that if the column of the SCESet is a factor it is coerced into a one-hot encoding. The intercept term is then removed as this is taken care of by the \(\lambda\) coefficient automatically, and the columns centre-scaled.
The posterior distribution is naturally multi-modal and the use of variational inference means we essentially dive straight into the nearest local maximum. Therefore, correct initialisation of the latent space is important. This is controlled through the z_init
argument.
PhenoPath allows for three different initialisations:
For our example these three look like
fit <- phenopath(sim$y, sim$x, z_init = 1) # 1, initialise to first principal component
fit <- phenopath(sim$y, sim$x, z_init = sim$z) # 2, initialise to true values
fit <- phenopath(sim$y, sim$x, z_init = "random") # 3, random initialisation
There are several parameters that tune the coordinate ascent variational inference (CAVI):
maxiter
maximum number of iterations to run CAVI forelbo_tol
the percentage change in the ELBO below which the model is considered convergedthin
Computing the ELBO is expensive, so only compute the ELBO every thin
iterationsFor example:
fit <- phenopath(sim$y, sim$x,
maxiter = 1000, # 1000 iterations max
elbo_tol = 1e-2, # consider model converged when change in ELBO < 0.02%
thin = 20 # calculate ELBO every 20 iterations
)
sessionInfo()
## R version 3.4.0 (2017-04-21)
## Platform: x86_64-apple-darwin15.6.0 (64-bit)
## Running under: OS X El Capitan 10.11.6
##
## Matrix products: default
## BLAS: /Library/Frameworks/R.framework/Versions/3.4/Resources/lib/libRblas.0.dylib
## LAPACK: /Library/Frameworks/R.framework/Versions/3.4/Resources/lib/libRlapack.dylib
##
## locale:
## [1] en_GB.UTF-8/en_GB.UTF-8/en_GB.UTF-8/C/en_GB.UTF-8/en_GB.UTF-8
##
## attached base packages:
## [1] parallel stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] scater_1.5.0 Biobase_2.37.2 BiocGenerics_0.23.1
## [4] bindrcpp_0.2 phenopath_0.99.4 forcats_0.2.0
## [7] tidyr_0.7.1 ggplot2_2.2.1 dplyr_0.7.4
## [10] BiocStyle_2.5.39
##
## loaded via a namespace (and not attached):
## [1] viridis_0.4.0 edgeR_3.19.1 viridisLite_0.2.0
## [4] shiny_1.0.3 assertthat_0.2.0 stats4_3.4.0
## [7] vipor_0.4.5 yaml_2.1.14 progress_1.1.2
## [10] RSQLite_1.1-2 backports_1.1.0 lattice_0.20-35
## [13] glue_1.1.1 limma_3.33.9 digest_0.6.12
## [16] RColorBrewer_1.1-2 colorspace_1.3-2 htmltools_0.3.6
## [19] httpuv_1.3.3 Matrix_1.2-10 plyr_1.8.4
## [22] XML_3.98-1.9 pkgconfig_2.0.1 biomaRt_2.33.3
## [25] zlibbioc_1.23.0 bookdown_0.4 purrr_0.2.3
## [28] xtable_1.8-2 scales_0.4.1 tibble_1.3.4
## [31] IRanges_2.11.12 lazyeval_0.2.0 magrittr_1.5
## [34] mime_0.5 memoise_1.1.0 evaluate_0.10
## [37] beeswarm_0.2.3 shinydashboard_0.6.1 tools_3.4.0
## [40] data.table_1.10.4 prettyunits_1.0.2 matrixStats_0.52.2
## [43] stringr_1.2.0 S4Vectors_0.15.5 munsell_0.4.3
## [46] locfit_1.5-9.1 AnnotationDbi_1.39.1 compiler_3.4.0
## [49] rlang_0.1.2 rhdf5_2.21.2 grid_3.4.0
## [52] RCurl_1.95-4.8 tximport_1.5.0 rjson_0.2.15
## [55] bitops_1.0-6 labeling_0.3 rmarkdown_1.6
## [58] gtable_0.2.0 DBI_0.7 reshape2_1.4.2
## [61] R6_2.2.2 gridExtra_2.2.1 knitr_1.16
## [64] bindr_0.1 rprojroot_1.2 stringi_1.1.5
## [67] ggbeeswarm_0.5.3 Rcpp_0.12.13 tidyselect_0.2.0
Blei, David M, Alp Kucukelbir, and Jon D McAuliffe. 2016. “Variational Inference: A Review for Statisticians.” ArXiv Preprint ArXiv:1601.00670.
McCarthy, Davis J, Kieran R Campbell, Aaron TL Lun, and Quin F Wills. 2017. “Scater: Pre-Processing, Quality Control, Normalization and Visualization of Single-Cell Rna-Seq Data in R.” Bioinformatics 33 (8). Oxford University Press: 1179–86.