API Reference

Core Functions

Factorization

RISE.factorization.pf2(X: anndata.AnnData, rank: int, random_state=1, doEmbedding: bool = True, tolerance=1e-9, max_iter: int = 500)

Perform PARAFAC2 tensor decomposition on single-cell RNA-seq data.

This is the main function for running RISE analysis. It decomposes the multi-condition single-cell data into condition factors, eigen-state factors, and gene factors, revealing patterns across experimental conditions.

Parameters:
  • X (anndata.AnnData) – Preprocessed AnnData object containing single-cell RNA-seq data. Must have X.obs[“condition_unique_idxs”] indicating which condition each cell belongs to (0-indexed).

  • rank (int) – Number of components to extract. Determines the complexity of the decomposition. Typically chosen based on variance explained and Factor Match Score analysis.

  • random_state (int, optional) – Random seed for reproducibility of the decomposition.

  • doEmbedding (bool, optional) – If True, automatically computes PaCMAP embedding of cell projections and stores in X.obsm[“X_pf2_PaCMAP”].

  • tolerance (float, optional) – Convergence threshold for the optimization algorithm.

  • max_iter (int, optional) – Maximum number of iterations for the optimization algorithm.

Returns:

The input AnnData object with added RISE decomposition results in X.uns[“Pf2_weights”], X.uns[“Pf2_A”], X.uns[“Pf2_B”], X.varm[“Pf2_C”], X.obsm[“projections”], X.obsm[“weighted_projections”], and X.obsm[“X_pf2_PaCMAP”]

Return type:

anndata.AnnData

Preprocessing

parafac2.normalize.prepare_dataset(X: anndata.AnnData, condition_name: str, geneThreshold: float, deviance: bool = False)

Preprocess single-cell RNA-seq data for RISE analysis.

This function performs gene filtering, normalization, scaling, and log transformation. It also creates the condition_unique_idxs column required for RISE decomposition.

Parameters:
  • X (anndata.AnnData) – AnnData object containing raw count data in sparse matrix format

  • condition_name (str) – Name of the column in X.obs that specifies experimental conditions for each cell

  • geneThreshold (float) – Minimum mean expression threshold for gene filtering

  • deviance (bool, optional) – If True, applies deviance transformation instead of log normalization

Returns:

Preprocessed AnnData object with condition_unique_idxs and gene means

Return type:

anndata.AnnData

Visualization Functions

General Plotting

RISE.figures.commonFuncs.plotGeneral.plot_r2x(data: anndata.AnnData, rank_vec, ax: matplotlib.axes.Axes)

Plot variance explained (R²X) for RISE and PCA across different ranks.

This visualization helps determine the optimal number of components by showing how variance explained increases with rank. The elbow point where the curve flattens indicates a good balance between model complexity and explanatory power.

Parameters:
  • data (anndata.AnnData) – Preprocessed AnnData object containing single-cell RNA-seq data.

  • rank_vec (array-like of int) – Array of rank values to test (e.g., [1, 5, 10, 15, 20, 25, 30]).

  • ax (matplotlib.axes.Axes) – Matplotlib axes object to plot on.

Factor Plotting

RISE.figures.commonFuncs.plotFactors.plot_condition_factors(data: anndata.AnnData, ax: matplotlib.axes.Axes, cond: str = 'Condition', log_transform: bool = True, cond_group_labels: pandas.Series = None, ThomsonNorm=False, color_key=None, group_cond=False)[source]

Plot condition factors as a heatmap showing how conditions contribute to components.

This visualization shows how each experimental condition (rows) contributes to each RISE component (columns). High values indicate strong association between a condition and a component’s pattern.

Parameters:
  • data (anndata.AnnData) – AnnData object with RISE decomposition results. Must contain data.uns[“Pf2_A”] and data.obs[cond].

  • ax (matplotlib.axes.Axes) – Matplotlib axes object to plot on.

  • cond (str, optional) – Name of column in data.obs containing condition labels.

  • log_transform (bool, optional) – If True, applies log10 transformation to condition factors before plotting.

  • cond_group_labels (pandas.Series, optional) – Series mapping conditions to group labels for colored row annotations.

  • ThomsonNorm (bool, optional) – If True, normalizes factors using only control conditions.

  • color_key (list, optional) – Custom colors for condition group labels.

  • group_cond (bool, optional) – If True and cond_group_labels provided, sorts conditions by group.

RISE.figures.commonFuncs.plotFactors.plot_eigenstate_factors(data: anndata.AnnData, ax: matplotlib.axes.Axes)[source]

Plot eigen-state factors as a heatmap showing cell state patterns.

Eigen-state factors represent the underlying cell state patterns across components. Each row represents an eigen-state and each column represents a component.

Parameters:
  • data (anndata.AnnData) – AnnData object with RISE decomposition results. Must contain data.uns[“Pf2_B”].

  • ax (matplotlib.axes.Axes) – Matplotlib axes object to plot on.

RISE.figures.commonFuncs.plotFactors.plot_gene_factors(data: anndata.AnnData, ax: matplotlib.axes.Axes, weight=0.08, trim=True)[source]

Plot gene factors as a heatmap showing which genes contribute to each component.

This visualization reveals coordinated gene modules by showing which genes (rows) are highly weighted in each component (columns).

Parameters:
  • data (anndata.AnnData) – AnnData object with RISE decomposition results. Must contain data.varm[“Pf2_C”].

  • ax (matplotlib.axes.Axes) – Matplotlib axes object to plot on.

  • weight (float, optional) – Minimum absolute weight threshold for including genes.

  • trim (bool, optional) – If True, filters genes based on the weight parameter.

PaCMAP Visualization

RISE.figures.commonFuncs.plotPaCMAP.plot_labels_pacmap(X: anndata.AnnData, labelType: str, ax: matplotlib.axes.Axes, condition=None, cmap: str = 'tab20', color_key=None)[source]

Plot PaCMAP embedding colored by categorical labels (cell type or condition).

This visualization shows the overall structure of the cell embedding, revealing how cells cluster by cell type or experimental condition.

Parameters:
  • X (anndata.AnnData) – AnnData object with RISE decomposition results. Must contain X.obsm[“X_pf2_PaCMAP”] and X.obs[labelType].

  • labelType (str) – Name of column in X.obs containing categorical labels to color by.

  • ax (matplotlib.axes.Axes) – Matplotlib axes object to plot on.

  • condition (list of str, optional) – If provided, only highlights cells from these specific conditions.

  • cmap (str, optional) – Matplotlib colormap name for coloring categories.

  • color_key (list, optional) – Custom list of colors for categories.

RISE.figures.commonFuncs.plotPaCMAP.plot_gene_pacmap(gene: str, X: anndata.AnnData, ax: matplotlib.axes.Axes, clip_outliers=0.9995)[source]

Plot PaCMAP embedding colored by gene expression levels.

This visualization overlays gene expression onto the PaCMAP embedding of cells, revealing which cell populations express specific genes.

Parameters:
  • gene (str) – Name of gene to visualize. Must be present in X.var_names.

  • X (anndata.AnnData) – AnnData object with RISE decomposition results. Must contain X.obsm[“X_pf2_PaCMAP”].

  • ax (matplotlib.axes.Axes) – Matplotlib axes object to plot on.

  • clip_outliers (float, optional) – Quantile threshold for clipping extreme expression values.

RISE.figures.commonFuncs.plotPaCMAP.plot_wp_pacmap(X: anndata.AnnData, cmp: int, ax: matplotlib.axes.Axes, cbarMax: float = 1.0)[source]

Plot PaCMAP embedding colored by weighted projections for a component.

This visualization shows which cells contribute most strongly to a specific component by coloring them according to their weighted projections.

Parameters:
  • X (anndata.AnnData) – AnnData object with RISE decomposition results. Must contain X.obsm[“X_pf2_PaCMAP”] and X.obsm[“weighted_projections”].

  • cmp (int) – Component number to visualize (1-indexed).

  • ax (matplotlib.axes.Axes) – Matplotlib axes object to plot on.

  • cbarMax (float, optional) – Maximum value for the color scale.

Factor Stability

RISE.figures.figureS4.plot_fms_diff_ranks(X: anndata.AnnData, ax: matplotlib.axes.Axes, ranksList: list[int], runs: int)

Plot Factor Match Score (FMS) across different ranks to assess stability.

FMS measures the reproducibility of PARAFAC2 decomposition results across multiple runs. Values above ~0.6 indicate stable, reproducible components.

Parameters:
  • X (anndata.AnnData) – Preprocessed AnnData object containing single-cell RNA-seq data.

  • ax (matplotlib.axes.Axes) – Matplotlib axes object to plot on.

  • ranksList (list of int) – List of rank values to test (e.g., [1, 5, 10, 15, 20, 25, 30]).

  • runs (int) – Number of independent runs per rank to use for FMS calculation.