Title: | Stan implementation of BASiCS |
---|---|
Description: | Provides an interface to infer the parameters of BASiCS using the variational inference (ADVI), Markov chain Monte Carlo (NUTS), and maximum a posteriori (BFGS) inference engines in the Stan programming language. BASiCS is a Bayesian hierarchical model that uses an adaptive Metropolis within Gibbs sampling scheme. Alternative inference methods provided by Stan may be preferable in some situations, for example for particularly large data or posterior distributions with difficult geometries. |
Authors: | Alan O'Callaghan [aut, cre], Catalina Vallejos [aut] |
Maintainer: | Alan O'Callaghan <[email protected]> |
License: | GPL-3 |
Version: | 1.9.0 |
Built: | 2024-10-30 04:21:50 UTC |
Source: | https://github.com/bioc/BASiCStan |
Provides an interface to infer the parameters of BASiCS using the variational inference (ADVI), Markov chain Monte Carlo (NUTS), and maximum a posteriori (BFGS) inference engines in the Stan programming language. BASiCS is a Bayesian hierarchical model that uses an adaptive Metropolis within Gibbs sampling scheme. Alternative inference methods provided by Stan may be preferable in some situations, for example for particularly large data or posterior distributions with difficult geometries. See also BASiCS_MCMC
Stan Development Team (2020). RStan: the R interface to Stan. R package version 2.21.2. https://mc-stan.org
Vallejos, Marioni and Richardson (2015). PLoS Computational Biology. https://doi.org/10.1371/journal.pcbi.1004333.
Vallejos, Richardson and Marioni (2016). Genome Biology. https://doi.org/10.1186/s13059-016-0930-3
Eling et al (2018). Cell Systems. https://doi.org/10.1016/j.cels.2018.06.011
The stan programming language enables the use of MAP, VB, and HMC inference. Only the regression mode featuring a joint prior between mean and overdispersion parameters is implemented
BASiCStan( Data, Method = c("vb", "sampling", "optimizing"), WithSpikes = length(altExpNames(Data)) > 0, Regression = TRUE, BatchInfo = Data$BatchInfo, L = 12, PriorMu = c("EmpiricalBayes", "uninformative"), NormFactorFun = scran::calculateSumFactors, ReturnBASiCS = TRUE, Verbose = TRUE, ... )
BASiCStan( Data, Method = c("vb", "sampling", "optimizing"), WithSpikes = length(altExpNames(Data)) > 0, Regression = TRUE, BatchInfo = Data$BatchInfo, L = 12, PriorMu = c("EmpiricalBayes", "uninformative"), NormFactorFun = scran::calculateSumFactors, ReturnBASiCS = TRUE, Verbose = TRUE, ... )
Data |
SingleCellExperiment object |
Method |
Inference method. One of: |
WithSpikes |
Do the data contain spike-in genes? See BASiCS for details.
When |
Regression |
Use joint prior for mean and overdispersion parameters?
Included for compatibility with |
BatchInfo |
Vector describing which batch each cell is from. |
L |
Number of regression terms (including slope and intercept) to use in joint prior for mu and delta. |
PriorMu |
Type of prior to use for mean expression. Default is "EmpiricalBayes", but "uninformative" is the prior used in Eling et al. and previous work. |
NormFactorFun |
Function that returns cell-specific scaling
normalisation factors. See |
ReturnBASiCS |
Should the object be converted into a BASiCS_Chain object? |
Verbose |
Should output of the stan commands be printed to the terminal? |
... |
Passed to vb or sampling. |
An object of class BASiCS_Chain
.
library("BASiCS") sce <- BASiCS_MockSCE(NGenes = 10, NCells = 10) fit_spikes <- BASiCStan(sce, tol_rel_obj = 1) ## uses fixed scaling normalisation factors fit_nospikes <- BASiCStan(sce, WithSpikes = FALSE, tol_rel_obj = 1)
library("BASiCS") sce <- BASiCS_MockSCE(NGenes = 10, NCells = 10) fit_spikes <- BASiCStan(sce, tol_rel_obj = 1) ## uses fixed scaling normalisation factors fit_nospikes <- BASiCStan(sce, WithSpikes = FALSE, tol_rel_obj = 1)
BASiCS_Chain
objects.Convert Stan fits to BASiCS_Chain
objects.
Stan2BASiCS( x, gene_names = attr(x, "gene_names"), cell_names = attr(x, "cell_names"), size_factors = attr(x, "size_factors") )
Stan2BASiCS( x, gene_names = attr(x, "gene_names"), cell_names = attr(x, "cell_names"), size_factors = attr(x, "size_factors") )
x |
A stan object |
gene_names , cell_names
|
Gene and cell names. The reason this argument
exists is that by default, stan fit parameters are not named.
NOTE: this must be the same order as the
data supplied to |
size_factors |
Cell-specific scaling normalisation factors, to be
stored as part of the chain object when |
A BASiCS_Chain
object.
library("BASiCS") sce <- BASiCS_MockSCE(NGenes = 10, NCells = 10) fit_spikes <- BASiCStan(sce, ReturnBASiCS = FALSE, tol_rel_obj = 1) summary(fit_spikes)
library("BASiCS") sce <- BASiCS_MockSCE(NGenes = 10, NCells = 10) fit_spikes <- BASiCStan(sce, ReturnBASiCS = FALSE, tol_rel_obj = 1) summary(fit_spikes)