class TRAMModel¶
- class deeptime.markov.msm.TRAMModel(count_models, transition_matrices, biased_conf_energies, lagrangian_mult_log, modified_state_counts_log, therm_state_energies=None, markov_state_energies=None)¶
The TRAM model containing the estimated parameters, free energies, and the underlying Markov models for each thermodynamic state. The TRAM Model is the result of a TRAM estimation and can be used to calculate observables and recover an (unbiased) potential of mean force (PMF). TRAM is described in [1].
- Parameters:
count_models (list(TransitionCountModel)) – The transition count models for all thermodynamic states.
transition_matrices (ndarray(n, m, m), float64) – The estimated transition matrices for each thermodynamic state. The transition matrices and count models are combined into a MarkovStateModelCollection that holds the Markov models for each thermodynamic state.
biased_conf_energies (ndarray(n, m), float64) – The estimated free energies of each Markov state for all thermodynamic states. biased_conf_energies[k,i] contains the bias energy of Markov state in thermodynamic state .
lagrangian_mult_log (ndarray(n, m), float64) – The estimated logarithm of the lagrange multipliers of each Markov state for all thermodynamic states. lagrangian_mult_log[k,i] contains the lagrange multiplier of Markov state in thermodynamic state .
modified_state_counts_log (ndarray(n, m), float64) – The logarithm of the modified state counts of each Markov state for all thermodynamic states. modified_state_counts_log[k,i] contains the state counts of Markov state in thermodynamic state .
therm_state_energies (ndarray(n), float64) – The estimated free energy of each thermodynamic state, .
markov_state_energies (ndarray(m), float64) – The estimated free energy of each Markov state, .
See also
References
Attributes
The estimated free energy per thermodynamic state and Markov state, , where is the thermodynamic state index, and the Markov state index.
The estimated logarithm of the lagrange multipliers of each Markov state for all thermodynamic states.
The estimated free energy per Markov state, , where is the Markov state index.
The logarithm of the modified state counts of each Markov state for all thermodynamic states.
The underlying
MarkovStateModelCollection
.The estimated free energy per thermodynamic state, , where is the thermodynamic state index.
Methods
compute_PMF
(dtrajs, bias_matrices, bin_indices)Compute the potential of mean force over a number of bins.
compute_log_likelihood
(dtrajs, bias_matrices)The (parameter-dependent part of the) likelihood to observe the given data.
compute_observable
(observable_values, ...[, ...])Compute an observable value.
compute_sample_weights_log
(dtrajs, bias_matrices)Compute the log of the sample weight for all samples .
copy
()Makes a deep copy of this model.
get_params
([deep])Get the parameters.
set_params
(**params)Set the parameters of this estimator.
- compute_PMF(dtrajs, bias_matrices, bin_indices, n_bins=None, therm_state=-1)¶
Compute the potential of mean force over a number of bins.
- Parameters:
dtrajs (list(np.ndarray)) – The list of discrete trajectories. dtrajs[i][n] contains the Markov state index of the -th sample in the -th trajectory.
bias_matrices (list(np.ndarray)) – The bias energy matrices. bias_matrices[i][n, k] contains the bias energy of the -th sample from the -th trajectory, evaluated at thermodynamic state , . The bias energy matrices should have the same size as dtrajs in both the first and second dimension. The third dimension is of size n_therm_state, i.e. for each sample, the bias energy in every thermodynamic state is calculated and stored in the bias_matrices.
bin_indices (list(np.ndarray)) – The list of bin indices that the samples are binned into. The PMF is calculated as a distribution over these bins. binned_samples[i][n] contains the bin index for the -th sample in the -th trajectory.
n_bins (int, optional) – The number of bins. If None, n_bins is inferred from the maximum bin index. The PMF array
therm_state (int, optional, default=-1) – The index of the thermodynamic state in which the PMF need to be computed. If therm_state=-1, the PMF is computed for the unbiased (reference) state.
- Returns:
PMF – A ndarray of size n_bins containing the estimated PMF from the data. Is n_bins was None, the PMF is of size max(bin_indices) + 1.
- Return type:
np.ndarray
- compute_log_likelihood(dtrajs, bias_matrices)¶
The (parameter-dependent part of the) likelihood to observe the given data.
The definition can be found in [1], Equation (9).
- Parameters:
dtrajs (list(np.ndarray)) – The list of discrete trajectories. dtrajs[i][n] contains the Markov state index of the -th sample in the -th trajectory.
bias_matrices (list(np.ndarray)) – The bias energy matrices. bias_matrices[i][n, k] contains the bias energy of the -th sample from the -th trajectory, evaluated at thermodynamic state , . The bias energy matrices should have the same size as dtrajs in both the first and second dimension. The third dimension is of size n_therm_state, i.e. for each sample, the bias energy in every thermodynamic state is calculated and stored in the bias_matrices.
- Returns:
log_likelihood – The parameter-dependent part of the log-likelihood.
- Return type:
float
Notes
Parameter-dependent, i.e., the factor
does not occur in the log-likelihood as it is constant with respect to the parameters, leading to
- compute_observable(observable_values, dtrajs, bias_matrices, therm_state=-1)¶
Compute an observable value.
- Parameters:
observable_values (list(np.ndarray)) – The list of observable values. observable_values[i][n] contains the observable value for the -th sample in the -th trajectory.
dtrajs (list(np.ndarray)) – The list of discrete trajectories. dtrajs[i][n] contains the Markov state index of the -th sample in the -th trajectory.
bias_matrices (list(np.ndarray)) – The bias energy matrices. bias_matrices[i][n, k] contains the bias energy of the -th sample from the -th trajectory, evaluated at thermodynamic state , . The bias energy matrices should have the same size as dtrajs in both the first and second dimension. The third dimension is of size n_therm_state, i.e. for each sample, the bias energy in every thermodynamic state is calculated and stored in the bias_matrices.
therm_state (int, optional, default=-1) – The index of the thermodynamic state in which the observable need to be computed. If therm_state=-1, the observable is computed for the unbiased (reference) state.
- Returns:
observable – The expected value of the observable in the selected thermodynamic state therm_state.
- Return type:
float
- compute_sample_weights_log(dtrajs, bias_matrices, therm_state=-1)¶
Compute the log of the sample weight for all samples . If the thermodynamic state index is >= 0, the sample weights for that thermodynamic state will be computed, i.e. . Otherwise, this gives the unbiased sample weights .
- Parameters:
dtrajs (list(np.ndarray)) – The list of discrete trajectories. dtrajs[i][n] contains the Markov state index of the -th sample in the -th trajectory. Sample weights are computed for each of the samples in the dtrajs.
bias_matrices (list(np.ndarray)) – The bias energy matrices. bias_matrices[i][n, k] contains the bias energy of the -th sample from the -th trajectory, evaluated at thermodynamic state , i.e. .
therm_state (int, optional) – The index of the thermodynamic state in which the sample weights need to be computed. If therm_state=-1, the unbiased sample weights are computed.
- Returns:
sample_weights_log – The log of the statistical weight of each sample (i.e., the probability distribution over all samples: the sum over all sample weights equals one.)
- Return type:
np.ndarray
Notes
The statistical distribution is given by
- copy() Model ¶
Makes a deep copy of this model.
- Returns:
A new copy of this model.
- Return type:
copy
- get_params(deep=False)¶
Get the parameters.
- Returns:
params – Parameter names mapped to their values.
- Return type:
mapping of string to any
- 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 latter have parameters of the form
<component>__<parameter>
so that it’s possible to update each component of a nested object.- Parameters:
**params (dict) – Estimator parameters.
- Returns:
self – Estimator instance.
- Return type:
object
- property biased_conf_energies: ndarray¶
The estimated free energy per thermodynamic state and Markov state, , where is the thermodynamic state index, and the Markov state index.
- property lagrangian_mult_log: ndarray¶
The estimated logarithm of the lagrange multipliers of each Markov state for all thermodynamic states. lagrangian_mult_log[k,i] contains the lagrange multiplier of Markov state in thermodynamic state .
- property markov_state_energies: ndarray¶
The estimated free energy per Markov state, , where is the Markov state index.
- property modified_state_counts_log¶
The logarithm of the modified state counts of each Markov state for all thermodynamic states. modified_state_counts_log[k,i] contains the state counts of Markov state in thermodynamic state .
- property msm_collection¶
The underlying
MarkovStateModelCollection
. Contains one Markov state model for each sampled thermodynamic state.- Getter:
The collection of markov state models containing one model for each thermodynamic state.
- Type:
- property therm_state_energies: ndarray¶
The estimated free energy per thermodynamic state, , where is the thermodynamic state index.