Assign leaflets — lipyphilic.lib.assign_leaflets

Author

Paul Smith

Year

2021

Copyright

GNU Public License v2

This module provides methods for assigning lipids to leaflets in a bilayer.

Assigning leaflets in planar bilayers

The class lipyphilic.lib.assign_leaflets.AssignLeaflets assigns each lipid to a leaflet based on the distance in z to the midpoint of the bilayer. Lipids may be assigned to the upper leaflet (indicated by 1), the lower leaflet (-1) or the bilayer midplane (0).

Input

Required:
  • universe : an MDAnalysis Universe object

  • lipid_sel : atom selection for all lipids in the bilayer, including e.g. sterols

Options:
  • midplane_sel : atom selection for lipid that may occupy the midplane

  • midplane_cutoff : atoms within this distance from the midpoint are considered to be the midplane

  • n_bins : split the membrane into n_bins * n_bins patches, and calculate local membrane midpoints for each patch

Output

  • leaflets : leaflet to which each lipid is assigned at each frame

Leaflet data are returned in a numpy.ndarray, where each row corresponds to an individual lipid and each column corresponds to an individual frame, i.e. leaflets[i, j] refers to the leaflet of lipid i at frame j. The results are accessible via the AssignLeaflets.leaflets attribute.

Example usage of AssignLeaflets

An MDAnalysis Universe must first be created before using AssignLeaflets:

import MDAnalysis as mda
from lipyphilic.lib.assign_leaflets import AssignLeaflets

u = mda.Universe(tpr, trajectory)

If we have used the MARTINI forcefield to study a phospholipid/cholesterol mixture, we can assign lipids and cholesterol to the upper and lower as follows:

leaflets = AssignLeaflets(
  universe=u,
  lipid_sel="name GL1 GL2 ROH"
)

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

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

The results are then available in the leaflets.leaflets attribute as a numpy.ndarray. Each row corresponds to an individual lipid and each column to an individual frame, i.e leaflets.leaflets[i, j] contains the leaflet membership of lipid i at frame j. Lipid i, at frame j, is in the upper leaflet if leaflets.leaflets[i, j]==1 and in the lower leaflet if leaflets.leaflets[i, j]==-1.

Allowing lipids in the midplane

The above example will assign every lipid (including sterols) to either the upper or lower leaflet. To allow cholesterol to be in the midplane, we can provide a midplane_sel and midplane_cutoff to AssignLeaflets:

leaflets = AssignLeaflets(
  universe=u,
  lipid_sel="name GL1 GL2 ROH",
  midplane_sel="resname CHOL and name ROH C2",
  midplane_cutoff=12.0
)

A cholesterol molecule that has both its ROH and C2 atoms within 12 Å of membrane midpoint will be assigned to the midplane, i.e. for cholesterol i at frame j that is in the midplane, leaflets.leaflets[i, j]==0.

Changing the resolution of the membrane grid

The first two examples compute a global membrane midpoint based on all the atoms of the lipids in the membrane. Lipids are then assigned a leaflet based on their distance in \(z\) to this midpoint. This is okay for planar bilayers, but can lead to incorrect leaflet classification in membranes with undulations. If your bilayer has undulations, AssignLeaflets can account for this by creating a grid in \(xy\) of your membrane, calculating the local membrane midpoint in each patch, then assigning leaflet membership based on distance in \(z\) to the local membrane midpoint. This is done through use of n_bins:

leaflets = AssignLeaflets(
  universe=u,
  lipid_sel="name GL1 GL2 ROH",
  midplane_sel="resname CHOL and name ROH C2",
  midplane_cutoff=12.0,
  n_bins=10
)

In this example, the membrane will be split into a 10 x 10 grid and a lipid assigned a leaflet based on the distance to the midpoint of the patch the lipid is in.

Assigning leaflets in membranes with high curvature

If your membrane is a vesicle or bilayer with very large undulations, such as in a buckled membrane, lipyphilic.lib.assign_leaflets.AssignLeaflets will assign lipids to the wrong leaflet

The class lipyphilic.lib.assign_leaflets.AssignCurvedLeaflets can be used in these scenaries to assign each lipid to a leaflet using MDAnalysis’ Leaflet Finder. Lipids may still be assigned to the upper/outer leaflet (indicated by 1), the lower/inner leaflet (-1) or the membrane midplane (0).

Input

Required:
  • universe : an MDAnalysis Universe object

  • lipid_sel : atom selection for all lipids in the bilayer, including e.g. sterols

  • lf_cutoff : distance cutoff below which two neighbouring atoms will be considered to be in the same leaflet.

Options:
  • midplane_sel : atom selection for lipid that may occupy the midplane

  • midplane_cutoff : atoms further than this distance from the either leaflet are considered to be the midplane

  • pbc : bool, specifying whether or not to take periodic boundaries into account

Output

  • leaflets : leaflet to which each lipid is assigned at each frame

Leaflet data are returned in a numpy.ndarray, where each row corresponds to an individual lipid and each column corresponds to an individual frame, i.e. leaflets[i, j] refers to the leaflet of lipid i at frame j. The results are accessible via the AssignLeaflets.leaflets attribute.

Example usage of AssignCurvedLeaflets

An MDAnalysis Universe must first be created before using AssignCurvedLeaflets:

import MDAnalysis as mda
from lipyphilic.lib.assign_leaflets import AssignLeaflets

u = mda.Universe(tpr, trajectory)

If we have used the MARTINI forcefield to study a phospholipid/cholesterol mixture, we can assign lipids and cholesterol to the upper and lower as follows:

leaflets = AssignCurvedLeaflets(
  universe=u,
  lipid_sel="name GL1 GL2 ROH",
  lf_cutoff=12.0,
  midplane_sel="name ROH",
  midplane_cutoff=10.0
)

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

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

This will first use MDAnalysis’ Leaflet Finder to assign all lipids, excluding those in midplane_sel, to either the upper or lower leaflet. The LeafletFinder will consider two lipids to be in the same leaflet if they have GL1 or GL2 atoms within \(12\) Å of one another. From this, we find the two largest leaflets, then assign the remaining phospholipids to a leaflet based on whichever leaflet they are closest to.

The phospholipids do not change leaflets throughtout the trajectory, only cholesterol — as specified with midplane_sel and midplane_cutoff. Thus, at each frame, each cholesterol is assinged a leaflet based on it’s minimum distance to the leaflet. In the above example, if a cholesterol is within \(10\) Å of one leaflet it is assigned to that leaflet. If it is within \(10\) Å of neither or both leaflets then it is assigned to the midplane.

The results are then available in the leaflets.leaflets attribute as a numpy.ndarray. Each row corresponds to an individual lipid and each column to an individual frame, i.e leaflets.leaflets[i, j] contains the leaflet membership of lipid i at frame j. Lipid i, at frame j, is in the upper leaflet if leaflets.leaflets[i, j]==1 and in the lower leaflet if leaflets.leaflets[i, j]==-1.

The classes and their methods

class lipyphilic.lib.assign_leaflets.AssignLeaflets(universe, lipid_sel, midplane_sel=None, midplane_cutoff=None, n_bins=1)[source]

Assign lipids in a bilayer to the upper leaflet, lower leaflet, or midplane.

Set up parameters for assigning lipids to a leaflet.

Parameters
  • universe (Universe) – MDAnalysis Universe object

  • lipid_sel (str) – Selection string for the lipids in a membrane. The selection should cover all residues in the membrane, including cholesterol.

  • midplane_sel (str, optional) – Selection string for residues that may be midplane. Any residues not in this selection will be assigned to a leaflet regardless of its proximity to the midplane. The default is None, in which case all lipids will be assigned to either the upper or lower leaflet.

  • midplane_cutoff (float, optional) – Minimum distance in z an atom must be from the midplane to be assigned to a leaflet rather than the midplane. The default is 0, in which case all lipids will be assigned to either the upper or lower leaflet. Must be non-negative.

  • n_bins (int, optional) – Number of bins in x and y to use to create a grid of membrane patches. Local membrane midpoints are computed for each patch, and lipids assigned a leaflet based on the distance to their local membrane midpoint. The default is 1, which is equivalent to computing a single global midpoint.

Note

Typically, midplane_sel should select only sterols. Other lipids have flip-flop rates that are currently unaccessible with MD simulations, and thus should always occupy either the upper or lower leaflet.

class lipyphilic.lib.assign_leaflets.AssignCurvedLeaflets(universe, lipid_sel, lf_cutoff=15, midplane_sel=None, midplane_cutoff=None, pbc=True)[source]

Assign lipids in a membrane to the upper leaflet, lower leaflet, or midplane.

Set up parameters for assigning lipids to a leaflet.

Parameters
  • universe (Universe) – MDAnalysis Universe object

  • lipid_sel (str) – Selection string for the lipids in a membrane. The selection should cover all residues in the membrane, including cholesterol.

  • lf_cutoff (float, optional) – Cutoff to pass to MDAnalysis.analysis.leaflet.LeafletFinder. Lipids closer than this cutoff distance apart will be considered to be in the same leaflet. The default is 15.0

  • midplane_sel (str, optional) – Selection string for residues that may be midplane. Any residues not in this selection will be assigned to a leaflet at ever frame. The default is None, in which case no molecules will be considered to be in the midplane.

  • midplane_cutoff (float, optional) – Lipids with atoms selected in midplane_sel that are within this distance of a leaflet will be to that leaflet. If a molecule is within this distance of neither or both leaflets, it will be assigned to the midplane. The default is None.

  • pbc (bool, optional) – Take periodic boundary conditions into account. The default is True.

Note

Typically, midplane_sel should select only sterols. Other lipids have flip-flop rates that are currently unaccessible with MD simulations, and thus should always occupy either the upper or lower leaflet.

filter_leaflets(lipid_sel=None, start=None, stop=None, step=None)

Create a subset of the leaflets results array.

Filter either by lipid species or by the trajectory frames, or both.

Parameters
  • lipid_sel (str, optional) – MDAnalysis selection string that will be used to select a subset of lipids present in the leaflets results array. The default is None, in which case data for all lipids will be returned.

  • start (int, optional) – Start frame for filtering. The default is None, in which case the first frame is used as the start.

  • stop (int, optional) – Stop frame for filtering. The default is None, in which case the final frame is used as the stop.

  • step (int, optional) – Number of frames to skip when filtering frames. The deafult is None, in which case all frames between start and stop are used.