4.1.1. Analysis building blocks — MDAnalysis.analysis.base

A collection of useful building blocks for creating Analysis classes.

class MDAnalysis.analysis.base.AnalysisBase(trajectory, verbose=False, **kwargs)[source]

Base class for defining multi frame analysis

The class it is designed as a template for creating multiframe analyses. This class will automatically take care of setting up the trajectory reader for iterating, and it offers to show a progress meter.

To define a new Analysis, AnalysisBase needs to be subclassed _single_frame must be defined. It is also possible to define _prepare and _conclude for pre and post processing. See the example below.

class NewAnalysis(AnalysisBase):
    def __init__(self, atomgroup, parameter, **kwargs):
        super(NewAnalysis, self).__init__(atomgroup.universe.trajectory,
                                          **kwargs)
        self._parameter = parameter
        self._ag = atomgroup

    def _prepare(self):
        # OPTIONAL
        # Called before iteration on the trajectory has begun.
        # Data structures can be set up at this time
        self.result = []

    def _single_frame(self):
        # REQUIRED
        # Called after the trajectory is moved onto each new frame.
        # store result of `some_function` for a single frame
        self.result.append(some_function(self._ag, self._parameter))

    def _conclude(self):
        # OPTIONAL
        # Called once iteration on the trajectory is finished.
        # Apply normalisation and averaging to results here.
        self.result = np.asarray(self.result) / np.sum(self.result)

Afterwards the new analysis can be run like this.

na = NewAnalysis(u.select_atoms('name CA'), 35).run(start=10, stop=20)
print(na.result)
times

array of Timestep times. Only exists after calling run()

Type:np.ndarray
frames

array of Timestep frame indices. Only exists after calling run()

Type:np.ndarray
Parameters:
  • trajectory (mda.Reader) – A trajectory Reader
  • verbose (bool, optional) – Turn on more logging and debugging, default False

Changed in version 1.0.0: Support for setting start, stop, and step has been removed. These should now be directly passed to AnalysisBase.run().

run(start=None, stop=None, step=None, verbose=None)[source]

Perform the calculation

Parameters:
  • start (int, optional) – start frame of analysis
  • stop (int, optional) – stop frame of analysis
  • step (int, optional) – number of frames to skip between each analysed frame
  • verbose (bool, optional) – Turn on verbosity
class MDAnalysis.analysis.base.AnalysisFromFunction(function, trajectory=None, *args, **kwargs)[source]

Create an analysis from a function working on AtomGroups

results

results of calculation are stored after call to run

Deprecated since version 1.1.0: The structure of the results array will change in MDAnalysis 2.0.

Type:ndarray

Example

>>> def rotation_matrix(mobile, ref):
>>>     return mda.analysis.align.rotation_matrix(mobile, ref)[0]
>>> rot = AnalysisFromFunction(rotation_matrix, trajectory, mobile, ref).run()
>>> print(rot.results)
Raises:

ValueError : if function has the same kwargs as BaseAnalysis

Parameters:
  • function (callable) – function to evaluate at each frame
  • trajectory (mda.coordinates.Reader (optional)) – trajectory to iterate over. If None the first AtomGroup found in args and kwargs is used as a source for the trajectory.
  • *args (list) – arguments for function
  • **kwargs (dict) – arguments for function and AnalysisBase
  • versionchanged: (.) – 1.0.0: Support for directly passing the start, stop, and step arguments has been removed. These should instead be passed to AnalysisFromFunction.run().
MDAnalysis.analysis.base.analysis_class(function)[source]

Transform a function operating on a single frame to an analysis class

For an usage in a library we recommend the following style:

>>> def rotation_matrix(mobile, ref):
>>>     return mda.analysis.align.rotation_matrix(mobile, ref)[0]
>>> RotationMatrix = analysis_class(rotation_matrix)

It can also be used as a decorator:

>>> @analysis_class
>>> def RotationMatrix(mobile, ref):
>>>     return mda.analysis.align.rotation_matrix(mobile, ref)[0]
>>> rot = RotationMatrix(u.trajectory, mobile, ref).run(step=2)
>>> print(rot.results)