msmbuilder.hmm.VonMisesHMM¶
-
class
msmbuilder.hmm.
VonMisesHMM
¶ Hidden Markov Model with von Mises Emissions
The von Mises distribution, (also known as the circular normal distribution or Tikhonov distribution) is a continuous probability distribution on the circle. For multivariate signals, the emmissions distribution implemented by this model is a product of univariate von Mises distributuons – analogous to the multivariate Gaussian distribution with a diagonal covariance matrix.
This class allows for easy evaluation of, sampling from, and maximum-likelihood estimation of the parameters of a HMM.
Parameters: n_states : int
The number of components (states) in the model
n_init : int
Number of time the EM algorithm will be run with different random seeds. The final results will be the best output of n_init consecutive runs in terms of log likelihood.
n_iter : int
The maximum number of iterations of expectation-maximization to run during each fitting round.
thresh : float
Convergence threshold for the log-likelihood during expectation maximization. When the increase in the log-likelihood is less than thresh between subsequent rounds of E-M, fitting will finish.
reversible_type : str
Method by which the reversibility of the transition matrix is enforced. ‘mle’ uses a maximum likelihood method that is solved by numerical optimization (BFGS), and ‘transpose’ uses a more restrictive (but less computationally complex) direct symmetrization of the expected number of counts.
random_state : int, optional
Random state, used during sampling.
Notes
The formulas for the maximization step of the E-M algorithim are adapted from [R29], especially equations (11) and (13).
References
[R29] (1, 2) Prati, Andrea, Simone Calderara, and Rita Cucchiara. “Using circular statistics for trajectory shape analysis.” Computer Vision and Pattern Recognition, 2008. CVPR 2008. IEEE Conference on. IEEE, 2008. .. [R30] Murray, Richard F., and Yaniv Morgenstern. “Cue combination on the circle and the sphere.” Journal of vision 10.11 (2010).
Attributes
means_ (array, shape (n_states, n_features)) Mean parameters for each state. kappas_ (array, shape (n_states, n_features)) Concentration parameter for each state. If kappa is zero, the distriution is uniform. If large, the distribution is very concentrated around the mean. transmat_ (array, shape (n_states, n_states)) Matrix of transition probabilities between states. populations_ (array, shape(n_states)) Population of each state. fit_logprob_ (array) The log probability of the model after each iteration of fitting. Methods
fit
Estimate model parameters. fit_predict
Find most likely hidden-state sequence corresponding to each data timeseries. predict
Find most likely hidden-state sequence corresponding to each data timeseries. score
Log-likelihood of sequences under the model summarize
Get a string summarizing the model. -
fit
()¶ Estimate model parameters.
Parameters: sequences : list
List of 2-dimensional array observation sequences, each of which has shape (n_samples_i, n_features), where n_samples_i is the length of the i_th observation.
-
fit_predict
()¶ Find most likely hidden-state sequence corresponding to each data timeseries.
Uses the Viterbi algorithm.
Parameters: sequences : list
List of 2-dimensional array observation sequences, each of which has shape (n_samples_i, n_features), where n_samples_i is the length of the i_th observation.
Returns: viterbi_logprob : float
Log probability of the maximum likelihood path through the HMM.
hidden_sequences : list of np.ndarrays[dtype=int, shape=n_samples_i]
Index of the most likely states for each observation.
-
overlap_
¶ Compute the matrix of normalized log overlap integrals between the hidden state distributions
Returns: noverlap : array, shape=(n_states, n_states)
noverlap[i,j] gives the log of normalized overlap integral between states i and j. The normalized overlap integral is `integral(f(i)*f(j)) / sqrt[integral(f(i)*f(i)) * integral(f(j)*f(j))]
Notes
The analytic formula used here follows from equation (A4) in 2
-
predict
()¶ Find most likely hidden-state sequence corresponding to each data timeseries.
Uses the Viterbi algorithm.
Parameters: sequences : list
List of 2-dimensional array observation sequences, each of which has shape (n_samples_i, n_features), where n_samples_i is the length of the i_th observation.
Returns: viterbi_logprob : float
Log probability of the maximum likelihood path through the HMM.
hidden_sequences : list of np.ndarrays[dtype=int, shape=n_samples_i]
Index of the most likely states for each observation.
-
score
()¶ Log-likelihood of sequences under the model
Parameters: sequences : list
List of 2-dimensional array observation sequences, each of which has shape (n_samples_i, n_features), where n_samples_i is the length of the i_th observation.
-
summarize
()¶ Get a string summarizing the model.
-
timescales_
¶ The implied relaxation timescales of the hidden Markov transition matrix
By diagonalizing the transition matrix, its propagation of an arbitrary initial probability vector can be written as a sum of the eigenvectors of the transition weighted by per-eigenvector term that decays exponentially with time. Each of these eigenvectors describes a “dynamical mode” of the transition matrix and has a characteristic timescale, which gives the timescale on which that mode decays towards equilibrium. These timescales are given by \(-1/log(u_i)\) where \(u_i\) are the eigenvalues of the transition matrix. In a reversible HMM with N states, the number of timescales is at most N-1. (The -1 comes from the fact that the stationary distribution of the chain is associated with an eigenvalue of 1, and an infinite characteristic timescale). The number of timescales can be less than N-1 for every eigenvalue of the transition matrix that is negative (which is allowable by detailed balance).
Returns: timescales : array, shape=[n_timescales]
The characteristic timescales of the transition matrix. If the model has not been fit or does not have a transition matrix, the return value will be None.
-