Trajectory transformations — lipyphilic.transformations
¶
This module contains methods for applying on-the-fly trajectory transformations with MDAnalysis.
Prevent atoms from jumping across periodic boundaries¶
lipyphilic.transformations.nojump
can be used to prevent atoms from jumping across
periodic boundaries. It is equivalent to using the
GROMACS command
trjconv with the flag
-pbc nojump.
The on-the-fly transformation can be added to your trajectory after loading it with MDAnalysis:
import MDAnalysis as mda
from lipyphilic.transformations import nojump
u = mda.Universe("production.tpr", "production.xtc")
ag = u.select_atoms("name GL1 GL2 ROH")
u.trajectory.add_transformations(nojump(ag))
Upon adding this transformation to your trajectory, lipyphilic will determine at which frames each atom crosses a boundary, keeping a record of the net movement across each boundary. Then, every time a new frame is loaded into memory by MDAnalysis — such as when you iterate over the trajectory — the transformation is applied.
This transformation is required when calculating the lateral diffusion of lipids in a membrane
using, for example, lipyphilic.lib.lateral_diffusion.MSD
. It can be used to remove the
need to create an unwrapped trajectory using GROMACS.
Fix membranes broken across periodic boundaries¶
The callable class lipyphilic.transformations.center_membrane
can be used to fix a membrane
split across periodic boundaries and then center it in the unit cell. The membrane is iteratively
shifted along a dimension until it is no longer split across periodic boundaries. It is then
moved it to the center of the box in this dimension.
The on-the-fly transformation can be added to your trajectory after loading it with MDAnalysis:
import MDAnalysis as mda
from lipyphilic.transformations import center_membrane
u = mda.Universe("production.tpr", "production.xtc")
ag = u.select_atoms("resname DPPC DOPC CHOL")
u.trajectory.add_transformations(center_membrane(ag))
This will center a DPPC/DOPC/cholesterol membrane in \(z\) every time a new frame is loaded into memory by MDAnalysis, such as when you iterate over the trajectory:
for ts in u.trajectory:
# do some nice analysis with your centered membrane
Note
ag should be an AtomGroup that contains all atoms in the membrane.
Transform triclinic coordinates to their orthorhombic representation¶
lipyphilic.transformations.triclinic_to_orthorhombic
can be used to transform
triclinic coordinates to their orthorhombic representation. It is equivalent to using the
GROMACS command
trjconv with the flag
-ur rect.
The on-the-fly transformation can be added to your trajectory after loading it with MDAnalysis:
import MDAnalysis as mda
from lipyphilic.transformations import triclinic_to_orthorhombic
u = mda.Universe("production.tpr", "production.xtc")
ag = u.select_atoms("resname DPPC DOPC CHOL")
u.trajectory.add_transformations(triclinic_to_orthorhombic(ag=ag))
After adding this transformation, upon load a new frame into memory the coordinates of the selected atoms will be transformed, and the dimensions of your system will be modified so that the angles are all 90°. Further analysis may then be performed using the orthorhombic coordinate system.
Some analyses in lipyphilic create a surface of the membrane plane using a two-dimensional rectangular grid. This includes
lipyphilic.lib.assign_leaflet.AssignLeaflets
lipyphilic.lib.memb_thickness.MembThicnkess
These analyses will fail with triclinic boxes - the triclinic_to_orthorhombic transformation must be applied to triclinic systems before these tools can be used.
Another case that will fail with triclinic systems is the lipyphilic.transformations.nojump
transformation - this transformation can currently only unwrap coordinates for orthorhombic
systems.
See lipyphilic.transformations.triclinic_to_orthorhombic
for the full list.
- class lipyphilic.transformations.nojump(ag, nojump_x=True, nojump_y=True, nojump_z=False, filename=None)[source]¶
Prevent atoms jumping across periodic boundaries.
This is useful if you would like to calculate the diffusion coefficient of lipids in your membrane.
This transformation does an initial pass over the trajectory to determine at which frames each atom crosses a boundary, keeping a record of the net movement across each boundary. Then, as a frame is loaded into memory, atom positions are translated according to their total displacement, taking into account crossing of boundaries as well box fluctuations in the box volume.
By default, atoms are only unwrapped in \(xy\), as it is assumed the membrane is a bilayer. To unwrap in all dimensions,
center_z
must also be set to True.- Parameters
ag (AtomGroup) – MDAnalysis AtomGroup to which to apply the transformation
nojump_x (bool, optional) – If true, atoms will be prevented from jumping across periodic boundaries in the x dimension.
nojump_y (bool, optional) – If true, atoms will be prevented from jumping across periodic boundaries in the y dimension.
nojump_z (bool, optional) – If true, atoms will be prevented from jumping across periodic boundaries in the z dimension.
filename (str, optional) – File in which to write the unwrapped, nojump trajectory. The default is None, in which case the transformation will be applied on-the-fly.py
- Returns
MDAnalysis.coordinates.base.Timestep
object, or None if a filename is provided.
Note
The nojump transformation is memory intensive to perform on-the-fly. If you have a long trajectory or a large number of atoms to be unwrapped, you can write the unwrapped coordinates to a new file by providing a
filename
tonojump
.Warning
The current implementation of nojump can only unwrap coordinates in orthorhombic systems.
- class lipyphilic.transformations.center_membrane(ag, shift=20, center_x=False, center_y=False, center_z=True, min_diff=10)[source]¶
Fix a membrane split across periodic boundaries and center it in the primary unit cell.
If, for example, the bilayer is split across \(z\), it will be iteratively translated in \(z\) until it is no longer broken. Then it will be moved to the center of the box.
A membrane with a maximum extent almost the same size as the box length in a given dimension will be considered to be split across that dimension.
By default, the membrane is only centered in \(z\), as it is assumed the membrane is a bilayer. To center a micelle,
center_x
andcenter_y
must also be set to True.- Parameters
ag (AtomGroup) – MDAnalysis AtomGroup containing all atoms in the membrane.
shift (float, optional) – The distance by which a bilayer will be iteratively translated. This must be smaller than the thickness of your bilayer or the diameter of your micelle.
min_diff (float, optional) – Minimum difference between the box size and the maximum extent of the membrane in order for the membrane to be considered unwrapped.
center_x (bool, optional) – If true, the membrane will be iteratively shifted in x until it is not longer split across periodic boundaries.
center_y (bool, optional) – If true, the membrane will be iteratively shifted in y until it is not longer split across periodic boundaries.
center_z (bool, optional) – If true, the membrane will be iteratively shifted in z until it is not longer split across periodic boundaries.
- Returns
MDAnalysis.coordinates.base.Timestep
object
- class lipyphilic.transformations.triclinic_to_orthorhombic(ag)[source]¶
Transform triclinic coordinates to their orthorhombic representation.
If you have a triclinic system, it is essential to apply this transformation before using the following analyses:
lipyphilic.lib.assign_leaflet.AssignLeaflets
lipyphilic.lib.area_per_lipid.AreaPerLipid
lipyphilic.lib.memb_thickness.MembThicnkess
lipyphilic.lib.registration.Registration
as well as before the following on-the-fly transformations:
lipyphilic.transformations.nojump
lipyphilic.transformations.center_membrane
The above tools will fail unless provided with an orthorhombic system.
This transformation is equivalent to using the GROMACS command trjconv with the flag -ur rect.
Note
triclinic_to_rectangular will put all selected atoms into the primary (orthorhombic) unit cell - molecules will not be kept whole or unwrapped.
Warning
If you wish to apply the triclinic_to_orthorhombic transformation along with other on-the-fly transformations, triclinic_to_orthorhombic must be the first one applied.
- Parameters
ag (AtomGroup) – MDAnalysis AtomGroup to which to apply the transformation