dynamo.tl.cell_velocities

dynamo.tl.cell_velocities(adata, ekey=None, vkey=None, X=None, V=None, X_embedding=None, transition_matrix=None, use_mnn=False, n_pca_components=None, transition_genes=None, min_r2=None, min_alpha=None, min_gamma=None, min_delta=None, basis='umap', neighbor_key_prefix='', adj_key='distances', add_transition_key=None, add_velocity_key=None, n_neighbors=30, method='pearson', neg_cells_trick=True, calc_rnd_vel=False, xy_grid_nums=(50, 50), correct_density=True, scale=True, sample_fraction=None, random_seed=19491001, enforce=False, preserve_len=False, **kernel_kwargs)[source]
Project high dimensional velocity vectors onto given low dimensional embeddings, and/or compute cell transition

probabilities.

Parameters:
  • adata (AnnData) – An AnnData object.

  • ekey (Optional[str]) – The dictionary key that corresponds to the gene expression in the layer attribute. If set to be None, ekey will be automatically detected from the adata object. Defaults to None.

  • vkey (Optional[str]) – The dictionary key that corresponds to the estimated velocity values in the layers attribute. If set to be None, vkey will be automatically detected from the adata object. Defaults to None.

  • X (Union[array, csr_matrix, None]) – The expression states of single cells (or expression states in reduced dimension, like pca, of single cells). Defaults to None.

  • V (Union[array, csr_matrix, None]) – The RNA velocity of single cells (or velocity estimates projected to reduced dimension, like pca, of single cells). Note that X, V need to have the exact dimensionalities. Defaults to None.

  • X_embedding (Optional[ndarray]) – The low expression reduced space (pca, umap, tsne, etc.) of single cells that RNA velocity will be projected onto. Note X_embedding, X and V has to have the same cell/sample dimension and X_embedding should have less feature dimension comparing that of X or V. Defaults to None.

  • transition_matrix (Union[ndarray, csr_matrix, None]) – The set of genes used for projection of high dimensional velocity vectors. If None, transition genes are determined based on the R2 of linear regression on phase planes. The argument can be either a dictionary key of .var, a list of gene names, or a list of booleans of length .n_vars. Defaults to None.

  • use_mnn (bool) – Whether to use mutual nearest neighbors for projecting the high dimensional velocity vectors. By default, we don’t use the mutual nearest neighbors. Mutual nearest neighbors are calculated from nearest neighbors across different layers, which accounts for cases where, for example, the cells from spliced expression may be nearest neighbors but far from nearest neighbors on unspliced data. Using mnn assumes your data from different layers are reliable (otherwise it will destroy real signals). Defaults to False.

  • n_pca_components (Optional[int]) – The number of pca components to project the high dimensional X, V before calculating transition matrix for velocity visualization. By default, it is None and if method is kmc, n_pca_components will be reset to 30; otherwise use all high dimensional data for velocity projection. Defaults to None.

  • transition_genes (Union[str, List[str], List[bool], None]) – The set of genes used for projection of high dimensional velocity vectors. If None, transition genes are determined based on the R2 of linear regression on phase planes. The argument can be either a dictionary key of .var, a list of gene names, or a list of booleans of length .n_vars. Defaults to None.

  • min_r2 (Optional[float]) – The minimal value of r-squared of the parameter fits for selecting transition genes. Defaults to None.

  • min_alpha (Optional[float]) – The minimal value of alpha kinetic parameter for selecting transition genes. Defaults to None.

  • min_gamma (Optional[float]) – The minimal value of gamma kinetic parameter for selecting transition genes. Defaults to None.

  • min_delta (Optional[float]) – The minimal value of delta kinetic parameter for selecting transition genes. Defaults to None.

  • basis (str) – The dictionary key that corresponds to the reduced dimension in .obsm attribute. Can be X_spliced_umap or X_total_umap, etc. Defaults to “umap”.

  • neighbor_key_prefix (str) – The dictionary key prefix in .uns. Connectivity and distance matrix keys are also generate with this prefix in adata.obsp. Defaults to “”.

  • adj_key (str) – The dictionary key for the adjacency matrix of the nearest neighbor graph in .obsp. Defaults to “distances”.

  • add_transition_key (Optional[str]) – The dictionary key that will be used for storing the transition matrix in .obsp. Defaults to None.

  • add_velocity_key (Optional[str]) – The dictionary key that will be used for storing the low dimensional velocity projection matrix in .obsm. Defaults to None.

  • n_neighbors (int) – The number of neighbors to be used to calculate velocity projection. Defaults to 30.

  • method (Literal['kmc', 'fp', 'cosine', 'pearson', 'transform']) – The method to calculate the transition matrix and project high dimensional vector to low dimension, either kmc, fp, cosine, pearson, or transform. “kmc” is our new approach to learn the transition matrix via diffusion approximation or an Itô kernel. “cosine” or “pearson” are the methods used in the original RNA velocity paper or the scvelo paper (Note that scVelo implementation actually centers both dX and V, so its cosine kernel is equivalent to pearson correlation kernel, but we also provide the raw cosine kernel). “kmc” option is arguable better than “correlation” or “cosine” as it not only considers the correlation but also the distance of the nearest neighbors to the high dimensional velocity vector. Finally, the “transform” method uses umap’s transform method to transform new data points to the UMAP space. “transform” method is NOT recommended. Kernels that are based on the reconstructed vector field in high dimension is also possible. Defaults to “pearson”.

  • neg_cells_trick (bool) – Whether to handle cells having negative correlations in gene expression difference with high dimensional velocity vector separately. This option was borrowed from scVelo package (https://github.com/theislab/scvelo) and use in conjunction with “pearson” and “cosine” kernel. Not required if method is set to be “kmc”. Defaults to True.

  • calc_rnd_vel (bool) – Whether to calculate the random velocity vectors which can be plotted downstream as a negative control and used to adjust the quiver scale of the velocity field. Defaults to False.

  • xy_grid_nums (Tuple[int]) – A tuple of number of grids on each dimension. Defaults to (50, 50).

  • correct_density (bool) – Whether to correct density when calculating the markov transition matrix. Defaults to True.

  • scale (bool) – whether to scale velocity when calculating the markov transition matrix, applicable to the kmc kernel. Defaults to True.

  • sample_fraction (Optional[float]) – The downsampled fraction of kNN for the purpose of acceleration, applicable to the kmc kernel. Defaults to None.

  • random_seed (int) – The random seed for numba to ensure consistency of the random velocity vectors. Default value 19491001 is a special day for those who care. Defaults to 19491001.

  • enforce (bool) –

    Whether to enforce 1) redefining use_for_transition column in obs attribute; However this is NOT executed if the argument

    ’transition_genes’ is not None.

    1. recalculation of the transition matrix. Defaults to False.

  • preserve_len (bool) – Whether to preserve the length of high dimension vector length. When set to be True, the length of low dimension projected vector will be proportionally scaled to that of the high dimensional vector. Defaults to False.

  • kernel_kwargs – Kwargs that would be passed to the kernel for constructing the transition matrix.

Raises:
  • Exception – Neighborhood info is invalid.

  • TypeError – Transition gene list is invalid.

  • ValueError – Provided transition genes do not have velocity data.

  • ExceptionX and V have different dimensions.

  • ExceptionX and X_embedding has different number of samples.

  • Exception – Number of dimension of X is smaller than the one of X_embedding.

  • Exception – Most calculated velocity is theoretically invalid.

  • NotImplementedError – The mode provided in kernel_kwargs is invalid.

Return type:

AnnData

Returns:

An updated AnnData object with projected velocity vectors, and a cell transition matrix calculated using either the Itô kernel method or similar methods from (La Manno et al. 2018).