molecular_simulations.analysis.cov_ppi module

Covariance-based protein-protein interaction analysis.

This module provides tools for analyzing protein-protein interactions based on covariance analysis of molecular dynamics trajectories. Adapted from https://www.biorxiv.org/content/10.1101/2025.03.24.644990v1.full.pdf

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)[source]

Bases: object

Analyze protein-protein interactions using covariance analysis.

Takes an input topology and trajectory file, computes the covariance matrix between two selections, filters interactions by distance (11Å for positive covariance, 13Å for negative covariance), and evaluates each based on distance and angle cutoffs for various interaction types.

u

MDAnalysis Universe object.

n_frames

Number of frames in the trajectory.

out

Output path for results.

mapping

Residue index to resID mapping for both selections.

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

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

  • out (Path | str) – Path to output results file.

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

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

  • cov_cutoff (tuple[float, float]) – Distance cutoffs for positive and negative covariance filtering respectively. Defaults to (11.0, 13.0) Angstroms.

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

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

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

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

  • plot (bool) – Whether to generate and save plots. Defaults to True.

Example

>>> ppi = PPInteractions("complex.prmtop", "traj.dcd", "results.json")
>>> ppi.run()

Initialize the protein-protein interaction analyzer.

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

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

  • out (Path | str) – Path to output results file.

  • sel1 (str) – MDAnalysis selection string for first selection.

  • sel2 (str) – MDAnalysis selection string for second selection.

  • cov_cutoff (tuple[float, float]) – Tuple of distance cutoffs for (positive, negative) covariance.

  • sb_cutoff (float) – Salt bridge distance cutoff in Angstroms.

  • hbond_cutoff (float) – Hydrogen bond distance cutoff in Angstroms.

  • hbond_angle (float) – Hydrogen bond angle cutoff in degrees.

  • hydrophobic_cutoff (float) – Hydrophobic interaction cutoff in Angstroms.

  • plot (bool) – Whether to generate plots.

__init__(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)[source]

Initialize the protein-protein interaction analyzer.

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

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

  • out (Path | str) – Path to output results file.

  • sel1 (str) – MDAnalysis selection string for first selection.

  • sel2 (str) – MDAnalysis selection string for second selection.

  • cov_cutoff (tuple[float, float]) – Tuple of distance cutoffs for (positive, negative) covariance.

  • sb_cutoff (float) – Salt bridge distance cutoff in Angstroms.

  • hbond_cutoff (float) – Hydrogen bond distance cutoff in Angstroms.

  • hbond_angle (float) – Hydrogen bond angle cutoff in degrees.

  • hydrophobic_cutoff (float) – Hydrophobic interaction cutoff in Angstroms.

  • plot (bool) – Whether to generate plots.

run()[source]

Execute the full interaction analysis workflow.

Obtains a covariance matrix, screens for close interactions, evaluates each pairwise interaction, and reports contact probabilities. Optionally generates plots.

Return type:

None

compute_interactions(res1, res2)[source]

Compute interaction probabilities between two residues.

Generates MDAnalysis AtomGroups for each residue, identifies relevant non-bonded interactions (hydrogen bonds, salt bridges, hydrophobic), and computes the fraction of simulation time each interaction is engaged.

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

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

Return type:

dict[str, dict[str, float]]

Returns:

Nested dictionary containing the results of each interaction type. Keys are residue pair names, values are dictionaries mapping interaction type to probability.

get_covariance()[source]

Compute the positional covariance matrix between selections.

Loops over all C-alpha atoms and computes 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>

The magnitude indicates correlation strength and the sign indicates positive or negative correlation.

Return type:

ndarray

Returns:

Covariance matrix with shape (N_residues_sel1, N_residues_sel2).

res_map(ag1, ag2)[source]

Create mapping from covariance matrix indices to resIDs.

Maps covariance matrix indices to AtomGroup resIDs to ensure correct residue pairs are examined.

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

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

Return type:

None

interpret_covariance(cov_mat)[source]

Identify residue pairs with positive or negative correlations.

Parameters:

cov_mat (ndarray) – Covariance matrix from get_covariance().

Returns:

(positive_pairs, negative_pairs). Each pair is a tuple of (resID_sel1, resID_sel2).

Return type:

tuple[list[tuple[int, int]], list[tuple[int, int]]]

identify_interaction_type(res1, res2)[source]

Determine which analyses to compute for a residue pair.

Identifies what analyses to compute based on residue types (hydrophobic interactions, hydrogen bonds, salt bridges).

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.

Return type:

tuple[list[Callable[..., int | float]], list[str]]

Returns:

Tuple containing (list of function calls, list of labels).

analyze_saltbridge(res1, res2)[source]

Calculate salt bridge occupancy between two residues.

Uses a simple distance cutoff to determine salt bridge formation.

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

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

Return type:

float

Returns:

Proportion of simulation time spent in salt bridge contact.

analyze_hbond(res1, res2)[source]

Calculate hydrogen bond occupancy between two residues.

Identifies all potential donor/acceptor atoms, filters by distance, then evaluates each pair over the trajectory using distance and angle cutoffs.

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

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

Return type:

float

Returns:

Proportion of simulation time spent in hydrogen bond contact.

analyze_hydrophobic(res1, res2)[source]

Calculate hydrophobic interaction occupancy between residues.

Uses a simple distance cutoff between carbon atoms to determine hydrophobic contact.

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

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

Return type:

float

Returns:

Proportion of simulation time spent in hydrophobic contact.

survey_donors_acceptors(res1, res2)[source]

Identify potential hydrogen bond donors and acceptors.

First-pass distance threshold to identify potential hydrogen bonds. Should be followed by querying H-bond angles but this serves to reduce the search space.

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

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

Return type:

tuple[AtomGroup, AtomGroup]

Returns:

Tuple of (donors, acceptors) AtomGroups containing atoms that pass the crude distance cutoff.

evaluate_hbond(donor, acceptor)[source]

Evaluate hydrogen bond formation in the current frame.

Checks whether there is a defined hydrogen bond between any donor and acceptor atoms using distance and angle criteria. Returns early when a valid H-bond is detected.

Parameters:
  • donor (AtomGroup) – AtomGroup of potential H-bond donors.

  • acceptor (AtomGroup) – AtomGroup of potential H-bond acceptors.

Return type:

int

Returns:

1 if a valid hydrogen bond is found, else 0.

save(results)[source]

Save results to a JSON file.

Parameters:

results (dict[str, dict[str, dict[str, float]]]) – Dictionary of results to be saved.

Return type:

None

plot_results(results)[source]

Generate and save plots of the results.

Creates bar plots for each combination of covariance type (positive/negative) and interaction type (hydrophobic, hydrogen bond, salt bridge).

Parameters:

results (dict[str, dict[str, dict[str, float]]]) – Dictionary of results to be plotted.

Return type:

None

parse_results(results)[source]

Prepare results for plotting.

Removes entries with all-zero interactions and converts to a Polars DataFrame for easier plotting.

Parameters:

results (dict[str, dict[str, dict[str, float]]]) – Dictionary of results to be prepped.

Return type:

DataFrame

Returns:

Polars DataFrame with columns for residue pair, interaction probabilities, and covariance type.

make_plot(data, column, name, fs=15)[source]

Generate a bar plot for a specified interaction type.

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

  • column (str) – Column name for the interaction type to plot.

  • name (Path | str) – Path to save the plot.

  • fs (int) – Font size for plot labels. Defaults to 15.

Return type:

None