pyemma.coordinates.tica¶
-
pyemma.coordinates.
tica
(data=None, lag=10, dim=-1, var_cutoff=0.95, kinetic_map=True, stride=1, force_eigenvalues_le_one=False, mean=None)¶ Time-lagged independent component analysis (TICA).
TICA is a linear transformation method. In contrast to PCA, which finds coordinates of maximal variance, TICA finds coordinates of maximal autocorrelation at the given lag time. Therefore, TICA is useful in order to find the slow components in a dataset and thus an excellent choice to transform molecular dynamics data before clustering data for the construction of a Markov model. When the input data is the result of a Markov process (such as thermostatted molecular dynamics), TICA finds in fact an approximation to the eigenfunctions and eigenvalues of the underlying Markov operator [1].
It estimates a TICA transformation from data. When input data is given as an argument, the estimation will be carried out straight away, and the resulting object can be used to obtain eigenvalues, eigenvectors or project input data onto the slowest TICA components. If no data is given, this object is an empty estimator and can be put into a
pipeline()
in order to use TICA in the streaming mode.Parameters: - data (ndarray (T, d) or list of ndarray (T_i, d) or a reader created by) – source function array with the data, if available. When given, the TICA transformation is immediately computed and can be used to transform data.
- lag (int, optional, default = 10) – the lag time, in multiples of the input time step
- dim (int, optional, default -1) – the number of dimensions (independent components) to project onto. A
call to the
map
function reduces the d-dimensional input to only dim dimensions such that the data preserves the maximum possible autocorrelation amongst dim-dimensional linear projections. -1 means all numerically available dimensions will be used unless reduced by var_cutoff. Setting dim to a positive value is exclusive with var_cutoff. - var_cutoff (float in the range [0,1], optional, default 0.95) – Determines the number of output dimensions by including dimensions until their cumulative kinetic variance exceeds the fraction subspace_variance. var_cutoff=1.0 means all numerically available dimensions (see epsilon) will be used, unless set by dim. Setting var_cutoff smaller than 1.0 is exclusive with dim
- kinetic_map (bool, optional, default False) – Eigenvectors will be scaled by eigenvalues. As a result, Euclidean distances in the transformed data approximate kinetic distances [4]. This is a good choice when the data is further processed by clustering.
- stride (int, optional, default = 1) – If set to 1, all input data will be used for estimation. Note that this could cause this calculation to be very slow for large data sets. Since molecular dynamics data is usually correlated at short timescales, it is often sufficient to estimate transformations at a longer stride. Note that the stride option in the get_output() function of the returned object is independent, so you can parametrize at a long stride, and still map all frames through the transformer.
- force_eigenvalues_le_one (boolean) – Compute covariance matrix and time-lagged covariance matrix such that the generalized eigenvalues are always guaranteed to be <= 1.
- mean (ndarray, optional, default None) – Optionally pass pre-calculated means to avoid their re-computation. The shape has to match the input dimension.
Returns: tica – Object for time-lagged independent component (TICA) analysis. it contains TICA eigenvalues and eigenvectors, and the projection of input data to the dominant TICA
Return type: a
TICA
transformation objectNotes
Given a sequence of multivariate data \(X_t\), it computes the mean-free covariance and time-lagged covariance matrix:
\[\begin{split}C_0 &= (X_t - \mu)^T (X_t - \mu) \\ C_{\tau} &= (X_t - \mu)^T (X_t + \tau - \mu)\end{split}\]and solves the eigenvalue problem
\[C_{\tau} r_i = C_0 \lambda_i r_i,\]where \(r_i\) are the independent components and \(\lambda_i\) are their respective normalized time-autocorrelations. The eigenvalues are related to the relaxation timescale by
\[t_i = -\frac{\tau}{\ln |\lambda_i|}.\]When used as a dimension reduction method, the input data is projected onto the dominant independent components.
TICA was originally introduced for signal processing in [2]. It was introduced to molecular dynamics and as a method for the construction of Markov models in [1] and [3]. It was shown in [1] that when applied to molecular dynamics data, TICA is an approximation to the eigenvalues and eigenvectors of the true underlying dynamics.
Examples
Invoke TICA transformation with a given lag time and output dimension:
>>> import numpy as np >>> from pyemma.coordinates import tica >>> data = np.random.random((100,3)) >>> projected_data = tica(data, lag=2, dim=1).get_output()[0]
For a brief explaination why TICA outperforms PCA to extract a good reaction coordinate have a look here.
-
class
pyemma.coordinates.transform.tica.
TICA
(lag, dim=-1, var_cutoff=0.95, kinetic_map=True, epsilon=1e-06, force_eigenvalues_le_one=False, mean=None)¶ Time-lagged independent component analysis (TICA)
Methods
describe
(*args, **kwargs)Get a descriptive string representation of this class. dimension
()output dimension fit
(X, **kwargs)For compatibility with sklearn fit_transform
(X, **kwargs)For compatibility with sklearn get_output
([dimensions, stride])Maps all input data of this transformer and returns it as an array or list of arrays. iterator
([stride, lag])Returns an iterator that allows to access the transformed data. map
(X)Deprecated: use transform(X) n_frames_total
([stride])Returns total number of frames. number_of_trajectories
()Returns the number of trajectories. output_type
()By default transformers return single precision floats. parametrize
([stride])Parametrize this Transformer register_progress_callback
(call_back[, stage])Registers the progress reporter. trajectory_length
(itraj[, stride])Returns the length of trajectory of the requested index. trajectory_lengths
([stride])Returns the length of each trajectory. transform
(X)Maps the input data through the transformer to correspondingly shaped output data array/list. Attributes
chunksize
chunksize defines how much data is being processed at once. cumvar
Cumulative sum of the the TICA eigenvalues data_producer
where the transformer obtains its data. eigenvalues
Eigenvalues of the TICA problem (usually denoted \(\lambda\) eigenvectors
Eigenvectors of the TICA problem, columnwise feature_TIC_correlation
Instantaneous correlation matrix between input features and TICs in_memory
are results stored in memory? lag
lag time of correlation matrix \(C_{ au}\) logger
The logger for this class instance mean
mean of input features name
The name of this instance ntraj
timescales
Implied timescales of the TICA transformation -
chunksize
¶ chunksize defines how much data is being processed at once.
-
cumvar
¶ Cumulative sum of the the TICA eigenvalues
Returns: cumvar Return type: 1D np.array
-
data_producer
¶ where the transformer obtains its data.
-
describe
(*args, **kwargs)¶ Get a descriptive string representation of this class.
-
dimension
()¶ output dimension
-
eigenvalues
¶ Eigenvalues of the TICA problem (usually denoted \(\lambda\)
Returns: eigenvalues Return type: 1D np.array
-
eigenvectors
¶ Eigenvectors of the TICA problem, columnwise
Returns: eigenvectors Return type: (N,M) ndarray
-
feature_TIC_correlation
¶ Instantaneous correlation matrix between input features and TICs
Denoting the input features as \(X_i\) and the TICs as \(\theta_j\), the instantaneous, linear correlation between them can be written as
\[\mathbf{Corr}(X_i, \mathbf{\theta}_j) = \frac{1}{\sigma_{X_i}}\sum_l \sigma_{X_iX_l} \mathbf{U}_{li}\]The matrix \(\mathbf{U}\) is the matrix containing, as column vectors, the eigenvectors of the TICA generalized eigenvalue problem .
Returns: feature_TIC_correlation – correlation matrix between input features and TICs. There is a row for each feature and a column for each TIC. Return type: ndarray(n,m)
-
fit
(X, **kwargs)¶ For compatibility with sklearn
-
fit_transform
(X, **kwargs)¶ For compatibility with sklearn
-
get_output
(dimensions=slice(0, None, None), stride=1)¶ Maps all input data of this transformer and returns it as an array or list of arrays.
Parameters: - dimensions (list-like of indexes or slice) – indices of dimensions you like to keep, default = all
- stride (int) – only take every n’th frame, default = 1
Returns: output – the mapped data, where T is the number of time steps of the input data, or if stride > 1, floor(T_in / stride). d is the output dimension of this transformer. If the input consists of a list of trajectories, Y will also be a corresponding list of trajectories
Return type: ndarray(T, d) or list of ndarray(T_i, d)
Notes
- This function may be RAM intensive if stride is too large or too many dimensions are selected.
- if in_memory attribute is True, then results of this methods are cached.
Example
plotting trajectories
>>> import pyemma.coordinates as coor >>> import matplotlib.pyplot as plt
Fill with some actual data!
>>> tica = coor.tica() >>> trajs = tica.get_output(dimensions=(0,), stride=100) >>> for traj in trajs: ... plt.figure() ... plt.plot(traj[:, 0])
-
in_memory
¶ are results stored in memory?
-
iterator
(stride=1, lag=0)¶ Returns an iterator that allows to access the transformed data.
Parameters: - stride (int) – Only transform every N’th frame, default = 1
- lag (int) – Configure the iterator such that it will return time-lagged data with a lag time of lag. If lag is used together with stride the operation will work as if the striding operation is applied before the time-lagged trajectory is shifted by lag steps. Therefore the effective lag time will be stride*lag.
Returns: iterator – If lag = 0, a call to the .next() method of this iterator will return the pair (itraj, X) : (int, ndarray(n, m)), where itraj corresponds to input sequence number (eg. trajectory index) and X is the transformed data, n = chunksize or n < chunksize at end of input.
If lag > 0, a call to the .next() method of this iterator will return the tuple (itraj, X, Y) : (int, ndarray(n, m), ndarray(p, m)) where itraj and X are the same as above and Y contain the time-lagged data.
Return type: a
TransformerIterator
-
lag
¶ lag time of correlation matrix \(C_{ au}\)
-
logger
¶ The logger for this class instance
-
map
(X)¶ Deprecated: use transform(X)
Maps the input data through the transformer to correspondingly shaped output data array/list.
-
mean
¶ mean of input features
-
n_frames_total
(stride=1)¶ Returns total number of frames.
Parameters: stride (int) – return value is the number of frames in trajectories when running through them with a step size of stride. Returns: int Return type: n_frames_total
-
name
¶ The name of this instance
-
ntraj
¶
-
number_of_trajectories
()¶ Returns the number of trajectories.
Returns: int Return type: number of trajectories
-
output_type
()¶ By default transformers return single precision floats.
-
parametrize
(stride=1)¶ Parametrize this Transformer
-
register_progress_callback
(call_back, stage=0)¶ Registers the progress reporter.
Parameters: - call_back (function) –
This function will be called with the following arguments:
- stage (int)
- instance of pyemma.utils.progressbar.ProgressBar
- 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.
- call_back (function) –
-
timescales
¶ Implied timescales of the TICA transformation
For each \(i\)-th eigenvalue, this returns
\[t_i = -\frac{\tau}{\log(|\lambda_i|)}\]where \(\tau\) is the
lag
of the TICA object and \(\lambda_i\) is the i-theigenvalue
of the TICA object.Returns: timescales – numpy array with the implied timescales. In principle, one should expect as many timescales as input coordinates were available. However, less eigenvalues will be returned if the TICA matrices were not full rank or var_cutoff
was parsedReturn type: 1D np.array
-
trajectory_length
(itraj, stride=1)¶ Returns the length of trajectory of the requested index.
Parameters: - itraj (int) – trajectory index
- stride (int) – return value is the number of frames in the trajectory when running through it with a step size of stride.
Returns: int
Return type: length of trajectory
-
trajectory_lengths
(stride=1)¶ Returns the length of each trajectory.
Parameters: stride (int) – return value is the number of frames of the trajectories when running through them with a step size of stride. Returns: array(dtype=int) Return type: containing length of each trajectory
-
transform
(X)¶ Maps the input data through the transformer to correspondingly shaped output data array/list.
Parameters: X (ndarray(T, n) or list of ndarray(T_i, n)) – The input data, where T is the number of time steps and n is the number of dimensions. If a list is provided, the number of time steps is allowed to vary, but the number of dimensions are required to be to be consistent. Returns: Y – The mapped data, where T is the number of time steps of the input data and d is the output dimension of this transformer. If called with a list of trajectories, Y will also be a corresponding list of trajectories Return type: ndarray(T, d) or list of ndarray(T_i, d)
-
References
[1] (1, 2, 3) Perez-Hernandez G, F Paul, T Giorgino, G De Fabritiis and F Noe. 2013. Identification of slow molecular order parameters for Markov model construction J. Chem. Phys. 139, 015102. doi:10.1063/1.4811489 [2] L. Molgedey and H. G. Schuster. 1994. Separation of a mixture of independent signals using time delayed correlations Phys. Rev. Lett. 72, 3634. [3] Schwantes C, V S Pande. 2013. Improvements in Markov State Model Construction Reveal Many Non-Native Interactions in the Folding of NTL9 J. Chem. Theory. Comput. 9, 2000-2009. doi:10.1021/ct300878a [4] Noe, F. and C. Clementi. 2015. Kinetic distance and kinetic maps from molecular dynamics simulation (in preparation).