pyemma.msm.estimate_hidden_markov_model¶
-
pyemma.msm.
estimate_hidden_markov_model
(dtrajs, nstates, lag, reversible=True, stationary=False, connectivity=None, mincount_connectivity='1/n', separate=None, observe_nonempty=True, stride=1, dt_traj='1 step', accuracy=0.001, maxit=1000)¶ Estimates a Hidden Markov state model from discrete trajectories
Returns a
MaximumLikelihoodHMSM
that contains a transition matrix between a few (hidden) metastable states. Each metastable state has a probability distribution of visiting the discrete ‘microstates’ contained in the input trajectories. The resulting object is a hidden Markov model that allows to compute a large number of quantities.- 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) – lagtime for the MSM estimation in multiples of trajectory steps
nstates (int) – the number of metastable states in the resulting HMM
reversible (bool, optional, default = True) – If true compute reversible MSM, else non-reversible MSM
stationary (bool, optional, default=False) – If True, the initial distribution of hidden states is self-consistently computed as the stationary distribution of the transition matrix. If False, it will be estimated from the starting states. Only set this to true if you’re sure that the observation trajectories are initiated from a global equilibrium distribution.
connectivity (str, optional, default = None) –
Defines if the resulting HMM will be defined on all hidden states or on a connected subset. Connectivity is defined by counting only transitions with at least mincount_connectivity counts. If a subset of states is used, all estimated quantities (transition matrix, stationary distribution, etc) are only defined on this subset and are correspondingly smaller than nstates. Following modes are available: * None or ‘all’ : The active set is the full set of states.
Estimation is done on all weakly connected subsets separately. The resulting transition matrix may be disconnected.
’largest’ : The active set is the largest reversibly connected set.
- ’populous’The active set is the reversibly connected set with
most counts.
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.
separate (None or iterable of int) – Force the given set of observed states to stay in a separate hidden state. The remaining nstates-1 states will be assigned by a metastable decomposition.
observe_nonempty (bool) – If True, will restricted the observed states to the states that have at least one observation in the lagged input trajectories.
dt_traj (str, optional, default='1 step') –
Description of the physical time corresponding to the trajectory 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*’accuracy (float) – convergence threshold for EM iteration. When two the likelihood does not increase by more than accuracy, the iteration is stopped successfully.
maxit (int) – stopping criterion for EM iteration. When so many iterations are performed without reaching the requested accuracy, the iteration is stopped without convergence (a warning is given)
- Returns
hmsm – Estimator object containing the HMSM and estimation information.
- Return type
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_hidden_markov_model(dtrajs, 2, 2)
We have estimated a 2x2 hidden transition matrix between the metastable states:
>>> print(mm.transition_matrix) [[ 0.684 0.316] [ 0.242 0.758]]
With the equilibrium distribution:
>>> print(mm.stationary_distribution) [ 0.43... 0.56...]
The observed states are the three discrete clusters that we have in our discrete trajectory:
>>> print(mm.observable_set) [0 1 2]
The metastable distributions (mm.metastable_distributions), or equivalently the observation probabilities are the probability to be in a given cluster (‘microstate’) if we are in one of the hidden metastable states. So it’s a 2 x 3 matrix:
>>> print(mm.observation_probabilities) [[ 0.9620883 0.0379117 0. ] [ 0. 0.28014352 0.71985648]]
The first metastable state ist mostly in cluster 0, and a little bit in the transition state cluster 1. The second metastable state is less well defined, but mostly in cluster 2 and less prominently in the transition state cluster 1.
We can print the lifetimes of the metastable states:
>>> print(mm.lifetimes) [ 5... 7...]
And the timescale of the hidden transition matrix - now we only have one relaxation timescale:
>>> print(mm.timescales()) [ 2.4...]
The mean first passage times can also be computed between metastable states:
>>> print(mm.mfpt(0, 1)) 6.3...
-
class
pyemma.msm.estimators.maximum_likelihood_hmsm.
MaximumLikelihoodHMSM
(*args, **kwargs)¶ Maximum likelihood estimator for a Hidden MSM given a MSM
Methods
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
save
(file_name[, model_name, overwrite, …])saves the current state of this object to given file and name.
set_model_params
([P, pobs, pi, reversible, …])- param P
coarse-grained or hidden transition matrix
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
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
-
property
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
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)
-
property
discrete_trajectories_full
¶ A list of integer arrays with the original trajectories.
-
property
discrete_trajectories_lagged
¶ Transformed original trajectories that are used as an input into the HMM estimation
-
property
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
-
property
dt_traj
¶
-
property
dtrajs_full
¶ A list of integer arrays with the original trajectories.
-
property
dtrajs_lagged
¶ Transformed original trajectories that are used as an input into the HMM estimation
-
property
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
-
property
lagtime
¶ The lag time in steps
-
property
msm_init
¶
-
property
nstates
¶ Number of active states on which all computations and estimations are done
-
property
nstates_obs
¶ Number of states in discrete trajectories
-
property
observable_set
¶ The active set of states on which all computations and estimations will be done
-
property
observable_state_indexes
¶ Ensures that the observable states are indexed and returns the indices
-
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) )
-
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
-
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})
References
[1]_ is an excellent review of estimation algorithms for discrete Hidden Markov Models. This function estimates a discrete HMM on the discrete input states using the Baum-Welch algorithm [2]_. We use a maximum-likelihood Markov state model to initialize the HMM estimation as described in 3.
- 1
L. R. Rabiner: A Tutorial on Hidden Markov Models and Selected Applications in Speech Recognition. Proc. IEEE 77, 257-286 (1989)
- 2
L. Baum, T. Petrie, G. Soules and N. Weiss N: A maximization technique occurring in the statistical analysis of probabilistic functions of Markov chains. Ann. Math. Statist. 41, 164-171 (1970)
- 3
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)