pygpcca.GPCCA¶
-
class
pygpcca.
GPCCA
(P, eta=None, z='LM', method='brandts')[source]¶ G-PCCA [Reuter18] spectral clustering method with optimized memberships.
Clusters the dominant m Schur vectors of a transition matrix.
This algorithm generates a fuzzy clustering such that the resulting membership functions are as crisp (characteristic) as possible.
- Parameters
P (
Union
[ndarray
,spmatrix
]) – Transition matrix (row-stochastic).Input probability distribution of the (micro)states. In theory eta can be an arbitrary distribution as long as it is a valid probability distribution (i.e., sums up to 1). A neutral and valid choice would be the uniform distribution (default).
In case of a reversible transition matrix, the stationary distribution \(\pi\) can (but don’t has to) be used here. In case of a non-reversible P, some initial or average distribution of the states might be chosen instead of the uniform distribution. If None, uniform distribution is used.
z (
str
) –Specifies which portion of the eigenvalue spectrum of P is to be sought. The returned invariant subspace of P will be associated with this part of the spectrum. Valid options are:
’LM’: largest magnitude (default).
’LR’: largest real part.
method (
str
) –Which method to use to determine the invariant subspace. Valid options are:
’brandts’: perform a full Schur decomposition of P utilizing scipy.schur (without the intrinsic sorting option, since it is flawed) and sort the returned Schur form R and Schur vector matrix Q afterwards using a routine published by Brandts [Brandts02]. This is well tested and thus the default method, although it is also the slowest choice.
’krylov’: calculate an orthonormal basis of the subspace associated with the m dominant eigenvalues of P using the Krylov-Schur method as implemented in SLEPc. This is the fastest choice and especially suitable for very large P, but it is still experimental.
See the installation instructions for more information.
References
If you use this code or parts of it, please cite [Reuter18].
Methods
minChi
(m_min, m_max)Calculate the minChi indicator (see [Reuter18]) for every \(m \in [m_{min},m_{max}]\).
optimize
(m)Full G-PCCA [Reuter18] spectral clustering method with optimized memberships.
Attributes
Crispness \(\xi \in [0,1]\) quantifying the optimality of the solution (higher is better).
Coarse grained input distribution of shape (m,).
Coarse grained stationary distribution of shape (m,).
Coarse grained transition matrix of shape (m, m).
Input probability distribution of the (micro)states.
Array of shape (n, m) containing the membership (or probability) of each state (to be assigned) to each cluster.
Crisp clustering using G-PCCA.
Crisp clustering using G-PCCA.
Number of clusters (macrostates) to group the n microstates into.
Optimized rotation matrix.
Ordered top left part of shape (m, m) of the real Schur matrix of P.
Array of shape (n, m) with m sorted Schur vectors in the columns.
Stationary probability distribution of the (micro)states.
Array of shape (m,) containing the m dominant eigenvalues of P.
Transition matrix (row-stochastic).