Interactions and Stability (InSty)¶
This module called InSty defines functions for calculating different types of interactions in protein structure, between proteins or between protein and ligand. The following interactions are available for protein interactions:
- Hydrogen bonds
- Salt Bridges
- Repulsive Ionic Bonding
- Pi stacking interactions
- Pi-cation interactions
- Hydrophobic interactions
- Disulfide Bonds
-
class
Interactions
(title='Unknown')[source]¶ Class for Interaction analysis of proteins.
-
buildInteractionMatrix
(**kwargs)[source]¶ - Build matrix with protein interactions. Interactions are counted as follows:
- Hydrogen bonds (HBs) +1
- Salt Bridges (SBs) +1 (Salt bridges might be included in hydrogen bonds)
- Repulsive Ionic Bonding (RIB) +1
- Pi stacking interactions (PiStack) +1
- Pi-cation interactions (PiCat) +1
- Hydrophobic interactions (HPh) +1
- Disulfide bonds (DiBs) +1
Some type of interactions could be removed from the analysis by replacing value 1 by 0. Example: >>> interactions = Interactions() >>> atoms = parsePDB(PDBfile).select(‘protein’) >>> interactions.calcProteinInteractions(atoms) >>> interactions.buildInteractionMatrix(RIB=0, HBs=0, HPh=0)
Parameters: - HBs (int, float) – score per single hydrogen bond
- SBs (int, float) – score per single salt bridge
- RIB (int, float) – score per single repulsive ionic bonding
- PiStack (int, float) – score per pi-stacking interaction
- PiCat (int, float) – score per pi-cation interaction
- HPh (int, float) – score per hydrophobic interaction
- DiBs (int, float) – score per disulfide bond
-
buildInteractionMatrixEnergy
(**kwargs)[source]¶ Build matrix with interaction energy comming from energy of pairs of specific residues.
Parameters: energy_list_type (str) – name of the list with energies Default is ‘IB_solv’ acceptable values are ‘IB_nosolv’, ‘IB_solv’, ‘CS’ ‘IB_solv’ and ‘IB_nosolv’ are derived from empirical potentials from O Keskin, I Bahar and colleagues from [OK98] and have RT units.
‘CS’ is from MD simulations of amino acid pairs from Carlos Simmerling and Gary Wu in the InSty paper (under preparation) and have units of kcal/mol.
-
calcProteinInteractions
(atoms, **kwargs)[source]¶ - Compute all protein interactions (shown below) using default parameters.
- Hydrogen bonds
- Salt Bridges
- RepulsiveIonicBonding
- Pi stacking interactions
- Pi-cation interactions
- Hydrophobic interactions
- Disulfide Bonds
Parameters: atoms ( Atomic
) – an Atomic object from which residues are selected
-
getDisulfideBonds
(**kwargs)[source]¶ Returns the list of disulfide bonds.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getFrequentInteractors
(contacts_min=2)[source]¶ Provide a list of residues with the most frequent interactions based on the following interactions:
- Hydrogen bonds (hb)
- Salt Bridges (sb)
- Repulsive Ionic Bonding (rb)
- Pi stacking interactions (ps)
- Pi-cation interactions (pc)
- Hydrophobic interactions (hp)
- Disulfide bonds (disb)
Parameters: contacts_min (int) – Minimal number of contacts which residue may form with other residues, Default is 2.
-
getHydrogenBonds
(**kwargs)[source]¶ Returns the list of hydrogen bonds.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getHydrophobic
(**kwargs)[source]¶ Returns the list of hydrophobic interactions.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getInteractions
(**kwargs)[source]¶ Returns the list of all interactions.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getInteractors
(residue_name)[source]¶ Provide information about interactions for a particular residue
Parameters: residue_name (str) – name of a resiude example: LEU234A, where A is a chain name
-
getPiCation
(**kwargs)[source]¶ Returns the list of Pi-cation interactions.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getPiStacking
(**kwargs)[source]¶ Returns the list of Pi-stacking interactions.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getRepulsiveIonicBonding
(**kwargs)[source]¶ Returns the list of repulsive ionic bonding.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getSaltBridges
(**kwargs)[source]¶ Returns the list of salt bridges.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
saveInteractionsPDB
(**kwargs)[source]¶ Save the number of potential interactions to PDB file in occupancy column.
Parameters:
-
setNewDisulfideBonds
(interaction)[source]¶ Replace default calculation of hydrophobic interactions by the one provided by user.
-
setNewHydrogenBonds
(interaction)[source]¶ Replace default calculation of hydrogen bonds by the one provided by user.
-
setNewHydrophobic
(interaction)[source]¶ Replace default calculation of hydrophobic interactions by the one provided by user.
-
setNewPiCation
(interaction)[source]¶ Replace default calculation of pi-cation interactions by the one provided by user.
-
setNewPiStacking
(interaction)[source]¶ Replace default calculation of pi-stacking interactions by the one provided by user.
-
setNewRepulsiveIonicBonding
(interaction)[source]¶ Replace default calculation of repulsive ionic bonding by the one provided by user.
-
setNewSaltBridges
(interaction)[source]¶ Replace default calculation of salt bridges by the one provided by user.
-
showCumulativeInteractionTypes
(**kwargs)[source]¶ Bar plot with the number of potential inetractions per residue. Particular type of interactions can be excluded by using keywords HBs=0, RIB=0, etc.
Parameters: - HBs (int, float) – score per single hydrogen bond
- SBs (int, float) – score per single salt bridge
- RIB (int, float) – score per single repulsive ionic bonding
- PiStack (int, float) – score per pi-stacking interaction
- PiCat (int, float) – score per pi-cation interaction
- HPh (int, float) – score per hydrophobic interaction
- DiBs (int, float) – score per disulfide bond
- selstr – selection string for focusing the plot
- energy (bool) – sum of the energy between residues default is False
- energy_list_type (str) – name of the list with energies default is ‘IB_solv’ acceptable values are ‘IB_nosolv’, ‘IB_solv’, ‘CS’
- overwrite_energies (bool) – whether to overwrite energies default is False
- percentile (float) – a percentile threshold to remove outliers, i.e. only showing data within p-th to 100-p-th percentile. Default is None, so no axis limits.
- vmin (float) – a minimum value threshold to remove outliers, i.e. only showing data greater than vmin This overrides percentile. Default is None, so no axis limits and little padding at the bottom when energy=True.
- vmax (float) – a maximum value threshold to remove outliers, i.e. only showing data less than vmax This overrides percentile. Default is None, so no axis limits and a little padding for interaction type labels.
‘IB_solv’ and ‘IB_nosolv’ are derived from empirical potentials from O Keskin, I Bahar and colleagues from [OK98] and have RT units.
‘CS’ is from MD simulations of amino acid pairs from Carlos Simmerling and Gary Wu for the InSty paper (under preparation) and have units kcal/mol.
-
showFrequentInteractors
(cutoff=4, **kwargs)[source]¶ Plots regions with the most frequent interactions.
Parameters: cutoff (int, float) – minimal score per residue which will be displayed. If cutoff value is too big, top 30% with the higest values will be returned. Default is 4. Nonstandard resiudes can be updated in a following way: d = {‘CYX’: ‘X’, ‘CEA’: ‘Z’} >>> name.showFrequentInteractors(d)
-
-
class
InteractionsTrajectory
(name='Unknown')[source]¶ Class for Interaction analysis of DCD trajectory or multi-model PDB (Ensemble PDB).
-
calcProteinInteractionsTrajectory
(atoms, trajectory=None, filename=None, **kwargs)[source]¶ - Compute all protein interactions (shown below) for DCD trajectory or multi-model PDB
- using default parameters. (1) Hydrogen bonds (2) Salt Bridges (3) RepulsiveIonicBonding (4) Pi stacking interactions (5) Pi-cation interactions (6) Hydrophobic interactions (7) Disulfide Bonds
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - trajectory (class:.Trajectory) – trajectory file
- filename (pkl) – Name of pkl filename in which interactions will be storage
- selection (str) – selection string
- selection2 (str) – selection string
- start_frame (int) – index of first frame to read
- stop_frame (int) – index of last frame to read
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getDisulfideBonds
(**kwargs)[source]¶ Return the list of disulfide bonds computed from DCD trajectory.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getHydrogenBonds
(**kwargs)[source]¶ Return the list of hydrogen bonds computed from DCD trajectory.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getHydrophobic
(**kwargs)[source]¶ Return the list of hydrophobic interactions computed from DCD trajectory.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getInteractions
(**kwargs)[source]¶ Return the list of all interactions.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getPiCation
(**kwargs)[source]¶ Return the list of Pi-cation interactions computed from DCD trajectory.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getPiStacking
(**kwargs)[source]¶ Return the list of Pi-stacking interactions computed from DCD trajectory.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getRepulsiveIonicBonding
(**kwargs)[source]¶ Return the list of repulsive ionic bonding computed from DCD trajectory.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getSaltBridges
(**kwargs)[source]¶ Return the list of salt bridges computed from DCD trajectory.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
getTimeInteractions
(filename=None, **kwargs)[source]¶ Return a bar plots with the number of interactions per each frame.
Parameters: filename (str) – PNG file name
-
parseInteractions
(filename)[source]¶ Import interactions from analysis of trajectory which was saved via calcProteinInteractionsTrajectory().
Parameters: filename (pkl) – Name of pkl file in which interactions will be storage
-
setNewDisulfideBondsTrajectory
(interaction)[source]¶ Replace default calculation of disulfide bonds by the one provided by user.
-
setNewHydrogenBondsTrajectory
(interaction)[source]¶ Replace default calculation of hydrogen bonds by the one provided by user.
-
setNewHydrophobicTrajectory
(interaction)[source]¶ Replace default calculation of hydrophobic interactions by the one provided by user.
-
setNewPiCationTrajectory
(interaction)[source]¶ Replace default calculation of pi-cation interactions by the one provided by user.
-
setNewPiStackingTrajectory
(interaction)[source]¶ Replace default calculation of pi-stacking interactions by the one provided by user.
-
-
class
LigandInteractionsTrajectory
(name='Unknown')[source]¶ Class for protein-ligand interaction analysis of DCD trajectory or multi-model PDB (Ensemble PDB). This class is using PLIP to provide the interactions. Install PLIP before using it.
## Instalation of PLIP using conda: >>> conda install -c conda-forge plip ## https://anaconda.org/conda-forge/plip # https://github.com/pharmai/plip/blob/master/DOCUMENTATION.md
## Instalation using PIP: >>> pip install plip
[SS15] Salentin S., Schreiber S., Haupt V. J., Adasme M. F., Schroeder M. PLIP: fully automated protein–ligand interaction profiler Nucl. Acids Res. 2015 43:W443-W447.
-
calcFrequentInteractors
(**kwargs)[source]¶ Returns a dictonary with residues involved in the interaction with ligand and their number of counts.
Parameters: selection (str) – selection string of ligand with chain ID e.g. “MESA” where MES is ligand resname and A is chain ID. Selection pointed as None will return all interactions together without ligands separation.
-
calcLigandInteractionsTrajectory
(atoms, trajectory=None, **kwargs)[source]¶ - Compute protein-ligand interactions for DCD trajectory or multi-model PDB
- using PLIP library.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - trajectory (class:.Trajectory) – trajectory file
- filename (pkl) – Name of pkl filename in which interactions will be storage
- start_frame (int) – index of first frame to read
- stop_frame (int) – index of last frame to read
- output (bool) – parameter to print the interactions on the screen while analyzing the structure, use ‘info’.
-
getLigandInteractions
(**kwargs)[source]¶ Return the list of protein-ligand interactions.
Parameters: - filters (str) – selection string of ligand with chain ID or interaction type e.g. ‘SBs’ (HBs, SBs, HPh, PiStack, PiCat, HPh, watBridge)
- include_frames (bool) – used with filters, it will leave selected keyword in orginal lists, if False it will collect selected interactions in one list, Use True to assign new selection using setLigandInteractions. Default is True
-
getLigandInteractionsNumber
(**kwargs)[source]¶ Return the number of interactions per each frame. Number of interactions can be a total number of interactions or it can be divided into interaction types.
Parameters: types (bool) – Interaction types can be included (True) or not (False). Default is True.
-
parseLigandInteractions
(filename)[source]¶ Import interactions from analysis of trajectory which was saved via calcLigandInteractionsTrajectory().
Parameters: filename (pkl) – Name of pkl file in which interactions will be storage
-
-
calcHydrogenBonds
(atoms, **kwargs)[source]¶ Compute hydrogen bonds for proteins and other molecules.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - distDA (int, float) – non-zero value, maximal distance between donor and acceptor. default is 3.5 distA also works
- angleDHA (int, float) – non-zero value, maximal (180 - D-H-A angle) (donor, hydrogen, acceptor). default is 40. angle and angleDA also work
- seq_cutoff_HB (int) – non-zero value, interactions will be found between atoms with index differences that are higher than seq_cutoff_HB. default is 25 atoms. seq_cutoff also works
- donors (list) – which atoms to count as donors default is [‘N’, ‘O’, ‘S’, ‘F’]
- acceptors (list) – which atoms to count as acceptors default is [‘N’, ‘O’, ‘S’, ‘F’]
- selection (str) – selection string
- selection2 (str) – selection string
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
Structure should contain hydrogens. If not they can be added using addMissingAtoms(pdb_name) function available in ProDy after Openbabel installation. conda install -c conda-forge openbabel
Note that the angle which it is considering is 180-defined angle D-H-A (in a good agreement with VMD) Results can be displayed in VMD by using showVMDinteraction()
- atoms (
-
calcChHydrogenBonds
(atoms, **kwargs)[source]¶ Finds hydrogen bonds between different chains. See more details in calcHydrogenBonds().
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - distA (int, float) – non-zero value, maximal distance between donor and acceptor. default is 3.5
- angle (int, float) – non-zero value, D-H-A angle (donor, hydrogen, acceptor). default is 40.
- seq_cutoff (int) – non-zero value, interactions will be found between atoms with index differences that are higher than seq_cutoff. default is 25 atoms.
Structure should contain hydrogens. If not they can be added using addMissingAtoms(pdb_name) function available in ProDy after Openbabel installation. conda install -c conda-forge openbabel
Note that the angle which it is considering is 180-defined angle D-H-A (in a good agreement with VMD) Results can be displayed in VMD.
- atoms (
-
calcSaltBridges
(atoms, **kwargs)[source]¶ Finds salt bridges in protein structure. Histidines are considered as charge residues.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - distSB (int, float) – non-zero value, maximal distance between center of masses of N and O atoms of negatively and positevely charged residues. default is 5. distA also works
- selection (str) – selection string
- selection2 (str) – selection string
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
Results can be displayed in VMD.
- atoms (
-
calcRepulsiveIonicBonding
(atoms, **kwargs)[source]¶ Finds repulsive ionic bonding in protein structure i.e. between positive-positive or negative-negative residues. Histidine is not considered as a charged residue.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
Results can be displayed in VMD.
-
calcPiStacking
(atoms, **kwargs)[source]¶ Finds π–π stacking interactions (between aromatic rings).
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - distPS (int, float) – non-zero value, maximal distance between center of masses of residues aromatic rings. default is 5. distA works too
- angle_min_PS (int, float) – minimal angle between aromatic rings. default is 0. angle_min works too
- angle_max_PS (int, float) – maximal angle between rings. default is 360. angle_max works too
- selection (str) – selection string
- selection2 (str) – selection string
- non_standard_PS (dict) – dictionary of non-standard residue in the protein structure that need to be included in calculations non_standard works too
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
Results can be displayed in VMD. By default three residues are included TRP, PHE, TYR and HIS. Additional selection can be added:
>>> non_standard = {"HSE": "noh and not backbone and not name CB", "HSD": "noh and not backbone and not name CB"} >>> calcPiStacking(atoms, non_standard)
Predictions for proteins only. To compute protein-ligand interactions use calcLigandInteractions() or define **kwargs.
- atoms (
-
calcPiCation
(atoms, **kwargs)[source]¶ Finds cation-Pi interaction i.e. between aromatic ring and positively charged residue (ARG and LYS).
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - distPC (int, float) – non-zero value, maximal distance between center of masses of aromatic ring and positively charge group. default is 5. distA works too
- selection (str) – selection string
- selection2 (str) – selection string
- non_standard_PC (dict) – dictionary of non-standard residue in the protein structure that need to be included in calculations non_standard also works
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
By default three residues are included TRP, PHE, TYR and HIS. Additional selection can be added:
>>> calcPiCation(atoms, 'HSE'='noh and not backbone and not name CB') or >>> non_standard = {"HSE": "noh and not backbone and not name CB", "HSD": "noh and not backbone and not name CB"} >>> calcPiCation(atoms, non_standard)
Results can be displayed in VMD. Predictions for proteins only. To compute protein-ligand interactions use calcLigandInteractions() or define **kwargs
- atoms (
-
calcHydrophobic
(atoms, **kwargs)[source]¶ Prediction of hydrophobic interactions between hydrophobic residues (ALA, ILE, LEU, MET, PHE, TRP, VAL, CG of ARG, and CG and CD of LYS).
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - distHPh (int, float) – non-zero value, maximal distance between atoms of hydrophobic residues. default is 4.5. distA works too
- non_standard_Hph (dict) – dictionary of non-standard residue in the protein structure that need to be included in calculations non_standard works too
- zerosHPh (bool) – zero values of hydrophobic overlapping areas included default is False
Last value in the output corresponds to the total hydrophobic overlapping area for two residues not only for the atoms that are included in the list. Atoms that which are listed are the closest between two residues and they will be inluded to draw the line in VMD.
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
Additional selection can be added as shown below (with selection that includes only hydrophobic part):
>>> calcHydrophobic(atoms, non_standard_Hph={'XLE'='noh and not backbone', 'XLI'='noh and not backbone'})
Predictions for proteins only. To compute protein-ligand interactions use calcLigandInteractions(). Results can be displayed in VMD by using showVMDinteraction()
Note that interactions between aromatic residues are omitted becasue they are provided by calcPiStacking().
- atoms (
-
calcDisulfideBonds
(atoms, **kwargs)[source]¶ Prediction of disulfide bonds.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - distDB (int, float) – non-zero value, maximal distance between atoms of cysteine residues. default is 3. distA works too
- atoms (
-
calcMetalInteractions
(atoms, distA=3.0, extraIons=['FE'], excluded_ions=['SOD', 'CLA'])[source]¶ Interactions with metal ions (includes water, ligands and other ions).
Parameters:
-
calcHydrogenBondsTrajectory
(atoms, trajectory=None, **kwargs)[source]¶ Compute hydrogen bonds for DCD trajectory or multi-model PDB using default parameters.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - trajectory (class:.Trajectory) – trajectory file
- distA (int, float) – non-zero value, maximal distance between donor and acceptor. default is 3.5 distDA also works
- angle (int, float) – non-zero value, maximal (180 - D-H-A angle) (donor, hydrogen, acceptor). default is 40. angleDHA also works
- seq_cutoff (int) – non-zero value, interactions will be found between atoms with index differences that are higher than seq_cutoff. default is 20 atoms.
- selection (str) – selection string
- selection2 (str) – selection string
- start_frame (int) – index of first frame to read
- stop_frame (int) – index of last frame to read
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
- atoms (
-
calcSaltBridgesTrajectory
(atoms, trajectory=None, **kwargs)[source]¶ Compute salt bridges for DCD trajectory or multi-model PDB using default parameters.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - trajectory (class:.Trajectory) – trajectory file
- distA (int, float) – non-zero value, maximal distance between center of masses of N and O atoms of negatively and positevely charged residues. default is 5. distSB also works
- selection (str) – selection string
- selection2 (str) – selection string
- start_frame (int) – index of first frame to read
- stop_frame (int) – index of last frame to read
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
- atoms (
-
calcRepulsiveIonicBondingTrajectory
(atoms, trajectory=None, **kwargs)[source]¶ Compute repulsive ionic bonding for DCD trajectory or multi-model PDB using default parameters.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - trajectory (class:.Trajectory) – trajectory file
- distA (int, float) – non-zero value, maximal distance between center of masses between N-N or O-O atoms of residues. default is 4.5. distRB also works
- selection (str) – selection string
- selection2 (str) – selection string
- start_frame (int) – index of first frame to read
- stop_frame (int) – index of last frame to read
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
- atoms (
-
calcPiStackingTrajectory
(atoms, trajectory=None, **kwargs)[source]¶ Compute Pi-stacking interactions for DCD trajectory or multi-model PDB using default parameters.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - trajectory (class:.Trajectory) – trajectory file
- distA (int, float) – non-zero value, maximal distance between center of masses of residues aromatic rings. default is 5. distPS also works
- angle_min (int) – minimal angle between aromatic rings. default is 0.
- angle_max (int, float) – maximal angle between rings. default is 360.
- selection (str) – selection string
- selection2 (str) – selection string
- start_frame (int) – index of first frame to read
- stop_frame (int) – index of last frame to read
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
- atoms (
-
calcPiCationTrajectory
(atoms, trajectory=None, **kwargs)[source]¶ Compute Pi-cation interactions for DCD trajectory or multi-model PDB using default parameters.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - trajectory (class:.Trajectory) – trajectory file
- distA (int, float) – non-zero value, maximal distance between center of masses of aromatic ring and positively charge group. default is 5. distPC also works
- selection (str) – selection string
- selection2 (str) – selection string
- start_frame (int) – index of first frame to read
- stop_frame (int) – index of last frame to read
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
- atoms (
-
calcHydrophobicTrajectory
(atoms, trajectory=None, **kwargs)[source]¶ Compute hydrophobic interactions for DCD trajectory or multi-model PDB using default parameters.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - trajectory (class:.Trajectory) – trajectory file
- distA (int, float) – non-zero value, maximal distance between atoms of hydrophobic residues. default is 4.5. distHPh also works
- selection (str) – selection string
- selection2 (str) – selection string
- start_frame (int) – index of first frame to read
- stop_frame (int) – index of last frame to read
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
- atoms (
-
calcDisulfideBondsTrajectory
(atoms, trajectory=None, **kwargs)[source]¶ Compute disulfide bonds for DCD trajectory or multi-model PDB using default parameters.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - trajectory (class:.Trajectory) – trajectory file
- distA (int, float) – non-zero value, maximal distance between atoms of cysteine residues. default is 2.5. distDB also works
- start_frame (int) – index of first frame to read
- stop_frame (int) – index of last frame to read
- atoms (
-
calcProteinInteractions
(atoms, **kwargs)[source]¶ - Compute all protein interactions (shown below).
- Hydrogen bonds
- Salt Bridges
- RepulsiveIonicBonding
- Pi stacking interactions
- Pi-cation interactions
- Hydrophobic interactions
- Disulfide Bonds
kwargs can be passed on to the underlying functions as described in their documentation. For example, distDA and angleDHA can be used to control hydrogen bonds, or distA and angle can be used across types.
Parameters: Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
-
calcStatisticsInteractions
(data, **kwargs)[source]¶ Return the statistics of interactions from PDB Ensemble or trajectory including: (1) the weight for each residue pair: corresponds to the number of counts divided by the number of frames (values >1 are obtained when the residue pair creates multiple contacts); (2) average distance of interactions for each pair [in Ang], (3) standard deviation [Ang.], (4) Energy [in kcal/mol] that is not distance dependent. Energy by default is solvent-mediated from [OK98] (‘IB_solv’) in RT units. To use non-solvent-mediated (residue-mediated) entries (‘IB_nosolv’) from [OK98] in RT units or solvent-mediated values obtained from MD for InSty paper (‘CS’, under preparation) in kcal/mol, change energy_list_type parameter. If energy information is not available, please check whether the pair of residues is listed in the “tabulated_energies.txt” file, which is localized in the ProDy directory.
Parameters: - data (list) – list with interactions from calcHydrogenBondsTrajectory() or other types
- weight_cutoff (int, float) – value above which results will be displayed 1 or more means that residue contact is present in all conformations/frames default value is 0.2 (in 20% of conformations contact appeared)
- energy_list_type ('IB_nosolv', 'IB_solv', 'CS') – name of the list with energies default is ‘IB_solv’
Example of usage: >>> atoms = parsePDB(‘PDBfile.pdb’) >>> dcd = Trajectory(‘DCDfile.dcd’) >>> dcd.link(atoms) >>> dcd.setCoords(atoms)
>>> data = calcPiCationTrajectory(atoms, dcd, distA=5) or >>> interactionsTrajectory = InteractionsTrajectory() >>> data = interactionsTrajectory.getPiCation()
-
calcDistribution
(interactions, residue1, residue2=None, **kwargs)[source]¶ Distributions/histograms of pairs of amino acids. Histograms are normalized.
Parameters:
-
calcSASA
(atoms, **kwargs)[source]¶ Provide information about solvent accessible surface area (SASA) based on the sasa program. To use this function compiled hpb.so is needed.
Parameters:
-
calcVolume
(atoms, **kwargs)[source]¶ Provide information about volume for each residue/molecule/chain or other selection”. To use this function compiled hpb.so is needed.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - selection (str) – selection string Default all the protein structure is used
- volume_cutoff – cutoff for volume Default is 0.0 to include all the results
- split_residues (bool) – it will provide values for each residue Default is False
- resnames (bool) – residues name included Default is False True - will give residue names and values for each residue False - will give only the values for each residues
- atoms (
-
compareInteractions
(data1, data2, **kwargs)[source]¶ Comparison of two outputs from interactions. It will provide information about the disappearance and appearance of some interactions as well as the similarities in the interactions for the same system. Two conformations can be compared.
Parameters: type filename: str
Example of usage: >>> atoms1 = parsePDB(‘PDBfile1.pdb’).select(‘protein’) >>> atoms2 = parsePDB(‘PDBfile2.pdb’).select(‘protein’) >>> compareInteractions(calcHydrogenBonds(atoms1), calcHydrogenBonds(atoms2))
-
showInteractionsGraph
(statistics, **kwargs)[source]¶ Return residue-residue interactions as graph/network.
Parameters: - statistics (list) – Results obtained from calcStatisticsInteractions analysis
- cutoff (int, float) – Minimal number of counts per residue in the trajectory Default is 0.1.
- code (str) – representation of the residues, 3-letter or 1-letter Default is 3-letter
- multiple_chains (bool) – display chain name for structure with many chains Default is False
- edge_cmap (str) – color of the residue connection Default is plt.cm.Blues (blue color)
- node_size (int) – size of the nodes which describes residues Default is 300
- node_distance (int) – value which will scale residue-residue interactions Default is 5
- font_size (int) – size of the font Default is 14
- seed (int) – random number which affect the distribution of residues Default is 42
-
showInteractionsHist
(statistics, **kwargs)[source]¶ Return information about residue-residue interactions as a bar plot.
Parameters:
-
calcLigandInteractions
(atoms, **kwargs)[source]¶ Provide ligand interactions with other elements of the system including protein, water and ions. Results are computed by PLIP [SS15] which should be installed. Note that PLIP will not recognize ligand unless it will be HETATM in the PDB file.
Parameters: To display results as a list of interactions use listLigandInteractions() and for visualization in VMD please use showLigandInteraction_VMD()
Requirements of usage: ## Instalation of Openbabel: >>> conda install -c conda-forge openbabel ## https://anaconda.org/conda-forge/openbabel
## Instalation of PLIP >>> conda install -c conda-forge plip ## https://anaconda.org/conda-forge/plip # https://github.com/pharmai/plip/blob/master/DOCUMENTATION.md
[SS15] Salentin S., Schreiber S., Haupt V. J., Adasme M. F., Schroeder M. PLIP: fully automated protein–ligand interaction profiler Nucl. Acids Res. 2015 43:W443-W447.
-
listLigandInteractions
(PLIP_output, **kwargs)[source]¶ Create a list of interactions from PLIP output created using calcLigandInteractions(). Results can be displayed in VMD.
Parameters: - PLIP_output (PLIP object obtained from calcLigandInteractions()) – Results from PLIP for protein-ligand interactions.
- output (bool) – parameter to print the interactions on the screen while analyzing the structure (True | False) Default is False
Note that five types of interactions are considered: hydrogen bonds, salt bridges, pi-stacking, cation-pi, hydrophobic and water bridges.
-
showProteinInteractions_VMD
(atoms, interactions, color='red', **kwargs)[source]¶ Save information about protein interactions to a TCL file (filename) which can be further use in VMD to display all intercations in a graphical interface (in TKConsole: play script_name.tcl). Different types of interactions can be saved separately (color can be selected) or all at once for all types of interactions (hydrogen bonds - blue, salt bridges - yellow, pi stacking - green, cation-pi - orangem, hydrophobic - silver, and disulfide bonds - black).
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - interactions (List of lists) – List of interactions for protein interactions.
- color (str) – color to draw interactions in VMD, not used only for single interaction type. default “red”
- filename (str) – name of TCL file where interactions will be saved.
Example (hydrogen bonds for protein only): >>> interactions = calcHydrogenBonds(atoms.protein, distA=3.2, angle=30) or all interactions at once: >>> interactions = calcProteinInteractions(atoms.protein)
- atoms (
-
showLigandInteraction_VMD
(atoms, interactions, **kwargs)[source]¶ Save information from PLIP for ligand-protein interactions in a TCL file which can be further used in VMD to display all intercations in a graphical interface (hydrogen bonds - blue, salt bridges - yellow, pi stacking - green, cation-pi - orange, hydrophobic - silver and water bridges - cyan).
Parameters: To obtain protein-ligand interactions: >>> calculations = calcLigandInteractions(atoms) >>> interactions = listLigandInteractions(calculations)
-
calcHydrogenBondsTrajectory
(atoms, trajectory=None, **kwargs)[source] Compute hydrogen bonds for DCD trajectory or multi-model PDB using default parameters.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - trajectory (class:.Trajectory) – trajectory file
- distA (int, float) – non-zero value, maximal distance between donor and acceptor. default is 3.5 distDA also works
- angle (int, float) – non-zero value, maximal (180 - D-H-A angle) (donor, hydrogen, acceptor). default is 40. angleDHA also works
- seq_cutoff (int) – non-zero value, interactions will be found between atoms with index differences that are higher than seq_cutoff. default is 20 atoms.
- selection (str) – selection string
- selection2 (str) – selection string
- start_frame (int) – index of first frame to read
- stop_frame (int) – index of last frame to read
Selection: If we want to select interactions for the particular residue or group of residues:
selection=’chain A and resid 1 to 50’- If we want to study chain-chain interactions:
- selection=’chain A’, selection2=’chain B’
- atoms (
-
calcHydrophobicOverlapingAreas
(atoms, **kwargs)[source]¶ Provide information about hydrophobic contacts between pairs of residues based on the regsurf program. To use this function compiled hpb.so is needed.
Parameters: - atoms (
Atomic
) – an Atomic object from which residues are selected - selection (str) – selection string of hydrophobic residues
- hpb_cutoff (float, int) – cutoff for hydrophobic overlapping area values default is 0.0
- cumulative_values ('pairs' or 'single_residues') – sum of results for pairs of residues or for single residues without distinguishing pairs, default is None
- atoms (
-
calcSminaBindingAffinity
(atoms, trajectory=None, **kwargs)[source]¶ Computing binding affinity of ligand toward protein structure using SMINA package [DRK13].
Parameters: - atoms (
Atomic
,LigandInteractionsTrajectory
) – an Atomic object from which residues are selected - protein_selection (str) – selection string for the protein and other compoment of the system that should be included, e.g. “protein and chain A”, Default is “protein”
- ligand_selection (str) – selection string for ligand, e.g. “resname ADP”, Default is “all not protein”
- ligand_selection – scoring function (vina or vinardo) Default is “vina”
- atom_terms (bool) – write per-atom interaction term values. It will provide more information as dictionary for each frame/model, and details for atom terms will be saved in terms_*frame_number*.txt, Default is False
SMINA installation is required to compute ligand binding affinity: >> conda install -c conda-forge smina (for Anaconda)
For more information on SMINA see https://sourceforge.net/projects/smina/. If you benefited from SMINA, please consider citing [DRK13].
[DRK13] (1, 2, 3) Koes D. R., Baumgartner M. P., Camacho C. J., Lessons Learned in Empirical Scoring with smina from the CSAR 2011 Benchmarking Exercise, J. Chem. Inf. Model. 2013 53: 1893–1904.
- atoms (
-
calcSminaPerAtomInteractions
(atoms, list_terms)[source]¶ Computing the summary of per-atom interaction term values using SMINA package [DRK13]. It will return dictionary with per-atom interaction term values.
Parameters: - atoms (
Atomic
,LigandInteractionsTrajectory
) – an Atomic object from which residues are selected - list_terms (list) – List of terms.txt files obtained from meth:.calcSminaBindingAffinity using atom_terms = True
Important: First text file in the list should be reference structure which correspond to the provided coordinates as atoms.
- atoms (
-
calcSminaTermValues
(data)[source]¶ Computing weights multiplied by term values, before weighting for each Term. As a results will are obtaining a dictionary.
Parameters: data (list) – List of results provided by Smina using meth:.calcSminaBindingAffinity with atom_terms = True
-
showSminaTermValues
(data)[source]¶ Display a histogram of weights multiplied by term values, before weighting for each Term. As a results will are obtaining a dictionary.
Parameters: data (list) – List of results provided by Smina using meth:.calcSminaBindingAffinity with atom_terms = True
-
showPairEnergy
(data, **kwargs)[source]¶ Return energies when a list of interactions is given. Energies will be added to each pair of residues at the last position in the list. Energy is based on the residue types and not on the distances. The unit of energy is kcal/mol. The energies defined as ‘IB_nosolv’ (non-solvent-mediated) and ‘IB_solv’ (solvent-mediated) are taken from [OK98] and ‘CS’ from InSty paper (under preparation). Protonation of residues is not distinguished. The protonation of residues is not distinguished. ‘IB_solv’ and ‘IB_nosolv’ have RT units and ‘CS’ has units of kcal/mol.
Known residues such as HSD, HSE, HIE, and HID (used in MD simulations) are treated as HIS.
Parameters: - data (list) – list with interactions from calcHydrogenBonds() or other types
- energy_list_type ('IB_nosolv', 'IB_solv', 'CS') – name of the list with energies default is ‘IB_solv’
[OK98] (1, 2, 3, 4, 5) Keskin O., Bahar I., Badretdinov A.Y., Ptitsyn O.B., Jernigan R.L., Empirical solvet-mediated potentials hold for both intra-molecular and inter-molecular inter-residues interactions, Protein Science 1998 7: 2578–2586.
-
checkNonstandardResidues
(atoms)[source]¶ Check whether the atomic structure contains non-standard residues and inform to replace the name to the standard one so that non-standard residues are treated in a correct way while computing interactions.
Parameters: atoms ( Atomic
) – an Atomic object from which residues are selected
-
saveInteractionsAsDummyAtoms
(atoms, interactions, filename, **kwargs)[source]¶ Creates a PDB file which will contain protein structure and dummy atoms that will be placed between pairs of interacting residues.
Parameters:
-
createFoldseekAlignment
(prot_seq, prot_foldseek, **kwargs)[source]¶ Aligns sequences from prot_seq with homologous sequences identified in prot_foldseek, generating a multiple sequence alignment.
Parameters:
-
runFoldseek
(pdb_file, chain, **kwargs)[source]¶ This script processes a PDB file to extract a specified chain’s sequence, searches for homologous structures using foldseek, and prepares alignment outputs for further analysis.
Before using the function, install Foldseek: >>> conda install bioconda::foldseek More information can be found: https://github.com/steineggerlab/foldseek?tab=readme-ov-file#databasesand
This function will not work under Windows. Example usage: runFoldseek(‘5kqm.pdb’, ‘A’, database_folder=’~/Downloads/foldseek/pdb’) where previous a folder called ‘foldseek’ were created and PDB database was uploaded using: >>> foldseek databases PDB pdb tmp (Linux console)
Parameters: - pdb_file (str) – A PDB file path
- chain (str) – Chain identifier
- coverage_threshold (float) – Coverage threshold Default is 0.3
- tm_threshold (float) – TM-score threshold Default is 0.5
- database_folder (str) – Folder with the database Default is ‘~/Downloads/foldseek/pdb’ To download the database use: ‘foldseek databases PDB pdb tmp’ in the console
- fixer ('pdbfixer' or 'openbabel') – The method for fixing lack of hydrogen bonds Default is ‘pdbfixer’
- subset (str) – subsets of atoms: ‘ca’, ‘bb’, ‘heavy’, ‘noh’, ‘all’ (see matchChains()) Default is ‘ca’
- seqid (float) – Minimum value of the sequence identity (see matchChains()) Default is 5
- overlap (float) – percent overlap (see matchChains()) Default is 50
- folder_name (str) – Folder where the results will be collected Default is ‘struc_homologs’
-
runDali
(pdb, chain, **kwargs)[source]¶ This function calls searchDali() and downloads all the PDB files, separate by chains, add hydrogens and missing side chains, and finally align them and put into the newly created folder.
Parameters: - pdb (str) – A PDB code
- chain (str) – chain identifier
- cutoff_len (float) – Length of aligned residues < cutoff_len (must be an integer or a float between 0 and 1) See searchDali for more details Default is 0.5
- cutoff_rmsd (float) – RMSD cutoff (see searchDali) Default is 1.0
- subset_Dali (str) – fullPDB, PDB25, PDB50, PDB90 Default is ‘fullPDB’
- fixer ('pdbfixer' or 'openbabel') – The method for fixing lack of hydrogen bonds Default is ‘pdbfixer’
- subset (str) – subsets of atoms: ‘ca’, ‘bb’, ‘heavy’, ‘noh’, ‘all’ (see matchChains()) Default is ‘ca’
- seqid (float) – Minimum value of the sequence identity (see matchChains()) Default is 5
- overlap (float) – percent overlap (see matchChains()) Default is 50
- folder_name (str) – Folder where the results will be collected Default is ‘struc_homologs’
-
runBLAST
(pdb, chain, **kwargs)[source]¶ This function calls blastPDB to find homologs and downloads all of them in PDB format to the local directory, separate chains that were identified by BLAST, add hydrogens and missing side chains, and finally align them and put into the newly created folder.
Parameters: - pdb (str) – A PDB code
- chain (str) – chain identifier
- fixer ('pdbfixer' or 'openbabel') – The method for fixing lack of hydrogen bonds Default is ‘pdbfixer’
- subset (str) – subsets of atoms: ‘ca’, ‘bb’, ‘heavy’, ‘noh’, ‘all’ (see matchChains()) Default is ‘bb’
- seqid (float) – Minimum value of the sequence identity (see matchChains()) Default is 90
- overlap (float) – percent overlap (see matchChains()) Default is 50
- folder_name (str) – Folder where the results will be collected Default is ‘struc_homologs’
-
extractMultiModelPDB
(multimodelPDB, **kwargs)[source]¶ Extracts individual PDB models from multimodel PDB and places them into the pointed directory. If used for calculating calcSignatureInteractions align the models.
Parameters:
-
calcSignatureInteractions
(PDB_folder, **kwargs)[source]¶ Analyzes protein structures to identify various interactions using InSty. Processes data from the MSA file and folder with selected models.
Example usage: >>> calcSignatureInteractions(‘./struc_homologs’) >>> calcSignatureInteractions(‘./struc_homologs’, mapping_file=’shortlisted_resind.msa’)
Parameters: