VelocityAutocorr

VelocityAutocorr — transport_analysis.analysis.VelocityAutocorr

This module offers a class to efficiently calculate a velocity autocorrelation function (VACF). Averaging over all atoms in the atom group ag, regardless of type, it will calculate

\[C(j \Delta t) = {1 \over N - j} \sum_{i=0}^{N-1-j} v(i \Delta t)v((i+j)\Delta t)\]

where \(N\) is the number of time frames and \(\Delta t\) are discrete time intervals between data points.

Basic usage

This example uses the files PRM_NCBOX and TRJ_NCBOX from the MDAnalysis test suite. To get started, execute

>>> import transport_analysis as ta
>>> from MDAnalysis.tests.datafiles import PRM_NCBOX, TRJ_NCBOX

We will calculate the VACF of an atom group of all water atoms in residues 1-5. To select these atoms:

>>> u = mda.Universe(PRM_NCBOX, TRJ_NCBOX)
>>> ag = u.select_atoms("resname WAT and resid 1-5")

We can run the calculation using any variable of choice such as wat_vacf and access our results using wat_vacf.results.timeseries:

>>> wat_vacf = ta.VelocityAutocorr(ag).run()
>>> wat_vacf.results.timeseries
[275.62075467 -18.42008255 -23.94383428  41.41415381  -2.3164344
-35.66393559 -22.66874897  -3.97575003   6.57888933  -5.29065096]

Notice that this example data is insufficient to provide a well-defined VACF. When working with real data, ensure that the frames are captured frequently enough to obtain a VACF suitable for your needs.

class transport_analysis.analysis.velocityautocorr.VelocityAutocorr(atomgroup: AtomGroup, dim_type: str = 'xyz', fft: bool = True, **kwargs)[source]

Class to calculate a velocity autocorrelation function (VACF).

Parameters
  • atomgroup (AtomGroup) – An MDAnalysis AtomGroup. Note that UpdatingAtomGroup instances are not accepted.

  • dim_type ({‘xyz’, ‘xy’, ‘yz’, ‘xz’, ‘x’, ‘y’, ‘z’}) – Desired dimensions to be included in the VACF. Defaults to ‘xyz’.

  • fft (bool) – If True, uses a fast FFT based algorithm for computation of the VACF. Otherwise, use the simple “windowed” algorithm. The tidynamics package is required for fft=True. Defaults to True.

Variables
  • atomgroup (AtomGroup) – The atoms to which this analysis is applied

  • dim_fac (int) – Dimensionality \(d\) of the VACF.

  • results.timeseries (numpy.ndarray) – The averaged VACF over all the particles with respect to lag-time. Obtained after calling VelocityAutocorr.run()

  • results.vacf_by_particle (numpy.ndarray) – The VACF of each individual particle with respect to lag-time.

  • start (Optional[int]) – The first frame of the trajectory used to compute the analysis.

  • stop (Optional[int]) – The frame to stop at for the analysis.

  • step (Optional[int]) – Number of frames to skip between each analyzed frame.

  • n_frames (int) – Number of frames analysed in the trajectory.

  • n_particles (int) – Number of particles VACF was calculated over.