molecular_simulations.analysis.cov_ppi module

class molecular_simulations.analysis.cov_ppi.PPInteractions(top, traj, out, sel1='chainID A', sel2='chainID B', cov_cutoff=(11.0, 13.0), sb_cutoff=6.0, hbond_cutoff=3.5, hbond_angle=30.0, hydrophobic_cutoff=8.0, plot=True)

Bases: object

Code herein adapted from:

https://www.biorxiv.org/content/10.1101/2025.03.24.644990v1.full.pdf

Takes an input topology file and trajectory file, and highlights relevant interactions between two selections. To this end we first compute the covariance matrix between the two selections, filter out all interactions which occur too far apart (11Å for positive covariance, 13Å for negative covariance), and examines each based on a variety of distance and angle cutoffs defined in the literature.

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

  • traj (PathLike) – Path to trajectory file.

  • out (PathLike) – Path to outputs.

  • sel1 (str) – Defaults to ‘chainID A’. MDAnalysis selection string for the first selection.

  • sel2 (str) – Defaults to ‘chainID B’. MDAnalysis selection string for the second selection.

  • cov_cutoff (tuple[float]) – Defaults to (11., 13.). Tuple of the distance cutoffs to use for positive and negative covariance respectively.

  • sb_cutoff (float) – Defaults to 6.0Å. Distance cutoff for salt bridges.

  • hbond_cutoff (float) – Defaults to 3.5Å. Distance cutoff for hydrogen bonds.

  • hbond_angle (float) – Defaults to 30.0 degrees. Angle cutoff for hydrogen bonds.

  • hydrophobic_cutoff (float) – Defaults to 8.0Å. Distance cutoff for hydrophobic interactions.

  • plot (bool) – Defaults to True. Whether or not to plot results. Saves plots at the output directory.

analyze_hbond(res1, res2)

Identifies all potential donor/acceptor atoms between two residues. Culls this list based on distance array across simulation and then evaluates each pair over the trajectory utilizing a distance and angle cutoff.

Parameters:
  • res1 (mda.AtomGroup) – AtomGroup for a residue from selection 1.

  • res2 (mda.AtomGroup) – AtomGroup for a residue from selection 2.

Returns:

Proportion of simulation time spent in hydrogen bond.

Return type:

(float)

analyze_hydrophobic(res1, res2)

Uses a simple distance cutoff to highlight the occupancy of hydrophobic interaction between two residues. Returns the fraction of simulation time spent engaged in interaction.

Parameters:
  • res1 (mda.AtomGroup) – AtomGroup for a residue from selection 1.

  • res2 (mda.AtomGroup) – AtomGroup for a residue from selection 2.

Returns:

Proportion of simulation time spent in interaction.

Return type:

(float)

analyze_saltbridge(res1, res2)

Uses a simple distance cutoff to highlight the occupancy of saltbridge between two residues. Returns the fraction of simulation time spent engaged in saltbridge.

Parameters:
  • res1 (mda.AtomGroup) – AtomGroup for a residue from selection 1.

  • res2 (mda.AtomGroup) – AtomGroup for a residue from selection 2.

Returns:

Proportion of simulation time spent in salt bridge.

Return type:

(float)

compute_interactions(res1, res2)

Ingests two resIDs, generates MDAnalysis AtomGroups for each, identifies relevant non-bonded interactions (HBonds, saltbridge, hydrophobic) and computes each. Returns a dict containing the proportion of simulation time that each interaction is engaged.

Parameters:
  • res1 (int) – ResID for a residue in sel1.

  • res2 (int) – ResID for a residue in sel2.

Returns:

A nested dictionary containing the results of each interaction

type.

Return type:

(Results)

evaluate_hbond(donor, acceptor)

Evaluates whether there is a defined hydrogen bond between any donor and acceptor atoms in a given frame. Must pass a distance cutoff as well as an angle cutoff. Returns early when a legal HBond is detected.

Parameters:
  • donor (mda.AtomGroup) – AtomGroup of HBond donor.

  • acceptor (mda.AtomGroup) – AtomGroup of HBond acceptor.

Returns:

1 if legal hbond found, else 0

Return type:

(int)

get_covariance()

Loop over all C-alpha atoms and compute the positional covariance using the functional form:

C = <(R1 - <R1>)(R2 - <R2>)T>

where each element corresponds to the ensemble average movement

C_ij = <deltaR_i * deltaR_j>

with the magnitude being the strength of correlation and the sign corresponding to positive and negative correlation respectively.

Returns:

Covariance matrix.

Return type:

(np.ndarray)

identify_interaction_type(res1, res2)

Identifies what analyses to compute for a given pair of protein residues (i.e. hydrophobic interactions, hydrogen bonds, saltbridges).

Parameters:
  • res1 (str) – 3-letter code resname for a residue from selection 1.

  • res2 (str) – 3-letter code resname for a residue from selection 2.

Returns:

Tuple containing list of function calls and list of labels.

Return type:

(TaskTree)

interpret_covariance(cov_mat)

Identify pairs of residues with positive or negative correlations. Returns a tuple comprised of pairs for each.

Parameters:

cov_mat (np.ndarray) – Covariance matrix.

Returns:

Tuple of positively and negatively correlated

pairs of residues coming from each selection.

Return type:

(tuple[tuple[int, int]])

make_plot(data, column, name, fs=15)

Generates a seaborn barplot from a dataframe for a specified column.

Parameters:
  • data (pl.DataFrame) – Polars dataframe of data.

  • column (str) – Label for desired column.

  • name (PathLike) – Path to file to save plot to.

  • fs (int) – Defaults to 15. Size of font for plot.

Return type:

None

Returns:

None

parse_results(results)

Prepares results for plotting. Removes any entries which are all 0. and returns as a pandas DataFrame for easier plotting.

Parameters:

results (Results) – Dictionary of results to be prepped.

Returns:

Polars dataframe of results.

Return type:

(pl.DataFrame)

plot_results(results)

Plot results.

Parameters:

results (Results) – Dictionary of results to be plotted.

Return type:

None

Returns:

None

res_map(ag1, ag2)

Map covariance matrix indices to AtomGroup resIDs so that we are examining the correct pairs of residues.

Parameters:
  • ag1 (mda.AtomGroup) – AtomGroup of the first selection.

  • ag2 (mda.AtomGroup) – AtomGroup of the second selection.

Return type:

None

Returns:

None

run()

Main function that runs the workflow. Obtains a covariance matrix, screens for close interactions, evaluates each pairwise interaction for each amino acid and report the contact probability of each.

Return type:

None

Returns:

None

save(results)

Save results as a json file.

Parameters:

results (Results) – Dictionary of results to be saved.

Return type:

None

Returns:

None

survey_donors_acceptors(res1, res2)

First pass distance threshhold to identify potential Hydrogen bonds. Should be followed by querying HBond angles but this serves to reduce our search space and time complexity. Only returns donors/acceptors which are within the distance cutoff in at least a single frame.

Parameters:
  • res1 (mda.AtomGroup) – AtomGroup for a residue from selection 1.

  • res2 (mda.AtomGroup) – AtomGroup for a residue from selection 2.

Returns:

Tuple of AtomGroups for residues which pass

crude distance cutoff for hydrogen bond donors/acceptors.

Return type:

(tuple[mda.AtomGroup])