pyemma.msm.estimate_markov_model

pyemma.msm.estimate_markov_model(dtrajs, lag, reversible=True, statdist=None, count_mode='sliding', sparse=False, connectivity='largest', dt_traj='1 step', maxiter=1000000, maxerr=1e-08)

Estimates a Markov model from discrete trajectories

Returns a MaximumLikelihoodMSM that contains the estimated transition matrix and allows to compute a large number of quantities related to Markov models.

Parameters:
  • dtrajs (list containing ndarrays(dtype=int) or ndarray(n, dtype=int)) – discrete trajectories, stored as integer ndarrays (arbitrary size) or a single ndarray for only one trajectory.
  • lag (int) – lag time at which transitions are counted and the transition matrix is estimated.
  • reversible (bool, optional) – If true compute reversible MSM, else non-reversible MSM
  • statdist ((M,) ndarray, optional) – Stationary vector on the full state-space. Transition matrix will be estimated such that statdist is its equilibrium distribution.
  • count_mode (str, optional, default='sliding') –

    mode to obtain count matrices from discrete trajectories. Should be one of:

    • ‘sliding’ : A trajectory of length T will have \(T-tau\) counts at time indexes
      \[(0 \rightarrow \tau), (1 \rightarrow \tau+1), ..., (T-\tau-1 \rightarrow T-1)\]
    • ‘effective’ : Uses an estimate of the transition counts that are statistically uncorrelated. Recommended when used with a Bayesian MSM.
    • ‘sample’ : A trajectory of length T will have \(T/tau\) counts at time indexes
      \[(0 \rightarrow \tau), (\tau \rightarrow 2 \tau), ..., (((T/\tau)-1) \tau \rightarrow T)\]
  • sparse (bool, optional) – If true compute count matrix, transition matrix and all derived quantities using sparse matrix algebra. In this case python sparse matrices will be returned by the corresponding functions instead of numpy arrays. This behavior is suggested for very large numbers of states (e.g. > 4000) because it is likely to be much more efficient.
  • connectivity (str, optional) –

    Connectivity mode. Three methods are intended (currently only ‘largest’ is implemented)

    • ‘largest’ : The active set is the largest reversibly connected set. All estimation will be done on this subset and all quantities (transition matrix, stationary distribution, etc) are only defined on this subset and are correspondingly smaller than the full set of states
    • ‘all’ : The active set is the full set of states. Estimation will be conducted on each reversibly connected set separately. That means the transition matrix will decompose into disconnected submatrices, the stationary vector is only defined within subsets, etc. Currently not implemented.
    • ‘none’ : The active set is the full set of states. Estimation will be conducted on the full set of states without ensuring connectivity. This only permits nonreversible estimation. Currently not implemented.
  • estimate (bool, optional) – If true estimate the MSM when creating the MSM object.
  • dt (str, optional) –

    Description of the physical time corresponding to the lag. May be used by analysis algorithms such as plotting tools to pretty-print the axes. By default ‘1 step’, i.e. there is no physical time unit. Specify by a number, whitespace and unit. Permitted units are (* is an arbitrary string):

    ‘fs’, ‘femtosecond*’
    ‘ps’, ‘picosecond*’
    ‘ns’, ‘nanosecond*’
    ‘us’, ‘microsecond*’
    ‘ms’, ‘millisecond*’
    ‘s’, ‘second*’
  • maxiter (int, optional) – Optional parameter with reversible = True. maximum number of iterations before the transition matrix estimation method exits
  • maxerr (float, optional) – Optional parameter with reversible = True. convergence tolerance for transition matrix estimation. This specifies the maximum change of the Euclidean norm of relative stationary probabilities (\(x_i = \sum_k x_{ik}\)). The relative stationary probability changes \(e_i = (x_i^{(1)} - x_i^{(2)})/(x_i^{(1)} + x_i^{(2)})\) are used in order to track changes in small probabilities. The Euclidean norm of the change vector, \(|e_i|_2\), is compared to maxerr.
Returns:

msm – Estimator object containing the MSM and estimation information.

Return type:

MaximumLikelihoodMSM

See also

MaximumLikelihoodMSM()
An MSM object that has been estimated from data
class pyemma.msm.estimators.maximum_likelihood_msm.MaximumLikelihoodMSM(lag=1, reversible=True, statdist_constraint=None, count_mode='sliding', sparse=False, connectivity='largest', dt_traj='1 step', maxiter=1000000, maxerr=1e-08)

Maximum likelihood estimator for MSMs given discrete trajectory statistics

Methods

cktest(nsets[, memberships, mlags, conf, ...]) Conducts a Chapman-Kolmogorow test.
coarse_grain(ncoarse[, method]) Returns a coarse-grained Markov model.
committor_backward(A, B) Backward committor from set A to set B
committor_forward(A, B) Forward committor (also known as p_fold or splitting probability) from set A to set B
correlation(a[, b, maxtime, k, ncv]) Time-correlation for equilibrium experiment.
eigenvalues([k]) Compute the transition matrix eigenvalues
eigenvectors_left([k]) Compute the left transition matrix eigenvectors
eigenvectors_right([k]) Compute the right transition matrix eigenvectors
estimate(dtrajs, **parms)
param dtrajs:discrete trajectories, stored as integer ndarrays (arbitrary size)
expectation(a) Equilibrium expectation value of a given observable.
fingerprint_correlation(a[, b, k, ncv]) Dynamical fingerprint for equilibrium time-correlation experiment.
fingerprint_relaxation(p0, a[, k, ncv]) Dynamical fingerprint for perturbation/relaxation experiment.
fit(X) Estimates parameters - for compatibility with sklearn.
generate_traj(N[, start, stop, stride]) Generates a synthetic discrete trajectory of length N and simulation time stride * lag time * N
get_model_params([deep]) Get parameters for this model.
get_params([deep]) Get parameters for this estimator.
hmm(nhidden) Estimates a hidden Markov state model as described in [1]_
mfpt(A, B) Mean first passage times from set A to set B, in units of the input trajectory time step
pcca(m) Runs PCCA++ [1]_ to compute a metastable decomposition of MSM states
propagate(p0, k) Propagates the initial distribution p0 k times
relaxation(p0, a[, maxtime, k, ncv]) Simulates a perturbation-relaxation experiment.
sample_by_distributions(distributions, nsample) Generates samples according to given probability distributions
sample_by_state(nsample[, subset, replace]) Generates samples of the connected states.
set_model_params(P[, pi, reversible, ...]) Call to set all basic model parameters.
set_params(**params) Set the parameters of this estimator.
simulate(N[, start, stop, dt]) Generates a realization of the Markov Model
timescales([k]) The relaxation timescales corresponding to the eigenvalues
trajectory_weights() Uses the MSM to assign a probability weight to each trajectory frame.
update_model_params(**params) Update given model parameter if they are set to specific values

Attributes

P The transition matrix on the active set.
active_count_fraction The fraction of counts in the largest connected set.
active_set The active set of states on which all computations and estimations will be done
active_state_fraction The fraction of states in the largest connected set.
active_state_indexes Ensures that the connected states are indexed and returns the indices
connected_sets The reversible connected sets of states, sorted by size (descending)
connectivity Returns the connectivity mode of the MSM
count_matrix_active The count matrix on the active set given the connectivity mode used.
count_matrix_full The count matrix on full set of discrete states, irrespective as to whether they are connected or not.
discrete_trajectories_active A list of integer arrays with the discrete trajectories mapped to the connectivity mode used.
discrete_trajectories_full A list of integer arrays with the original (unmapped) discrete trajectories
dt_model
dtrajs_active A list of integer arrays with the discrete trajectories mapped to the connectivity mode used.
dtrajs_full A list of integer arrays with the original (unmapped) discrete trajectories
effective_count_matrix Statistically uncorrelated transition counts within the active set of states
is_reversible Returns whether the MSM is reversible
is_sparse Returns whether the MSM is sparse
lagtime The lag time at which the Markov model was estimated
largest_connected_set The largest reversible connected set of states
logger The logger for this class instance
metastable_assignments Assignment of states to metastable sets using PCCA++
metastable_distributions Probability of metastable states to visit an MSM state by PCCA++
metastable_memberships Probabilities of MSM states to belong to a metastable state by PCCA++
metastable_sets Metastable sets using PCCA++
model The model estimated by this Estimator
name The name of this instance
neig number of eigenvalues to compute.
nstates Number of active states on which all computations and estimations are done
nstates_full Number of states in discrete trajectories
pi
reversible Returns whether the MSM is reversible
sparse Returns whether the MSM is sparse
stationary_distribution
timestep_model Physical time corresponding to one transition matrix step, e.g.
transition_matrix The transition matrix on the active set.
P

The transition matrix on the active set.

active_count_fraction

The fraction of counts in the largest connected set.

active_set

The active set of states on which all computations and estimations will be done

active_state_fraction

The fraction of states in the largest connected set.

active_state_indexes

Ensures that the connected states are indexed and returns the indices

cktest(nsets, memberships=None, mlags=10, conf=0.95, err_est=False, n_jobs=1, show_progress=True)

Conducts a Chapman-Kolmogorow test.

Parameters:
  • nsets (int) – number of sets to test on
  • memberships (ndarray(nstates, nsets), optional) – optional state memberships. By default (None) will conduct a cktest on PCCA (metastable) sets.
  • mlags (int or int-array, optional) – multiples of lag times for testing the Model, e.g. range(10). A single int will trigger a range, i.e. mlags=10 maps to mlags=range(10). The setting None will choose mlags automatically according to the longest available trajectory
  • conf (float, optional) – confidence interval
  • err_est (bool, optional) – compute errors also for all estimations (computationally expensive) If False, only the prediction will get error bars, which is often sufficient to validate a model.
  • n_jobs (int, default=1) – how many jobs to use during calculation
  • show_progress (bool, optional) – Show progress bars for calculation?
Returns:

cktest

Return type:

ChapmanKolmogorovValidator

References

This test was suggested in [1]_ and described in detail in [2]_.

[1]F. Noe, Ch. 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)
[2]Prinz, J H, H Wu, M Sarich, B Keller, M Senne, M Held, J D Chodera, C Schuette and F Noe. 2011. Markov models of molecular kinetics: Generation and validation. J Chem Phys 134: 174105
coarse_grain(ncoarse, method='hmm')

Returns a coarse-grained Markov model.

Currently only the HMM method described in [1]_ is available for coarse-graining MSMs.

Parameters:ncoarse (int) – number of coarse states
Returns:hmsm
Return type:MaximumLikelihoodHMSM

References

[1]F. Noe, H. Wu, J.-H. Prinz and N. Plattner: Projected and hidden Markov models for calculating kinetics and metastable states of complex molecules J. Chem. Phys. 139, 184114 (2013)
committor_backward(A, B)

Backward committor from set A to set B

Parameters:
  • A (int or int array) – set of starting states
  • B (int or int array) – set of target states
committor_forward(A, B)

Forward committor (also known as p_fold or splitting probability) from set A to set B

Parameters:
  • A (int or int array) – set of starting states
  • B (int or int array) – set of target states
connected_sets

The reversible connected sets of states, sorted by size (descending)

connectivity

Returns the connectivity mode of the MSM

correlation(a, b=None, maxtime=None, k=None, ncv=None)

Time-correlation for equilibrium experiment.

In order to simulate a time-correlation experiment (e.g. fluorescence correlation spectroscopy [NDD11], dynamical neutron scattering [LYP13], ...), first compute the mean values of your experimental observable \(a\) by MSM state:

\[a_i = \frac{1}{N_i} \sum_{x_t \in S_i} f(x_t)\]

where \(S_i\) is the set of configurations belonging to MSM state \(i\) and \(f()\) is a function that computes the experimental observable of interest for configuration \(x_t\). If a cross-correlation function is wanted, also apply the above computation to a second experimental observable \(b\).

Then the accurate (i.e. without statistical error) autocorrelation function of \(f(x_t)\) given the Markov model is computed by correlation(a), and the accurate cross-correlation function is computed by correlation(a,b). This is done by evaluating the equation

\[\begin{split}acf_a(k\tau) &= \mathbf{a}^\top \mathrm{diag}(\boldsymbol{\pi}) \mathbf{P(\tau)}^k \mathbf{a} \\ ccf_{a,b}(k\tau) &= \mathbf{a}^\top \mathrm{diag}(\boldsymbol{\pi}) \mathbf{P(\tau)}^k \mathbf{b}\end{split}\]

where \(acf\) stands for autocorrelation function and \(ccf\) stands for cross-correlation function, \(\mathbf{P(\tau)}\) is the transition matrix at lag time \(\tau\), \(\boldsymbol{\pi}\) is the equilibrium distribution of \(\mathbf{P}\), and \(k\) is the time index.

Note that instead of using this method you could generate a long synthetic trajectory from the MSM and then estimating the time-correlation of your observable(s) directly from this trajectory. However, there is no reason to do this because the present method does that calculation without any sampling, and only in the limit of an infinitely long synthetic trajectory the two results will agree exactly. The correlation function computed by the present method still has statistical uncertainty from the fact that the underlying MSM transition matrix has statistical uncertainty when being estimated from data, but there is no additional (and unnecessary) uncertainty due to synthetic trajectory generation.

Parameters:
  • a ((n,) ndarray) – Observable, represented as vector on state space
  • maxtime (int or float) – Maximum time (in units of the input trajectory time step) until which the correlation function will be evaluated. Internally, the correlation function can only be computed in integer multiples of the Markov model lag time, and therefore the actual last time point will be computed at \(\mathrm{ceil}(\mathrm{maxtime} / \tau)\) By default (None), the maxtime will be set equal to the 5 times the slowest relaxation time of the MSM, because after this time the signal is almost constant.
  • b ((n,) ndarray (optional)) – Second observable, for cross-correlations
  • k (int (optional)) – Number of eigenvalues and eigenvectors to use for computation. This option is only relevant for sparse matrices and long times for which an eigenvalue decomposition will be done instead of using the matrix power.
  • ncv (int (optional)) – Only relevant for sparse matrices and large lag times where the relaxation will be computed using an eigenvalue decomposition. The number of Lanczos vectors generated, ncv must be greater than k; it is recommended that ncv > 2*k.
Returns:

  • times (ndarray (N)) – Time points (in units of the input trajectory time step) at which the correlation has been computed
  • correlations (ndarray (N)) – Correlation values at given times

Examples

This example computes the autocorrelation function of a simple observable on a three-state Markov model and plots the result using matplotlib:

>>> import numpy as np
>>> import pyemma.msm as msm
>>>
>>> P = np.array([[0.99, 0.01, 0], [0.01, 0.9, 0.09], [0, 0.1, 0.9]])
>>> a = np.array([0.0, 0.5, 1.0])
>>> M = msm.markov_model(P)
>>> times, acf = M.correlation(a)
>>>
>>> import matplotlib.pylab as plt
>>> plt.plot(times, acf)  

References

[NDD11]Noe, F., S. Doose, I. Daidone, M. Loellmann, J. D. Chodera, M. Sauer and J. C. Smith. 2011 Dynamical fingerprints for probing individual relaxation processes in biomolecular dynamics with simulations and kinetic experiments. Proc. Natl. Acad. Sci. USA 108, 4822-4827.
[LYP13]Lindner, B., Z. Yi, J.-H. Prinz, J. C. Smith and F. Noe. 2013. Dynamic Neutron Scattering from Conformational Dynamics I: Theory and Markov models. J. Chem. Phys. 139, 175101.
count_matrix_active

The count matrix on the active set given the connectivity mode used.

For example, for connectivity=’largest’, the count matrix is given only on the largest reversibly connected set. Attention: This count matrix has been obtained by sliding a window of length tau across the data. It contains a factor of tau more counts than are statistically uncorrelated. It’s fine to use this matrix for maximum likelihood estimated, but it will give far too small errors if you use it for uncertainty calculations. In order to do uncertainty calculations, use the effective count matrix, see: effective_count_matrix

See also

effective_count_matrix
For a count matrix with effective (statistically uncorrelated) counts.
count_matrix_full

The count matrix on full set of discrete states, irrespective as to whether they are connected or not. Attention: This count matrix has been obtained by sliding a window of length tau across the data. It contains a factor of tau more counts than are statistically uncorrelated. It’s fine to use this matrix for maximum likelihood estimated, but it will give far too small errors if you use it for uncertainty calculations. In order to do uncertainty calculations, use the effective count matrix, see: effective_count_matrix (only implemented on the active set), or divide this count matrix by tau.

See also

effective_count_matrix
For a active-set count matrix with effective (statistically uncorrelated) counts.
discrete_trajectories_active

A list of integer arrays with the discrete trajectories mapped to the connectivity mode used. For example, for connectivity=’largest’, the indexes will be given within the connected set. Frames that are not in the connected set will be -1.

discrete_trajectories_full

A list of integer arrays with the original (unmapped) discrete trajectories

dt_model
dtrajs_active

A list of integer arrays with the discrete trajectories mapped to the connectivity mode used. For example, for connectivity=’largest’, the indexes will be given within the connected set. Frames that are not in the connected set will be -1.

dtrajs_full

A list of integer arrays with the original (unmapped) discrete trajectories

effective_count_matrix

Statistically uncorrelated transition counts within the active set of states

You can use this count matrix for Bayesian estimation or error perturbation.

References

[1] Noe, F. (2015) Statistical inefficiency of Markov model count matrices
http://publications.mi.fu-berlin.de/1699/1/autocorrelation_counts.pdf
eigenvalues(k=None)

Compute the transition matrix eigenvalues

Parameters:k (int) – number of eigenvalues to be returned. By default will return all available eigenvalues
Returns:ts – transition matrix eigenvalues \(\lambda_i, i = 1, ..., k\)., sorted by descending norm.
Return type:ndarray(k,)
eigenvectors_left(k=None)

Compute the left transition matrix eigenvectors

Parameters:k (int) – number of eigenvectors to be returned. By default all available eigenvectors.
Returns:L – left eigenvectors in a row matrix. l_ij is the j’th component of the i’th left eigenvector
Return type:ndarray(k,n)
eigenvectors_right(k=None)

Compute the right transition matrix eigenvectors

Parameters:k (int) – number of eigenvectors to be computed. By default all available eigenvectors.
Returns:R – right eigenvectors in a column matrix. r_ij is the i’th component of the j’th right eigenvector
Return type:ndarray(n,k)
estimate(dtrajs, **parms)
Parameters:
  • dtrajs (list containing ndarrays(dtype=int) or ndarray(n, dtype=int) or pyemma.msm.util.dtraj_states.DiscreteTrajectoryStats) – discrete trajectories, stored as integer ndarrays (arbitrary size) or a single ndarray for only one trajectory.
  • **params – Other keyword parameters if different from the settings when this estimator was constructed
Returns:

MSM

Return type:

pyemma.msm.MaximumlikelihoodMSM

expectation(a)

Equilibrium expectation value of a given observable.

Parameters:a ((n,) ndarray) – Observable vector on the MSM state space
Returns:val – Equilibrium expectation value fo the given observable
Return type:float

Notes

The equilibrium expectation value of an observable \(a\) is defined as follows

\[\mathbb{E}_{\mu}[a] = \sum_i \pi_i a_i\]

\(\pi=(\pi_i)\) is the stationary vector of the transition matrix \(P\).

fingerprint_correlation(a, b=None, k=None, ncv=None)

Dynamical fingerprint for equilibrium time-correlation experiment.

Parameters:
  • a ((n,) ndarray) – Observable, represented as vector on MSM state space
  • b ((n,) ndarray, optional) – Second observable, for cross-correlations
  • k (int, optional) – Number of eigenvalues and eigenvectors to use for computation. This option is only relevant for sparse matrices and long times for which an eigenvalue decomposition will be done instead of using the matrix power
  • ncv (int, optional) – Only relevant for sparse matrices and large lag times, where the relaxation will be computed using an eigenvalue decomposition. The number of Lanczos vectors generated, ncv must be greater than k; it is recommended that ncv > 2*k
Returns:

  • timescales ((k,) ndarray) – Time-scales (in units of the input trajectory time step) of the transition matrix
  • amplitudes ((k,) ndarray) – Amplitudes for the correlation experiment

References

Spectral densities are commonly used in spectroscopy. Dynamical fingerprints are a useful representation for computational spectroscopy results and have been introduced in [NDD11].

[NDD11]Noe, F, S Doose, I Daidone, M Loellmann, M Sauer, J D Chodera and J Smith. 2010. Dynamical fingerprints for probing individual relaxation processes in biomolecular dynamics with simulations and kinetic experiments. PNAS 108 (12): 4822-4827.
fingerprint_relaxation(p0, a, k=None, ncv=None)

Dynamical fingerprint for perturbation/relaxation experiment.

Parameters:
  • p0 ((n,) ndarray) – Initial distribution for a relaxation experiment
  • a ((n,) ndarray) – Observable, represented as vector on state space
  • k (int, optional) – Number of eigenvalues and eigenvectors to use for computation
  • ncv (int, optional) – Only relevant for sparse matrices and large lag times, where the relaxation will be computes using an eigenvalue decomposition. The number of Lanczos vectors generated, ncv must be greater than k; it is recommended that ncv > 2*k
Returns:

  • timescales ((k,) ndarray) – Time-scales (in units of the input trajectory time step) of the transition matrix
  • amplitudes ((k,) ndarray) – Amplitudes for the relaxation experiment

References

Spectral densities are commonly used in spectroscopy. Dynamical fingerprints are a useful representation for computational spectroscopy results and have been introduced in [NDD11].

[NDD11]Noe, F, S Doose, I Daidone, M Loellmann, M Sauer, J D Chodera and J Smith. 2010. Dynamical fingerprints for probing individual relaxation processes in biomolecular dynamics with simulations and kinetic experiments. PNAS 108 (12): 4822-4827.
fit(X)

Estimates parameters - for compatibility with sklearn.

Parameters:X (object) – A reference to the data from which the model will be estimated
Returns:estimator – The estimator (self) with estimated model.
Return type:object
generate_traj(N, start=None, stop=None, stride=1)

Generates a synthetic discrete trajectory of length N and simulation time stride * lag time * N

This information can be used in order to generate a synthetic molecular dynamics trajectory - see pyemma.coordinates.save_traj()

Note that the time different between two samples is the Markov model lag time tau. When comparing quantities computing from this synthetic trajectory and from the input trajectories, the time points of this trajectory must be scaled by the lag time in order to have them on the same time scale.

Parameters:
  • N (int) – Number of time steps in the output trajectory. The total simulation time is stride * lag time * N
  • start (int, optional, default = None) – starting state. If not given, will sample from the stationary distribution of P
  • stop (int or int-array-like, optional, default = None) – stopping set. If given, the trajectory will be stopped before N steps once a state of the stop set is reached
  • stride (int, optional, default = 1) – Multiple of lag time used as a time step. By default, the time step is equal to the lag time
Returns:

indexes – trajectory and time indexes of the simulated trajectory. Each row consist of a tuple (i, t), where i is the index of the trajectory and t is the time index within the trajectory. Note that the time different between two samples is the Markov model lag time tau

Return type:

ndarray( (N, 2) )

See also

pyemma.coordinates.save_traj()
in order to save this synthetic trajectory as a trajectory file with molecular structures
get_model_params(deep=True)

Get parameters for this model.

Parameters:deep (boolean, optional) – If True, will return the parameters for this estimator and contained subobjects that are estimators.
Returns:params – Parameter names mapped to their values.
Return type:mapping of string to any
get_params(deep=True)

Get parameters for this estimator.

Parameters:deep (boolean, optional) – If True, will return the parameters for this estimator and contained subobjects that are estimators.
Returns:params – Parameter names mapped to their values.
Return type:mapping of string to any
hmm(nhidden)

Estimates a hidden Markov state model as described in [1]_

Parameters:nhidden (int) – number of hidden (metastable) states
Returns:hmsm
Return type:MaximumLikelihoodHMSM

References

[1]F. Noe, H. Wu, J.-H. Prinz and N. Plattner: Projected and hidden Markov models for calculating kinetics and metastable states of complex molecules J. Chem. Phys. 139, 184114 (2013)
is_reversible

Returns whether the MSM is reversible

is_sparse

Returns whether the MSM is sparse

lagtime

The lag time at which the Markov model was estimated

largest_connected_set

The largest reversible connected set of states

logger

The logger for this class instance

metastable_assignments

Assignment of states to metastable sets using PCCA++

Computes the assignment to metastable sets for active set states using the PCCA++ method [1]_. pcca() needs to be called first to make this attribute available.

This is only recommended for visualization purposes. You cannot compute any actual quantity of the coarse-grained kinetics without employing the fuzzy memberships!

Returns:assignments – For each MSM state, the metastable state it is located in.
Return type:ndarray (n,)

See also

pcca
to compute the metastable decomposition

References

[1]Roeblitz, S and M Weber. 2013. Fuzzy spectral clustering by PCCA+: application to Markov state models and data classification. Advances in Data Analysis and Classification 7 (2): 147-179
metastable_distributions

Probability of metastable states to visit an MSM state by PCCA++

Computes the probability distributions of active set states within each metastable set by combining the the PCCA++ method [1]_ with Bayesian inversion as described in [2]_.

pcca() needs to be called first to make this attribute available.

Returns:p_out – A matrix containing the probability distribution of each active set state, given that we are in one of the m metastable sets, i.e. p(state | metastable). The row sums of p_out are 1.
Return type:ndarray (m,n)

See also

pcca
to compute the metastable decomposition

References

[1]Roeblitz, S and M Weber. 2013. Fuzzy spectral clustering by PCCA+: application to Markov state models and data classification. Advances in Data Analysis and Classification 7, 147-179.
[2]F. Noe, H. Wu, J.-H. Prinz and N. Plattner. 2013. Projected and hidden Markov models for calculating kinetics and metastable states of complex molecules J. Chem. Phys. 139, 184114.
metastable_memberships

Probabilities of MSM states to belong to a metastable state by PCCA++

Computes the memberships of active set states to metastable sets with the PCCA++ method [1]_.

pcca() needs to be called first to make this attribute available.

Returns:M – A matrix containing the probability or membership of each state to be assigned to each metastable set, i.e. p(metastable | state). The row sums of M are 1.
Return type:ndarray((n,m))

See also

pcca
to compute the metastable decomposition

References

[1]Roeblitz, S and M Weber. 2013. Fuzzy spectral clustering by PCCA+: application to Markov state models and data classification. Advances in Data Analysis and Classification 7 (2): 147-179
metastable_sets

Metastable sets using PCCA++

Computes the metastable sets of active set states within each metastable set using the PCCA++ method [1]_. pcca() needs to be called first to make this attribute available.

This is only recommended for visualization purposes. You cannot compute any actual quantity of the coarse-grained kinetics without employing the fuzzy memberships!

Returns:sets – A list of length equal to metastable states. Each element is an array with microstate indexes contained in it
Return type:list of ndarray

See also

pcca
to compute the metastable decomposition

References

[1]Roeblitz, S and M Weber. 2013. Fuzzy spectral clustering by PCCA+: application to Markov state models and data classification. Advances in Data Analysis and Classification 7 (2): 147-179
mfpt(A, B)

Mean first passage times from set A to set B, in units of the input trajectory time step

Parameters:
  • A (int or int array) – set of starting states
  • B (int or int array) – set of target states
model

The model estimated by this Estimator

name

The name of this instance

neig

number of eigenvalues to compute.

nstates

Number of active states on which all computations and estimations are done

nstates_full

Number of states in discrete trajectories

pcca(m)

Runs PCCA++ [1]_ to compute a metastable decomposition of MSM states

After calling this method you can access metastable_memberships(), metastable_distributions(), metastable_sets() and metastable_assignments().

Parameters:m (int) – Number of metastable sets
Returns:pcca_obj – An object containing all PCCA quantities. However, you can also ignore this return value and instead retrieve the quantities of your interest with the following MSM functions: metastable_memberships(), metastable_distributions(), metastable_sets() and metastable_assignments().
Return type:PCCA

Notes

If you coarse grain with PCCA++, the order of the obtained memberships might not be preserved. This also applies for metastable_memberships(), metastable_distributions(), metastable_sets(), metastable_assignments()

References

[1]Roeblitz, S and M Weber. 2013. Fuzzy spectral clustering by PCCA+: application to Markov state models and data classification. Advances in Data Analysis and Classification 7 (2): 147-179
pi
propagate(p0, k)

Propagates the initial distribution p0 k times

Computes the product

\[p_k = p_0^T P^k\]

If the lag time of transition matrix \(P\) is \(\tau\), this will provide the probability distribution at time \(k \tau\).

Parameters:
  • p0 (ndarray(n,)) – Initial distribution. Vector of size of the active set.
  • k (int) – Number of time steps
Returns:

pk – Distribution after k steps. Vector of size of the active set.

Return type:

ndarray(n,)

relaxation(p0, a, maxtime=None, k=None, ncv=None)

Simulates a perturbation-relaxation experiment.

In perturbation-relaxation experiments such as temperature-jump, pH-jump, pressure jump or rapid mixing experiments, an ensemble of molecules is initially prepared in an off-equilibrium distribution and the expectation value of some experimental observable is then followed over time as the ensemble relaxes towards equilibrium.

In order to simulate such an experiment, first determine the distribution of states at which the experiment is started, \(p_0\) and compute the mean values of your experimental observable \(a\) by MSM state:

\[a_i = \frac{1}{N_i} \sum_{x_t \in S_i} f(x_t)\]

where \(S_i\) is the set of configurations belonging to MSM state \(i\) and \(f()\) is a function that computes the experimental observable of interest for configuration \(x_t\).

Then the accurate (i.e. without statistical error) time-dependent expectation value of \(f(x_t)\) given the Markov model is computed by relaxation(p0, a). This is done by evaluating the equation

\[E_a(k\tau) = \mathbf{p_0}^{\top} \mathbf{P(\tau)}^k \mathbf{a}\]

where \(E\) stands for the expectation value that relaxes to its equilibrium value that is identical to expectation(a), \(\mathbf{P(\tau)}\) is the transition matrix at lag time \(\tau\), \(\boldsymbol{\pi}\) is the equilibrium distribution of \(\mathbf{P}\), and \(k\) is the time index.

Note that instead of using this method you could generate many synthetic trajectory from the MSM with starting points drawn from the initial distribution and then estimating the time-dependent expectation value by an ensemble average. However, there is no reason to do this because the present method does that calculation without any sampling, and only in the limit of an infinitely many trajectories the two results will agree exactly. The relaxation function computed by the present method still has statistical uncertainty from the fact that the underlying MSM transition matrix has statistical uncertainty when being estimated from data, but there is no additional (and unnecessary) uncertainty due to synthetic trajectory generation.

Parameters:
  • p0 ((n,) ndarray) – Initial distribution for a relaxation experiment
  • a ((n,) ndarray) – Observable, represented as vector on state space
  • maxtime (int or float, optional) – Maximum time (in units of the input trajectory time step) until which the correlation function will be evaluated. Internally, the correlation function can only be computed in integer multiples of the Markov model lag time, and therefore the actual last time point will be computed at \(\mathrm{ceil}(\mathrm{maxtime} / \tau)\). By default (None), the maxtime will be set equal to the 5 times the slowest relaxation time of the MSM, because after this time the signal is constant.
  • k (int (optional)) – Number of eigenvalues and eigenvectors to use for computation
  • ncv (int (optional)) – Only relevant for sparse matrices and large lag times, where the relaxation will be computes using an eigenvalue decomposition. The number of Lanczos vectors generated, ncv must be greater than k; it is recommended that ncv > 2*k
Returns:

  • times (ndarray (N)) – Time points (in units of the input trajectory time step) at which the relaxation has been computed
  • res (ndarray) – Array of expectation value at given times

reversible

Returns whether the MSM is reversible

sample_by_distributions(distributions, nsample)

Generates samples according to given probability distributions

Parameters:
  • distributions (list or array of ndarray ( (n) )) – m distributions over states. Each distribution must be of length n and must sum up to 1.0
  • nsample (int) – Number of samples per distribution. If replace = False, the number of returned samples per state could be smaller if less than nsample indexes are available for a state.
Returns:

indexes – List of the sampled indices by distribution. Each element is an index array with a number of rows equal to nsample, with rows consisting of a tuple (i, t), where i is the index of the trajectory and t is the time index within the trajectory.

Return type:

length m list of ndarray( (nsample, 2) )

sample_by_state(nsample, subset=None, replace=True)

Generates samples of the connected states.

For each state in the active set of states, generates nsample samples with trajectory/time indexes. This information can be used in order to generate a trajectory of length nsample * nconnected using pyemma.coordinates.save_traj() or nconnected trajectories of length nsample each using pyemma.coordinates.save_traj()

Parameters:
  • nsample (int) – Number of samples per state. If replace = False, the number of returned samples per state could be smaller if less than nsample indexes are available for a state.
  • subset (ndarray((n)), optional, default = None) – array of states to be indexed. By default all states in the connected set will be used
  • replace (boolean, optional) – Whether the sample is with or without replacement
Returns:

indexes – list of trajectory/time index arrays with an array for each state. Within each index array, each row consist of a tuple (i, t), where i is the index of the trajectory and t is the time index within the trajectory.

Return type:

list of ndarray( (N, 2) )

See also

pyemma.coordinates.save_traj()
in order to save the sampled frames sequentially in a trajectory file with molecular structures
pyemma.coordinates.save_trajs()
in order to save the sampled frames in nconnected trajectory files with molecular structures
set_model_params(P, pi=None, reversible=None, dt_model='1 step', neig=None)

Call to set all basic model parameters.

Sets or updates given model parameters. This argument list of this method must contain the full list of essential, or independent model parameters. It can additionally contain derived parameters, e.g. in order to save computational costs of re-computing them.

Parameters:
  • P (ndarray(n,n)) – transition matrix
  • pi (ndarray(n), optional, default=None) – stationary distribution. Can be optionally given in case if it was already computed, e.g. by the estimator.
  • reversible (bool, optional, default=None) – whether P is reversible with respect to its stationary distribution. If None (default), will be determined from P
  • dt_model (str, optional, default='1 step') –

    Description of the physical time corresponding to the model time step. May be used by analysis algorithms such as plotting tools to pretty-print the axes. By default ‘1 step’, i.e. there is no physical time unit. Specify by a number, whitespace and unit. Permitted units are (* is an arbitrary string):

    ‘fs’, ‘femtosecond*’
    ‘ps’, ‘picosecond*’
    ‘ns’, ‘nanosecond*’
    ‘us’, ‘microsecond*’
    ‘ms’, ‘millisecond*’
    ‘s’, ‘second*’
  • neig (int or None) – The number of eigenvalues / eigenvectors to be kept. If set to None, defaults will be used. For a dense MSM the default is all eigenvalues. For a sparse MSM the default is 10.

Notes

Explicitly define all independent model parameters in the argument list of this function (by mandatory or keyword arguments)

set_params(**params)

Set the parameters of this estimator. The method works on simple estimators as well as on nested objects (such as pipelines). The former have parameters of the form <component>__<parameter> so that it’s possible to update each component of a nested object. :returns: :rtype: self

simulate(N, start=None, stop=None, dt=1)

Generates a realization of the Markov Model

Parameters:
  • N (int) – trajectory length in steps of the lag time
  • start (int, optional, default = None) – starting hidden state. If not given, will sample from the stationary distribution of the hidden transition matrix.
  • stop (int or int-array-like, optional, default = None) – stopping hidden set. If given, the trajectory will be stopped before N steps once a hidden state of the stop set is reached
  • dt (int) – trajectory will be saved every dt time steps. Internally, the dt’th power of P is taken to ensure a more efficient simulation.
Returns:

htraj – The state trajectory with length N/dt

Return type:

(N/dt, ) ndarray

sparse

Returns whether the MSM is sparse

stationary_distribution
timescales(k=None)

The relaxation timescales corresponding to the eigenvalues

Parameters:k (int) – number of timescales to be returned. By default all available eigenvalues, minus 1.
Returns:ts – relaxation timescales in units of the input trajectory time step, defined by \(-\tau / ln | \lambda_i |, i = 2,...,k+1\).
Return type:ndarray(m)
timestep_model

Physical time corresponding to one transition matrix step, e.g. ‘10 ps’

trajectory_weights()

Uses the MSM to assign a probability weight to each trajectory frame.

This is a powerful function for the calculation of arbitrary observables in the trajectories one has started the analysis with. The stationary probability of the MSM will be used to reweigh all states. Returns a list of weight arrays, one for each trajectory, and with a number of elements equal to trajectory frames. Given \(N\) trajectories of lengths \(T_1\) to \(T_N\), this function returns corresponding weights:

\[(w_{1,1}, ..., w_{1,T_1}), (w_{N,1}, ..., w_{N,T_N})\]

that are normalized to one:

\[\sum_{i=1}^N \sum_{t=1}^{T_i} w_{i,t} = 1\]

Suppose you are interested in computing the expectation value of a function \(a(x)\), where \(x\) are your input configurations. Use this function to compute the weights of all input configurations and obtain the estimated expectation by:

\[\langle a \rangle = \sum_{i=1}^N \sum_{t=1}^{T_i} w_{i,t} a(x_{i,t})\]

Or if you are interested in computing the time-lagged correlation between functions \(a(x)\) and \(b(x)\) you could do:

\[\langle a(t) b(t+\tau) \rangle_t = \sum_{i=1}^N \sum_{t=1}^{T_i} w_{i,t} a(x_{i,t}) a(x_{i,t+\tau})\]
Returns:weights – The normalized trajectory weights. Given \(N\) trajectories of lengths \(T_1\) to \(T_N\), returns the corresponding weights:
\[(w_{1,1}, ..., w_{1,T_1}), (w_{N,1}, ..., w_{N,T_N})\]
Return type:list of ndarray
transition_matrix

The transition matrix on the active set.

update_model_params(**params)

Update given model parameter if they are set to specific values

References

The mathematical theory of Markov (state) model estimation was introduced in [1]_. Further theoretical developments were made in [2]_. The term Markov state model was coined in [3]. Continuous-time Markov models (Master equation models) were suggested in [4]. Reversible Markov model estimation was introduced in [5], and further developed in [6],[7]_,[9]_. It was shown in [8] that the quality of Markov state models does in fact not depend on memory loss, but rather on where the discretization is suitable to approximate the eigenfunctions of the Markov operator (the ‘reaction coordinates’). With a suitable choice of discretization and lag time, MSMs can thus become very accurate. [9] introduced a number of methodological improvements and gives a good overview of the methodological basics of Markov state modeling today. [10] is a more extensive review book of theory, methods and applications.

[1]Schuette, C. , A. Fischer, W. Huisinga and P. Deuflhard: A Direct Approach to Conformational Dynamics based on Hybrid Monte Carlo. J. Comput. Phys., 151, 146-168 (1999)
[2]Swope, W. C., J. W. Pitera and F. Suits: Describing protein folding kinetics by molecular dynamics simulations: 1. Theory J. Phys. Chem. B 108, 6571-6581 (2004)
[3]Singhal, N., C. D. Snow, V. S. Pande: Using path sampling to build better Markovian state models: Predicting the folding rate and mechanism of a tryptophan zipper beta hairpin. J. Chem. Phys. 121, 415 (2004).
[4]Sriraman, S., I. G. Kevrekidis and G. Hummer, G. J. Phys. Chem. B 109, 6479-6484 (2005)
[5]Noe, F.: Probability Distributions of Molecular Observables computed from Markov Models. J. Chem. Phys. 128, 244103 (2008)
[6]Buchete, N.-V. and Hummer, G.: Coarse master equations for peptide folding dynamics. J. Phys. Chem. B 112, 6057–6069 (2008)
[7]Bowman, G. R., K. A. Beauchamp, G. Boxer and V. S. Pande: Progress and challenges in the automated construction of Markov state models for full protein systems. J. Chem. Phys. 131, 124101 (2009)
[8]Sarich, M., F. Noe and C. Schuette: On the approximation quality of Markov state models. SIAM Multiscale Model. Simul. 8, 1154-1177 (2010)
[9]Prinz, J.-H., H. Wu, M. Sarich, B. Keller, M. Senne, M. Held, J. D. Chodera, C. Schuette and F. Noe: Markov models of molecular kinetics: Generation and Validation J. Chem. Phys. 134, 174105 (2011)
[10]Bowman, G. R., V. S. Pande and F. Noe: An Introduction to Markov State Models and Their Application to Long Timescale Molecular Simulation. Advances in Experimental Medicine and Biology 797, Springer, Heidelberg (2014)

Example

>>> from pyemma import msm
>>> import numpy as np
>>> np.set_printoptions(precision=3)
>>> dtrajs = [[0,1,2,2,2,2,1,2,2,2,1,0,0,0,0,0,0,0], [0,0,0,0,1,1,2,2,2,2,2,2,2,1,0,0]]  # two trajectories
>>> mm = msm.estimate_markov_model(dtrajs, 2)

Which is the active set of states we are working on?

>>> print(mm.active_set)
[0 1 2]

Show the count matrix

>>> print(mm.count_matrix_active)
[[ 7.  2.  1.]
 [ 2.  0.  4.]
 [ 2.  3.  9.]]

Show the estimated transition matrix

>>> print(mm.transition_matrix)
[[ 0.7    0.167  0.133]
 [ 0.388  0.     0.612]
 [ 0.119  0.238  0.643]]

Is this model reversible (i.e. does it fulfill detailed balance)?

>>> print(mm.is_reversible)
True

What is the equilibrium distribution of states?

>>> print(mm.stationary_distribution)
[ 0.393  0.17   0.437]

Relaxation timescales?

>>> print(mm.timescales())
[ 3.415  1.297]

Mean first passage time from state 0 to 2:

>>> print(mm.mfpt(0, 2))  
9.929...