pyemma.thermo.TRAM¶
-
class
pyemma.thermo.
TRAM
(lag, count_mode='sliding', connectivity='post_hoc_RE', nstates_full=None, equilibrium=None, maxiter=10000, maxerr=1e-15, save_convergence_info=0, dt_traj='1 step', nn=None, connectivity_factor=1.0, direct_space=False, N_dtram_accelerations=0, callback=None, init='mbar', init_maxiter=5000, init_maxerr=1e-08, overcounting_factor=1.0)¶ Transition(-based) Reweighting Analysis Method.
-
__init__
(lag, count_mode='sliding', connectivity='post_hoc_RE', nstates_full=None, equilibrium=None, maxiter=10000, maxerr=1e-15, save_convergence_info=0, dt_traj='1 step', nn=None, connectivity_factor=1.0, direct_space=False, N_dtram_accelerations=0, callback=None, init='mbar', init_maxiter=5000, init_maxerr=1e-08, overcounting_factor=1.0)¶ Transition(-based) Reweighting Analysis Method
Parameters: - lag (int) – Integer lag time at which transitions are counted.
- count_mode (str, optional, default='sliding') –
mode to obtain count matrices from discrete trajectories. Should be one of: * ‘sliding’ : A trajectory of length T will have \(T-\tau\) counts at time indexes
\[(0 \rightarrow \tau), (1 \rightarrow \tau+1), ..., (T-\tau-1 \rightarrow T-1)\]- ’sample’ : A trajectory of length T will have \(T/\tau\) counts
at time indexes\[(0 \rightarrow \tau), (\tau \rightarrow 2 \tau), ..., ((T/\tau-1) \tau \rightarrow T)\]
Currently only ‘sliding’ is supported.
- ’sample’ : A trajectory of length T will have \(T/\tau\) counts
at time indexes
- connectivity (str, optional, default='post_hoc_RE') –
One of ‘post_hoc_RE’, ‘BAR_variance’, ‘reversible_pathways’ or ‘summed_count_matrix’. Defines what should be considered a connected set in the joint (product) space of conformations and thermodynamic ensembles. * ‘reversible_pathways’ : requires that every state in the connected set
can be reached by following a pathway of reversible transitions. A reversible transition between two Markov states (within the same thermodynamic state k) is a pair of Markov states that belong to the same strongly connected component of the count matrix (from thermodynamic state k). A pathway of reversible transitions is a list of reversible transitions [(i_1, i_2), (i_2, i_3),…, (i_(N-2), i_(N-1)), (i_(N-1), i_N)]. The thermodynamic state where the reversible transitions happen, is ignored in constructing the reversible pathways. This is equivalent to assuming that two ensembles overlap at some Markov state whenever there exist frames from both ensembles in that Markov state.- ’post_hoc_RE’ : similar to ‘reversible_pathways’ but with a more strict requirement for the overlap between thermodynamic states. It is required that every state in the connected set can be reached by following a pathway of reversible transitions or jumping between overlapping thermodynamic states while staying in the same Markov state. A reversible transition between two Markov states (within the same thermodynamic state k) is a pair of Markov states that belong to the same strongly connected component of the count matrix (from thermodynamic state k). Two thermodynamic states k and l are defined to overlap at Markov state n if a replica exchange simulation [2]_ restricted to state n would show at least one transition from k to l or one transition from from l to k. The expected number of replica exchanges is estimated from the simulation data. The minimal number required of replica exchanges per Markov state can be increased by decreasing connectivity_factor.
- ’BAR_variance’ : like ‘post_hoc_RE’ but with a different condition to define the thermodynamic overlap based on the variance of the BAR estimator [3]_. Two thermodynamic states k and l are defined to overlap at Markov state n if the variance of the free energy difference Delta f_{kl} computed with BAR (and restricted to conformations form Markov state n) is less or equal than one. The minimally required variance can be controlled with connectivity_factor.
- ’summed_count_matrix’ : all thermodynamic states are assumed to overlap. The connected set is then computed by summing the count matrices over all thermodynamic states and taking it’s largest strongly connected set. Not recommended!
For more details see
pyemma.thermo.extensions.cset.compute_csets_TRAM()
. - nstates_full (int, optional, default=None) – Number of cluster centers, i.e., the size of the full set of states.
- equilibrium (list of booleans, optional) – For every trajectory triple (ttraj[i], dtraj[i], btraj[i]), indicates whether to assume global equilibrium. If true, the triple is not used for computing kinetic quantities (but only thermodynamic quantities). By default, no trajectory is assumed to be in global equilibrium. This is the TRAMMBAR extension.
- maxiter (int, optional, default=10000) – The maximum number of self-consistent 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 log-likelihood; 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*’ - connectivity_factor (float, optional, default=1.0) – Only needed if connectivity=’post_hoc_RE’ or ‘BAR_variance’. Values greater than 1.0 weaken the connectivity conditions. For ‘post_hoc_RE’ this multiplies the number of hypothetically observed transitions. For ‘BAR_variance’ this scales the threshold for the minimal allowed variance of free energy differences.
- direct_space (bool, optional, default=False) – Whether to perform the self-consistent iteration with Boltzmann factors (direct space) or free energies (log-space). When analyzing data from multi-temperature simulations, direct-space is not recommended.
- N_dtram_accelerations (int, optional, default=0) – Convergence of TRAM can be speeded up by interleaving the updates in the self-consistent iteration with a dTRAM-like update step. N_dtram_accelerations says how many times the dTRAM-like update step should be applied in every iteration of the TRAM equations. Currently this is only effective if direct_space=True.
- init (str, optional, default=None) –
Use a specific initialization for self-consistent iteration:
None: use a hard-coded guess for free energies and Lagrangian multipliers’mbar’: perform a short MBAR estimate to initialize the free energies - init_maxiter (int, optional, default=5000) – The maximum number of self-consistent iterations during the initialization.
- init_maxerr (float, optional, default=1.0E-8) – Convergence criterion for the initialization.
- overcounting_factor (double, default = 1.0) – Only needed if equilibrium contains True (TRAMMBAR). Sets the relative statistical weight of equilibrium and non-equilibrium frames. An overcounting_factor of value n means that every non-equilibrium frame is counted n times. Values larger than 1 increase the relative weight of the non-equilibrium data. Values less than 1 increase the relative weight of the equilibrium data.
References
[1] Wu, H. et al 2016 Multiensemble Markov models of molecular thermodynamics and kinetics Proc. Natl. Acad. Sci. USA 113 E3221–E3230
Methods
__init__
(lag[, count_mode, connectivity, …])Transition(-based) Reweighting Analysis Method estimate
(X, **params)param X: Simulation trajectories. ttrajs contain the indices of the thermodynamic state, dtrajs expectation
(a)Equilibrium expectation value of a given observable. 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. log_likelihood
()Returns the value of the log-likelihood of the converged TRAM estimate. mbar_pointwise_free_energies
([therm_state])meval
(f, *args, **kw)Evaluates the given function call for all models Returns the results of the calls in a list pointwise_free_energies
([therm_state])Computes the pointwise free energies \(-\log(\mu^k(x))\) for all points x. 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. set_params
(**params)Set the parameters of this estimator. update_model_params
(**params)Update given model parameter if they are set to specific values Attributes
active_set
The active set of states on which all computations and estimations will be done. dt_traj
f
The free energies (in units of kT) on the configuration states. f_full_state
force_constants
The individual force matrices labelled accordingly to ttrajs. free_energies
The free energies (in units of kT) on the configuration states. free_energies_full_state
label
Human-readable description for the thermodynamic state of this model. logger
The logger for this class instance model
The model estimated by this Estimator msm
MSM of the unbiased thermodynamic state; only present when unbiased data available. name
The name of this instance nstates
Number of active states on which all computations and estimations are done. nstates_full
Size of the full set of states. pi
The stationary distribution on the configuration states. pi_full_state
stationary_distribution
The stationary distribution on the configuration states. stationary_distribution_full_state
temperatures
The individual temperatures labelled accordingly to ttrajs. umbrella_centers
The individual umbrella centers labelled accordingly to ttrajs. unbiased_state
Index of the unbiased thermodynamic state. -
active_set
¶ The active set of states on which all computations and estimations will be done.
-
estimate
(X, **params)¶ Parameters: X (tuple of (ttrajs, dtrajs, btrajs)) – Simulation trajectories. ttrajs contain the indices of the thermodynamic state, dtrajs contains the indices of the configurational states and btrajs contain the biases.
- ttrajs : list of numpy.ndarray(X_i, dtype=int)
- Every element is a trajectory (time series). ttrajs[i][t] is the index of the thermodynamic state visited in trajectory i at time step t.
- dtrajs : list of numpy.ndarray(X_i, dtype=int)
- dtrajs[i][t] is the index of the configurational state (Markov state) visited in trajectory i at time step t.
- btrajs : list of numpy.ndarray((X_i, T), dtype=numpy.float64)
- For every simulation frame seen in trajectory i and time step t, btrajs[i][t,k] is the bias energy of that frame evaluated in the k’th thermodynamic state (i.e. at the k’th Umbrella/Hamiltonian/temperature).
-
expectation
(a)¶ Equilibrium expectation value of a given observable.
Parameters: a ((M,) ndarray) – Observable vector Returns: val – Equilibrium expectation value of the given observable Return type: float Notes
The equilibrium expectation value of an observable a is defined as follows
\[\mathbb{E}_{\mu}[a] = \sum_i \mu_i a_i\]\(\mu=(\mu_i)\) is the stationary vector of the transition matrix \(T\).
-
f
¶ The free energies (in units of kT) on the configuration states.
-
fit
(X, y=None)¶ Estimates parameters - for compatibility with sklearn.
Parameters: X (object) – A reference to the data from which the model will be estimated Returns: estimator – The estimator (self) with estimated model. Return type: object
-
force_constants
¶ The individual force matrices labelled accordingly to ttrajs. (only set, when estimated from umbrella data).
-
free_energies
¶ The free energies (in units of kT) on the configuration states.
-
get_model_params
(deep=True)¶ Get parameters for this model.
Parameters: deep (boolean, optional) – If True, will return the parameters for this estimator and contained subobjects that are estimators. Returns: params – Parameter names mapped to their values. Return type: mapping of string to any
-
get_params
(deep=True)¶ Get parameters for this estimator.
Parameters: deep (boolean, optional) – If True, will return the parameters for this estimator and contained subobjects that are estimators. Returns: params – Parameter names mapped to their values. Return type: mapping of string to any
-
label
¶ Human-readable description for the thermodynamic state of this model.
-
classmethod
load
(file_name, model_name='default')¶ Loads a previously saved PyEMMA object from disk.
Parameters: - file_name (str or file like object (has to provide read method)) – The file like object tried to be read for a serialized object.
- model_name (str, default='default') – if multiple models are contained in the file, these can be accessed by
their name. Use
pyemma.list_models()
to get a representation of all stored models.
Returns: obj
Return type: the de-serialized object
-
log_likelihood
()¶ Returns the value of the log-likelihood of the converged TRAM estimate.
-
logger
¶ The logger for this class instance
-
meval
(f, *args, **kw)¶ Evaluates the given function call for all models Returns the results of the calls in a list
-
model
¶ The model estimated by this Estimator
-
msm
¶ MSM of the unbiased thermodynamic state; only present when unbiased data available.
-
name
¶ The name of this instance
-
nstates
¶ Number of active states on which all computations and estimations are done.
-
nstates_full
¶ Size of the full set of states.
-
pi
¶ The stationary distribution on the configuration states.
-
pointwise_free_energies
(therm_state=None)¶ Computes the pointwise free energies \(-\log(\mu^k(x))\) for all points x.
\(\mu^k(x)\) is the optimal estimate of the Boltzmann distribution of the k’th ensemble defined on the set of all samples.
Parameters: therm_state (int or None, default=None) – Selects the thermodynamic state k for which to compute the pointwise free energies. None selects the “unbiased” state which is defined by having zero bias energy. Returns: mu_k – list of the same layout as dtrajs (or ttrajs). mu_k[i][t] contains the pointwise free energy of the frame seen in trajectory i and time step t. Frames that are not in the connected sets get assiged an infinite pointwise free energy. Return type: list of numpy.ndarray(X_i, dtype=numpy.float64)
-
save
(file_name, model_name='default', overwrite=False, save_streaming_chain=False)¶ saves the current state of this object to given file and name.
Parameters: - file_name (str) – path to desired output file
- model_name (str, default='default') – creates a group named ‘model_name’ in the given file, which will contain all of the data. If the name already exists, and overwrite is False (default) will raise a RuntimeError.
- overwrite (bool, default=False) – Should overwrite existing model names?
- save_streaming_chain (boolean, default=False) – if True, the data_producer(s) of this object will also be saved in the given file.
Examples
>>> import pyemma, numpy as np >>> from pyemma.util.contexts import named_temporary_file >>> m = pyemma.msm.MSM(P=np.array([[0.1, 0.9], [0.9, 0.1]]))
>>> with named_temporary_file() as file: # doctest: +SKIP ... m.save(file, 'simple') # doctest: +SKIP ... inst_restored = pyemma.load(file, 'simple') # doctest: +SKIP >>> np.testing.assert_equal(m.P, inst_restored.P) # doctest: +SKIP
-
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’.
-
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 former have parameters of the form
<component>__<parameter>
so that it’s possible to update each component of a nested object. :returns: :rtype: self
-
stationary_distribution
¶ The stationary distribution on the configuration states.
-
temperatures
¶ The individual temperatures labelled accordingly to ttrajs. (only set, when estimated from multi-temperature data).
-
umbrella_centers
¶ The individual umbrella centers labelled accordingly to ttrajs. (only set, when estimated from umbrella data).
-
unbiased_state
¶ Index of the unbiased thermodynamic state.
-
update_model_params
(**params)¶ Update given model parameter if they are set to specific values
-