Lipid z thickness — lipyphilic.lib.z_thickness
¶
- Author
Paul Smith
- Year
2021
- Copyright
GNU Public License v2
This module provides methods for calculating the thickness in \(z\) of lipids or lipid tails.
The thickness of lipid tails is a useful input feature for creating Hidden Markov Models (HMM) to detect phase separation in lipid bilayers. See Park and Im (2019) for a description of using HMMs in lipid membrane analysis.
Input¶
- Required:
universe : an MDAnalysis Universe object
lipid_sel : atom selection for the atoms to be used in calculating the thickness of a lipid
Output¶
z_thickness : thickness in \(z\) of each lipid in the bilayer
The \(z\) thickness 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 ZThickness
¶
An MDAnalysis Universe must first be created before using ZThickness:
import MDAnalysis as mda
from lipyphilic.lib.z_thickness import ZThickness
u = mda.Universe(tpr, trajectory)
If we have used the MARTINI forcefield to study a phospholipid/cholesterol mixture, we can calculate the thickness of cholesterol and sn1 tails in the bilayer as follows:
z_thickness_sn1 = ZThickness(
universe=u,
lipid_sel="(name ??1 ??A) or (resname CHOL and not name ROH)"
)
Above, our lipid_sel selection will select sn1 beads and cholesterol beads in the MARTINI forcefield, making use of the powerful MDAnalysis atom selection language.
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_thickness_sn1.run(
start=None,
stop=None,
step=None,
verbose=True
)
The results are then available in the z_thickness_sn1.z_thickness
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.
Averaging the thickness of two tails¶
Above we saw how to calculate the thickness of the sn1 tail of lipids along with cholesterol. Similarly, we can calculate the thickness of the sn2 tails:
z_thickness_sn2 = ZThickness(
universe=u,
lipid_sel="(name ??1 ??A) or (resname CHOL and not name ROH)"
)
z_thickness_sn2.run(verbose=True)
Now, if we would like to know the mean thickness of acyl tails across both sn1 and sn2 tails,
we can use the average()
method of ZThickness
:
z_thickness = ZThickness.average(
z_thickness_sn1,
z_thickness_sn2
)
This will average the thickness of the two tails, leaving the cholesterol thicknesses (from
z_thickness_sn1) unchanged, and return a new ZThickness
object containing the averaged data
in its z_thickness
attribute.
The class and its methods¶
- class lipyphilic.lib.z_thickness.ZThickness(universe, lipid_sel)[source]¶
Calculate the thickness in z of lipids in a bilayer.
Set up parameters for calculating lipid thicknesses.
- Parameters
universe (Universe) – MDAnalysis Universe object
lipid_sel (str) – Selection string for atoms to use in calculating lipid thickneses
- static average(sn1_thickness, sn2_thickness)[source]¶
Calculate the average thickness of two tails.
Given two ZThickness objects, typically each representing either the sn1 or sn2 tails of the lipids, an averagte thickness of each lipid is calculated.
- Parameters
sn1_thickness (ZThickness) – A ZThickness object for which the thicknesses have been calculated.
sn2_thickness (ZThickness) – A ZThickness object for which the thicknesses have been calculated.
- Returns
z_thickness (ZThickness) – A new ZThickness object containing the averaged data in its z_thickness attribute.
Warning
The frames used in analysing ‘sn1_thickness’ and ‘sn2_thickness’ must be the same - i.e. the ‘start’, ‘stop’, and ‘step’ parameters passed to the ‘.run()’ methods must be identical.