spateo.alignment.methods.paste
#
Module Contents#
Functions#
|
Calculates and returns optimal alignment of two slices. |
|
|
|
Computes center alignment of slices. |
|
Finds and applies optimal rotation between spatial coordinates of two layers (may also do a reflection). |
- spateo.alignment.methods.paste.paste_pairwise_align(sampleA: anndata.AnnData, sampleB: anndata.AnnData, layer: str = 'X', genes: list | numpy.ndarray | None = None, spatial_key: str = 'spatial', alpha: float = 0.1, dissimilarity: str = 'kl', G_init=None, a_distribution=None, b_distribution=None, norm: bool = False, numItermax: int = 200, numItermaxEmd: int = 100000, dtype: str = 'float32', device: str = 'cpu', verbose: bool = True) Tuple[numpy.ndarray, int | None] [source]#
Calculates and returns optimal alignment of two slices.
- Parameters:
- sampleA
Sample A to align.
- sampleB
Sample B to align.
- layer
If ‘X’, uses
sample.X
to calculate dissimilarity between spots, otherwise uses the representation given bysample.layers[layer]
.- genes
Genes used for calculation. If None, use all common genes for calculation.
- spatial_key
The key in .obsm that corresponds to the raw spatial coordinates.
- alpha
Alignment tuning parameter. Note: 0 <= alpha <= 1. When α = 0 only the gene expression data is taken into account, while when α =1 only the spatial coordinates are taken into account.
- dissimilarity
Expression dissimilarity measure:
'kl'
or'euclidean'
.- G_init
Initial mapping to be used in FGW-OT, otherwise default is uniform mapping.
- a_distribution
Distribution of sampleA spots, otherwise default is uniform.
- b_distribution
Distribution of sampleB spots, otherwise default is uniform.
- norm
If
True
, scales spatial distances such that neighboring spots are at distance 1. Otherwise, spatial distances remain unchanged.- numItermax
Max number of iterations for cg during FGW-OT.
- numItermaxEmd
Max number of iterations for emd during FGW-OT.
- dtype
The floating-point number type. Only float32 and float64.
- device
Equipment used to run the program. You can also set the specified GPU for running. E.g.: ‘0’.
- verbose
If
True
, print progress updates.
- Returns:
Alignment of spots. obj: Objective function output of FGW-OT.
- Return type:
pi
- spateo.alignment.methods.paste.paste_center_align(init_center_sample: anndata.AnnData, samples: List[anndata.AnnData], layer: str = 'X', genes: list | numpy.ndarray | None = None, spatial_key: str = 'spatial', lmbda: numpy.ndarray | None = None, alpha: float = 0.1, n_components: int = 15, threshold: float = 0.001, max_iter: int = 10, numItermax: int = 200, numItermaxEmd: int = 100000, dissimilarity: str = 'kl', norm: bool = False, random_seed: int | None = None, pis_init: List[numpy.ndarray] | None = None, distributions: List[numpy.ndarray] | None = None, dtype: str = 'float32', device: str = 'cpu', verbose: bool = True) Tuple[anndata.AnnData, List[numpy.ndarray]] [source]#
Computes center alignment of slices.
- Parameters:
- init_center_sample
Sample to use as the initialization for center alignment; Make sure to include gene expression and spatial information.
- samples
List of samples to use in the center alignment.
- layer
If ‘X’, uses
sample.X
to calculate dissimilarity between spots, otherwise uses the representation given bysample.layers[layer]
.- genes
Genes used for calculation. If None, use all common genes for calculation.
- spatial_key
The key in .obsm that corresponds to the raw spatial coordinates.
- lmbda
List of probability weights assigned to each slice; If
None
, use uniform weights.- alpha
Alignment tuning parameter. Note: 0 <= alpha <= 1. When α = 0 only the gene expression data is taken into account, while when α =1 only the spatial coordinates are taken into account.
- n_components
Number of components in NMF decomposition.
- threshold
Threshold for convergence of W and H during NMF decomposition.
- max_iter
Maximum number of iterations for our center alignment algorithm.
- numItermax
Max number of iterations for cg during FGW-OT.
- numItermaxEmd
Max number of iterations for emd during FGW-OT.
- dissimilarity
Expression dissimilarity measure:
'kl'
or'euclidean'
.- norm
If
True
, scales spatial distances such that neighboring spots are at distance 1. Otherwise, spatial distances remain unchanged.- random_seed
Set random seed for reproducibility.
- pis_init
Initial list of mappings between ‘A’ and ‘slices’ to solver. Otherwise, default will automatically calculate mappings.
- distributions
Distributions of spots for each slice. Otherwise, default is uniform.
- dtype
The floating-point number type. Only float32 and float64.
- device
Equipment used to run the program. You can also set the specified GPU for running. E.g.: ‘0’.
- verbose
If
True
, print progress updates.
- Returns:
Inferred center sample with full and low dimensional representations (W, H) of the gene expression matrix.
List of pairwise alignment mappings of the center sample (rows) to each input sample (columns).
- spateo.alignment.methods.paste.generalized_procrustes_analysis(X, Y, pi)[source]#
Finds and applies optimal rotation between spatial coordinates of two layers (may also do a reflection).
- Parameters:
- X
np array of spatial coordinates.
- Y
np array of spatial coordinates.
- pi
mapping between the two layers output by PASTE.
- Returns:
Aligned spatial coordinates of X, Y and the mapping relations.