alg
– Algorithms for Structures¶
Topics¶
Submodules¶
chain_mapping
– Chain Mappingcontact_score
– Contact-Based Scoresdockq
– DockQ Implementationhelix_kinks
– Algorithms to Calculate Helix Kinksligand_scoring
– Ligand scoring functionsqsscore
– New QS Score Implementationscoring
– Specialized Scoring Functionsscoring_base
– Helper functions for scoringstereochemistry
– Stereochemistry Checksstructure_analysis
– Functions to Analyze Structurestrajectory_analysis
– DRMSD, Pairwise Distances and More
Algorithms on Structures¶
- Accessibility(ent, probe_radius=1.4, include_hydrogens=False, include_hetatm=False, include_water=False, oligo_mode=False, selection='', asa_abs='asaAbs', asa_rel='asaRel', asa_atom='asaAtom', algorithm=NACCESS)¶
Calculates the accesssible surface area for ever atom in ent. The algorithm mimics the behaviour of the bindings available for the NACCESS and DSSP tools and has been tested to reproduce the numbers accordingly.
- Parameters:
ent (
EntityView
/EntityHandle
) – Entity on which to calculate the surfaceprobe_radius (
float
) – Radius of probe to determine accessible surface areainclude_hydrogens (
bool
) – Whether to include hydrogens in the solvent accessibility calculations. By default, every atom with ele=H,D is simply neglected.include_hetatms (
bool
) – Whether to include atoms flagged as hetatms , i.e. ligands, in the solvent accessibility calculations. They are neglected by default.include_water (
bool
) – Whether to include water in the solvent accessibility calculations. By default, every residue with name “HOH” is neglected.oligo_mode (
bool
) – A typical used case of accessibility calculations is to determine the solvent accessibility of a full complex and then the accessibility of each chain individually. Lots of calculations can be cached because only the values of the atoms close to an interface change. This is exactly what happens when you activate the oligo mode. It returns exactly the same value but adds, additionally to the values estimated in full complex, the values from each individual chain as float properties on every residue and atom. Example for atom accessible surface if the according property name is set to “asaAtom”: Accessibility in the full complex is stored as “asaAtom”, the accessibility when only considering that particular chain is stored as “asaAtom_single_chain”. The other properties follow the same naming scheme.selection (
str
) – Selection statement, that gets applied on ent before doing anything. Everything that is not selected is neglected. The default value of “” results in no selection at all.asa_abs (
str
) – Float property name to assign the summed solvent accessible surface from each atom to a residue.asa_rel (
str
) –Float property name to assign the relative solvent accessibility to a residue. This is the absolute accessibility divided by the maximum solvent accessibility of that particular residue. This maximum solvent accessibility is dependent on the chosen
AccessibilityAlgorithm
. Only residues of the 20 standarad amino acids can be handled.In case of the NACCESS algorithm you can expect a value in range [0.0, 100.0] and a value of -99.9 for non standard residues.
In case of the DSSP algorithm you can expect a value in range [0.0, 1.0], no float property is assigned in case of a non standard residue.
asa_atom (
str
) – Float property name to assign the solvent accessible area to each atom.algorithm (
AccessibilityAlgorithm
) – Specifies the used algorithm for solvent accessibility calculations
- Returns:
The summed solvent accessibilty of each atom in ent.
- class AccessibilityAlgorithm¶
The accessibility algorithm enum specifies the algorithm used by the respective tools. Available are:
NACCESS, DSSP
- AssignSecStruct(ent)¶
Assigns secondary structures to all residues based on hydrogen bond patterns as described by DSSP.
- Parameters:
ent (
EntityView
/EntityHandle
) – Entity on which to assign secondary structures
- class FindMemParam¶
Result object for the membrane detection algorithm described below
- axis¶
initial search axis from which optimal membrane slab could be found
- tilt_axis¶
Axis around which we tilt the membrane starting from the initial axis
- tilt¶
Angle to tilt around tilt axis
- angle¶
After the tilt operation we perform a rotation around the initial axis with this angle to get the final membrane axis
- membrane_axis¶
The result of applying the tilt and rotation procedure described above. The membrane_axis is orthogonal to the membrane plane and has unit length.
- pos¶
Real number that describes the membrane center point. To get the actual position you can do: pos * membrane_axis
- width¶
Total width of the membrane in A
- energy¶
Pseudo energy of the implicit solvation model
- membrane_asa¶
Membrane accessible surface area
- membrane_representation¶
Dummy atoms that represent the membrane. This entity is only valid if the according flag has been set to True when calling FindMembrane.
- FindMembrane(ent, assign_membrane_representation=True, fast=False)¶
Estimates the optimal membrane position of a protein by using an implicit solvation model. The original algorithm and the used energy function are described in: Lomize AL, Pogozheva ID, Lomize MA, Mosberg HI (2006) Positioning of proteins in membranes: A computational approach.
There are some modifications in this implementation and the procedure is as follows:
Initial axis are constructed that build the starting point for initial parameter grid searches.
For every axis, the protein is rotated so that the axis builds the z-axis
In order to exclude internal hydrophilic pores, only the outermost atoms with respect the the z-axis enter an initial grid search
The width and position of the membrane is optimized for different combinations of tilt and rotation angles (further described in
FindMemParam
). The top 20 parametrizations (only top parametrization if fast is True) are stored for further processing.
The 20 best membrane parametrizations from the initial grid search (only the best if fast is set to True) enter a final minimization step using a Levenberg-Marquardt minimizer.
- Parameters:
ent (
ost.mol.EntityHandle
/ost.mol.EntityView
) –Entity of a transmembrane protein, you’ll get weird results if this is not the case. The energy term of the result is typically a good indicator whether ent is an actual transmembrane protein. The following float properties will be set on the atoms:
’asaAtom’ on all atoms that are selected with ent.Select(‘peptide=true and ele!=H’) as a result of envoking
Accessibility()
.’membrane_e’ the contribution of the potentially membrane facing atoms to the energy function.
assign_membrane_representation (
bool
) – Whether to construct a membrane representation using dummy atomsfast – If set to false, the 20 best results of the initial grid search undergo a Levenberg-Marquardt minimization and the parametrization with optimal minimized energy is returned. If set to yes, only the best result of the initial grid search is selected and returned after Levenberg-Marquardt minimization.
- Returns:
The results object
- Return type:
Trajectory Analysis¶
This is a set of functions used for basic trajectory analysis such as extracting
positions, distances, angles and RMSDs. The organization is such that most
functions have their counterpart at the individual frame level
so that they can also be called on one frame instead of
the whole trajectory.
All these functions have a “stride” argument that defaults to stride=1, which is used to skip frames in the analysis.
Additional functions are available in the ost.mol.alg.trajectory_analysis
submodule.
- AnalyzeAngle(traj, atom1, atom2, atom3, stride=1)¶
This function extracts the angle between three atoms from a trajectory and returns it as a vector. The second atom is taken as being the central atom, so that the angle is between the vectors (atom1.pos-atom2.pos) and (atom3.pos-atom2.pos).
- Parameters:
traj (
CoordGroupHandle
) – The trajectory to be analyzed.atom1 – The first
AtomHandle
.atom2 – The second
AtomHandle
.atom3 – The third
AtomHandle
.stride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- AnalyzeAromaticRingInteraction(traj, view_ring1, view_ring2, stride=1)¶
This function is a crude analysis of aromatic ring interactions. For each frame in a trajectory, it calculates the minimal distance between the atoms in one view and the center of mass of the other and vice versa, and returns the minimum between these two minimal distances. Concretely, if the two views are the heavy atoms of two rings, then it returns the minimal center of mass - heavy atom distance betweent he two rings
- Parameters:
traj (
CoordGroupHandle
) – The trajectory to be analyzed.view_ring1 (
EntityView
.) – First group of atomsview_ring2 (
EntityView
.) – Second group of atomsstride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- AnalyzeAtomPos(traj, atom1, stride=1)¶
This function extracts the position of an atom from a trajectory. It returns a vector containing the position of the atom for each analyzed frame.
- Parameters:
traj (
CoordGroupHandle
) – The trajectory to be analyzed.atom1 – The
AtomHandle
.stride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- AnalyzeCenterOfMassPos(traj, sele, stride=1)¶
This function extracts the position of the center-of-mass of a selection (
EntityView
) from a trajectory and returns it as a vector.- Parameters:
traj (
CoordGroupHandle
) – The trajectory to be analyzed.sele (
EntityView
.) – The selection from which the center of mass is computedstride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- AnalyzeDihedralAngle(traj, atom1, atom2, atom3, atom4, stride=1)¶
This function extracts the dihedral angle between four atoms from a trajectory and returns it as a vector. The angle is between the planes containing the first three and the last three atoms.
- Parameters:
traj (
CoordGroupHandle
) – The trajectory to be analyzed.atom1 – The first
AtomHandle
.atom2 – The second
AtomHandle
.atom3 – The third
AtomHandle
.atom4 – The fourth
AtomHandle
.stride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- AnalyzeDistanceBetwAtoms(traj, atom1, atom2, stride=1)¶
This function extracts the distance between two atoms from a trajectory and returns it as a vector.
- Parameters:
traj (
CoordGroupHandle
) – The trajectory to be analyzed.atom1 – The first
AtomHandle
.atom2 – The second
AtomHandle
.stride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- AnalyzeDistanceBetwCenterOfMass(traj, sele1, sele2, stride=1)¶
This function extracts the distance between the center-of-mass of two selections (
EntityView
) from a trajectory and returns it as a vector.- Parameters:
traj (
CoordGroupHandle
) – The trajectory to be analyzed.sele1 (
EntityView
.) – The selection from which the first center of mass is computedsele2 (
EntityView
.) – The selection from which the second center of mass is computedstride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- AnalyzeMinDistance(traj, view1, view2, stride=1)¶
This function extracts the minimal distance between two sets of atoms (view1 and view2) for each frame in a trajectory and returns it as a vector.
- Parameters:
traj (
CoordGroupHandle
) – The trajectory to be analyzed.view1 (
EntityView
.) – The first group of atomsview2 (
EntityView
.) – The second group of atomsstride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- AnalyzeMinDistanceBetwCenterOfMassAndView(traj, view_cm, view_atoms, stride=1)¶
This function extracts the minimal distance between a set of atoms (view_atoms) and the center of mass of a second set of atoms (view_cm) for each frame in a trajectory and returns it as a vector.
- Parameters:
traj (
CoordGroupHandle
) – The trajectory to be analyzed.view_cm (
EntityView
.) – The group of atoms from which the center of mass is takenview_atoms (
EntityView
.) – The second group of atomsstride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- AnalyzeRMSD(traj, reference_view, sele_view, stride=1)¶
This function extracts the rmsd between two
EntityView
and returns it as a vector. The views don’t have to be from the same entity. The reference positions are taken directly from the reference_view, evaluated only once. The positions from the sele_view are evaluated for each frame. If you want to compare to frame i of the trajectory t, first use t.CopyFrame(i) for example:eh = io.LoadPDB(...) t = io.LoadCHARMMTraj(eh, ...) sele = eh.Select(...) t.CopyFrame(0) mol.alg.AnalyzeRMSD(t, sele, sele)
- Parameters:
traj (
CoordGroupHandle
) – The trajectory to be analyzed.reference_view (
EntityView
.) – The selection used as reference structuresele_view (
EntityView
.) – The selection compared to the reference_viewstride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- SuperposeFrames(frames, sel, from=0, to=-1, ref=-1)¶
This function superposes the frames of the given coord group and returns them as a new coord group.
- Parameters:
frames (
CoordGroupHandle
) – The source coord group.sel (
ost.mol.EntityView
) – An entity view containing the selection of atoms to be used for superposition. If set to an invalid view, all atoms in the coord group are used.from – index of the first frame
to – index of the last frame plus one. If set to -1, the value is set to the number of frames in the coord group
ref – The index of the reference frame to use for superposition. If set to -1, the each frame is superposed to the previous frame.
- Returns:
A newly created coord group containing the superposed frames.
- SuperposeFrames(frames, sel, ref_view, from=0, to=-1)
Same as SuperposeFrames above, but the superposition is done on a reference view and not on another frame of the trajectory.
- Parameters:
frames (
CoordGroupHandle
) – The source coord group.sel (
ost.mol.EntityView
) – An entity view containing the selection of atoms of the frames to be used for superposition.ref_view (
ost.mol.EntityView
) – The reference view on which the frames will be superposed. The number of atoms in this reference view should be equal to the number of atoms in sel.from – index of the first frame
to – index of the last frame plus one. If set to -1, the value is set to the number of frames in the coord group
- Returns:
A newly created coord group containing the superposed frames.