def initial_model_gaussian1d(observations, nstates, reversible=True): """Generate an initial model with 1D-Gaussian output densities Parameters ---------- observations : list of ndarray((T_i), dtype=float) list of arrays of length T_i with observation data nstates : int The number of states. Examples -------- Generate initial model for a gaussian output model. >>> from bhmm import testsystems >>> [model, observations, states] = testsystems.generate_synthetic_observations(output_model_type='gaussian') >>> initial_model = initial_model_gaussian1d(observations, model.nstates) """ ntrajectories = len(observations) # Concatenate all observations. collected_observations = np.array([], dtype=config.dtype) for o_t in observations: collected_observations = np.append(collected_observations, o_t) # Fit a Gaussian mixture model to obtain emission distributions and state stationary probabilities. from bhmm._external.sklearn import mixture gmm = mixture.GMM(n_components=nstates) gmm.fit(collected_observations[:, None]) from bhmm import GaussianOutputModel output_model = GaussianOutputModel(nstates, means=gmm.means_[:, 0], sigmas=np.sqrt(gmm.covars_[:, 0])) logger().info("Gaussian output model:\n" + str(output_model)) # Extract stationary distributions. Pi = np.zeros([nstates], np.float64) Pi[:] = gmm.weights_[:] logger().info("GMM weights: %s" % str(gmm.weights_)) # Compute fractional state memberships. Nij = np.zeros([nstates, nstates], np.float64) for o_t in observations: # length of trajectory T = o_t.shape[0] # output probability pobs = output_model.p_obs(o_t) # normalize pobs /= pobs.sum(axis=1)[:, None] # Accumulate fractional transition counts from this trajectory. for t in range(T - 1): Nij[:, :] = Nij[:, :] + np.outer(pobs[t, :], pobs[t + 1, :]) logger().info("Nij\n" + str(Nij)) # Compute transition matrix maximum likelihood estimate. import msmtools.estimation as msmest Tij = msmest.transition_matrix(Nij, reversible=reversible) # Update model. model = HMM(Tij, output_model, reversible=reversible) return model
def initial_model_discrete(observations, nstates, lag=1, reversible=True): """Generate an initial model with discrete output densities Parameters ---------- observations : list of ndarray((T_i), dtype=int) list of arrays of length T_i with observation data nstates : int The number of states. lag : int, optional, default=1 The lag time to use for initializing the model. TODO ---- * Why do we have a `lag` option? Isn't the HMM model, by definition, lag=1 everywhere? Why would this be useful instead of just having the user subsample the data? Examples -------- Generate initial model for a discrete output model. >>> from bhmm import testsystems >>> [model, observations, states] = testsystems.generate_synthetic_observations(output_model_type='discrete') >>> initial_model = initial_model_discrete(observations, model.nstates) """ # check input if not reversible: warnings.warn( "nonreversible initialization of discrete HMM currently not supported. Using a reversible matrix for initialization." ) reversible = True # estimate Markov model C_full = msmtools.estimation.count_matrix(observations, lag) lcc = msmtools.estimation.largest_connected_set(C_full) Clcc = msmtools.estimation.largest_connected_submatrix(C_full, lcc=lcc) T = msmtools.estimation.transition_matrix(Clcc, reversible=True).toarray() # pcca pcca = msmtools.analysis.dense.pcca.PCCA(T, nstates) # default output probability, in order to avoid zero columns nstates_full = C_full.shape[0] eps = 0.01 * (1.0 / nstates_full) # Use PCCA distributions, but avoid 100% assignment to any state (prevents convergence) B_conn = np.maximum(pcca.output_probabilities, eps) # full state space output matrix B = eps * np.ones((nstates, nstates_full), dtype=np.float64) # expand B_conn to full state space B[:, lcc] = B_conn[:, :] # renormalize B to make it row-stochastic B /= B.sum(axis=1)[:, None] # coarse-grained transition matrix M = pcca.memberships W = np.linalg.inv(np.dot(M.T, M)) A = np.dot(np.dot(M.T, T), M) P_coarse = np.dot(W, A) # symmetrize and renormalize to eliminate numerical errors X = np.dot(np.diag(pcca.coarse_grained_stationary_probability), P_coarse) X = 0.5 * (X + X.T) # if there are values < 0, set to eps X = np.maximum(X, eps) # turn into coarse-grained transition matrix A = X / X.sum(axis=1)[:, None] logger().info('Initial model: ') logger().info('transition matrix = \n' + str(A)) logger().info('output matrix = \n' + str(B.T)) # initialize HMM # -------------- output_model = DiscreteOutputModel(B) model = HMM(A, output_model) return model