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).
- molecular_simulations.analysis.fingerprinter.dist_mat(xyz1, xyz2)¶
Compute distance matrix between two coordinate sets.
- 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:
- Return type:
- 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.
- 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:
- 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:
- 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:
- 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:
objectCalculate 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:
- __init__(topology, trajectory=None, target_selection='segid A', binder_selection=None, out_path=None, out_name=None)[source]¶
Initialize the Fingerprinter.
- 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:
- 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:
- 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:
- iterate_frames()[source]¶
Run fingerprint calculations over all trajectory frames.
Initializes fingerprint arrays and iterates through all frames, calculating fingerprints for each.
- Return type: