pyemma.thermo.DTRAM

class pyemma.thermo.DTRAM(bias_energies_full, lag, count_mode='sliding', connectivity='reversible_pathways', maxiter=10000, maxerr=1e-15, save_convergence_info=0, dt_traj='1 step', init=None, init_maxiter=10000, init_maxerr=1e-08)

Discrete Transition(-based) Reweighting Analysis Method.

__init__(bias_energies_full, lag, count_mode='sliding', connectivity='reversible_pathways', maxiter=10000, maxerr=1e-15, save_convergence_info=0, dt_traj='1 step', init=None, init_maxiter=10000, init_maxerr=1e-08)

Discrete Transition(-based) Reweighting Analysis Method

Parameters:
  • bias_energies_full (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.
  • 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.

  • connectivity (str, optional, default='reversible_pathways') –

    One of ‘reversible_pathways’, ‘summed_count_matrix’ or None. 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.
    • ‘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!
    • None : assume that everything is connected. For debugging.

    For more details see thermotools.cset.compute_csets_dTRAM().

  • maxiter (int, optional, default=10000) – The maximum number of self-consistent iterations before the estimator exits unsuccessfully.
  • maxerr (float, optional, default=1.0E-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*’
  • 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
    ‘wham’: perform a short WHAM estimate to initialize the free energies
  • init_maxiter (int, optional, default=10000) – The maximum number of self-consistent iterations during the initialization.
  • init_maxerr (float, optional, default=1.0E-8) – Convergence criterion for the initialization.

Example

>>> from pyemma.thermo import DTRAM
>>> import numpy as np
>>> B = np.array([[0, 0],[0.5, 1.0]])
>>> dtram = DTRAM(B, 1)
>>> 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])]
>>> dtram = dtram.estimate((ttrajs, dtrajs))
>>> dtram.log_likelihood() 
-9.805...
>>> dtram.count_matrices 
array([[[5, 1],
        [1, 2]],
[[1, 4],
[3, 1]]], dtype=int32)
>>> dtram.stationary_distribution 
array([ 0.38...,  0.61...])
>>> dtram.meval('stationary_distribution') 
[array([ 0.38...,  0.61...]), array([ 0.50...,  0.49...])]

References

[1]Wu, H. et al 2014 Statistically optimal analysis of state-discretized trajectory data from multiple thermodynamic states J. Chem. Phys. 141, 214106

Methods

__init__(bias_energies_full, lag[, ...]) Discrete Transition(-based) Reweighting Analysis Method
estimate(trajs)
param X:Simulation trajectories. ttrajs contain the indices of the thermodynamic state and
expectation(a) Equilibrium expectation value of a given observable.
fit(X) Estimates parameters - for compatibility with sklearn.
get_model_params([deep]) Get parameters for this model.
get_params([deep]) Get parameters for this estimator.
log_likelihood()
meval(f, *args, **kw) Evaluates the given function call for all models
register_progress_callback(call_back[, stage]) Registers the progress reporter.
set_model_params([models, f_therm, pi, f, label])
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.
f The free energies (in units of kT) on the configuration states.
f_full_state
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
model_active_set
msm MSM of the unbiased thermodynamic state; only present when unbiased data available.
msm_active_set
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
show_progress whether to show the progress of heavy calculations on this object.
stationary_distribution The stationary distribution on the configuration states.
stationary_distribution_full_state
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(trajs)
Parameters:X (tuple of (ttrajs, dtrajs)) –

Simulation trajectories. ttrajs contain the indices of the thermodynamic state and dtrajs contains the indices of the configurational states.

ttrajs : list of numpy.ndarray(X_i, dtype=int)
Every elements 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.
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)

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
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.

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.

register_progress_callback(call_back, stage=0)

Registers the progress reporter.

Parameters:
  • call_back (function) –

    This function will be called with the following arguments:

    1. stage (int)
    2. instance of pyemma.utils.progressbar.ProgressBar
    3. optional *args and named keywords (**kw), for future changes
  • stage (int, optional, default=0) – The stage you want the given call back function to be fired.
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

show_progress

whether to show the progress of heavy calculations on this object.

stationary_distribution

The stationary distribution on the configuration states.

unbiased_state

Index of the unbiased thermodynamic state.

update_model_params(**params)

Update given model parameter if they are set to specific values