molecular_simulations.analysis.interaction_energy module

Interaction energy calculation module.

This module provides tools for computing linear interaction energies between specified chains and other simulation components using OpenMM. Supports both static models and dynamic trajectory analysis.

class molecular_simulations.analysis.interaction_energy.InteractionEnergy[source]

Bases: ABC

Abstract base class for interaction energy calculations.

Defines the interface for all interaction energy calculation classes.

Initialize the InteractionEnergy base class.

abstractmethod __init__()[source]

Initialize the InteractionEnergy base class.

abstractmethod compute()[source]

Compute the interaction energy.

Must be implemented by subclasses.

abstractmethod energy(*args, **kwargs)[source]

Get the computed energy.

Must be implemented by subclasses.

abstractmethod get_selection(*args, **kwargs)[source]

Get the atom selection for energy calculation.

Must be implemented by subclasses.

class molecular_simulations.analysis.interaction_energy.StaticInteractionEnergy(pdb, chain='A', platform='CUDA', first_residue=None, last_residue=None)[source]

Bases: InteractionEnergy

Compute linear interaction energy for a static model.

Computes the linear interaction energy between a specified chain and other simulation components. Can specify a range of residues to limit the calculation. Works on a static model but can be adapted for trajectory data.

pdb

Path to the input PDB file.

chain

Chain identifier for interaction calculations.

platform

OpenMM platform for computation.

lj

Lennard-Jones energy after compute().

coulomb

Coulombic energy after compute().

selection

Atom indices for the selected chain.

Parameters:
  • pdb (str) – Path to input PDB file.

  • chain (str) – Chain identifier for energy calculation. Computes energy between this chain and all other components. Use whitespace if there are no chains. Defaults to ‘A’.

  • platform (str) – OpenMM platform name. Defaults to ‘CUDA’.

  • first_residue (int | None) – If set, restricts calculation to residues starting from this resid. Defaults to None.

  • last_residue (int | None) – If set, restricts calculation to residues ending at this resid. Defaults to None.

Example

>>> ie = StaticInteractionEnergy("complex.pdb", chain="B")
>>> ie.compute()
>>> print(f"LJ: {ie.lj}, Coulomb: {ie.coulomb}")

Initialize the StaticInteractionEnergy calculator.

Parameters:
  • pdb (str) – Path to input PDB file.

  • chain (str) – Chain identifier for calculation.

  • platform (str) – OpenMM platform name.

  • first_residue (int | None) – Starting residue for calculation.

  • last_residue (int | None) – Ending residue for calculation.

__init__(pdb, chain='A', platform='CUDA', first_residue=None, last_residue=None)[source]

Initialize the StaticInteractionEnergy calculator.

Parameters:
  • pdb (str) – Path to input PDB file.

  • chain (str) – Chain identifier for calculation.

  • platform (str) – OpenMM platform name.

  • first_residue (int | None) – Starting residue for calculation.

  • last_residue (int | None) – Ending residue for calculation.

get_system()[source]

Build an implicit solvent OpenMM system.

Loads the PDB file and creates an OpenMM system with GB/n2 implicit solvent model. Automatically fixes the PDB if needed.

Return type:

System

Returns:

OpenMM System object configured with implicit solvent.

compute(positions=None)[source]

Compute interaction energy of the system.

Computes both Lennard-Jones and Coulombic interaction energies between the selected chain and all other atoms.

Parameters:

positions (ndarray | None) – Optional atomic positions to use instead of those from the PDB file. Useful for trajectory analysis.

Return type:

None

get_selection(topology)[source]

Get indices of atoms for pairwise interaction calculation.

Uses OpenMM’s selection capabilities to identify atoms in the specified chain and residue range.

Parameters:

topology (Topology) – OpenMM Topology object.

Return type:

None

fix_pdb()[source]

Repair the input PDB using PDBFixer.

Adds missing residues, atoms, and hydrogens to create a complete structure suitable for OpenMM.

Return type:

tuple

Returns:

Tuple of (positions, topology) after fixing.

property interactions: ndarray

Get LJ and Coulombic energies as an array.

Returns:

Array of shape (2, 1) containing [lj, coulomb] energies.

static energy(context, solute_coulomb_scale=0, solute_lj_scale=0, solvent_coulomb_scale=0, solvent_lj_scale=0)[source]

Compute potential energy for the given context.

Parameters:
  • context (Context) – OpenMM Context object.

  • solute_coulomb_scale (int) – Scale for solute Coulombic energy (0 or 1).

  • solute_lj_scale (int) – Scale for solute LJ energy (0 or 1).

  • solvent_coulomb_scale (int) – Scale for solvent Coulombic energy (0 or 1).

  • solvent_lj_scale (int) – Scale for solvent LJ energy (0 or 1).

Returns:

Computed energy term with units.

class molecular_simulations.analysis.interaction_energy.InteractionEnergyFrame(system, top, chain='A', platform='CUDA', first_residue=None, last_residue=None)[source]

Bases: StaticInteractionEnergy

Interaction energy calculator for trajectory frames.

Inherits from StaticInteractionEnergy and overloads get_system to allow for easier trajectory analysis. Requires the OpenMM system and topology to be built externally.

Parameters:
  • system (System) – Pre-built OpenMM System object.

  • top (Topology) – OpenMM Topology object.

  • chain (str) – Chain identifier for calculation. Defaults to ‘A’.

  • platform (str) – OpenMM platform name. Defaults to ‘CUDA’.

  • first_residue (int | None) – Starting residue for calculation. Defaults to None.

  • last_residue (int | None) – Ending residue for calculation. Defaults to None.

Example

>>> system = build_system(topology)
>>> ie = InteractionEnergyFrame(system, topology, chain="A")
>>> ie.compute(positions)

Initialize the InteractionEnergyFrame calculator.

Parameters:
  • system (System) – Pre-built OpenMM System object.

  • top (Topology) – OpenMM Topology object.

  • chain (str) – Chain identifier for calculation.

  • platform (str) – OpenMM platform name.

  • first_residue (int | None) – Starting residue for calculation.

  • last_residue (int | None) – Ending residue for calculation.

__init__(system, top, chain='A', platform='CUDA', first_residue=None, last_residue=None)[source]

Initialize the InteractionEnergyFrame calculator.

Parameters:
  • system (System) – Pre-built OpenMM System object.

  • top (Topology) – OpenMM Topology object.

  • chain (str) – Chain identifier for calculation.

  • platform (str) – OpenMM platform name.

  • first_residue (int | None) – Starting residue for calculation.

  • last_residue (int | None) – Ending residue for calculation.

get_system()[source]

Return the pre-built OpenMM system.

Sets self.selection via get_selection and returns the existing system object.

Return type:

System

Returns:

The pre-built OpenMM System object.

class molecular_simulations.analysis.interaction_energy.DynamicInteractionEnergy(top, traj, stride=1, chain='A', platform='CUDA', first_residue=None, last_residue=None, progress_bar=False)[source]

Bases: object

Compute interaction energies over a trajectory.

Uses InteractionEnergyFrame to run per-frame energy calculations and orchestrates trajectory operations.

system

OpenMM System object.

coordinates

Trajectory coordinate array.

stride

Frame stride for calculations.

energies

Energy array after compute_energies().

IE

InteractionEnergyFrame instance.

Parameters:
  • top (Path | str) – Path to prmtop topology file.

  • traj (Path | str) – Path to DCD trajectory file.

  • stride (int) – Stride for moving through trajectory. Defaults to 1.

  • chain (str) – Chain identifier for calculation. Defaults to ‘A’.

  • platform (str) – OpenMM platform name. Defaults to ‘CUDA’.

  • first_residue (int | None) – Starting residue for calculation. Defaults to None.

  • last_residue (int | None) – Ending residue for calculation. Defaults to None.

  • progress_bar (bool) – Whether to display a tqdm progress bar. Defaults to False.

Example

>>> die = DynamicInteractionEnergy("system.prmtop", "traj.dcd")
>>> die.compute_energies()
>>> print(die.energies.shape)  # (n_frames, 2)

Initialize the DynamicInteractionEnergy calculator.

Parameters:
  • top (Path | str) – Path to topology file.

  • traj (Path | str) – Path to trajectory file.

  • stride (int) – Frame stride.

  • chain (str) – Chain identifier.

  • platform (str) – OpenMM platform name.

  • first_residue (int | None) – Starting residue.

  • last_residue (int | None) – Ending residue.

  • progress_bar (bool) – Whether to show progress.

__init__(top, traj, stride=1, chain='A', platform='CUDA', first_residue=None, last_residue=None, progress_bar=False)[source]

Initialize the DynamicInteractionEnergy calculator.

Parameters:
  • top (Path | str) – Path to topology file.

  • traj (Path | str) – Path to trajectory file.

  • stride (int) – Frame stride.

  • chain (str) – Chain identifier.

  • platform (str) – OpenMM platform name.

  • first_residue (int | None) – Starting residue.

  • last_residue (int | None) – Ending residue.

  • progress_bar (bool) – Whether to show progress.

compute_energies()[source]

Compute energies for each frame in the trajectory.

Stores results in self.energies with shape (n_frames, 2) where columns are [LJ, Coulomb].

Return type:

None

build_system(top)[source]

Build an OpenMM system from the topology file.

Handles both PDB and prmtop topology files.

Parameters:

top (Path | str) – Path to topology file.

Return type:

System

Returns:

OpenMM System object.

Raises:

NotImplementedError – If topology file type is not supported.

load_traj(top, traj)[source]

Load trajectory into mdtraj and extract coordinates.

Parameters:
  • top (Path | str) – Path to topology file.

  • traj (Path | str) – Path to trajectory file.

Return type:

ndarray

Returns:

Coordinate array with shape (n_frames, n_atoms, 3).

setup_pbar()[source]

Build a tqdm progress bar for trajectory iteration.

Return type:

None