4.1.1. Analysis building blocks — MDAnalysis.analysis.base¶
The 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 is designed as a template for creating multi-frame analyses. This class will automatically take care of setting up the trajectory reader for iterating, and it offers to show a progress meter. Computed results are stored inside the
resultsattribute.To define a new Analysis, AnalysisBase needs to be subclassed and
_single_frame()must be defined. It is also possible to define_prepare()and_conclude()for pre- and post-processing. All results should be stored as attributes of theResultscontainer.Parameters: - trajectory (MDAnalysis.coordinates.base.ReaderBase) – A trajectory Reader
- verbose (bool, optional) – Turn on more logging and debugging
-
times¶ array of Timestep times. Only exists after calling
AnalysisBase.run()Type: numpy.ndarray
-
frames¶ array of Timestep frame indices. Only exists after calling
AnalysisBase.run()Type: numpy.ndarray
-
results¶ results of calculation are stored after call to
AnalysisBase.run()Type: Results
Example
from MDAnalysis.analysis.base import AnalysisBase 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.results.example_result = [] def _single_frame(self): # REQUIRED # Called after the trajectory is moved onto each new frame. # store an example_result of `some_function` for a single frame self.results.example_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.results.example_result = np.asarray(self.example_result) self.results.example_result /= np.sum(self.result)
Afterwards the new analysis can be run like this
import MDAnalysis as mda from MDAnalysisTests.datafiles import PSF, DCD u = mda.Universe(PSF, DCD) na = NewAnalysis(u.select_atoms('name CA'), 35) na.run(start=10, stop=20) print(na.results.example_result) # results can also be accessed by key print(na.results["example_result"])
Changed in version 1.0.0: Support for setting
start,stop, andstephas been removed. These should now be directly passed toAnalysisBase.run().Changed in version 2.0.0: Added
results
-
class
MDAnalysis.analysis.base.AnalysisFromFunction(function, trajectory=None, *args, **kwargs)[source]¶ Create an
AnalysisBasefrom a function working on AtomGroupsParameters: - function (callable) – function to evaluate at each frame
- trajectory (mda.coordinates.Reader, optional) – trajectory to iterate over. If
Nonethe first AtomGroup found in args and kwargs is used as a source for the trajectory. - *args (list) – arguments for
function - **kwargs (dict) – arguments for
functionandAnalysisBase
-
results.frames¶ simulation frames used in analysis
Type: numpy.ndarray
-
results.times¶ simulation times used in analysis
Type: numpy.ndarray
-
results.timeseries¶ Results for each frame of the wrapped function, stored after call to
AnalysisFromFunction.run().Type: numpy.ndarray
Raises: ValueError– iffunctionhas the samekwargsasAnalysisBaseExample
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.timeseries)
Changed in version 1.0.0: Support for directly passing the
start,stop, andsteparguments has been removed. These should instead be passed toAnalysisFromFunction.run().Changed in version 2.0.0: Former
resultsare now stored asresults.timeseries
-
class
MDAnalysis.analysis.base.Results(*args, **kwargs)[source]¶ Container object for storing results.
Resultsare dictionaries that provide two ways by which values can be accessed: by dictionary keyresults["value_key"]or by object attribute,results.value_key.Resultsstores all results obtained from an analysis after callingrun().The implementation is similar to the
sklearn.utils.Bunchclass in scikit-learn.Raises: AttributeError– If an assigned attribute has the same name as a default attribute.ValueError– If a key is not of typestrand therefore is not able to be accessed by attribute.
Examples
>>> from MDAnalysis.analysis.base import Results >>> results = Results(a=1, b=2) >>> results['b'] 2 >>> results.b 2 >>> results.a = 3 >>> results['a'] 3 >>> results.c = [1, 2, 3, 4] >>> results['c'] [1, 2, 3, 4]
New in version 2.0.0.
-
MDAnalysis.analysis.base.analysis_class(function)[source]¶ Transform a function operating on a single frame to an
AnalysisBaseclass.Parameters: function (callable) – function to evaluate at each frame -
results.frames¶ simulation frames used in analysis
Type: numpy.ndarray
-
results.times¶ simulation times used in analysis
Type: numpy.ndarray
-
results.timeseries¶ Results for each frame of the wrapped function, stored after call to
AnalysisFromFunction.run().Type: numpy.ndarray
Raises: ValueError– iffunctionhas the samekwargsasAnalysisBaseExamples
For use 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.timeseries)
Changed in version 2.0.0: Former
resultsare now stored asresults.timeseries-