pyrovelocity.utils#

pyrovelocity.utils.anndata_counts_to_df(adata)[source]#
pyrovelocity.utils.attributes(obj)[source]#

get object attributes

pyrovelocity.utils.debug(x)[source]#
pyrovelocity.utils.ensure_numpy_array(obj)[source]#
pyrovelocity.utils.filter_startswith_dict(dictionary_with_underscore_keys)[source]#

Remove entries from a dictionary whose keys start with an underscore.

Parameters:

dictionary_with_underscore_keys (dict) – Dictionary to be filtered.

Returns:

Filtered dictionary.

Return type:

dict

pyrovelocity.utils.generate_sample_data(n_obs=100, n_vars=12, alpha=5, beta=0.5, gamma=0.3, alpha_=0, noise_model='gillespie', random_seed=0)[source]#

Generate synthetic single-cell RNA sequencing data with spliced and unspliced layers. If using the “iid” noise model, the data will be generated with scvi.data.synthetic_iid. If using the “normal” or “gillespie” noise model, the data will be generated with scvelo.datasets.simulation accounting for the given expression dynamics parameters.

Parameters:
  • n_obs (int, optional) – Number of observations (cells). Default is 100.

  • n_vars (int, optional) – Number of variables (genes). Default is 12.

  • alpha (float, optional) – Transcription rate. Default is 5.

  • beta (float, optional) – Splicing rate. Default is 0.5.

  • gamma (float, optional) – Degradation rate. Default is 0.3.

  • alpha – Additional transcription rate. Default is 0.

  • noise_model (str, optional) – Noise model to be used. Must be one of ‘iid’, ‘gillespie’, or ‘normal’. Default is ‘gillespie’.

  • random_seed (int, optional) – Random seed for reproducibility. Default is 0.

Returns:

An AnnData object containing the generated synthetic data.

Return type:

anndata.AnnData

Raises:

ValueError – If noise_model is not one of ‘iid’, ‘gillespie’, or ‘normal’.

Examples

>>> from pyrovelocity.utils import generate_sample_data, print_anndata
>>> adata = generate_sample_data(random_seed=99)
>>> print_anndata(adata)
>>> adata = generate_sample_data(n_obs=50, n_vars=10, noise_model="normal", random_seed=99)
>>> print_anndata(adata)
>>> adata = generate_sample_data(n_obs=50, n_vars=10, noise_model="iid", random_seed=99)
>>> print_anndata(adata)
>>> adata = generate_sample_data(noise_model="wishful thinking")
Traceback (most recent call last):
    ...
ValueError: noise_model must be one of 'iid', 'gillespie', 'normal'
pyrovelocity.utils.get_pylogger(name='pyrovelocity.utils', log_level='DEBUG')[source]#

Initializes multi-GPU-friendly python command line logger.

Return type:

Logger

pyrovelocity.utils.get_velocity_samples(posterior_samples, model)[source]#

Computes the velocity samples from the posterior samples.

Parameters:
  • posterior_samples (dict) – Dictionary containing posterior samples.

  • model – Model used for predictions.

Returns:

Velocity samples.

Return type:

torch.Tensor

Examples

>>> import torch
>>> posterior_samples = {
...     "beta": torch.tensor([[0.4]]),
...     "gamma": torch.tensor([[0.3]]),
...     "u_scale": torch.tensor([[[1.0, 2.0]]]),
...     "s_scale": torch.tensor([[[2.0, 4.0]]]),
...     "u": torch.tensor([[1.0, 2.0]]),
...     "s": torch.tensor([[0.5, 1.0]])
... }
>>> model = None  # Model is not used in the function
>>> get_velocity_samples(posterior_samples, model)
tensor([[0.6500, 1.3000]])
pyrovelocity.utils.init_with_all_cells(adata, input_type='knn', shared_time=True, latent_factor='linear', latent_factor_size=10, plate_size=2, num_aux_cells=200, init_smooth=True)[source]#
pyrovelocity.utils.inv(x)[source]#

Computes the element-wise reciprocal of a tensor.

Parameters:

x (torch.Tensor) – Input tensor.

Returns:

Tensor with element-wise reciprocal of x.

Return type:

torch.Tensor

Examples

>>> import torch
>>> x = torch.tensor([2., 4., 0.5])
>>> inv(x)
tensor([0.5000, 0.2500, 2.0000])
pyrovelocity.utils.log(x)[source]#

Computes the element-wise natural logarithm of a tensor, while clipping its values to avoid numerical instability.

Parameters:

x (torch.Tensor) – Input tensor.

Returns:

Tensor with element-wise natural logarithm of x.

Return type:

torch.Tensor

Examples

>>> import torch
>>> x = torch.tensor([0.0001, 0.5, 0.9999])
>>> log(x)
tensor([-9.2103e+00, -6.9315e-01, -1.0002e-04])
pyrovelocity.utils.mRNA(tau, u0, s0, alpha, beta, gamma)[source]#

Computes the mRNA dynamics given the parameters and initial conditions.

Parameters:
  • tau (torch.Tensor) – Time points.

  • u0 (torch.Tensor) – Initial value of u.

  • s0 (torch.Tensor) – Initial value of s.

  • alpha (torch.Tensor) – Alpha parameter.

  • beta (torch.Tensor) – Beta parameter.

  • gamma (torch.Tensor) – Gamma parameter.

Returns:

Tuple containing the final values of u and s.

Return type:

Tuple[torch.Tensor, torch.Tensor]

Examples

>>> import torch
>>> tau = torch.tensor(2.0)
>>> u0 = torch.tensor(1.0)
>>> s0 = torch.tensor(0.5)
>>> alpha = torch.tensor(0.5)
>>> beta = torch.tensor(0.4)
>>> gamma = torch.tensor(0.3)
>>> mRNA(tau, u0, s0, alpha, beta, gamma)
(tensor(1.1377), tensor(0.9269))
pyrovelocity.utils.mae(pred_counts, true_counts)[source]#

Computes the mean average error between predicted counts and true counts.

Parameters:
  • pred_counts (np.ndarray) – Predicted counts.

  • true_counts (np.ndarray) – True counts.

Returns:

Mean average error.

Return type:

float

Examples

>>> import numpy as np
>>> pred_counts = np.array([[1, 2], [3, 4]])
>>> true_counts = np.array([[2, 3], [4, 5]])
>>> mae(pred_counts, true_counts)
1.0
pyrovelocity.utils.mae_evaluate(posterior_samples, adata)[source]#
pyrovelocity.utils.mse_loss_sum(u_model, s_model, u_data, s_data)[source]#

Computes the mean squared error loss sum between the model and data.

Parameters:
  • u_model (torch.Tensor) – Predicted values of u from the model.

  • s_model (torch.Tensor) – Predicted values of s from the model.

  • u_data (torch.Tensor) – True values of u from the data.

  • s_data (torch.Tensor) – True values of s from the data.

Returns:

Mean squared error loss sum.

Return type:

torch.Tensor

Examples

>>> import torch
>>> u_model = torch.tensor([0.5, 0.6])
>>> s_model = torch.tensor([0.7, 0.8])
>>> u_data = torch.tensor([0.4, 0.5])
>>> s_data = torch.tensor([0.6, 0.7])
>>> mse_loss_sum(u_model, s_model, u_data, s_data)
tensor(0.0200)
pyrovelocity.utils.ode_mRNA(tau, u0, s0, alpha, beta, gamma)[source]#

Solves the ODE system for mRNA dynamics.

Parameters:
  • tau (torch.Tensor) – Time points.

  • u0 (torch.Tensor) – Initial value of u.

  • s0 (torch.Tensor) – Initial value of s.

  • alpha (torch.Tensor) – Alpha parameter.

  • beta (torch.Tensor) – Beta parameter.

  • gamma (torch.Tensor) – Gamma parameter.

Returns:

Tuple containing the final values of u and s.

Return type:

Tuple[torch.Tensor, torch.Tensor]

pyrovelocity.utils.pretty_print_dict(d)[source]#
pyrovelocity.utils.print_anndata(anndata_obj)[source]#

Print a formatted representation of an AnnData object.

This function produces a custom output for the AnnData object with each element of obs, var, uns, obsm, varm, layers, obsp, varp indented and displayed on a new line.

Parameters:

anndata_obj (anndata.AnnData) – The AnnData object to be printed.

Raises:

AssertionError – If the input object is not an instance of anndata.AnnData.

Examples

>>> import anndata
>>> import numpy as np
>>> import pandas as pd
>>> np.random.seed(42)
>>> X = np.random.randn(10, 5)
>>> obs = pd.DataFrame({"clusters_coarse": np.random.randint(0, 2, 10),
...                     "clusters": np.random.randint(0, 2, 10),
...                     "S_score": np.random.rand(10),
...                     "G2M_score": np.random.rand(10)})
>>> var = pd.DataFrame({"gene_name": [f"gene_{i}" for i in range(5)]})
>>> adata = anndata.AnnData(X, obs=obs, var=var)
>>> print_anndata(adata)  
AnnData object with n_obs × n_vars = 10 × 5
    obs:
        clusters_coarse,
        clusters,
        S_score,
        G2M_score,
    var:
        gene_name,
pyrovelocity.utils.print_attributes(obj)[source]#

print object attributes

pyrovelocity.utils.rescale_time(dx_dt, t_start, t_end)[source]#

Converts an ODE to be solved on a batch of different time intervals into an equivalent system of ODEs to be solved on [0, 1].

Parameters:
  • dx_dt (Callable) – Function representing the ODE system.

  • t_start (torch.Tensor) – Start time of the time interval.

  • t_end (torch.Tensor) – End time of the time interval.

Returns:

Function representing the rescaled ODE system.

Return type:

Callable

Examples

>>> import torch
>>> def dx_dt(t, x):
...     return -x
>>> t_start = torch.tensor(0.0)
>>> t_end = torch.tensor(1.0)
>>> rescaled_dx_dt = rescale_time(dx_dt, t_start, t_end)
>>> rescaled_dx_dt(torch.tensor(0.5), torch.tensor(2.0))
tensor(-2.)
pyrovelocity.utils.site_is_discrete(site)[source]#
Return type:

bool

pyrovelocity.utils.tau_inv(u=None, s=None, u0=None, s0=None, alpha=None, beta=None, gamma=None)[source]#

Computes the inverse tau given the parameters and initial conditions.

Parameters:
  • u (torch.Tensor) – Value of u.

  • s (torch.Tensor) – Value of s.

  • u0 (torch.Tensor) – Initial value of u.

  • s0 (torch.Tensor) – Initial value of s.

  • alpha (torch.Tensor) – Alpha parameter.

  • beta (torch.Tensor) – Beta parameter.

  • gamma (torch.Tensor) – Gamma parameter.

Returns:

Inverse tau.

Return type:

torch.Tensor

Examples

>>> import torch
>>> u = torch.tensor(0.6703)
>>> s = torch.tensor(0.4596)
>>> u0 = torch.tensor(1.0)
>>> s0 = torch.tensor(0.5)
>>> alpha = torch.tensor(0.5)
>>> beta = torch.tensor(0.4)
>>> gamma = torch.tensor(0.3)
>>> tau_inv(u, s, u0, s0, alpha, beta, gamma)
tensor(3.9736e-07)
pyrovelocity.utils.trace(func)[source]#
pyrovelocity.utils.velocity_dus_dt(alpha, beta, gamma, tau, x)[source]#

Computes the velocity du/dt and ds/dt.

Parameters:
  • alpha (torch.Tensor) – Alpha parameter.

  • beta (torch.Tensor) – Beta parameter.

  • gamma (torch.Tensor) – Gamma parameter.

  • tau (torch.Tensor) – Time points.

  • x (Tuple[torch.Tensor, torch.Tensor]) – Tuple containing u and s.

Returns:

Tuple containing du/dt and ds/dt.

Return type:

Tuple[torch.Tensor, torch.Tensor]

Examples

>>> import torch
>>> alpha = torch.tensor(0.5)
>>> beta = torch.tensor(0.4)
>>> gamma = torch.tensor(0.3)
>>> tau = torch.tensor(2.0)
>>> x = (torch.tensor(1.0), torch.tensor(0.5))
>>> velocity_dus_dt(alpha, beta, gamma, tau, x)
(tensor(0.1000), tensor(0.2500))