molecular_simulations.analysis.interaction_energy module
- class molecular_simulations.analysis.interaction_energy.DynamicInteractionEnergy(top, traj, stride=1, chain='A', platform='CUDA', first_residue=None, last_residue=None, progress_bar=False)
Bases:
objectClass for obtaining interaction energies of a trajectory. Utilizes the InteractionEnergyFrame child class to run per-frame energy calculations and orchestrates the trajectory operations.
- Parameters:
top (PathLike) – Path to prmtop topology file.
traj (PathLike) – Path to DCD trajectory file.
stride (int) – Defaults to 1. The stride with which to move through the trajectory.
chain (str) – Defaults to A. The chain for which to compute the energy between. Computes energy between this chain and all other components in PDB file. Set to a whitespace if there are no chains in your PDB.
platform (str) – Defaults to CUDA. Supports running on GPU for speed.
first_residue (int | None) – Defaults to None. If set, will restrict the calculation to residues beginning with resid first_residue.
last_residue (int | None) – Defaults to None. If set, will restrict the calculation to residues ending with resid last_residue.
progress_bar (bool) – Defaults to False. If True a tqdm progress bar will display progress.
- build_system(top)
Builds OpenMM system for both pdb and prmtop topology files.
- Parameters:
top (PathLike) – Path to topology file.
- Returns:
OpenMM system object.
- Return type:
(System)
- compute_energies()
Computes the energy for each frame in trajectory, storing internally.
- Return type:
- Returns:
None
- load_traj(top, traj)
Loads trajectory into mdtraj and extracts full coordinate array.
- Parameters:
top (PathLike) – Path to topology file.
traj (PathLike) – Path to trajectory file.
- Returns:
Coordinate array of shape (n_frames, n_atoms, 3)
- Return type:
(np.ndarray)
- class molecular_simulations.analysis.interaction_energy.InteractionEnergy
Bases:
ABC- abstractmethod compute()
- abstractmethod energy()
- abstractmethod get_selection()
- class molecular_simulations.analysis.interaction_energy.InteractionEnergyFrame(system, top, chain='A', platform='CUDA', first_residue=None, last_residue=None)
Bases:
StaticInteractionEnergyInherits from StaticInteractionEnergy and overloads get_system to allow for more easily running this analysis on a trajectory of frames. Requires the OpenMM system and topology to be built externally and passed in rather than beginning from a PDB file.
- Parameters:
system (System) – OpenMM system object.
top (Topology) – OpenMM topology object.
chain (str) – Defaults to A. The chain for which to compute the energy between. Computes energy between this chain and all other components in PDB file. Set to a whitespace if there are no chains in your PDB.
platform (str) – Defaults to CUDA. Supports running on GPU for speed.
first_residue (int | None) – Defaults to None. If set, will restrict the calculation to residues beginning with resid first_residue.
last_residue (int | None) – Defaults to None. If set, will restrict the calculation to residues ending with resid last_residue.
- get_system()
Sets self.selection via self.get_selection and returns existing OpenMM system object.
- Returns:
OpenMM system object.
- Return type:
(System)
- class molecular_simulations.analysis.interaction_energy.StaticInteractionEnergy(pdb, chain='A', platform='CUDA', first_residue=None, last_residue=None)
Bases:
InteractionEnergyComputes the linear interaction energy between specified chain and other simulation components. Can specify a range of residues in chain to limit calculation to. Works on a static model but can be adapted to run on dynamics data.
- Inputs:
pdb (str): Path to input PDB file chain (str): Defaults to A. The chain for which to compute the energy between.
Computes energy between this chain and all other components in PDB file. Set to a whitespace if there are no chains in your PDB.
platform (str): Defaults to CUDA. Supports running on GPU for speed. first_residue (int | None): Defaults to None. If set, will restrict the
calculation to residues beginning with resid first_residue.
- last_residue (int | None): Defaults to None. If set, will restrict the
calculation to residues ending with resid last_residue.
- Parameters:
- compute(positions=None)
Compute interaction energy of system. Can optionally provide atomic positions such that this operation can be scaled onto a trajectory of frames rather than a static model.
- Parameters:
positions (np.ndarray | None) – Defaults to None. If provided, inject the positions into the OpenMM context.
- Return type:
- Returns:
None
- static energy(context, solute_coulomb_scale=0, solute_lj_scale=0, solvent_coulomb_scale=0, solvent_lj_scale=0)
Computes the potential energy for provided context object.
- Parameters:
context (Context) – OpenMM context object.
solute_coulomb_scale (int) – Defaults to 0. If 1 we will consider solute contributions to coulombic non-bonded energy.
solute_lj_scale (int) – Defaults to 0. If 1 we will consider solute contributions to LJ non-bonded energy.
solvent_coulomb_scale (int) – Defaults to 0. If 1 we will consider solvent contributions to coulombic non-bonded energy.
solvent_lj_scale (int) – Defaults to 0. If 1 we will consider solvent contributions to LJ non-bonded energy.
- Returns:
Computed energy term.
- Return type:
(float)
- fix_pdb()
Using the OpenMM adjacent tool, PDBFixer, repair input PDB by adding hydrogens, and missing atoms such that we can actually construct an OpenMM system.
- Return type:
- Returns:
None
- get_selection(topology)
Using the poorly documented OpenMM selection language, get indices of atoms that we want to isolate for pairwise interaction energy calculation.
- Parameters:
topology (Topology) – OpenMM topology object.
- Return type:
- Returns:
None
- get_system()
Builds implicit solvent OpenMM system.
- Returns:
OpenMM system object.
- Return type:
(System)