phenopath: PhenoPath - genomic trajectories with heterogeneous...

Description Usage Arguments Details Value See Also Examples

View source: R/phenopath.R

Description

PhenoPath learns genomic trajectories in the presence of heterogenous environmental and genetic backgrounds. It takes input gene expression measurements that are modelled by a single unobserved factor (the "trajectory"). The regulation of genes along the trajectory is perturbed by an additional set of covariates (such as genetic or environmental status) allowing for the identification of covariate-trajectory interactions. The model is fitted using mean-field co-ordinate ascent variational inference.

Usage

1
2
phenopath(exprs_obj, x, sce_assay = "exprs", elbo_tol = 1e-05, z_init = 1,
  ...)

Arguments

exprs_obj

Input gene expression, either

  1. An SummarizedExperiment object, or

  2. A cell-by-gene matrix of normalised expression values in log form.

x

The covariate vector, either

  1. The name of a column of colData(exprs_obj) if exprs_obj is an SummarizedExperiment, or

  2. A numeric of factor vector of length equal to the number of cells, or

  3. A formula from which to build a model matrix from colData(exprs_obj), if exprs_obj is a SummarizedExperiment

sce_assay

The assay from exprs_obj to use as expression if exprs_object is a SummarizedExperiment

elbo_tol

The relative pct change in the ELBO below which is considered converged. See convergence section in details below.

z_init

The initialisation of the latent trajectory. Should be one of

  1. A positive integer describing which principal component of the data should be used for initialisation (default 1), or

  2. A numeric vector of length number of samples to be used directly for initialisation, or

  3. The text character "random", for random initialisation from a standard normal distribution.

...

Additional arguments to be passed to clvm. See description below for more details or call ?clvm.

Details

Input expression

If an SummarizedExperiment is provided, assay(exprs_obj, sce_assay) is used. This is assumed to be in a form that is suitably normalised and approximately normal, such as the log of TPM values (plus a suitable offset) or similar.

Encoding covariates

See the vignette.

Convergence of variational inference

It is strongly recommended to call plot_elbo(phenopath_fit) after the fitting procedure to ensure the ELBO has approximately converged (though convergence metrics are printed during the fitting process). For a good introduction to variational inference see Blei, D.M., Kucukelbir, A. & McAuliffe, J.D., 2017. Variational Inference: A Review for Statisticians. Journal of the American Statistical Association.

Additional arguments

Addition arguments to clvm are passed via .... For full documentation, call ?clvm. Some notable options:

Value

An S3 structure with the following entries:

See Also

clvm for the underlying CAVI function, trajectory to extract the latent trajectory, interaction_effects for the interaction effect sizes, significant_interactions for the results of Bayesian significance testing.

Examples

1
2
3
4
5
sim <- simulate_phenopath() # returns a list with gene expression in y and covariates in x
fit <- phenopath(sim$y, sim$x, elbo_tol = 1e-2)

# Extract the trajectory
z <- trajectory(fit)

phenopath documentation built on Nov. 8, 2020, 6:53 p.m.