pyemma.msm.BayesianHMSM

class pyemma.msm.BayesianHMSM(nstates=2, lag=1, stride='effective', p0_prior='mixed', transition_matrix_prior='mixed', nsamples=100, init_hmsm=None, reversible=True, stationary=False, connectivity='largest', mincount_connectivity='1/n', separate=None, observe_nonempty=True, dt_traj='1 step', conf=0.95, store_hidden=False, show_progress=True)

Estimator for a Bayesian Hidden Markov state model

__init__(nstates=2, lag=1, stride='effective', p0_prior='mixed', transition_matrix_prior='mixed', nsamples=100, init_hmsm=None, reversible=True, stationary=False, connectivity='largest', mincount_connectivity='1/n', separate=None, observe_nonempty=True, dt_traj='1 step', conf=0.95, store_hidden=False, show_progress=True)

Estimator for a Bayesian HMSM

Parameters
  • nstates (int, optional, default=2) – number of hidden states

  • lag (int, optional, default=1) – lagtime to estimate the HMSM at

  • stride (str or int, default=1) –

    stride between two lagged trajectories extracted from the input trajectories. Given trajectory s[t], stride and lag will result in trajectories

    s[0], s[tau], s[2 tau], … s[stride], s[stride + tau], s[stride + 2 tau], …

    Setting stride = 1 will result in using all data (useful for maximum likelihood estimator), while a Bayesian estimator requires a longer stride in order to have statistically uncorrelated trajectories. Setting stride = None ‘effective’ uses the largest neglected timescale as an estimate for the correlation time and sets the stride accordingly.

  • p0_prior (None, str, float or ndarray(n)) –

    Prior for the initial distribution of the HMM. Will only be active if stationary=False (stationary=True means that p0 is identical to the stationary distribution of the transition matrix). Currently implements different versions of the Dirichlet prior that is conjugate to the Dirichlet distribution of p0. p0 is sampled from:

    where \(n_i\) are the number of times a hidden trajectory was in state \(i\) at time step 0 and \(a_i\) is the prior count. Following options are available:

    • ’mixed’ (default), \(a_i = p_{0,init}\), where \(p_{0,init}\) is the initial distribution of initial_model.

    • ndarray(n) or float, the given array will be used as A.

    • ’uniform’, \(a_i = 1\)

    • None, \(a_i = 0\). This option ensures coincidence between sample mean an MLE. Will sooner or later lead to sampling problems, because as soon as zero trajectories are drawn from a given state, the sampler cannot recover and that state will never serve as a starting state subsequently. Only recommended in the large data regime and when the probability to sample zero trajectories from any state is negligible.

  • transition_matrix_prior (str or ndarray(n, n)) –

    Prior for the HMM transition matrix. Currently implements Dirichlet priors if reversible=False and reversible transition matrix priors as described in 3 if reversible=True. For the nonreversible case the posterior of transition matrix \(P\) is:

    where \(c_{ij}\) are the number of transitions found for hidden trajectories and \(b_{ij}\) are prior counts.

    • ’mixed’ (default), \(b_{ij} = p_{ij,init}\), where \(p_{ij,init}\) is the transition matrix of initial_model. That means one prior count will be used per row.

    • ndarray(n, n) or broadcastable, the given array will be used as B.

    • ’uniform’, \(b_{ij} = 1\)

    • None, \(b_ij = 0\). This option ensures coincidence between sample mean an MLE. Will sooner or later lead to sampling problems, because as soon as a transition \(ij\) will not occur in a sample, the sampler cannot recover and that transition will never be sampled again. This option is not recommended unless you have a small HMM and a lot of data.

  • init_hmsm (HMSM, default=None) – Single-point estimate of HMSM object around which errors will be evaluated. If None is give an initial estimate will be automatically generated using the given parameters.

  • store_hidden (bool, optional, default=False) – store hidden trajectories in sampled HMMs

  • show_progress (bool, default=True) – Show progressbars for calculation?

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)

2

J. D. Chodera Et Al: Bayesian hidden Markov model analysis of single-molecule force spectroscopy: Characterizing kinetics under measurement uncertainty. arXiv:1108.1430 (2011)

3

Trendelkamp-Schroer, B., H. Wu, F. Paul and F. Noe: Estimation and uncertainty of reversible Markov models. J. Chem. Phys. 143, 174101 (2015).

Methods

__init__([nstates, lag, stride, p0_prior, …])

Estimator for a Bayesian HMSM

cktest([mlags, conf, err_est, n_jobs, …])

Conducts a Chapman-Kolmogorow test.

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(X, **params)

Estimates the model given the data X

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[, y])

Estimates parameters - for compatibility with sklearn.

get_model_params([deep])

Get parameters for this model.

get_params([deep])

Get parameters for this estimator.

load(file_name[, model_name])

Loads a previously saved PyEMMA object from disk.

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_observation_probabilities(nsample)

Generates samples according to the current observation probability distribution

sample_conf(f, *args, **kwargs)

Sample confidence interval of numerical method f over all samples

sample_f(f, *args, **kwargs)

Evaluated method f for all samples

sample_mean(f, *args, **kwargs)

Sample mean of numerical method f over all samples

sample_std(f, *args, **kwargs)

Sample standard deviation of numerical method f over all samples

save(file_name[, model_name, overwrite, …])

saves the current state of this object to given file and name.

set_model_params([samples, conf, P, pobs, …])

param samples

sampled MSMs

set_params(**params)

Set the parameters of this estimator.

simulate(N[, start, stop, dt])

Generates a realization of the Hidden Markov Model

submodel([states, obs, …])

Returns a HMM with restricted state space

submodel_disconnect([mincount_connectivity])

Disconnects sets of hidden states that are barely connected

submodel_largest([strong, mincount_connectivity])

Returns the largest connected sub-HMM (convenience function)

submodel_populous([strong, …])

Returns the most populous connected sub-HMM (convenience function)

timescales([k])

The relaxation timescales corresponding to the eigenvalues

trajectory_weights()

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

transition_matrix_obs([k])

Computes the transition matrix between observed states

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_set

The active set of hidden states on which all hidden state computations are done

discrete_trajectories_full

A list of integer arrays with the original trajectories.

discrete_trajectories_lagged

Transformed original trajectories that are used as an input into the HMM estimation

discrete_trajectories_obs

A list of integer arrays with the discrete trajectories mapped to the observation mode used.

dt_model

Description of the physical time corresponding to the lag.

dt_traj

dtrajs_full

A list of integer arrays with the original trajectories.

dtrajs_lagged

Transformed original trajectories that are used as an input into the HMM estimation

dtrajs_obs

A list of integer arrays with the discrete trajectories mapped to the observation mode used.

eigenvectors_left_obs

eigenvectors_right_obs

init_hmsm

is_reversible

Returns whether the MSM is reversible

is_sparse

Returns whether the MSM is sparse

lagtime

The lag time in steps

lifetimes

Lifetimes of states of the hidden transition matrix

logger

The logger for this class instance

metastable_assignments

Computes the assignment to metastable sets for observable states

metastable_distributions

Returns the output probability distributions. Identical to

metastable_memberships

Computes the memberships of observable states to metastable sets by

metastable_sets

Computes the metastable sets of observable states within each

model

The model estimated by this Estimator

msm_init

n_metastable

Number of states chosen for PCCA++ computation.

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_obs

Number of states in discrete trajectories

observable_set

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

observable_state_indexes

Ensures that the observable states are indexed and returns the indices

observation_probabilities

returns the output probability matrix

pi

The stationary distribution on the MSM states

reversible

Returns whether the MSM is reversible

samples

show_progress

whether to show the progress of heavy calculations on this object.

sparse

Returns whether the MSM is sparse

stationary_distribution

The stationary distribution on the MSM states

stationary_distribution_obs

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_set

The active set of hidden states on which all hidden state computations are done

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

Conducts a Chapman-Kolmogorow test.

Parameters
  • mlags (int or int-array, default=10) – 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, default = 0.95) – confidence interval

  • err_est (bool, default=False) – 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=None) – how many jobs to use during calculation

  • show_progress (bool, default=True) – Show progressbars for calculation?

Returns

cktest

Return type

ChapmanKolmogorovValidator

References

This is an adaption of the Chapman-Kolmogorov Test described in detail in [1]_ to Hidden MSMs as described in [2]_.

1

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

2

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

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)  # doctest: +SKIP

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.

discrete_trajectories_full

A list of integer arrays with the original trajectories.

discrete_trajectories_lagged

Transformed original trajectories that are used as an input into the HMM estimation

discrete_trajectories_obs

A list of integer arrays with the discrete trajectories mapped to the observation mode used. When using observe_active = True, the indexes will be given on the MSM active set. Frames that are not in the observation set will be -1. When observe_active = False, this attribute is identical to discrete_trajectories_full

dt_model

Description of the physical time corresponding to the lag.

dtrajs_full

A list of integer arrays with the original trajectories.

dtrajs_lagged

Transformed original trajectories that are used as an input into the HMM estimation

dtrajs_obs

A list of integer arrays with the discrete trajectories mapped to the observation mode used. When using observe_active = True, the indexes will be given on the MSM active set. Frames that are not in the observation set will be -1. When observe_active = False, this attribute is identical to discrete_trajectories_full

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(X, **params)

Estimates the model given the data X

Parameters
  • X (object) – A reference to the data from which the model will be estimated

  • params (dict) – New estimation parameter values. The parameters must that have been announced in the __init__ method of this estimator. The present settings will overwrite the settings of parameters given in the __init__ method, i.e. the parameter values after this call will be those that have been used for this estimation. Use this option if only one or a few parameters change with respect to the __init__ settings for this run, and if you don’t need to remember the original settings of these changed parameters.

Returns

estimator – The estimated estimator with the model being available.

Return type

object

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, y=None)

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

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

is_reversible

Returns whether the MSM is reversible

is_sparse

Returns whether the MSM is sparse

lagtime

The lag time in steps

lifetimes

Lifetimes of states of the hidden transition matrix

Returns

l – state lifetimes in units of the input trajectory time step, defined by \(-\tau / ln \mid p_{ii} \mid, i = 1,...,nstates\), where \(p_{ii}\) are the diagonal entries of the hidden transition matrix.

Return type

ndarray(nstates)

classmethod load(file_name, model_name='default')

Loads a previously saved PyEMMA object from disk.

Parameters
  • file_name (str or file like object (has to provide read method)) – The file like object tried to be read for a serialized object.

  • model_name (str, default='default') – if multiple models are contained in the file, these can be accessed by their name. Use pyemma.list_models() to get a representation of all stored models.

Returns

obj

Return type

the de-serialized object

logger

The logger for this class instance

metastable_assignments

Computes the assignment to metastable sets for observable states

Notes

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

Returns

For each observable state, the metastable state it is located in.

Return type

ndarray((n) ,dtype=int)

metastable_distributions
Returns the output probability distributions. Identical to

observation_probabilities()

Returns

Pout – output probability matrix from hidden to observable discrete states

Return type

ndarray (m,n)

metastable_memberships
Computes the memberships of observable states to metastable sets by

Bayesian inversion as described in [1]_.

Returns

M – A matrix containing the probability or membership of each observable state to be assigned to each metastable or hidden state. The row sums of M are 1.

Return type

ndarray((n,m))

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)

metastable_sets
Computes the metastable sets of observable states within each

metastable set

Notes

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 observable state indexes contained in it

Return type

list of int-arrays

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

n_metastable

Number of states chosen for PCCA++ computation.

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_obs

Number of states in discrete trajectories

observable_set

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

observable_state_indexes

Ensures that the observable states are indexed and returns the indices

observation_probabilities

returns the output probability matrix

Returns

Pout – output probability matrix from hidden to observable discrete states

Return type

ndarray (m,n)

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

The stationary distribution on the MSM states

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_observation_probabilities(nsample)

Generates samples according to the current observation probability distribution

Parameters

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_conf(f, *args, **kwargs)

Sample confidence interval of numerical method f over all samples

Calls f(*args, **kwargs) on all samples and computes the confidence interval. Size of confidence interval is given in the construction of the SampledModel. f must return a numerical value or an ndarray.

Parameters
  • f (method reference or name (str)) – Model method to be evaluated for each model sample

  • args (arguments) – Non-keyword arguments to be passed to the method in each call

  • kwargs (keyword-argments) – Keyword arguments to be passed to the method in each call

Returns

  • L (float or ndarray) – lower value or array of confidence interval

  • R (float or ndarray) – upper value or array of confidence interval

sample_f(f, *args, **kwargs)

Evaluated method f for all samples

Calls f(*args, **kwargs) on all samples.

Parameters
  • f (method reference or name (str)) – Model method to be evaluated for each model sample

  • args (arguments) – Non-keyword arguments to be passed to the method in each call

  • kwargs (keyword-argments) – Keyword arguments to be passed to the method in each call

Returns

vals – list of results of the method calls

Return type

list

sample_mean(f, *args, **kwargs)

Sample mean of numerical method f over all samples

Calls f(*args, **kwargs) on all samples and computes the mean. f must return a numerical value or an ndarray.

Parameters
  • f (method reference or name (str)) – Model method to be evaluated for each model sample

  • args (arguments) – Non-keyword arguments to be passed to the method in each call

  • kwargs (keyword-argments) – Keyword arguments to be passed to the method in each call

Returns

mean – mean value or mean array

Return type

float or ndarray

sample_std(f, *args, **kwargs)

Sample standard deviation of numerical method f over all samples

Calls f(*args, **kwargs) on all samples and computes the standard deviation. f must return a numerical value or an ndarray.

Parameters
  • f (method reference or name (str)) – Model method to be evaluated for each model sample

  • args (arguments) – Non-keyword arguments to be passed to the method in each call

  • kwargs (keyword-argments) – Keyword arguments to be passed to the method in each call

Returns

std – standard deviation or array of standard deviations

Return type

float or ndarray

save(file_name, model_name='default', overwrite=False, save_streaming_chain=False)

saves the current state of this object to given file and name.

Parameters
  • file_name (str) – path to desired output file

  • model_name (str, default='default') – creates a group named ‘model_name’ in the given file, which will contain all of the data. If the name already exists, and overwrite is False (default) will raise a RuntimeError.

  • overwrite (bool, default=False) – Should overwrite existing model names?

  • save_streaming_chain (boolean, default=False) – if True, the data_producer(s) of this object will also be saved in the given file.

Examples

>>> import pyemma, numpy as np
>>> from pyemma.util.contexts import named_temporary_file
>>> m = pyemma.msm.MSM(P=np.array([[0.1, 0.9], [0.9, 0.1]]))
>>> with named_temporary_file() as file: # doctest: +SKIP
...    m.save(file, 'simple') # doctest: +SKIP
...    inst_restored = pyemma.load(file, 'simple') # doctest: +SKIP
>>> np.testing.assert_equal(m.P, inst_restored.P) # doctest: +SKIP
set_model_params(samples=None, conf=0.95, P=None, pobs=None, pi=None, reversible=None, dt_model='1 step', neig=None)
Parameters
  • samples (list of MSM objects) – sampled MSMs

  • conf (float, optional, default=0.68) – Confidence interval. By default one-sigma (68.3%) is used. Use 95.4% for two sigma or 99.7% for three sigma.

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

show_progress

whether to show the progress of heavy calculations on this object.

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

Generates a realization of the Hidden 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 ((N/dt, ) ndarray) – The hidden state trajectory with length N/dt

  • otraj ((N/dt, ) ndarray) – The observable state discrete trajectory with length N/dt

sparse

Returns whether the MSM is sparse

stationary_distribution

The stationary distribution on the MSM states

submodel(states=None, obs=None, mincount_connectivity='1/n', inplace=False)

Returns a HMM with restricted state space

Parameters
  • states (None, str or int-array) – Hidden states to restrict the model to. In addition to specifying the subset, possible options are: * None : all states - don’t restrict * ‘populous-strong’ : strongly connected subset with maximum counts * ‘populous-weak’ : weakly connected subset with maximum counts * ‘largest-strong’ : strongly connected subset with maximum size * ‘largest-weak’ : weakly connected subset with maximum size

  • obs (None, str or int-array) – Observed states to restrict the model to. In addition to specifying an array with the state labels to be observed, possible options are: * None : all states - don’t restrict * ‘nonempty’ : all states with at least one observation in the estimator

  • mincount_connectivity (float or '1/n') – minimum number of counts to consider a connection between two states. Counts lower than that will count zero in the connectivity check and may thus separate the resulting transition matrix. Default value: 1/nstates.

  • inplace (Bool) – if True, submodel is estimated in-place, overwriting the original estimator and possibly discarding information. Default value: False

Returns

hmm – The restricted HMM.

Return type

HMM

submodel_disconnect(mincount_connectivity='1/n')

Disconnects sets of hidden states that are barely connected

Runs a connectivity check excluding all transition counts below mincount_connectivity. The transition matrix and stationary distribution will be re-estimated. Note that the resulting transition matrix may have both strongly and weakly connected subsets.

Parameters

mincount_connectivity (float or '1/n') – minimum number of counts to consider a connection between two states. Counts lower than that will count zero in the connectivity check and may thus separate the resulting transition matrix. The default evaluates to 1/nstates.

Returns

hmm – The restricted HMM.

Return type

HMM

submodel_largest(strong=True, mincount_connectivity='1/n')

Returns the largest connected sub-HMM (convenience function)

Returns

hmm – The restricted HMM.

Return type

HMM

submodel_populous(strong=True, mincount_connectivity='1/n')

Returns the most populous connected sub-HMM (convenience function)

Returns

hmm – The restricted HMM.

Return type

HMM

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 HMSM 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

  • The normalized trajectory weights. Given \(N\) trajectories of lengths \(T_1\) to \(T_N\),

  • returns the corresponding weights

  • .. math:: – (w_{1,1}, …, w_{1,T_1}), (w_{N,1}, …, w_{N,T_N})

transition_matrix

The transition matrix on the active set.

transition_matrix_obs(k=1)

Computes the transition matrix between observed states

Transition matrices for longer lag times than the one used to parametrize this HMSM can be obtained by setting the k option. Note that a HMSM is not Markovian, thus we cannot compute transition matrices at longer lag times using the Chapman-Kolmogorow equality. I.e.:

\[P (k \tau) \neq P^k (\tau)\]

This function computes the correct transition matrix using the metastable (coarse) transition matrix \(P_c\) as:

\[P (k \tau) = {\Pi}^-1 \chi^{\top} ({\Pi}_c) P_c^k (\tau) \chi\]

where \(\chi\) is the output probability matrix, \(\Pi_c\) is a diagonal matrix with the metastable-state (coarse) stationary distribution and \(\Pi\) is a diagonal matrix with the observable-state stationary distribution.

Parameters

k (int, optional, default=1) – Multiple of the lag time for which the By default (k=1), the transition matrix at the lag time used to construct this HMSM will be returned. If a higher power is given,

update_model_params(**params)

Update given model parameter if they are set to specific values