6.2.1. Ensemble Analysis base class

New in version 0.8.0.

The Analysis modules help in the implementation analyses of MDPOW simulations. To simplify the process of analyzing a collection of systems generated by a free energy simulation the objects in Ensemble Objects allow for a molecule directory’s systems to be loaded into MDAnalysis Universes and be analyzed as a group.

EnsembleAnalysis is a class inspired by the AnalysisBase from MDAnalysis which iterates over the systems in the ensemble or the frames in the systems. It sets up iterations between universes or universe frames allowing for analysis to be run on either whole systems or the frames of those systems. This allows for users to easily run analyses on MDPOW simulations.

NotImplementedError will detect whether _single_universe() or _single_frame() should be implemented, based on which is defined in the EnsembleAnalysis. Only one of the two methods should be defined for an EnsembleAnalysis. For verbose functionality, the analysis may show two iteration bars, where only one of which will actually be iterated, while the other will load to completion instantaneously, showing the system that is being worked on.

class mdpow.analysis.ensemble.EnsembleAnalysis(ensemble=None)[source]

Base class for running multi-system analyses

The class is designed based on the AnalysisBase class in MDAnalysis and is a template for creating multi-universe multi-frame analyses using the Ensemble object

Keywords:
ensemble
The Ensemble being analyzed in the class

Example Analysis

Dihedral Analysis Demonstration:

class DihedralAnalysis(mdpow.ensemble.EnsembleAnalysis):
    def __init__(self, DihedralEnsembleGroup):
        super(DihedralAnalysis, self).__init__(DihedralEnsembleGroup.ensemble)

        self._sel = DihedralEnsembleGroup

    def _prepare_ensemble(self):
        self.result_dict = {}
        for s in ['water', 'octanol']:
            self.result_dict[s] = {'Coulomb': {},
                                   'VDW': {}}
        for key in self._sel.group_keys():
            self.result_dict[key[0]][key[1]][key[2]] = None

    def _prepare_universe(self):
        self.angle_dict = {'angle': None,
                           'time': None}
        self.angles = []

    def _single_frame(self):
        angle = calc_dihedrals(self._sel[self._key].positions[0], self._sel[self._key].positions[1],
                               self._sel[self._key].positions[2], self._sel[self._key].positions[3])
        self.angles.append(angle)

    def _conclude_universe(self):
        self.angle_dict['time'] = self.times
        self.angle_dict['angle'] = self.angles
        self.result_dict[self._key[0]][self._key[1]][self._key[2]] = self.angle_dict

    def _conclude_ensemble(self):
        self.results = pd.DataFrame(data=self.result_dict)

D = DihedralAnalysis.run(start=0 stop=10, step=1)
_prepare_ensemble()[source]

For establishing data structures used in running analysis on the entire ensemble.

Data structures will not be overwritten upon moving to next system in ensemble.

_prepare_universe()[source]

For establishing data structures used in running analysis on each trajectory in ensemble

Data structures will be overwritten between upon after each trajectory has been run

_single_universe()[source]

Calculations on a single MDAnalysis.Universe object.

Run on each MDAnalysis.Universe in the Ensemble during when run() in called. NotImplementedError will detect whether _single_universe() or _single_frame() should be implemented, based on which is defined in the EnsembleAnalysis.

_single_frame()[source]

Calculate data from a single frame of trajectory.

Called on each frame for each MDAnalysis.Universe in the Ensemble.

NotImplementedError will detect whether _single_universe() or _single_frame() should be implemented, based on which is defined in the EnsembleAnalysis.

_conclude_ensemble()[source]

Run after all trajectories in ensemble are finished

_conclude_universe()[source]

Run after each trajectory is finished

static check_groups_from_common_ensemble(groups: List[mdpow.analysis.ensemble.EnsembleAtomGroup])[source]

Checks if inputted list of EnsembleAtomGroup originate from the same Ensemble

Checks every EnsembleAtomGroup in list to determine if their ensemble() references the same object in memory. If two EnsembleAtomGroup object don’t have a common Ensemble ValueError is raised.

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

Runs _single_universe() on each system or _single_frame() on each frame in the system.

First iterates through keys of ensemble, then runs _single_universe() or _single_frame().