Lipid order parameter — lipyphilic.lib.order_parameter

Author

Paul Smith

Year

2021

Copyright

GNU Public License v2

This module provides methods for calculating the orientational order parameter of lipid tails in a bilayer.

Coarse-grained order parameter

The class liyphilic.lib.order_parameter.SCC calculates the coarse-grained order parameter, as defined in Seo et al. (2020). The coarse-grained order parameter, \(S_{CC}\), is defined as:

\[S_{CC} = \displaystyle \frac{\big \langle 3 \cos^2 \theta - 1 \big \rangle}{2}\]

where \(\theta\) is the angle between the membrane normal and the vector connecting two consecutive tail beads. Angular brackets denote averages over all beads in an acyl tail.

See Piggot et al. (2017) for an excellent discussion on calculating acyl tali order parameters in molecular dynamics simulations.

Input

Required:
  • universe : an MDAnalysis Universe object

  • tail_sel : atom selection for beads in the acyl tail

Options:
  • normals : local membrane normals for each tail at each frame

Output

  • SCC : order parameter of each tail at each frame

The order parameter data are returned in a numpy.ndarray, where each row corresponds to an individual lipid and each column corresponds to an individual frame.

Warning

tail_sel should select beads in either the sn1 or sn2 tails, not both tails.

Example usage of Scc

An MDAnalysis Universe must first be created before using SCC:

import MDAnalysis as mda
from lipyphilic.lib.order_parameter import SCC

u = mda.Universe(tpr, trajectory)

If we have used the MARTINI forcefield to study a DPPC/DOPC/cholesterol mixture, we can calculate the order parameter of the sn1 of tails of DPPC and DOPC as follows:

scc_sn1 = SCC(
  universe=u,
  tail_sel="name ??A"  # selects C1A, C2A, D2A, C3A, and C4A
)

This will calculate \(S_{CC}\) of each DOPC and DPPC sn1 tail.

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

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

The results are then available in the scc_sn1.SCC attribute as a numpy.ndarray. The array has the shape (n_residues, n_frames). Each row corresponds to an individual lipid and each column to an individual frame.

Likewise, to calculate the \(S_{CC}\) of the sn2 tails, we can do:

scc_sn2 = SCC(
  universe=u,
  tail_sel="name ??B"  # selects C1B, C2B, D2B, C3B, and C4B
)
scc_sn2.run(verbose=True)

And then get a weighted-average \(S_{CC}\) we can do:

SCC.weighted_average(scc_sn1, scc_sn2)

which will take into account the number of beads in each tail and return a new SCC object whose SCC attribute contains the weighted-average \(S_{CC}\) for each lipid at each frame.

Local membrane normals

By default, the \(S_{CC}\) is calculated as the angle between the positive \(z\) axis and the vector between two consecutive beads in an acyl tail. However, it is also possible to pass to SCC local membrane normals to use instead of the positive \(z\) axis.

You can also calculate local membrane normals using, for example, MemSurfer. If you store the local membrane normals in a numpy.ndarray called normals, with shape (n_residues, n_frames, 3), then you can simply pass these normals to SCC:

scc_sn1 = SCC(
  universe=u,
  tail_sel="name ??A",
  normals=normals
)
scc_sn1.run(verbose=True)

\(S_{CC}\) projected onto the membrane plane

Once the \(S_{CC}\) has been calculated, it is possible to create a 2D plot of time-averaged \(S_{CC}\) values projected onto the membrane plane. This can be done using the liypphilic.lib.SCC.project_SCC() method, which is a wrapper around the more general liypphilic.lib.plotting.ProjectionPlot class.

If the lipids have been assigned to leaflets, and the weighted average of the sn1 and sn2 tails stored in an SCC object named scc, we can plot the projection of the coarse-grained order parameter onto the membrane plane as follows:

scc_projection = scc.project_SCC(
  lipid_sel="name ??A ??B",
  start=-100,
  stop=None,
  step=None,
  filter_by=leaflets.filter_by("name ??A ??B") == -1
)

The order parameter of each lipid parameter will be averaged over the final 100 frames, as specified by the start argument. The frame in the middle of the selected frames will be used for determining lipid positions. In the above case, the lipid positions at frame \(-50\) will be used. The lipid_sel specifies that the center of mass of the sn1 (“??A”) and sn2 (“??B”) atoms will be used for projecting lipid positions onto the membrane plane. And the filter_by argument is used here to specificy that only lipids in the lower (-1) leaflet should be used for plotting the projected \(S_{CC}\) values.

The class and its methods

class lipyphilic.lib.order_parameter.SCC(universe, tail_sel, normals=None)[source]

Calculate coarse-grained acyl tail order parameter.

Set up parameters for calculating the SCC.

Parameters
  • universe (Universe) – MDAnalysis Universe object

  • tail_sel (str) – Selection string for atoms in either the sn1 or sn2 tail of lipids in the membrane

  • normals (numpy.ndarray, optional) – Local membrane normals, a 3D array of shape (n_residues, n_frames, 3), containing x, y and z vector components of the local membrane normals.

project_SCC(lipid_sel=None, start=None, stop=None, step=None, filter_by=None, unwrap=True, bins=None, ax=None, cmap=None, vmin=None, vmax=None, cbar=True, cbar_kws={}, imshow_kws={})[source]

Project the SCC values onto the xy plane of the membrane.

The SCC values, averaged over a selected range of frames, are projected onto the xy plane based on the center of mass of each lipid. The atoms to be used in calculating the center of mass of the lipids can be specified using the lipid_sel arugment.

This method creates an instance of lipyphilic.lib.plotting.ProjectionPlot with the projected \(S_{CC}\) interpolated across periodic boundaries. The plot is returned so further modification can be performed if needed.

Note

The lipid positions are taken from the middle frame of the selected range.

Parameters
  • lipid_sel (MDAnalysis atom selection, optional) – The center of mass of each lipid will be determined via this selection. The default is None, in which case every atom of a lipid is used to determine its center of mass.

  • start (int, optional) – Start frame for averaging the SCC results.

  • stop (int, optional) – Final frame for averaging the SCC results.

  • step (int, optional) – Number of frames to skip

  • filter_by (array-like, optional) – A Boolean mask for selecting a subset of lipids. It may take the following shapes:

    (n_lipids) The mask is used to select a subset of lipids for projecting the SCC onto the membrane plane.

    (n_lipids, n_frames) This is the same shape as the NumPy array created by the lipyphilic.lib.SCC.run() method. Boolean values are used only from the column corresponding to the middle frame of the range selected by start, stop, and step.

    The default is None, in which case no filtering is applied.

  • unwrap (bool, optional) – If True, lipids will be unwrapped before computing their center of mass, which is.

  • bins (int or array_like or [int, int] or [array, array]) – The bin specification:

    int

    If int, the number of bins for the two dimensions (nx=ny=bins).

    array-like

    If array_like, the bin edges for the two dimensions (x_edges=y_edges=bins).

    [int, int]

    If [int, int], the number of bins in each dimension (nx, ny = bins).

    [array, array]

    If [array, array], the bin edges in each dimension (x_edges, y_edges = bins).

    combination

    A combination [int, array] or [array, int], where int is the number of bins and array is the bin edges.

    The default is None, in which case a grid with 1 x 1 Angstrom resolution is created.

  • ax (Axes, optional) – Matplotlib Axes on which to plot the projection. The default is None, in which case a new figure and axes will be created.

  • cmap (str or ~matplotlib.colors.Colormap, optional) – The Colormap instance or registered colormap name used to map scalar data to colors.

  • vmin, vmax (float, optional) – Define the data range that the colormap covers. By default, the colormap covers the complete value range of the supplied data.

  • cbar (bool, optional) – Whether or not to add a colorbar to the plot.

  • cbar_kws (dict, optional) – A dictionary of keyword options to pass to matplotlib.pyplot.colorbar.

  • imshow_kws (dict, optional) – A dictionary of keyword options to pass to matplotlib.pyplot.imshow, which is used to plot the 2D density map.

Returns

scc_projection (ProjectionPlot) – The ProjectionPlot object containing the SCC data and the matplotlob.pyplot.imshow plot of the projection.

static weighted_average(sn1_scc, sn2_scc)[source]

Calculate the weighted average Scc of two tails.

Given two SCC objects, a weighted average of the Scc of each lipid is calculated.

Parameters
  • sn1_scc (SCC) – An SCC object for which the order parameters have been calculated.

  • sn2_scc (SCC) – An SCC object for which the order parameters have been calculated.

Returns

scc (SCC) – An SCC object with the weighted average Scc of each lipid at each frame stored in the scc.SCC attirbute

Warning

The frames used in analysing ‘sn1_scc’ and ‘sn2_scc’ must be the same - i.e. the ‘start’, ‘stop’, and ‘step’ parameters passed to the ‘.run()’ methods must be identical.