% Generated by roxygen2: do not edit by hand % Please edit documentation in R/scvelo.R \name{scvelo} \alias{scvelo} \alias{scvelo,ANY-method} \alias{scvelo,SummarizedExperiment-method} \alias{scvelo,SingleCellExperiment-method} \title{RNA velocity with \pkg{scVelo}} \usage{ scvelo(x, ...) \S4method{scvelo}{ANY}( x, subset.row = NULL, sf.X = NULL, sf.spliced = NULL, sf.unspliced = NULL, use.theirs = FALSE, mode = c("steady_state", "deterministic", "stochastic", "dynamical"), scvelo.params = list(), dimred = NULL, ncomponents = 30, BPPARAM = SerialParam(), BSPARAM = bsparam() ) \S4method{scvelo}{SummarizedExperiment}( x, ..., assay.X = "counts", assay.spliced = "spliced", assay.unspliced = "unspliced" ) \S4method{scvelo}{SingleCellExperiment}(x, ..., sf.X = sizeFactors(x), dimred = NULL, use.dimred = NULL) } \arguments{ \item{x}{A named list of three matrices of the same dimensions where genes are in rows and cells are in columns. The list should contain \code{"spliced"} and \code{"unspliced"} entries containing spliced and unspliced counts, respectively. It should also contain an \code{"X"} entry containing the \dQuote{usual} count matrix, see details below. Alternatively, a \linkS4class{SummarizedExperiment} object containing three such matrices among its assays.} \item{...}{For the generic, further arguments to pass to specific methods. For the SummarizedExperiment and SingleCellExperiment methods, further arguments to pass to the ANY method.} \item{subset.row}{A character, integer or logical vector specifying the genes to use for the velocity calculations. Defaults to all genes.} \item{sf.X}{A numeric vector containing size factors for usual count matrix. Defaults to \code{\link{librarySizeFactors}} on the \code{"X"} matrix in \code{x}.} \item{sf.spliced}{A numeric vector containing size factors for the spliced counts for each cell. Defaults to \code{\link{librarySizeFactors}} on the \code{"spliced"} matrix in \code{x}.} \item{sf.unspliced}{A numeric vector containing size factors for the unspliced counts for each cell. Defaults to \code{\link{librarySizeFactors}} on the \code{"unspliced"} matrix in \code{x}.} \item{use.theirs}{Logical scalar indicating whether \pkg{scVelo}'s gene filtering and normalization should be used.} \item{mode}{String specifying the method to use to estimate the transcriptional dynamics.} \item{scvelo.params}{List of lists containing arguments for individual \pkg{scVelo} functions, see details below.} \item{dimred}{A low-dimensional representation of the cells with number of rows equal to the number of cells in \code{x}, used to find the nearest neighbors.} \item{ncomponents}{Numeric scalar indicating the number of principal components to obtain. Only used if \code{use.theirs=FALSE} and \code{dimred=NULL}.} \item{BPPARAM}{A \linkS4class{BiocParallelParam} object specifying whether the PCA calculations should be parallelized. Only used if \code{use.theirs=FALSE} and \code{dimred=NULL}.} \item{BSPARAM}{A \linkS4class{BiocSingularParam} object specifying which algorithm should be used to perform the PCA. Only used if \code{use.theirs=FALSE} and \code{dimred=NULL}.} \item{assay.X}{An integer scalar or string specifying the assay of \code{x} containing the usual count matrix.} \item{assay.spliced}{An integer scalar or string specifying the assay of \code{x} containing the spliced counts.} \item{assay.unspliced}{An integer scalar or string specifying the assay of \code{x} containing the unspliced counts.} \item{use.dimred}{String naming the entry of \code{\link{reducedDims}(x)} to use for nearest neighbor calculations. Ignored if \code{dimred} is supplied.} } \value{ A \linkS4class{SingleCellExperiment} is returned containing the output of the velocity calculations. Of particular interest are: \itemize{ \item the \code{velocity_pseudotime} field in the \code{\link{colData}}, containing the velocity pseudotime for each cell. \item the \code{velocity} entry of the \code{\link{assays}}, containing the velocity vectors for each cell. } The output will always have number of columns equal to the number of cells supplied in \code{x}, though the number of rows will depend on whether any subsetting (if \code{subset.row} is supplied) or feature selection (if \code{use.theirs=TRUE}) was performed. } \description{ Perform RNA velocity calculations with the \pkg{scVelo} package. } \details{ This function uses the \pkg{scVelo} Python package (\url{https://pypi.org/project/scvelo/}) for RNA velocity calculations. The main difference from the original \pkg{velocyto} approach is that the dynamical model of \pkg{scVelo} does not rely on the presence of observed steady-state populations, which should improve the reliability of the velocity calculations in general applications. For consistency with other Bioconductor workflows, we perform as many standard steps in R as we can before starting the velocity calculations with \pkg{scVelo}. This involves: \enumerate{ \item Size factor-based normalization with \code{sf.*} values and \code{\link{normalizeCounts}}. For \code{"X"}, log-transformation is performed as well, while for the others, only scaling normalization is performed. \item Subsetting all matrices to \code{subset.row}, most typically to a subset of interest, e.g., highly variable genes. Note that, if set, any subsetting is done \emph{after} normalization so that library sizes are correctly computed. \item If \code{dimred=NULL}, the PCA step on the log-expression values derived from the \code{"X"} matrix, using the specified \code{BSPARAM} to obtain the first \code{ncomponents} PCs. } This allows us to guarantee that, for example, the log-expression matrix of HVGs or the PCA coordinates are the same as that used in other applications like clustering or trajectory reconstruction. Nonetheless, one can set \code{use.theirs=TRUE} to directly use the entire \pkg{scVelo} normalization and filtering pipeline. This ignores all of the size factors arguments (\code{sf.*}), all of the PCA-related arguments (\code{ncomponents}, \code{BSPARAM}) and \code{subset.row}. However, if a low-dimensionality result is supplied via \code{dimred} or \code{use.dimred}, the \pkg{scVelo} PCA will always be omitted. Upon first use, this function will instantiate a conda environment containing the \pkg{scVelo} package. This is done via the \pkg{basilisk} package - see the documentation for that package for trouble-shooting. } \section{Comments on the three matrices}{ Strictly speaking, only the spliced and unspliced matrices are necessary for the velocity calculations. However, it is often the case that the spliced matrix is not actually the same as a \dQuote{usual} count matrix (e.g., generated by summing counts across all exons for all mapped genes). This is due to differences in the handling of ambiguous reads that map across exon-intron boundaries, or to genomic regions that can be either exonic or intronic depending on the isoform; the spliced count matrix is more likely to exclude such reads. We request the usual count matrix as the \code{"X"} entry of \code{x} to ensure that the PCA and nearest neighbor detection in \pkg{scVelo} are done on the same data as that used in other steps of the large analysis (e.g., clustering, visualization, trajectory reconstruction). In practice, if the usual count matrix is not available, one can often achieve satisfactory results by simply re-using the spliced count matrix as both the \code{"X"} and \code{"spliced"} entries of \code{x}. Note that if reduced dimensions are supplied in \code{dimred}, any \code{"X"} entry is only used to create the AnnData object and is not used in any actual calculations. } \section{Additional arguments to Python}{ Additional arguments to \pkg{scVelo} functions are provided via \code{scvelo.params}. This is a named list where each entry is named after a function and is itself a named list of arguments for that function. The following function names are currently recognized: \itemize{ \item \code{"filter_and_normalize"}, for gene selection and normalization. This is not used unless \code{use.theirs=TRUE}. \item \code{"moments"}, for PCA and nearest neighbor detection. The PCA is not performed if \code{dimred} or \code{use.dimred} is already supplied. \item \code{"recover_dynamics"} \item \code{"velocity"} \item \code{"velocity_graph"} \item \code{"velocity_pseudotime"} \item \code{"latent_time"} \item \code{"velocity_confidence"} } See the \pkg{scVelo} documentation for more details about the available arguments and the examples below for a syntax example. } \examples{ # Using mock data to demonstrate the process: library(scuttle) sce1 <- mockSCE() sce2 <- mockSCE() spliced <- counts(sce1) unspliced <- counts(sce2) out <- scvelo(list(X=spliced, spliced=spliced, unspliced=unspliced)) # make scvelo use 10 rather than the default 30 neighbors to compute moments for velocity estimation: out <- scvelo(list(X=spliced, spliced=spliced, unspliced=unspliced), scvelo.params=list(moments=list(n_neighbors=10L))) } \references{ Bergen, V., Lange, M., Peidli, S. et al. Generalizing RNA velocity to transient cell states through dynamical modeling. Nat Biotechnol 38, 1408–1414 (2020). \url{https://doi.org/10.1038/s41587-020-0591-3} } \author{ Aaron Lun, Charlotte Soneson }