pyemma.thermo.wham¶
-
pyemma.thermo.
wham
(ttrajs, dtrajs, bias, maxiter=100000, maxerr=1e-15, save_convergence_info=0, dt_traj='1 step')¶ Weighted histogram analysis method
- Parameters
ttrajs (numpy.ndarray(T) of int, or list of numpy.ndarray(T_i) of int) – A single discrete trajectory or a list of discrete trajectories. The integers are indexes in 0,…,num_therm_states-1 enumerating the thermodynamic states the trajectory is in at any time.
dtrajs (numpy.ndarray(T) of int, or list of numpy.ndarray(T_i) of int) – A single discrete trajectory or a list of discrete trajectories. The integers are indexes in 0,…,num_conf_states-1 enumerating the num_conf_states Markov states or the bins the trajectory is in at any time.
bias (numpy.ndarray(shape=(num_therm_states, num_conf_states)) object) – bias_energies_full[j, i] is the bias energy in units of kT for each discrete state i at thermodynamic state j.
maxiter (int, optional, default=10000) – The maximum number of dTRAM iterations before the estimator exits unsuccessfully.
maxerr (float, optional, default=1e-15) – Convergence criterion based on the maximal free energy change in a self-consistent iteration step.
save_convergence_info (int, optional, default=0) – Every save_convergence_info iteration steps, store the actual increment and the actual loglikelihood; 0 means no storage.
dt_traj (str, optional, default='1 step') –
Description of the physical time corresponding to the lag. 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*’
- Returns
A stationary model which consists of thermodynamic quantities at all temperatures/thermodynamic states.
- Return type
A
MultiThermModel
object
Example
Umbrella sampling: Suppose we simulate in K umbrellas, centered at positions \(y_0,...,y_{K-1}\) with bias energies
\[b_k(x) = \frac{c_k}{2 \textrm{kT}} \cdot (x - y_k)^2\]Suppose we have one simulation of length T in each umbrella, and they are ordered from 0 to K-1. We have discretized the x-coordinate into 100 bins. Then dtrajs and ttrajs should each be a list of \(K\) arrays. dtrajs would look for example like this:
[ (0, 0, 0, 0, 1, 1, 1, 0, 0, 0, ...), (0, 1, 0, 1, 0, 1, 1, 0, 0, 1, ...), ... ]
where each array has length T, and is the sequence of bins (in the range 0 to 99) visited along the trajectory. ttrajs would look like this:
[ (0, 0, 0, 0, 0, 0, 0, 0, 0, 0, ...), (1, 1, 1, 1, 1, 1, 1, 1, 1, 1, ...), ... ]
Because trajectory 1 stays in umbrella 1 (index 0), trajectory 2 stays in umbrella 2 (index 1), and so forth. bias is a \(K \times n\) matrix with all reduced bias energies evaluated at all centers:
\[\begin{split}\left(\begin{array}{cccc} b_0(y_0) & b_0(y_1) & ... & b_0(y_{n-1}) \\ b_1(y_0) & b_1(y_1) & ... & b_1(y_{n-1}) \\ ... \\ b_{K-1}(y_0) & b_{K-1}(y_1) & ... & b_{K-1}(y_{n-1}) \end{array}\right)\end{split}\]Let us try the above example:
>>> from pyemma.thermo import wham >>> import numpy as np >>> ttrajs = [np.array([0,0,0,0,0,0,0,0,0,0]), np.array([1,1,1,1,1,1,1,1,1,1])] >>> dtrajs = [np.array([0,0,0,0,1,1,1,0,0,0]), np.array([0,1,0,1,0,1,1,0,0,1])] >>> bias = np.array([[0.0, 0.0], [0.5, 1.0]]) >>> wham_obj = wham(ttrajs, dtrajs, bias) >>> wham_obj.log_likelihood() -6.6... >>> wham_obj.state_counts array([[7, 3], [5, 5]]) >>> wham_obj.stationary_distribution array([ 0.5..., 0.4...])
See
MultiThermModel
for a full documentation.-
class
pyemma.thermo.models.multi_therm.
MultiThermModel
(*args, **kwargs)¶ Coupled set of stationary models at multiple thermodynamic states
Methods
expectation
(a)Equilibrium expectation value of a given observable.
get_model_params
([deep])Get parameters for this model.
load
(file_name[, model_name])Loads a previously saved PyEMMA object from disk.
meval
(f, *args, **kw)Evaluates the given function call for all models Returns the results of the calls in a list
save
(file_name[, model_name, overwrite, …])saves the current state of this object to given file and name.
set_model_params
([models, f_therm, pi, f, label])Call to set all basic model parameters.
update_model_params
(**params)Update given model parameter if they are set to specific values
Attributes
-
meval
(f, *args, **kw)¶ Evaluates the given function call for all models Returns the results of the calls in a list
-
set_model_params
(models=None, f_therm=None, pi=None, f=None, label='ground state')¶ Call to set all basic model parameters.
- Parameters
pi (ndarray(n)) – Stationary distribution. If not already normalized, pi will be scaled to fulfill \(\sum_i \pi_i = 1\). The free energies f will then be computed from pi via \(f_i = - \log(\pi_i)\).
f (ndarray(n)) – Discrete-state free energies. If normalized_f = True, a constant will be added to normalize the stationary distribution. Otherwise f is left as given. Then, pi will be computed from f via \(\pi_i = \exp(-f_i)\) and, if necessary, scaled to fulfill \(\sum_i \pi_i = 1\). If both (pi and f) are given, f takes precedence over pi.
normalize_energy (bool, default=True) – If parametrized by free energy f, normalize them such that \(\sum_i \pi_i = 1\), which is achieved by \(\log \sum_i \exp(-f_i) = 0\).
label (str, default=None) – Human-readable description for the thermodynamic state of this model. May contain a temperature description, such as ‘300 K’ or a description of bias energy such as ‘unbiased’ or ‘Umbrella 1’.
-
property
unbiased_state
¶ Index of the unbiased thermodynamic state.
-
References
- 1
Ferrenberg, A.M. and Swensen, R.H. 1988. New Monte Carlo Technique for Studying Phase Transitions. Phys. Rev. Lett. 23, 2635–2638
- 2
Kumar, S. et al 1992. The Weighted Histogram Analysis Method for Free-Energy Calculations on Biomolecules. I. The Method. J. Comp. Chem. 13, 1011–1021