Table Of Contents

Previous topic

3.2.1. Hydrogen Bond analysis — MDAnalysis.analysis.hbonds.hbond_analysis

Next topic

3.3.1. Generation and Analysis of HOLE pore profiles — MDAnalysis.analysis.hole

This Page

3.2.2. Hydrogen bond autocorrelation — MDAnalysis.analysis.hbonds.hbond_autocorrel

Author:Richard J. Gowers
Year:2014
Copyright:GNU Public License v3

New in version 0.9.0.

3.2.2.1. Description

Calculates the time autocorrelation function, \(C_x(t)\), for the hydrogen bonds in the selections passed to it. The population of hydrogen bonds at a given startpoint, \(t_0\), is evaluated based on geometric criteria and then the lifetime of these bonds is monitored over time. Multiple passes through the trajectory are used to build an average of the behaviour.

\(C_x(t) = \left \langle \frac{h_{ij}(t_0) h_{ij}(t_0 + t)}{h_{ij}(t_0)^2} \right\rangle\)

The subscript \(x\) refers to the definition of lifetime being used, either continuous or intermittent. The continuous definition measures the time that a particular hydrogen bond remains continuously attached, whilst the intermittent definition allows a bond to break and then subsequently reform and be counted again. The relevent lifetime, \(\tau_x\), can then be found via integration of this function

\(\tau_x = \int_0^\infty C_x(t) dt\)

For this, the observed behaviour is fitted to a multi exponential function, using 2 exponents for the continuous lifetime and 3 for the intermittent lifetime.

\(C_x(t) = A_1 \exp( - t / \tau_1) + A_2 \exp( - t / \tau_2) [+ A_3 \exp( - t / \tau_3)]\)

Where the final pre expoential factor \(A_n\) is subject to the condition:

\(A_n = 1 - \sum\limits_{i=1}^{n-1} A_i\)

References

[notsure]Multiscale modelling of polymeric systems with hydrogen bonding: Selective removal of degrees of freedom

3.2.2.2. Input

Three AtomGroup selections representing the hydrogens, donors and acceptors that you wish to analyse. Note that the hydrogens and donors selections must be aligned, that is hydrogens[0] and donors[0] must represent a bonded pair. If a single donor therefore has two hydrogens, it must feature twice in the donors AtomGroup.

The keyword exclusions allows a tuple of array addresses to be provided, (Hidx, Aidx),these pairs of hydrogen-acceptor are then not permitted to be counted as part of the analysis. This could be used to exclude the consideration of hydrogen bonds within the same functional group, or to perform analysis on strictly intermolecular hydrogen bonding.

Hydrogen bonds are defined on the basis of geometric criteria; a Hydrogen-Acceptor distance of less then dist_crit and a Donor-Hydrogen-Acceptor angle of greater than angle_crit.

The length of trajectory to analyse in ps, sample_time, is used to choose what length to analyse.

Multiple passes, controlled by the keyword nruns, through the trajectory are performed and an average calculated. For each pass, nsamples number of points along the run are calculated.

3.2.2.3. Output

All results of the analysis are available through the solution attribute. This is a dictionary with the following keys

  • results The raw results of the time autocorrelation function.

  • time Time axis, in ps, for the results.

  • fit Results of the exponential curve fitting procedure. For the

    continuous lifetime these are (A1, tau1, tau2), for the intermittent lifetime these are (A1, A2, tau1, tau2, tau3).

  • tau Calculated time constant from the fit.

  • estimate Estimated values generated by the calculated fit.

The results and time values are only filled after the run() method, fit, tau and estimate are filled after the solve() method has been used.

3.2.2.4. Examples

from MDAnalysis.analysis import hbonds
import matplotlib.pyplot as plt
H = u.selectAtoms('name Hn')
O = u.selectAtoms('name O')
N = u.selectAtoms('name N')
hb_ac = hbonds.HydrogenBondAutoCorrel(u, acceptors = u.atoms.O,
            hydrogens = u.atoms.Hn, donors = u.atoms.N,bond_type='continuous',
            sample_time = 2, nruns = 20, nsamples = 1000)
hb_ac.run()
hb_ac.solve()
tau = hb_ac.solution['tau']
time = hb_ac.solution['time']
results = hb_ac.solution['results']
estimate = hb_ac.solution['estimate']
plt.plot(time, results, 'ro')
plt.plot(time, estimate)
plt.show()
class MDAnalysis.analysis.hbonds.hbond_autocorrel.HydrogenBondAutoCorrel(universe, hydrogens=None, acceptors=None, donors=None, bond_type=None, exclusions=None, angle_crit=130.0, dist_crit=3.0, sample_time=100, time_cut=None, nruns=1, nsamples=50, pbc=True)[source]

Perform a time autocorrelation of the hydrogen bonds in the system.

Arguments:
universe

The MDA universe

hydrogens

AtomGroup of Hydrogens which can form hydrogen bonds

acceptors

AtomGroup of all Acceptor atoms

donors

The atoms which are connected to the hydrogens. This group must be identical in length to the hydrogen group and matched, ie hydrogens[1] is bonded to donors[0]. For many cases, this will mean a donor appears twice in this group.

bond_type

Which definition of hydrogen bond lifetime to consider, either ‘continuous’ or ‘intermittent’.

Keywords:
exclusions

Indices of Hydrogen-Acceptor pairs to be excluded. With nH and nA Hydrogens and Acceptors, a (nH x nA) array of distances is calculated, exclusions is used as a mask on this array to exclude some pairs.

angle_crit

The angle (in degrees) which all bonds must be greater than [130.0]

dist_crit

The maximum distance (in Angstroms) for a hydrogen bond [3.0]

sample_time

The amount of time, in ps, that you wish to observe hydrogen bonds for [100]

nruns

The number of different start points within the trajectory to use [1]

nsamples

Within each run, the number of frames to analyse [50]

pbc

Whether to consider periodic boundaries in calculations [True]

run(force=False)[source]

Run all the required passes

Keywords:
force

Will overwrite previous results if they exist

solve(p_guess=None)[source]

Fit results to an multi exponential decay and integrate to find characteristic time

Keywords:
p_guess

Initial guess for the leastsq fit, must match the shape of the expected coefficients

Continuous defition results are fitted to a double exponential, intermittent definition are fit to a triple exponential.

The results of this fitting procedure are saved into the fit, tau and estimate keywords in the solution dict.

  • fit contains the coefficients, (A1, tau1, tau2) or (A1, A2, tau1, tau2, tau3)
  • tau contains the calculated lifetime in ps for the hydrogen bonding
  • estimate contains the estimate provided by the fit of the time autocorrelation function

In addition, the output of the leastsq function is saved into the solution dict

  • infodict
  • mesg
  • ier
save_results(filename='hbond_autocorrel')[source]

Saves the results to a numpy zipped array (.npz, see numpy.savez)

This can be loaded using numpy.load(filename)

Keywords:
filename

The desired filename [hbond_autocorrel]