Trajectory transformations — lipyphilic.transformations.transformations

This module contains methods for applying on-the-fly trajectory transformations with MDAnalysis.

Fix membranes broken across periodic boundaries

The callable class 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
import lipyphilic as lpp

u = mda.Universe("production.tpr", "production.xtc")

ag = u.select_atoms("resname DPPC DOPC CHOL")

u.trajectory.add_transformations(lpp.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.

class lipyphilic.transformations.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.

Deprecated since version 0.11.0: Please use MDAnalysis.transformations.NoJump instead.

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 to nojump.

Warning

The current implementation of nojump can only unwrap coordinates in orthorhombic systems.