msmtools.analysis.committor

msmtools.analysis.committor(T, A, B, forward=True, mu=None)

Compute the committor between sets of microstates.

The committor assigns to each microstate a probability that being at this state, the set B will be hit next, rather than set A (forward committor), or that the set A has been hit previously rather than set B (backward committor). See [1] for a detailed mathematical description. The present implementation uses the equations given in [2].

Parameters:
  • T ((M, M) ndarray or scipy.sparse matrix) – Transition matrix
  • A (array_like) – List of integer state labels for set A
  • B (array_like) – List of integer state labels for set B
  • forward (bool) – If True compute the forward committor, else compute the backward committor.
Returns:

q – Vector of comittor probabilities.

Return type:

(M,) ndarray

Notes

Committor functions are used to characterize microstates in terms of their probability to being visited during a reaction/transition between two disjoint regions of state space A, B.

Forward committor

The forward committor \(q^{(+)}_i\) is defined as the probability that the process starting in i will reach B first, rather than A.

Using the first hitting time of a set \(S\),

\[T_{S}=\inf\{t \geq 0 | X_t \in S \}\]

the forward committor \(q^{(+)}_i\) can be fromally defined as

\[q^{(+)}_i=\mathbb{P}_{i}(T_{A}<T_{B}).\]

The forward committor solves to the following boundary value problem

\[\begin{split}\begin{array}{rl} \sum_j L_{ij} q^{(+)}_{j}=0 & i \in X \setminus{(A \cup B)} \\ q_{i}^{(+)}=0 & i \in A \\ q_{i}^{(+)}=1 & i \in B \end{array}\end{split}\]

\(L=T-I\) denotes the generator matrix.

Backward committor

The backward committor is defined as the probability that the process starting in \(x\) came from \(A\) rather than from \(B\).

Using the last exit time of a set \(S\),

\[t_{S}=\sup\{t \geq 0 | X_t \notin S \}\]

the backward committor can be formally defined as

\[q^{(-)}_i=\mathbb{P}_{i}(t_{A}<t_{B}).\]

The backward comittor solves another boundary value problem

\[\begin{split}\begin{array}{rl} \sum_j K_{ij} q^{(-)}_{j}=0 & i \in X \setminus{(A \cup B)} \\ q_{i}^{(-)}=1 & i \in A \\ q_{i}^{(-)}=0 & i \in B \end{array}\end{split}\]

\(K=(D_{\pi}L)^{T}\) denotes the adjoint generator matrix.

References

[1]P. Metzner, C. Schuette and E. Vanden-Eijnden. Transition Path Theory for Markov Jump Processes. Multiscale Model Simul 7: 1192-1219 (2009).
[2]F. Noe, C. Schuette, E. Vanden-Eijnden, L. Reich and T.Weikl Constructing the Full Ensemble of Folding Pathways from Short Off-Equilibrium Simulations. Proc. Natl. Acad. Sci. USA, 106: 19011-19016 (2009).

Examples

>>> import numpy as np
>>> from msmtools.analysis import committor
>>> T = np.array([[0.89, 0.1, 0.01], [0.5, 0.0, 0.5], [0.0, 0.1, 0.9]])
>>> A = [0]
>>> B = [2]
>>> u_plus = committor(T, A, B)
>>> u_plus
array([ 0. ,  0.5,  1. ])
>>> u_minus = committor(T, A, B, forward=False)
>>> u_minus
array([ 1.        ,  0.45454545,  0.        ])