CodeEntropy.levels.mda module¶
MDAnalysis universe utilities.
This module contains helpers for creating reduced MDAnalysis Universe objects by sub-selecting frames and/or atoms, and for building a Universe that combines coordinates from one trajectory with forces sourced from a second trajectory.
- class CodeEntropy.levels.mda.UniverseOperations[source]¶
Bases:
objectFunctions to create and manipulate MDAnalysis Universe objects.
- This helper provides methods to:
Build reduced universes by selecting subsets of frames or atoms.
Extract a single fragment (molecule) into a standalone universe.
Merge coordinates from one trajectory with forces from another trajectory.
- convert_lammps(tprfile: str, trrfile, fileformat: str | None = None) MDAnalysis.Universe[source]¶
Update the units produced from the universe produced from LAMMPS format topology and trajectory files. MDA currently has a bug that results in forces not being converted to the correct units (see issue for more details: https://github.com/MDAnalysis/mdanalysis/issues/5115 ) The method currently expects the following additional columns in the lammps dump file: fx fy fz c_5 c_7 where c_5 and c_7 are the atom potential and kinetic energies respectively.
This method loads:
Coordinates and dimensions from the coordinate trajectory (
tprfile+trrfile).
- Parameters:
tprfile – Topology input file.
trrfile – Coordinate trajectory file(s). This can be a single path or a list, as accepted by MDAnalysis.
fileformat – Optional file format for the coordinate trajectory, as recognised by MDAnalysis.
- Returns:
A new Universe containing coordinates, forces and dimensions loaded into memory.
- Return type:
MDAnalysis.Universe
- Raises:
ValueError – If fileformat is not one of the supported values.
- extract_fragment(universe: MDAnalysis.Universe, molecule_id: int) MDAnalysis.Universe[source]¶
Extract a single molecule (fragment) as a standalone reduced universe.
- Parameters:
universe – The source universe.
molecule_id – Fragment index in universe.atoms.fragments.
- Returns:
A reduced universe containing only the atoms of the selected fragment.
- merge_forces(tprfile: str, trrfile, forcefile: str, fileformat: str | None = None, kcal: bool = False, *, force_format: str | None = None, fallback_to_positions_if_no_forces: bool = True) MDAnalysis.Universe[source]¶
Create a universe by merging coordinates and forces from different files.
This method loads:
Coordinates and dimensions from the coordinate trajectory (
tprfile+trrfile).Forces from the force trajectory (
tprfile+forcefile).
If the force trajectory does not expose forces in MDAnalysis (e.g., the file does not contain forces, or the reader does not provide them), then:
If
fallback_to_positions_if_no_forcesis True, positions from the force trajectory are used as the “forces” array (backwards-compatible behaviour with earlier implementations).Otherwise, the underlying
NoDataErroris raised.
- Parameters:
tprfile – Topology input file.
trrfile – Coordinate trajectory file(s). This can be a single path or a list, as accepted by MDAnalysis.
forcefile – Trajectory containing forces.
fileformat – Optional file format for the coordinate trajectory, as recognised by MDAnalysis.
kcal – If True, scale the force array by 4.184 to convert from kcal to kJ.
force_format – Optional file format for the force trajectory. If not provided, uses
fileformat.fallback_to_positions_if_no_forces – If True, and the force trajectory has no accessible forces, use positions from the force trajectory as a fallback (legacy behaviour).
- Returns:
A new Universe containing coordinates, forces and dimensions loaded into memory.
- Return type:
MDAnalysis.Universe
- select_atoms(u: MDAnalysis.Universe, select_string: str = 'all') MDAnalysis.Universe[source]¶
Create a reduced universe by dropping atoms according to user selection.
- Parameters:
u – A Universe object with topology, coordinates and (optionally) forces.
select_string – MDAnalysis select_atoms selection string.
- Returns:
A reduced universe containing only the selected atoms. Coordinates, forces (if present) and dimensions are loaded into memory.
- select_frames(u: MDAnalysis.Universe, start: int | None = None, end: int | None = None, step: int = 1) MDAnalysis.Universe[source]¶
Create a reduced universe by dropping frames according to user selection.
- Parameters:
u – A Universe object with topology, coordinates and (optionally) forces.
start – Frame index to start analysis. If None, defaults to 0.
end – Frame index to stop analysis (Python slicing semantics). If None, defaults to the full trajectory length.
step – Step size between frames.
- Returns:
A reduced universe containing the selected frames, with coordinates, forces (if present) and unit cell dimensions loaded into memory.