Lateral diffusion — lipyphilic.lib.lateral_diffusion

This module contains methods for calculating the lateral diffusion coefficient of lipids in a bilayer.

The class lipyphilic.lib.lateral_diffusion.MSD calculates the two-dimensional mean squared displacent (MSD) of lipids in a bilayer. The Fast Correlation Algorithm, implemented in tidynamics is used to calculate the MSD of each lipid, with optional removal of the center of mass motion of the bilayer.

lipyphilic.lib.lateral_diffusion.MSD also contains a method for calculating the lateral diffusion coefficient, \(D_{xy}\), via the Einstein relation:

\[D_{xy} = \frac{1}{4} \lim_{t\to\infty} \frac{d}{dt} \displaystyle \Bigg\langle \frac{1}{N} \sum_{i=1}^{N} \left | r_i(t_0 + \Delta t) - r_i(t_0) \right |^2 \displaystyle \Bigg\rangle_{t_0}\]

where \(N\) is the number of lipids, \(r_i(t0)\) is the center of mass in \(xy\) of lipids \(i\) at a time origin t_0, \(r_i(t0 + \Delta t)\) is the same lipid’s center of mass at a lagtime \(\Delta t\), and the angular brackets denote an average over all time origins, \(t_0\).

Typically, the MSD is averaged over all molecules. However, lipyphilic.lib.lateral_diffusion.MSD will return the MSD for each individual lipid. This makes it simple to later calculate the diffusion coefficient using a subset of the lipids, such as a specific lipid species or lipids near a protein.

Input

Required:
  • universe : an MDAnalysis Universe object

  • lipid_sel : atom selection for calculating the MSD

Optional:
  • com_removal_sel : atom selection for center of mass removal from the MSD

  • dt : time period betwen consecutive frames in the MSD analysis

Output

  • msd : the mean squared displacement of each lipid at each lagtime, \(\Delta t\), in \(nm^2\):

  • lagtimes : a NumPy array of lagtimes (in \(ns\))

The data are stored in the MSD.msd and MSD.lagtimes attributes.

Warning

Before using lipyphilic.lib.lateral_diffusion.MSD you must ensure that the coordinates have been unwrapped using, for example, lipyphilic.transformations.nojump.

Example usage of MSD

To calculate the MSD of each lipid in a bilayer we must first load a trajectory using MDAnalysis:

import MDAnalysis as mda
from lipyphilic.lib.lateral_diffusion import MSD

u = mda.Universe(tpr, trajectory)

If we have used the MARTINI forcefield to study a phospholipid/cholesterol mixture, we can calculate the MSD of each lipid as follows:

msd = MSD(
  universe=u,
  lipid_sel="resname DPPC DOPC CHOL"
)

We then select which frames of the trajectory to analyse (None will use every frame) and choose to display a progress bar (verbose=True):

msd.run(
  start=None,
  stop=None,
  step=None,
  verbose=True
)

The results are then available in the msd.MSD attribute as a numpy.ndarray. Each row corresponds to an individual lipid and each column to a different lagtime`.

Center of mass removal

During your simulation, it is likely that you removed the center of mass motion of your bilayer in the \(z\) direction. However, it is not possible to remove the \(x\) and \(y\) center of mass motions until you have unwrapped your lipid positions.

You may select which lipids to use for the center of mass motion removal using the com_removal_sel keyword:

msd = MSD(
  universe=u,
  lipid_sel="resname DPPC",
  com_removla_sel="resname DPPC DOPC CHOL"
)

In this case, the MSD of DPPC will be calculated with the center of mass motion of the bilayer will be subtracted from it.

Plotting the MSD of each species

If you have calculated the MSD of DPPC, DOPC and cholesterol as in the first example, you can plot the MSD of each species as follows:

for species in ["DPPC", "DOPC", "CHOL"]:

  plt.loglog(
    msd.lagtimes,
    np.mean(msd.msd[msd.membrane.residues.resnames == species], axis=0),
    label=species
  )

plt.legend()

The linear part of the log-log plot can be used for fitting a line and calculating the diffusion coefficient.

Calculating the lateral diffusion coefficient

After calculating the MSD and identifying the linear portion of the plot, the :func: diffusion_coefficient method of lipyphilic.lib.lateral_diffusion.MSD can be used to calculate \(D_{xy}\). We need to pass the time at which to start and stop the linear fit:

d, sem = msd.diffusion_coefficient(
  start_fit=400,
  end_fit=600
)

This will calculate a diffusion coefficient for each individual lipid and return the mean and standard error of the distribution of coefficients.

To calculate the diffusion coefficient of a subset of lipids we can use the :attr:lipid_sel keyword:

d, sem = msd.diffusion_coefficient(
  start_fit=400,
  end_fit=600,
  lipid_sel="resname CHOL"
)

which will calculate the lateral diffusion coefficient for cholesterol, using a fit to the MSD curve from lagtime \(\Delta t = 400\) to lagtime \(\Delta t = 600\).

class lipyphilic.lib.lateral_diffusion.MSD(universe, lipid_sel, com_removal_sel=None, dt=None)[source]

Calculate the mean-squared lateral displacement of lipids in a bilayer.

The MSD is returned in units of \(nm^2/ns\).

Parameters
  • universe (Universe) – MDAnalysis Universe object

  • lipid_sel (str) – Selection string for calculating the mean-squared displacemnt. IF multiple atoms per lipid are selected, the center-of-mass of these atoms will be used for calculating the MSD.

  • com_removal_sel (str, optional) – The MSD of the center of mass of atoms in this selection will be subtracted from all individual lipid MSDs. The default is None, in which case no center of mass motion removal is performed.

  • dt (float, optional) – The time, in nanoseconds, between consecutive frames in universe.trajectory. The defualt is None, in which case dt is taken to be universe.trajectory.dt divided by 1000.

diffusion_coefficient(start_fit=None, stop_fit=None, lipid_sel=None)[source]

Calculate the lateral diffusion coefficient via the Einstein relation.

A diffusion is calculated for each lipid through a linear fit to its MSD curve. The mean and standard error of the diffusion coefficient is returned.

Parameters
  • start_fit (float, optional) – The time at which to start the linear fit to the MSD curve. The default is None, in which case the fit will exclude the first 20% of the MSD data.

  • stop_fit (float, optional) – The time at which to stop the linear fit to the MSD curve. The default is None, in which case the fit will exclude the final 20% of the MSD data.

  • lipid_sel (str, optional) – Selection string for lipids to include in calculating the diffusion coefficient.

Returns

  • d (float) – The mean lateral diffusion coefficient, in \(cm^2/s\)., averaged over all lipids in lipid_sel.

  • sem (float) – The standard error of the diffusion coefficients.