8. Trajectory transformations (“on-the-fly” transformations)

In MDAnalysis, a transformation is a function/function-like class that modifies the data for the current Timestep and returns the Timestep. For instance, coordinate transformations, such as PBC corrections and molecule fitting are often required for some analyses and visualization. Transformation functions (transformation_1 and transformation_2 in the following example) can be called by the user for any given Timestep of the trajectory,

u = MDAnalysis.Universe(topology, trajectory)

for ts in u.trajectory:
   ts = transformation_2(transformation_1(ts))

where they change the coordinates of the timestep ts in place. There is nothing special about these transformations except that they have to be written in such a way that they change the Timestep in place.

As described under Workflows, multiple transformations can be grouped together and associated with a trajectory so that the trajectory is transformed on-the-fly, i.e., the data read from the trajectory file will be changed before it is made available in, say, the AtomGroup.positions attribute.

The submodule MDAnalysis.transformations contains a collection of transformations (see Transformations in MDAnalysis) that can be immediately used but one can always write custom transformations (see Creating transformations).

8.1. Workflows

Instead of manually applying transformations, it is much more convenient to associate a whole workflow of transformations with a trajectory and have the transformations be called automatically.

A workflow is a sequence (tuple or list) of transformation functions that will be applied in this order. For example,

workflow = [transformation_1, transformation_2]

would effectively result in

ts = transformation_2(transformation_1(ts))

for every time step in the trajectory.

One can add a workflow using the Universe.trajectory.add_transformations method of a trajectory (where the list workflow is taken from the example above),

u.trajectory.add_transformations(*workflow)

or upon Universe creation using the keyword argument transformations:

u = MDAnalysis.Universe(topology, trajectory, transformations=workflow)

Note that in these two cases, the workflow cannot be changed after having being added.

8.2. Creating transformations

A simple transformation can also be a function that takes a Timestep as input, modifies it, and returns it. If it takes no other arguments but a Timestep can be defined as the following example:

def up_by_2(ts):
    """
    Translate all coordinates by 2 angstroms up along the Z dimension.
    """
    ts.positions = ts.positions + np.array([0, 0, 2], dtype=np.float32)
    return ts

If the transformation requires other arguments besides the Timestep, the following two methods can be used to create such transformation:

8.2.1. Creating complex transformation classes

It is implemented by inheriting from MDAnalysis.transformations.base.TransformationBase, which defines __call__() for the transformation class and can be applied directly to a Timestep. _transform() has to be defined and include the operations on the MDAnalysis.coordinates.base.Timestep.

So, a transformation class can be roughly defined as follows:

from MDAnalysis.transformations import TransformationBase

class up_by_x_class(TransformationBase):
    def __init__(self, distance):
        self.distance = distance

    def _transform(self, ts):
        ts.positions = ts.positions + np.array([0, 0, self.distance], dtype=np.float32)
        return ts

It is the default construction method in MDAnalysis.transformations from release 2.0.0 onwards because it can be reliably serialized. See MDAnalysis.transformations.translate for a simple example.

8.2.2. Creating complex transformation closure functions

Transformation can also be a wrapped function takes the Timestep object as argument. So in this case, a transformation function (closure) can be roughly defined as follows:

def up_by_x_func(distance):
    """
    Creates a transformation that will translate all coordinates by a given amount along the Z dimension.
    """
    def wrapped(ts):
        ts.positions = ts.positions + np.array([0, 0, distance], dtype=np.float32)
        return ts
    return wrapped

An alternative to using a wrapped function is using partials from functools. The above function can be written as:

import functools

def up_by_x(ts, distance):
    ts.positions = ts.positions + np.array([0, 0, distance], dtype=np.float32)
    return ts

up_by_2 = functools.partial(up_by_x, distance=2)

Although functions (closures) work as transformations, they are not used in in MDAnalysis from release 2.0.0 onwards because they cannot be reliably serialized and thus a Universe with such transformations cannot be used with common parallelization schemes (e.g., ones based on multiprocessing). For detailed descriptions about how to write a closure-style transformation, please refer to MDAnalysis 1.x documentation.

8.3. Transformations in MDAnalysis

The module MDAnalysis.transformations contains transformations that can be immediately used in your own workflows. In order to use any of these transformations, the module must first be imported:

import MDAnalysis.transformations

A workflow can then be added to a trajectory as described above. Notably, the parameter max_threads can be defined when creating a transformation instance to limit the maximum threads. (See MDAnalysis.transformations.base.TransformationBase for more details) Whether a specific transformation can be used along with parallel analysis can be assessed by checking its parallelizable attribute.

See Currently implemented transformations for more on the existing transformations in MDAnalysis.transformations.

8.4. How to transformations

Translating the coordinates of a single frame (although one would normally add the transformation to a workflow, as shown in the subsequent examples):

u = MDAnalysis.Universe(topology, trajectory)
new_ts = MDAnalysis.transformations.translate([1,1,1])(u.trajectory.ts)

Create a workflow and add it to the trajectory:

u = MDAnalysis.Universe(topology, trajectory)
workflow = [MDAnalysis.transformations.translate([1,1,1]),
            MDAnalysis.transformations.translate([1,2,3])]
u.trajectory.add_transformations(*workflow)

Giving a workflow as a keyword argument when defining the universe:

workflow = [MDAnalysis.transformations.translate([1,1,1]),
            MDAnalysis.transformations.translate([1,2,3])]
u = MDAnalysis.Universe(topology, trajectory, transformations=workflow)

8.5. Building blocks for Transformation Classes

Transformations normally ultilize the power of NumPy to get better performance on array operations. However, when it comes to parallelism, NumPy will sometimes oversubscribe the threads, either by hyper threading (when it uses OpenBlas backend), or by working with other parallel engines (e.g. Dask).

In MDAnalysis, we use threadpoolctl inside TransformationBase to control the maximum threads for transformations.

It is also possible to apply a global thread limit by setting the external environmental varibale, e.g. OMP_NUM_THREADS=1 MKL_NUM_THREADS=1 OPENBLAS_NUM_THREADS=1 BLIS_NUM_THREADS=1 python script.py. Read more about parallelism and resource management in scikit-learn documentations.

Users are advised to benchmark code because interaction between different libraries can lead to sub-optimal performance with defaults.