Flip-flop — lipyphilic.lib.flip_flop

Author

Paul Smith

Year

2021

Copyright

GNU Public License v2

This module provides methods for finding flip-flop events in a lipid bilayer.

A flip-flop event occurs when a molecule - typically a sterol - moves from one leaflet of a bilayer into the opposing leaflet.

The class lipyphilic.lib.flip_flop.FlipFlop finds the frames at which a flip-flop event begins and ends, as well as the direction of travel (upper-to-lower or lower-to-upper). FlipFlop can also determine whether each event was successful (the molecule resides in the opposing leaflet for at least a given length of time), or not (the molecule went to the midplane but returned to its original leaflet).

See Baral et al. (2020) for further discussion on flip-flop in lipid bilayers, including the affect on the flip-flop rate of the buffer size used to assign molecules to the midplane of the bilayer.

Input

Required:
  • universe : an MDAnalysis Universe object.

  • lipid_sel : atom selection for atoms to use in detecting flip-flop

  • leaflets : leaflet membership (-1: lower leaflet, 0: midplane, 1: upper leaflet) of each lipid in the membrane at each frame

Output

  • resindex : residue index of a flip-flopping molecule

  • flip_flop_start_frame : final frame at which the molecule was present in its original leaflet

  • flip_flop_end_frame : first frame at which the molecule is present in the new leaflet

  • moves_to : direction of travel of the molecule: equal to 1 if the upper leaflet is the new lealet, equal to -1 if the lower leaflet is the new leaflet

Flip-flop data area returned in a numpy.ndarray, on a “one line, one observation” basis and can be accessed via FlipFlop.flip_flops:

flip_flops = [
    [
        <resindex (0-based)>,
        <end_frame (0-based)>,
        <start_frame (0-based)>,
        <moves_to>
    ],
    ...
]

moves_to is equal to 1 or -1 if the molecule flip-flops into the upper or the lower leaflet, respectively.

Additionaly, the success or failure of each flip-flop event is stored in the attribute FlipFlop.flip_flop_success.

Example usage of FlipFlop

An MDAnalysis Universe must first be created before using FlipFlop:

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

u = mda.Universe(tpr, trajectory)

Then we need to know which leaflet each lipid is in at each frame. This may be done using lipyphilic.lib.assign_leaflets.AssignLeaflets:

leaflets = AssignLeaflets(
  universe=u,
  lipid_sel="name GL1 GL2 ROH"  # assuming we are using the MARTINI forcefield
  midplane_sel="name ROH",       # only cholesterol is allowed to flip-flop
  midplane_cutoff=8.0,          # buffer size for assigning molecules to the midplane
)
leaflets.run()

The leaflet data are stored in the leaflets.leaflets attribute. We can now create our FlipFlop object:

flip_flop = FlipFlop(
  universe=u,
  lipid_sel="name ROH",
  leaflets=leaflets.filter_leaflets("name ROH")  # pass only the relevant leaflet data
)

We then select which frames of the trajectory to analyse (None will use every frame):

flip_flop.run(
  start=None,
  stop=None,
  step=None
)

The results are then available in the flipflop.flip_flop attribute as a numpy.ndarray. Each row corresponds to an individual flip-flop event, and the four columns correspond, respectively, to the molecule resindex, flip-flop start frame, flip-flop end frame, and the leaflet in which the molecule resides after the flip-flop.

Specify minimum residence time for successful flip-flops

We can also specify the minumum number of frames a molecule must reside in its new leaflet for the flip-flop to be considered successful. We do this using the frame_cutoff parameter:

flip_flop = FlipFlop(
  universe=u,
  lipid_sel="name ROH",
  leaflets=leaflets.filter_leaflets("name ROH")
  frame_cuotff=10,
)

With frame_cutoff=10, a molecule must remain in its new leaflet for at least 10 consecutive frames for the flip-flop to be considered successful. If this condition is not met, the flip-flop event is recorded as failing.

Calculating the flip-flop rate

The flip-flop rate can be calculatd directly from the number of successfull flip-flop evetns, which itself can be calculated as:

n_successful = sum(flip_flop.flip_flop_success == "Success")

The rate is then given by the total number of successful flip-flops divided by the total simulations time and the number of molecules of the translocating species.

The class and its methods

class lipyphilic.lib.flip_flop.FlipFlop(universe, lipid_sel, leaflets, frame_cutoff=1)[source]

Find flip-flop events in a lipid bilayer.

Set up parameters for finding flip-flop events.

Parameters
  • universe (Universe) – MDAnalysis Universe object

  • lipid_sel (str) – Selection string for atoms to use in detecting flip-flop.

  • leaflets (numpy.ndarray (n_lipids,, n_frames)) – An array of leaflet membership for each lipid as each frame, in which: -1 corresponds to the lower leaflet; 1 corresponds to the upper leaflet; and 0 corresponds to the midplane.

  • frame_cutoff (int, optional) – To be counted as a successful flip-flop, a molecule must reside in its new leaflet for at least ‘frame_cutoff’ consecutive frames. The default is 1, in which case the molecule only needs to move to the opposing leaflet for a single frame for the flip-flop to be successful.

Tip

Leaflet membership can be determined using lipyphilic.lib.assign_leaflets.AssignLeaflets.