Lipid z angles — lipyphilic.lib.z_angles

Author

Paul Smith

Year

2021

Copyright

GNU Public License v2

This module provides methods for calculating the angle lipids make with the positive \(z\) axis.

Two atoms must be selected per lipid, and the angle between the \(z\) axis and the vector joining the two atoms will be calculated for each lipid. The vector will always point from atom B to atom A, even for lipids in the lower leaflet. This means the angle \(\theta_{ABz}\) will be in the range \(-180° < \theta < 180°\).

Input

Required:
  • universe : an MDAnalysis Universe object

  • atom_A_sel : atom selection for atom A in each lipid

  • atom_B_sel : atom selection for atom B in each lipid

Options:
  • rad : boolean variable specifying whether to return the angle in radians

Output

  • z_angles : angle made between the \(z\)-axis and the vector from \(B\) to \(A\)

The \(z\) angles data are returned in a numpy.ndarray, where each row corresponds to an individual lipid and each column corresponds to an individual frame.

Example usage of ZAngles

An MDAnalysis Universe must first be created before using ZAngles:

import MDAnalysis as mda
from lipyphilic.lib.z_angles import ZAngles

u = mda.Universe(tpr, trajectory)

If we have used the MARTINI forcefield to study a phospholipid/cholesterol mixture, we can calculate the orientation of cholesterol in the bilayer as follows:

z_angles = ZAngles(
  universe=u,
  atom_A_sel="name ROH",
  atom_B_sel="name R5"
)

This will calculate the angle between the \(z\)-axis and the vector from the R5 bead to the ROH bead of cholesterol.

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

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

The results are then available in the z_angles.z_angles 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.

Calculate the angle in radians

By default, the results are returned in degrees. We can also specify that the results should be returned in radians:

z_angles = ZAngles(
  universe=u,
  atom_A_sel="name ROH",
  atom_B_sel="name R5",
  rad=True
)

The class and its methods

class lipyphilic.lib.z_angles.ZAngles(universe, atom_A_sel, atom_B_sel, rad=False)[source]

Calculate the orientation of lipids in a bilayer.

Set up parameters for calculating the orientations.

Parameters
  • universe (Universe) – MDAnalysis Universe object

  • atom_A_sel (str) – Selection string for atom A of lipids in the membrane.

  • atom_B_sel (str) – Selection string for atom B of lipids in the membrane.

  • rad (bool, optional) – Whether to return the angles in radians. The default is False, in which case the results are returned in degrees.

Note

The orientation is defined as the angle between \(z\) and the vector from atom B to atom A.