cospar.tmap.infer_Tmap_from_multitime_clones

cospar.tmap.infer_Tmap_from_multitime_clones(adata_orig, clonal_time_points=None, later_time_point=None, smooth_array=[15, 10, 5], CoSpar_KNN=20, sparsity_threshold=0.1, intraclone_threshold=0.05, normalization_mode=1, extend_Tmap_space=False, save_subset=True, use_full_Smatrix=True, trunca_threshold=[0.001, 0.01], compute_new=False, max_iter_N=5, epsilon_converge=0.05)

Infer transition map for clonal data with multiple time points.

It prepares adata object for cells of targeted time points by cospar.tmap._utils.select_time_points(), generates the similarity matrix via cospar.tmap._utils.generate_similarity_matrix(), and iteratively calls the core function refine_Tmap_through_cospar() to update the transition map.

  • If later_time_point=None, the inferred map allows transitions between neighboring time points. For example, if clonal_time_points=[‘day1’,’day2’,’day3’], then it computes transitions for pairs (day1, day2) and (day2, day3), but not (day1, day3).

  • If later_time_point is specified, the function produces a map between earlier time points and this later time point. For example, if later_time_point=’day3, the map allows transitions for pairs (day1, day3) and (day2, day3), but not (day1,day2).

Parameters
adata_orig : AnnData object

Should be prepared from our anadata initialization.

clonal_time_points : list of str, optional (default: all time points)

List of time points to be included for analysis. We assume that each selected time point has clonal measurements.

later_time_points : list, optional (default: None)

If specified, the function will produce a map T between these early time points among clonal_time_points and the later_time_point. If not specified, it produces a map T between neighboring clonal time points.

smooth_array : list, optional (default: [15,10,5])

List of smooth rounds at initial runs of iteration. Suppose that it has a length N. For iteration n<N, the n-th entry of smooth_array determines the kernel exponent to build the S matrix at the n-th iteration. When n>N, we use the last entry of smooth_array to compute the S matrix. We recommend starting with more smoothing depth and gradually reduce the depth, as inspired by simulated annealing. Data with higher clonal dispersion should start with higher smoothing depth. The final depth should depend on the manifold itself. For fewer cells, it results in a small KNN graph, and a small final depth should be used. We recommend to use a number at the multiple of 5 for computational efficiency i.e., smooth_array=[20, 15, 10, 5], or [20,15,10]

max_iter_N : int, optional (default: 5)

The maximum iterations used to compute the transition map, regardless of epsilon_converge.

epsilon_converge : float, optional (default: 0.05)

The convergence threshold for the change of map correlations between consecutive iterations. This convergence test is activated only when CoSpar has iterated for 3 times.

CoSpar_KNN : int, optional (default: 20)

The number of neighbors for KNN graph used for computing the similarity matrix.

trunca_threshold : list, optional (default: [0.001,0.01])

Threshold to reset entries of a matrix to zero. The first entry is for Similarity matrix; the second entry is for the Tmap. This is only for computational and storage efficiency.

sparsity_threshold : float, optional (default: 0.1)

The relative threshold to remove noises in the updated transition map, in the range [0,1].

intraclone_threshold : float, optional (default: 0.05)

The threshold to remove noises in the demultiplexed (un-smoothed) map, in the range [0,1]

normalization_mode : int, optional (default: 1)

Normalization method. Choice: [0,1]. 0, single-cell normalization; 1, Clone normalization. The clonal normalization suppresses the contribution of large clones, and is much more robust.

extend_Tmap_space : bool optional (default: False)

If true, the initial states for Tmap will include all states at initial time points, and the later states for Tmap will include all states at later time points. Otherwise, the initial and later state space of the Tmap will be restricted to cells with multi-time clonal information alone. The latter case speeds up the computation, which is recommended. This option is ignored when later_time_points is not None.

save_subset : bool, optional (default: True)

If true, save only Smatrix at smooth round [5,10,15,…]; Otherwise, save Smatrix at each round.

use_full_Smatrix : bool, optional (default: True)

If true, extract the relevant Smatrix from the full Smatrix defined by all cells. This tends to be more accurate. The package is optimized around this choice.

Compute_new : bool, optional (default: False)

If True, compute Smatrix from scratch, whether it was computed and saved before or not. This is activated only when use_full_Smatrix=False.

Returns

adata (AnnData object) – Store results at adata.uns[‘transition_map’] and adata.uns[‘intraclone_transition_map’]. This adata is different from the input adata_orig due to subsampling cells.