Table Of Contents

Previous topic

hydrophobicity module

Next topic

lipid_order_parameters module

This Page

lipid_analysis module

Code author: Niklaus Johner <niklaus.johner@a3.epfl.ch>

This module contains functions to determine lipid tilt and splay angles and Calculate the elastic properties of the membrane from there.

If you use this code please cite ref. [1], where the method and its implementation is described in details or ref. [2], where the method is applied to complex lipidic phases.

To ensure proper treatment of the periodic boundary conditions, the trajectory should first be extended to neighboring unit cells, and then aligned. The tilts and splays are then calculated only for the central cell, the others being only used to ascertain correct treatment of the PBCs.

References

[1]Niklaus Johner, D. Harries and G. Khelashvili, “Implementation of a methodology for determining elastic properties of lipid assemblies from molecular dynamics simulations”, submitted to BMC bioinformatics (2015).
[2](1, 2) Niklaus Johner, D. Harries and G. Khelashvili, “Curvature and lipid packing modulate the elastic properties of lipid assemblies: comparing the HII and lamellar phases.” The Journal of Physical Chemistry Letters 5, no. 23 (2014), 4201-6.
GetBoundaryBetweenViews(t, waters, lipids, outdir='', density_cutoff=None, stride=1, within_size_normals=5.0, filename_basis='')[source]

This function determines the interface between two views, typically the lipid-water interface and assigns normals to every point on the surface.

Parameters:
  • t (CoordGroupHandle) – the trajectory
  • waters (EntityView) – First view
  • lipids (EntityView) – Second view
  • outdir (str) – Path to output directory
  • density_cutoff (float) – Interface will not be calculated for regions where the density is lower than this cutoff.
  • stride (int) – stride used to calculate the average density from the trajectory.
  • within_size_normals – radius of the patch used to determine the normals on the water-lipid interface.
  • filename_basis (str) – used as first part in the name of all the files generated.
Returns:

A tuple (water_filtered,lipid_filtered,b_eh) containing the density for the first and second views and an entity for the boundary between the two views. Every atom in the boundary has an associated normal vector set as a Vec3 property ‘n’.

WARNING: Interface was changed. Taking 2 views instead of a list of lipid names and water name. The order of the two views was also inversed. I also removed the PBC, cell_center and cell_size parameters

AssignNormalsToLipids(t, eh, b_eh, lipid_names, head_group_dict)[source]
Parameters:
  • t (CoordGroupHandle) – The trajectory
  • eh (Entity) – The associated entity
  • b_eh (Entity) – The surface from which the normals will get assigned to the lipids.
  • lipid_names (list) – a list of the lipid residue names to which normals will be assigned.
  • head_group_dict (dict) – a dictionary containing the selections defining the headgroup for each lipid type. There should be one entry for each residue name in lipid_names
Returns:

A dictionary with one entry for each residue name in lipid_names. Each element in the dictionary is a list(Vec3List). Each element in the list corresponds to one residue and each Vec3 in the Vec3List is the normal for one frame of the trajectory.

AnalyzeLipidTilts(t, eh, lipid_names, lipid_normal_dict, head_group_dict, tail_dict, prot_cm=None, bool_prop='')[source]

This function calculates the lipid tilts from a trajectory.

Parameters:
  • t (CoordGroupHandle) – The trajectory
  • eh (EntityHandle) – The associated entity
  • lipid_names (str) – List of the residue names of the different lipids in the system
  • lipid_normal_dict (dict) – Dictionary of normal vectors. One entry for every lipid type (element in lipid_names) Every entry is a list(Vec3List) of normals for every frame for every lipid of that type (size of list:NLipidsx NFrames).
  • head_group_dict (dict) – Dictionary containing a selection string for each lipid type that is used to determine the position of the lipid headgroups (center of mass of the selection).
  • tail_dict (dict) – Dictionary containing a selection string for each lipid type that is used to determine the position of the lipid tails (center of mass of the selection).
  • prot_cm (Vec3List) – A list of position (one for each frame). If specified the each tilt the distance between this position and the lipid in question will also be returned. This is typically used to calculate local properties of the membrane around an insertion.
  • bool_prop (bool) – Boolean property assigned to lipids to determine whether they should be considered in the tilt calculations. This is typically used to treat the periodic boundary conditions, to differentiate lipids from the central unit cell, for which tilt and splay are calculated, from the lipids from neighboring unit cells, used only to ensure correct treatment of PBC.
Returns:

Dictionary of lipid tilts. One entry for every lipid type (element in lipid_names). Every entry is a list with two elements. The first one is a list(FloatList) of tilts for every frame for every lipid of that type (size of list:NLipidsx NFrames). The second element contains the list of distances to prot_cm if it was defined and is empty otherwise.

WARNING: Removed parameters PBC, cell_center, cell_size

AnalyzeLipidSplays(t, eh, lipid_names, head_group_dict, tail_dict, lipid_normal_dict, lipid_tilt_dict, distance_sele_dict, distance_cutoff=10, bool_prop='')[source]

This function calculates the lipid splays from a trajectory.

Parameters:
  • t – The trajectory
  • eh – The associated entity
  • lipid_names – A list of the residue names of the different lipids in the system
  • head_group_dict – Dictionary containing a selection string for each lipid type that is used to determine the position of the lipid headgroups (center of mass of the selection).
  • tail_dict – Dictionary containing a selection string for each lipid type that is used to determine the position of the lipid tails (center of mass of the selection).
  • lipid_normal_dict – Dictionary of normal vectors. One entry for every lipid type (element in lipid_names) Every entry is a list(Vec3List) of normals for every frame for every lipid of that type (size of list:NLipidsx NFrames).
  • lipid_tilt_dict – Dictionary of lipid tilts. One entry for every lipid type (element in lipid_names) Every entry is a list(FloatList) of tilts for every frame for every lipid of that type (size of list:NLipidsx NFrames).
  • distance_sele_dict – Dictionary containing a selection string for each lipid type that is used to calculate the distance between lipids (center of mass distance). The center of mass of these selections should lie on the neutral plane.
  • distance_cutoff – Lipid pairs further apart than this distance will not be considered for splay calculation.
  • bool_prop – Boolean property assigned to lipids to determine whether they should be considered in the splay calculations. This is typically used to treat the periodic boundary conditions, to differentiate lipids from the central unit cell, for which tilt and splay are calculated, from the lipids from neighboring unit cells, used only to ensure correct treatment of PBC.
Returns:

Dictionary of splays, containing one entry for each possible lipid pairs

WARNING: Changed the name of that function from AnalyzeLipidSplay to AnalyzeLipidSplays.

AnalyzeLipidTiltAndSplay(t, lipid_names, head_group_dict, tail_dict, distance_cutoff=10.0, within_size_normals=10.0, distance_sele_dict={}, water_name='TIP3', outdir='', density_cutoff=None, prot_sele=None, density_stride=10, tilt_bool_prop='', splay_bool_prop='', filename_basis='', sele_dict={})[source]

This function is a wrapper to determine the membrane elastic moduli from the lipid titls and splays. Periodic boundary conditions are not treated explicitely here and should be treated as suggested in the description of this module.

Parameters:
  • t (CoordGroupHandle) – The trajectory
  • lipid_names (str) – List of the residue names of the different lipids in the system
  • head_group_dict (dict) – Dictionary containing a selection string for each lipid type that is used to determine the position of the lipid headgroups (center of mass of the selection).
  • tail_dict (dict) – Dictionary containing a selection string for each lipid type that is used to determine the position of the lipid tails (center of mass of the selection).
  • distance_cutoff (float) – Lipid pairs further apart than this distance will not be considered for splay calculation.
  • within_size_normals (float) – radius of the patch used to determine the normals on the water-lipid interface.
  • distance_sele_dict (dict) – Dictionary containing a selection string for each lipid type that is used to calculate the distance between lipids (center of mass distance). The center of mass of these selections should lie on the neutral plane.
  • water_name – Residue name of the waters (used to calculate the water-lipid interface).
  • outdir (str) – Path to output directory. If none, no files will be written.
  • density_cutoff (float) – Interface will not be calculated for regions where the density is lower than this cutoff.
  • prot_sele (str) – Selection string used to determine the position of the protein. This is used to calculate the distance between lipids and the protein which will be returned together with splays and tilts and can be used to determine the vatiation in membrane properties around a protein
  • density_stride (int) – Stride to be used for the calculation of the average lipid and water densities used to determine the interface. Using every single frame can slow down the calculation.
  • tilt_bool_prop (bool) – Boolean property assigned to lipids to determine whether they should be considered in the titl calculations. This is typically used to treat the periodic boundary conditions, to differentiate lipids from the central unit cell, for which tilt and splay are calculated, from the lipids from neighboring unit cells, used only to ensure correct treatment of PBC.
  • splay_bool_prop (bool) – Same as tilt_bool_prop but for the calculation of splays
  • filename_basis (str) – used as first part in the name of all the files generated.
  • sele_dict (dict) – Dictionary containing selection strings used to separate the system in several parts for the calculation This can be used for example to make the tilt and splay calculations separately for each leaflet of a bilayer. In such a case sele_dict would be something like sele_dict={“upper”:”z>0”,”lower”:”z<0”}
Returns:

A tuple (lipid_tilt_dict,lipid_normal_dict,splay_dict,b_eh), where lipid_tilt_dict, lipid_normal_dict and lipid_splay_dict are dictionaries with keys corresponding to the elements in sele_dict. For more information about lipid_tilt_dict and lipid_splay_dict, refer to the the documentation for the AnalyzeLipidTilts and AnalyzeLipidSplays functions.

WARNING: Removed parameters PBC, cell_center, cell_size

WriteTiltDict(lipid_tilt_dict, outdir, filename_basis='')[source]

This function writes out the lipid_tilt_dict to the firectory outdir. It writes out one file for every key in lipid_tilt_dict, i.e. for every lipid type. Filenames are preceded by the filename_basis.

Parameters:
  • lipid_tilt_dict (dict) – the lipid tilt dictionary as produced by AnalyzeLipidTiltAndSplay().
  • outdir (str) – the directory to which the files will be written
  • filename_basis (str) – Will be prepended to all file names.
WriteSplayDict(splay_dict, outdir, filename_basis)[source]

This function writes out the splay_dict to the firectory outdir. It writes out one file for every key in splay_dict, i.e. for every lipid type pair. Filenames are preceded by the filename_basis.

Parameters:
  • splay_dict (dict) – the lipid splay dictionary as produced by AnalyzeLipidTiltAndSplay().
  • outdir (str) – the directory to which the files will be written
  • filename_basis (str) – Will be prepended to all file names.
FitTiltDistribution(tilt_list, nbins=90, x_range=None, outdir='', filename_basis='', title_complement='', degrees=False)[source]

This function extracts the tilt modulus from a list of lipid tilts. It will first fit a gaussian y=A exp[(x-mu)/sigma^2] to the distribution of tilts to determine the fitting range used to fit the potential of mean force (PMF). Different fitting ranges are used to estimate the error on the extracted tilt modulus

Parameters:
  • tilt_list (list) – A list of lipid splays.
  • nbins (int) – The number of bins used when determining the distribution of splays
  • x_range (tuple of 2 floats) – The range in which the distribution will be calculated. Defaults to [mean-3*std,mean+3*std].
  • outdir (str) – the directory to which output files will be written, i.e. plots of splay distribution and PMF. if it is not defined, plots will not be generated
  • filename_basis (str) – Will be prepended to all file names.
  • title_complement (str) – will be added in the title of the plots
  • degrees (bool) – Whether plots should be in degrees or radians.
Returns:

A tuple (Chi, DeltaChi,Chi_list), containing the tilt modulus Chi, the estimated uncertainty on Chi, and the list of Chi values obtained from the different fitting ranges.

Return type:

(float,:class:float,:class:list)

FitSplayDistribution(splay_list, lipid_area, nbins=100, x_range=None, outdir='', filename_basis='', title_complement='')[source]

This function extracts the bending modulus from a list of splays. It will first fit a gaussian y=A exp[(x-mu)/sigma^2] to the distribution of splays to determine the fitting range used to fit the potential of mean force (PMF). Different fitting ranges are used to estimate the error on the extracted bending rigidity

Parameters:
  • splay_list (list) – A list of lipid splays.
  • lipid_area (float) – The area per lipid.
  • nbins (int) – The number of bins used when determining the distribution of splays
  • x_range (tuple of 2 floats) – The range in which the distribution will be calculated. Defaults to [mean-3*std,mean+3*std].
  • outdir (str) – the directory to which output files will be written, i.e. plots of splay distribution and PMF. if it is not defined, plots will not be generated
  • filename_basis (str) – Will be prepended to all file names.
  • title_complement (str) – will be added in the title of the plots
Returns:

A tuple (K, DeltaK,K_list), containing the bending rigidity K, the estimated uncertainty on K, and the list of K values obtained from the different fitting ranges.

Return type:

(float,:class:float,:class:list)

AnalyzeAreaPerLipid(t, lipids)[source]

This function calculates the area per lipid simply from the number of lipids and the size of the simulation box. The area per lipid is calculated for each frame in the simulation and the average is returned. This is only suitable for bilayers

Parameters:
  • t (CoordGroupHandle) – The trajectory
  • lipids – Selection of the lipids in the system
Returns:

The area per lipid

Return type:

float

ExtractTiltAndSplayModuli(tilt_dict, splay_dict, lipid_area, outdir, nbins=100)[source]

This function extracts the tilt and splay moduli from the dictionaries of tilts and splays obtained from the AnalyzeLipidTiltAndSplay function. It will first fit a gaussian y=A exp[(x-mu)/sigma^2] to the distribution of tilts and the distribution of splays to determine the fitting range then used to fit the corresponding potential of mean force (PMF). Different fitting ranges are used to estimate the error on the extracted tilt and splay moduli. The function will calculate one tilt modulus for each lipid species and one splay modulus for each pair of lipid species. It will then combine these to calculate the overall tilt modulus and splay modulus (bending rigidity). More details about this procedure can be found in ref. [2]

Parameters:
  • tilt_dict (dict) – A dictionary of lipid tilts as returned by the AnalyzeLipidTiltAndSplay function.
  • splay_dict (dict) – A dictionary of lipid splays as returned by the AnalyzeLipidTiltAndSplay function.
  • lipid_area – The area per lipid.
  • outdir (str) – the directory to which output files will be written, i.e. plots of tilt and splay distributions and PMFs. and text files containing the elastic moduli
  • nbins (int) – The number of bins used when determining the distributions of tilts and splays