10.5. EDR auxiliary reader — MDAnalysis.auxiliary.EDR

New in version 2.4.0.

10.5.1. Background

EDR files are binary files following the XDR protocol. They are written by GROMACS during simulations and contain time-series non-trajectory data on the system, such as energies, temperature, or pressure.

pyedr is a Python package that reads EDR binary files and returns them human-readable form as a dictionary of NumPy arrays. It is used by the EDR auxiliary reader to parse EDR files. As such, a dictionary with string keys and numpy array values is loaded into the EDRReader. It is basically a Python-based version of the C++ code in GROMACS.

The classes in this module are based on the pyedr package. Pyedr is an optional dependency and must be installed to use this Reader. Use of the reader without pyedr installed will raise an ImportError. The variable HAS_PYEDR indicates whether this module has pyedr availble.

The EDR auxiliary reader takes the output from pyedr and makes this data available within MDAnalysis. The usual workflow starts with creating an EDRReader and passing it the file to read as such:

import MDAnalysis as mda
aux = mda.auxiliary.EDR.EDRReader(some_edr_file)

The newly created aux object contains all the data found in the EDR file. It is stored in the data_dict dictionary, which maps the names GROMACS gave each data entry to a NumPy array that holds the relevant data. These GROMACS-given names are stored in and available through the terms attribute. In addition to the numeric data, the new EDRReader also stores the units of each entry in the data_dict dictionary in its unit_dict dictionary.

Warning

Units are converted to MDAnalysis base units automatically unless otherwise specified. However, not all unit types have a defined base unit in MDAnalysis. (cf. MDAnalysis.units.MDANALYSIS_BASE_UNITS). Pressure units, for example, are not currently defined, and will not be converted. This might cause inconsistencies between units! Conversion can be switched off by passing convert_units=False when the EDRReader is created:

aux = mda.auxiliary.EDR.EDRReader(some_edr_file, convert_units=False)

10.5.2. Standalone Usage of the EDRReader

The EDRReader can be used to access the data stored in EDR files on its own, without association of data to trajectory frames. This is useful, for example, for plotting. The data for a single term, a list of terms, or for all terms can be returned in dictionary form. “Time” is always returned in this dictionary to make plotting easier:

temp = aux.get_data("Temperature")
plt.plot(temp["Time"], temp["Temperature"])

some_terms = aux.get_data(["Potential", "Kinetic En.", "Box-X"])
plt.plot(some_terms["Time"], some_terms["Potential"])

all_terms = aux.get_data()
plt.plot(all_terms["Time"], all_terms["Pressure"])

10.5.3. Adding EDR data to trajectories

Like other AuxReaders, the EDRReader can attach its data to a trajectory by associating it to the appropriate time steps. In general, to add EDR data to a trajectory, one needs to provide two arguments.

Note

The following will change with the release of MDAnalysis 3.0. From then on, the order of these two arguments will be reversed.

The first argument is the aux_spec dictionary. With it, users specify which entries from the EDR file they want to add, and they give it a more convenient name to be used in MDAnalysis (because GROMACS creates names like “#Surf*SurfTen” or “‘Constr. rmsd’” which may be inconvenient to use.) This dictionary might look like this:

aux_spec = {"epot": "Potential",
            "surf_tension": "#Surf*SurfTen"}

When provided as shown below, this would direct the EDRReader to take the data it finds under the “Potential” key in its data_dict dictionary and attach it to the trajectory time steps under u.trajectory.ts.aux.epot (and the same for the surface tension).

The second argument needed is the source of the EDR data itself. Here, either the path to the EDR file or a previously created EDRReader object can be provided.

Examples:

import MDAnalysis as mda
from MDAnalysisTests.datafiles import AUX_EDR, AUX_EDR_TPR, AUX_EDR_XTC
import matplotlib.pyplot as plt

A Universe and an EDRReader object are created and the data available in the EDR file is printed:

In [1]: u = mda.Universe(AUX_EDR_TPR, AUX_EDR_XTC)
In [2]: aux = mda.auxiliary.EDR.EDRReader(AUX_EDR)
In [3]: aux.terms
Out[3]: ['Time', 'Bond', 'Angle', ...]

Data is associated with the trajectory, using an aux_spec dictionary to specify which data to add under which name. Any number of terms can be added using this method. The data is then accessible in the ts.aux namespace via both attribute and dictionary syntax:

In [4]: u.trajectory.add_auxiliary({"epot": "Potential",
                                    "angle": "Angle"}, aux)
In [5]: u.trajectory.ts.aux.epot
Out[5]: -525164.0625
In [6]: u.trajectory.ts.aux.Angle
Out[6]: 3764.52734375
In [7]: u.trajectory.ts.aux["epot"]
Out[7]: -525164.0625

Note

Some GROMACS-provided terms have spaces. Unless an attribute name without a space is provided, these terms will not be accessible via the attribute syntax. Only the dictionary syntax will work in that case.

Further, it is possible to add all data from the EDR file to the trajectory. To do this, the aux_spec dictionary is omitted, and the data source (the second argument as explained above) is provided explicitly as auxdata. When adding data this way, the terms in terms become the names used in ts.aux:

In [7]: u.trajectory.add_auxiliary(auxdata=aux)
In [8]: u.trajectory.ts.aux["#Surf*SurfTen"]
Out[8]: -1857.519287109375

10.5.4. Classes

class MDAnalysis.auxiliary.EDR.EDRReader(filename: str, convert_units: bool = True, **kwargs)[source]

Auxiliary reader to read data from an .edr file.

EDR files are created by GROMACS during a simulation. They are binary files which contain time-series energy data and other data related to the simulation.

Default reader for .edr files. All data from the file will be read and stored on initialisation.

Parameters:
  • filename (str) – Location of the file containing the auxiliary data.

  • convert_units (bool, optional) – If True (default), units from the EDR file are automatically converted to MDAnalysis base units. If False, units are taken from the file as-is. Where no base unit is defined in MDAnalysis, no conversion takes place. Unit types in MDAnalysis.units.MDANALYSIS_BASE_UNITS will be converted automatically by default.

  • **kwargs – Other AuxReader options.

_auxdata

path at which the auxiliary data file is located

Type:

pathlib.PosixPath

data_dict

dictionary that contains the auxiliary data, mapping the names GROMACS gave the entries in the EDR file to a NumPy array containing this data

Type:

dict

unit_dict

dictionary that contains the units of the auxiliary data, mapping the data_selector of the Reader (i.e. the name of the dataset in the EDR file) to its unit.

Type:

dict

_n_steps

Number of steps for which auxdata is available

Type:

int

terms

Names of the auxiliary data entries available in data_dict. These are the names GROMACS set in the EDR file.

Type:

list

Note

The file is assumed to be of a size such that reading and storing the full contents is practical. A warning will be issued when memory usage exceeds 1 GB. This warning limit can be changed via the memory_limit kwarg.

get_data(data_selector: str | List[str] | None = None) Dict[str, ndarray][source]

Returns the auxiliary data contained in the EDRReader. Returns either all data or data specified as data_selector in form of a str or a list of any of EDRReader.terms. Time is always returned to allow easy plotting.

Parameters:

data_selector (str, List[str], None) – Keys to be extracted from the auxiliary reader’s data dictionary. If None, returns all data found in data_dict.

Returns:

data_dict – Dictionary mapping data_selector keys to NumPy arrays of the auxiliary data.

Return type:

dict

Raises:

KeyError – if an invalid data_selector key is passed.

read_all_times() ndarray[source]

Get list of time at each step.

Returns:

Time at each step.

Return type:

NumPy array of float

The actual data for each step is stored by instances of EDRStep.

class MDAnalysis.auxiliary.EDR.EDRStep(time_selector: str = 'Time', data_selector: str | None = None, **kwargs)[source]

AuxStep class for the .edr file format.

Extends the base AuxStep class to allow selection of time and data-of-interest fields (by dictionary key) from the full set of data read each step.

Parameters:
  • time_selector (str, optional) – Name of the dictionary key that links to the time values (assumed to be in ps). Default value is “Time”

  • data_selector (str | list of str | None, optional) – List of dictionary keys linking to data of interest in the EDR file to be stored in data. Default value is None.

  • **kwargs – Other AuxStep options.