spateo.tools.paste#

Module Contents#

Functions#

filter_common_genes(→ list)

Filters for the intersection of genes between all samples.

kl_divergence_backend(X, Y)

Returns pairwise KL divergence (over all pairs of samples) of two matrices X and Y.

check_backend([device, dtype])

Check the proper backend for the device.

check_spatial_coords(→ numpy.ndarray)

Check spatial coordinate information.

check_exp(→ numpy.ndarray)

Check expression matrix.

pairwise_align(→ Tuple[numpy.ndarray, Optional[int]])

Calculates and returns optimal alignment of two slices.

center_NMF(n_components, random_seed[, dissimilarity])

center_align(→ Tuple[anndata.AnnData, List[numpy.ndarray]])

Computes center alignment of slices.

generalized_procrustes_analysis(X, Y, pi)

Finds and applies optimal rotation between spatial coordinates of two layers (may also do a reflection).

_get_optimal_mapping_relationship(X, Y, pi[, keep_all])

mapping_aligned_coords(→ Tuple[dict, dict])

Optimal mapping coordinates between X and Y.

mapping_center_coords(→ dict)

Optimal mapping coordinates between X and Y based on intermediate coordinates.

Attributes#

spateo.tools.paste.intersect_lsts[source]#
spateo.tools.paste.to_dense_matrix[source]#
spateo.tools.paste.filter_common_genes(*genes) list[source]#

Filters for the intersection of genes between all samples.

Parameters
genes

List of genes.

spateo.tools.paste.kl_divergence_backend(X, Y)[source]#

Returns pairwise KL divergence (over all pairs of samples) of two matrices X and Y.

Takes advantage of POT backend to speed up computation.

Parameters
X

np array with dim (n_samples by n_features)

Y

np array with dim (m_samples by n_features)

Returns

np array with dim (n_samples by m_samples). Pairwise KL divergence matrix.

Return type

D

spateo.tools.paste.check_backend(device: str = 'cpu', dtype: str = 'float32')[source]#

Check the proper backend for the device.

Parameters
device

Equipment used to run the program. You can also set the specified GPU for running. E.g.: ‘0’.

dtype

The floating-point number type. Only float32 and float64.

Returns

The proper backend. type_as: The type_as.device is the device used to run the program and the type_as.dtype is the floating-point number type.

Return type

backend

spateo.tools.paste.check_spatial_coords(sample: anndata.AnnData, spatial_key: str = 'spatial') numpy.ndarray[source]#

Check spatial coordinate information.

Parameters
sample

An anndata object.

spatial_key

The key in .obsm that corresponds to the raw spatial coordinates.

Returns

The spatial coordinates.

spateo.tools.paste.check_exp(sample: anndata.AnnData, layer: str = 'X') numpy.ndarray[source]#

Check expression matrix.

Parameters
sample

An anndata object.

layer

The key in .layers that corresponds to the expression matrix.

Returns

The expression matrix.

spateo.tools.paste.pairwise_align(sampleA: anndata.AnnData, sampleB: anndata.AnnData, layer: str = 'X', genes: Optional[Union[list, numpy.ndarray]] = 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') Tuple[numpy.ndarray, Optional[int]][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 by sample.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 : array-like, optional

Initial mapping to be used in FGW-OT, otherwise default is uniform mapping.

a_distribution : array-like, optional

Distribution of sampleA spots, otherwise default is uniform.

b_distribution : array-like, optional

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’.

Returns

Alignment of spots. obj: Objective function output of FGW-OT.

Return type

pi

spateo.tools.paste.center_NMF(n_components, random_seed, dissimilarity='kl')[source]#
spateo.tools.paste.center_align(init_center_sample: anndata.AnnData, samples: List[anndata.AnnData], layer: str = 'X', genes: Optional[Union[list, numpy.ndarray]] = None, spatial_key: str = 'spatial', lmbda: Optional[numpy.ndarray] = 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: Optional[int] = None, pis_init: Optional[List[numpy.ndarray]] = None, distributions: Optional[List[numpy.ndarray]] = None, dtype: str = 'float32', device: str = 'cpu') 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 by sample.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’.

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.tools.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.

spateo.tools.paste._get_optimal_mapping_relationship(X: numpy.ndarray, Y: numpy.ndarray, pi: numpy.ndarray, keep_all: bool = False)[source]#
spateo.tools.paste.mapping_aligned_coords(X: numpy.ndarray, Y: numpy.ndarray, pi: numpy.ndarray, keep_all: bool = False) Tuple[dict, dict][source]#

Optimal mapping coordinates between X and Y.

Parameters
X

Aligned spatial coordinates.

Y

Aligned spatial coordinates.

pi

Mapping between the two layers output by PASTE.

keep_all

Whether to retain all the optimal relationships obtained only based on the pi matrix, If keep_all is False, the optimal relationships obtained based on the pi matrix and the nearest coordinates.

Returns

Two dicts of mapping_X, mapping_Y, pi_index, pi_value.

mapping_X is X coordinates aligned with Y coordinates. mapping_Y is the Y coordinate aligned with X coordinates. pi_index is index between optimal mapping points in the pi matrix. pi_value is the value of optimal mapping points.

spateo.tools.paste.mapping_center_coords(modelA: anndata.AnnData, modelB: anndata.AnnData, center_key: str) dict[source]#

Optimal mapping coordinates between X and Y based on intermediate coordinates.

Parameters
modelA

modelA aligned with center model.

modelB

modelB aligned with center model.

center_key

The key in .uns that corresponds to the alignment info between modelA/modelB and center model.

Returns

A dict of raw_X, raw_Y, mapping_X, mapping_Y, pi_value.

raw_X is the raw X coordinates. raw_Y is the raw Y coordinates. mapping_X is the Y coordinates aligned with X coordinates. mapping_Y is the X coordinates aligned with Y coordinates. pi_value is the value of optimal mapping points.