molecular_simulations.analysis.fingerprinter module

Interaction energy fingerprinting module.

This module calculates interaction energy fingerprints between target and binder chains in molecular structures, computing both electrostatic and Lennard-Jones contributions at the residue level.

molecular_simulations.analysis.fingerprinter.unravel_index(n1, n2)

Create unraveled indices for vectorized distance calculations.

Generates two arrays of indices that, when used together, represent all pairwise combinations of indices from ranges [0, n1) and [0, n2).

Parameters:
  • n1 (int) – Size of the first dimension.

  • n2 (int) – Size of the second dimension.

Return type:

tuple[ndarray, ndarray]

Returns:

Tuple of two arrays, each of shape (n1*n2,), containing the row and column indices respectively.

molecular_simulations.analysis.fingerprinter.dist_mat(xyz1, xyz2)

Compute distance matrix between two coordinate sets.

Parameters:
  • xyz1 (ndarray) – First coordinate array of shape (n1, ndim).

  • xyz2 (ndarray) – Second coordinate array of shape (n2, ndim).

Return type:

ndarray

Returns:

Distance matrix of shape (n1, n2) where element [i, j] is the Euclidean distance between xyz1[i] and xyz2[j].

molecular_simulations.analysis.fingerprinter.electrostatic(distance, charge_i, charge_j)

Calculate electrostatic energy between two particles.

Uses the reaction field method with a 10 Å cutoff and a solvent dielectric of 78.5 (water).

Parameters:
  • distance (float) – Distance between particles i and j in nm.

  • charge_i (float) – Charge of particle i in elementary charge units.

  • charge_j (float) – Charge of particle j in elementary charge units.

Return type:

float

Returns:

Electrostatic energy in kJ/mol.

Note

Conversion factors used: - Avogadro = 6.022e23 molecules/mol - e- to Coulomb = 1.602e-19 C/e- - nm to m = 1e-9 m/nm - 1/(4*pi*epsilon_0) = 8.988e9 J*m/C^2

molecular_simulations.analysis.fingerprinter.electrostatic_sum(distances, charge_is, charge_js)

Calculate sum of all electrostatic interactions between two groups.

Parameters:
  • distances (ndarray) – Distance matrix between particles with shape (len(charge_is), len(charge_js)).

  • charge_is (ndarray) – Array of charges for group i.

  • charge_js (ndarray) – Array of charges for group j.

Return type:

float

Returns:

Total electrostatic interaction energy in kJ/mol.

molecular_simulations.analysis.fingerprinter.lennard_jones(distance, sigma_i, sigma_j, epsilon_i, epsilon_j)

Calculate Lennard-Jones energy between two particles.

Uses standard Lorentz-Berthelot combining rules with a 12 Å cutoff.

Parameters:
  • distance (float) – Distance between particles i and j in nm.

  • sigma_i (float) – Sigma parameter for particle i in nm.

  • sigma_j (float) – Sigma parameter for particle j in nm.

  • epsilon_i (float) – Epsilon parameter for particle i in kJ/mol.

  • epsilon_j (float) – Epsilon parameter for particle j in kJ/mol.

Return type:

float

Returns:

Lennard-Jones interaction energy in kJ/mol.

molecular_simulations.analysis.fingerprinter.lennard_jones_sum(distances, sigma_is, sigma_js, epsilon_is, epsilon_js)

Calculate sum of all LJ interactions between two groups.

Parameters:
  • distances (ndarray) – Distance matrix between particles with shape (len(sigma_is), len(sigma_js)).

  • sigma_is (ndarray) – Array of sigma parameters for group i in nm.

  • sigma_js (ndarray) – Array of sigma parameters for group j in nm.

  • epsilon_is (ndarray) – Array of epsilon parameters for group i in kJ/mol.

  • epsilon_js (ndarray) – Array of epsilon parameters for group j in kJ/mol.

Return type:

float

Returns:

Total Lennard-Jones interaction energy in kJ/mol.

molecular_simulations.analysis.fingerprinter.fingerprints(xyzs, charges, sigmas, epsilons, target_resmap, binder_inds)

Calculate per-residue interaction fingerprints.

Computes electrostatic and Lennard-Jones energies between each target residue and all binder residues.

Parameters:
  • xyzs (ndarray) – Coordinate array for all atoms in nm.

  • charges (ndarray) – Array of partial charges for all atoms.

  • sigmas (ndarray) – Array of sigma parameters for all atoms in nm.

  • epsilons (ndarray) – Array of epsilon parameters for all atoms in kJ/mol.

  • target_resmap (list) – List of arrays, each containing atom indices for a target residue.

  • binder_inds (ndarray) – Array of atom indices for the binder.

Return type:

tuple[ndarray, ndarray]

Returns:

Tuple of (lj_fingerprint, es_fingerprint), each of shape (n_target_residues,) containing the per-residue interaction energies.

class molecular_simulations.analysis.fingerprinter.Fingerprinter(topology, trajectory=None, target_selection='segid A', binder_selection=None, out_path=None, out_name=None)[source]

Bases: object

Calculate interaction energy fingerprints between target and binder.

Computes per-residue Lennard-Jones and electrostatic interaction energies between a target selection and a binder selection over a molecular dynamics trajectory.

topology

Path to the topology file.

trajectory

Path to the trajectory or coordinate file.

target_selection

MDAnalysis selection string for the target.

binder_selection

MDAnalysis selection string for the binder.

out

Output path for the fingerprint data.

target_fingerprint

Target fingerprint array after run().

binder_fingerprint

Binder fingerprint array after run().

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

  • trajectory (Path | str | None) – Path to trajectory or coordinate file. If None, will look for inpcrd or rst7 with same stem as topology.

  • target_selection (str) – MDAnalysis selection string for target. Defaults to ‘segid A’.

  • binder_selection (str | None) – MDAnalysis selection string for binder. If None, binder is defined as everything not in target.

  • out_path (Path | str | None) – Output directory path. If None, uses topology parent.

  • out_name (str | None) – Output filename. If None, uses ‘fingerprint.npz’.

Example

>>> fp = Fingerprinter("complex.prmtop", "traj.dcd")
>>> fp.run()
>>> fp.save()

Initialize the Fingerprinter.

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

  • trajectory (Path | str | None) – Path to trajectory or coordinate file.

  • target_selection (str) – MDAnalysis selection for target.

  • binder_selection (str | None) – MDAnalysis selection for binder.

  • out_path (Path | str | None) – Output directory path.

  • out_name (str | None) – Output filename.

__init__(topology, trajectory=None, target_selection='segid A', binder_selection=None, out_path=None, out_name=None)[source]

Initialize the Fingerprinter.

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

  • trajectory (Path | str | None) – Path to trajectory or coordinate file.

  • target_selection (str) – MDAnalysis selection for target.

  • binder_selection (str | None) – MDAnalysis selection for binder.

  • out_path (Path | str | None) – Output directory path.

  • out_name (str | None) – Output filename.

assign_nonbonded_params()[source]

Extract nonbonded parameters from the topology.

Builds an OpenMM system from the topology and extracts charge, sigma, and epsilon parameters for all particles.

Return type:

None

load_pdb()[source]

Load the structure into an MDAnalysis Universe.

Can handle either PDB or AMBER prmtop files. For prmtop files, will look for inpcrd or rst7 coordinates if trajectory was not specified.

Return type:

None

assign_residue_mapping()[source]

Create mappings from residue indices to atom indices.

Creates residue-to-atom mappings for both target and binder selections.

Return type:

None

iterate_frames()[source]

Run fingerprint calculations over all trajectory frames.

Initializes fingerprint arrays and iterates through all frames, calculating fingerprints for each.

Return type:

None

calculate_fingerprints(frame_index)[source]

Calculate fingerprints for a single frame.

Parameters:

frame_index (int) – Index of the current frame (may differ from frame number if trajectory is discontinuous).

Return type:

None

run()[source]

Execute the complete fingerprinting workflow.

Obtains parameters, loads the structure, assigns residue mappings, and iterates through the trajectory to compute fingerprints.

Return type:

None

save()[source]

Save fingerprint data to an NPZ file.

Saves both target and binder fingerprints to the configured output path.

Return type:

None