mol.alg
– Algorithms for Structures¶
Local Distance Test scores (lDDT, DRMSD)¶
-
LocalDistDiffTest
(model, distance_list, tolerance_list, sequence_separation=0, local_lddt_property_string="")¶ This function counts the number of conserved local contacts between a model and a reference structure which is needed to compute the Local Distance Difference Test score.
The Local Distance Difference Test score is a number between zero and one, which measures the agreement of local contacts between a model and a reference structure. One means complete agreement, and zero means no agreement at all. The calculation of this score does not require any superposition between the model and the reference structures.
All distances between atoms in the reference structure that are shorter than a certain predefined length (inclusion radius) are compared with the corresponding distances in the model structure. If the difference between a reference distance and the corresponding model distance is smaller than a threshold value (tolerance), that distance is considered conserved. The final lDDT score is the fraction of conserved distances. Missing atoms in the model structure lead to non-conserved distances (and thus lower the final lDDT score).
This function takes as an input a list of distances to be checked for conservation. Any number of threshold values can be specified when the function is called. All thresholds are then applied in sequence and the return counts are averaged over all threshold values. A sequence separation parameter can be passed to the function. If this happens, only distances between residues whose separation in sequence is higher than the provided parameter are considered when the score is computed.
If a string is passed as the last parameter, residue-based counts and the value of the residue-based Local Distance Difference Test score are saved in each ResidueHandle as int and float properties. Specifically, the local residue-based lddt score is stored in a float property named as the provided string, while the residue-based number of conserved and total distances are saved in two int properties named <string>_conserved and <string>_total.
Parameters: - model (
EntityView
) – the model structure - distance_list (
GlobalRDMap
) – the list of distances to check for conservation - tolerance_list – a list of thresholds used to determine distance conservation
- sequence_separation – sequence separation parameter used when computing the score
- local_lddt_property_string – the base name for the ResidueHandle properties that store the local scores
Returns: a tuple containing the counts of the conserved distances in the model and of all the checked distances
- model (
-
LocalDistDiffTest
(model, target, cutoff, max_dist, local_lddt_property_string="") Wrapper around
LocalDistDiffTest()
above using: distance_list =CreateDistanceList()
with target and max_dist as parameters and tolerance_list = [cutoff].Parameters: - model (
EntityView
) – the model structure - target (
EntityView
) – the target structure from which distances are derived - cutoff (
float
) – single distance threshold to determine distance conservation - max_dist (
float
) – the inclusion radius in Angstroms (to determine which distances are checked for conservation) - local_lddt_property_string – the base name for the ResidueHandle properties that store the local scores
Returns: the Local Distance Difference Test score (conserved distances divided by all the checked distances)
Return type: float
- model (
-
LocalDistDiffTest
(alignment, tolerance, radius, ref_index=0, mdl_index=1) Calculates the Local Distance Difference Test score (see previous function) starting from an alignment between a reference structure and a model. The AlignmentHandle parameter used to provide the alignment to the function needs to have the two structures attached to it. By default the first structure in the alignment is considered to be the reference structure, and the second structure is taken as the model. This can however be changed by passing the indexes of the two structures in the AlignmentHandle as parameters to the function.
Note
This function uses the old implementation of the Local Distance Difference Test algorithm and will give slightly different results from the new one.
Parameters: - alignment (
AlignmentHandle
) – an alignment containing the sequences of the reference and of the model structures, with the structures themselves attached - tolerance – a list of thresholds used to determine distance conservation
- radius – the inclusion radius in Angstroms (to determine which distances are checked for conservation)
- ref_index – index of the reference structure in the alignment
- mdl_index – index of the model in the alignment
Returns: the Local Distance Difference Test score
- alignment (
-
LDDTHA
(model, distance_list, sequence_separation=0)¶ This function calculates the Local Distance Difference Test, using the same threshold values as the GDT-HA test (the default set of thresholds used for the lDTT score) (See previous functions). The thresholds are 0.5, 1, 2, and 4 Angstroms.
The function only compares the input distance list to the first chain of the model structure.
The local residue-based lDDT score values are stored in the ResidueHandles of the model passed to the function in a float property called “locallddt”.
A sequence separation parameter can be passed to the function. If this happens, only distances between residues whose separation is higher than the provided parameter are considered when computing the score.
Parameters: - model (
EntityView
) – the model structure - distance_list (
GlobalRDMap
) – the list of distances to check for conservation - sequence_separation – sequence separation parameter
Returns: the Local Distance Difference Test score
- model (
-
DistanceRMSDTest
(model, distance_list, cap_difference, sequence_separation=0, local_drmsd_property_string="")¶ This function performs a Distance RMSD Test on a provided model, and calculates the two values that are necessary to determine the Distance RMSD Score, namely the sum of squared distance deviations and the number of distances on which the sum was computed.
The Distance RMSD Test (or DRMSD Test) computes the deviation in the length of local contacts between a model and a reference structure and expresses it in the form of a score value. The score has an an RMSD-like form, with the deviations in the RMSD formula computed as contact distance differences. The score is open-ended, with a value of zero meaning complete agreement of local contact distances, and a positive value revealing a disagreement of magnitude proportional to the score value itself. This score does not require any superposition between the model and the reference.
This function processes a list of distances provided by the user, together with their length in the reference structure. For each distance that is found in the model, its difference with the reference length is computed and used as deviation term in the RMSD-like formula.When a distance is not present in the model because one or both the atoms are missing, a default deviation value provided by the user is used.
The function only processes distances between atoms that do not belong to the same residue, and considers only standard residues in the first chain of the model. For residues with symmetric sidechains (GLU, ASP, ARG, VAL, PHE, TYR), the naming of the atoms is ambiguous. For these residues, the function computes the Distance RMSD Test score that each naming convention would generate when considering all non-ambiguous surrounding atoms. The solution that gives the lower score is then picked to compute the final Distance RMSD Score for the whole model.
A sequence separation parameter can be passed to the function. If this happens, only distances between residues whose separation is higher than the provided parameter are considered when computing the score.
If a string is passed as last parameter to the function, the function computes the Distance RMSD Score for each residue and saves it as a float property in the ResidueHandle, with the passed string as property name. Additionally, the actual sum of squared deviations and the number of distances on which it was computed are stored as properties in the ResidueHandle. The property names are respectively <passed string>_sum (a float property) and <passed string>_count (an integer property).
Parameters: - model (
EntityView
) – the model structure - distance_list (
GlobalRDMap
) – the list of distances to check (here we only use the first of the two distance values stored, the second is ignored) - cap_difference – a default deviation value to be used when a distance is not found in the model
- sequence_separation – sequence separation parameter
- local_ldt_property_string – the base name for the ResidueHandle properties that store the local scores
Returns: a tuple containing the sum of squared distance deviations, and the number of distances on which it was computed.
- model (
-
DRMSD
(model, distance_list, cap_difference, sequence_separation=0)¶ This function calculates the Distance RMSD Test score (see
DistanceRMSDTest()
).The function only considers distances between atoms not belonging to the same residue, and only compares the input distance list to the first chain of the model structure. It requires, in addition to the model and the list themselves, a default deviation value to be used in the DRMSD Test when a distance is not found in the model.
The local Local Distance Difference Test score values are stored in the ResidueHandles of the model passed to the function in a float property called “localdrmsd”.
A sequence separation parameter can be passed to the function. If this happens, only distances between residues whose separation is higher than the provided parameter are considered when computing the score.
Parameters: - model (
EntityView
) – the model structure - distance_list (
GlobalRDMap
) – the list of distances as inDistanceRMSDTest()
- cap_difference – a default deviation value to be used when a distance is not found in the model
- sequence_separation – sequence separation parameter
Returns: the Distance RMSD Test score
- model (
-
CreateDistanceList
(reference, radius)¶ -
CreateDistanceListFromMultipleReferences
(reference_list, tolerance_list, sequence_separation, radius)¶ Both these functions create lists of distances to be checked during a Local Distance Difference Test (see description of the functions above).
Note
These functions process only standard residues present in the first chain of the reference structures.
The only difference between the two functions is that one takes a single reference structure and the other a list of reference structures. The structures in the list have to be properly prepared before being passed to the function. Corresponding residues in the structures must have the same residue number, the same chain name, etc. Gaps are allowed and automatically dealt with: if information about a distance is present in at least one of the structures, it will be considered.
If a distance between two atoms is shorter than the inclusion radius in all structures in which the two atoms are present, it is included in the list. However, if the distance is longer than the inclusion radius in at least one of the structures, it is not considered to be a local interaction and is excluded from the list.
The multiple-reference function takes care of residues with ambiguous symmetric sidechains. To decide which naming convention to use, the function computes a Local Distance Difference Test score foreach reference against the first reference structure in the list, using only non ambiguously-named atoms. It picks then the naming convention that gives the highest score, guaranteeing that all references are processed with the correct atom names.
The cutoff list that will later be used to compute the Local Distance Difference Test score and the sequence separation parameter must be passed to the multi-reference function. These parameters do not influence the output distance list, which always includes all distances within the provided radius (to make it consistent with the single-reference corresponding function). However, the parameters are used when dealing with the naming convention of residues with ambiguous nomenclature.
Parameters: - reference (
EntityView
) – a reference structure from which distances are derived - reference_list (list of
EntityView
) – a list of reference structures from which distances are derived - tolerance_list – a list of thresholds used to determine distance conservation when computing the lDDT score
- sequence_separation – sequence separation parameter used when computing the lDDT score
- radius – inclusion radius (in Angstroms) used to determine the distances included in the list
Returns: - reference (
-
class
UniqueAtomIdentifier
¶ Object containing enough information to uniquely identify an atom in a structure
-
UniqueAtomIdentifier
(chain, residue_number, residue_name, atom_name)¶ Creates an UniqueAtomIdentifier object starting from relevant atom information
Parameters: - chain – a string containing the name of the chain to which the atom belongs
- residue_number (
ResNum
) – the number of the residue to which the atom belongs - residue_name – a string containing the name of the residue to which the atom belongs
- atom_name – a string containing the name of the atom
-
GetChainName
()¶ Returns the name of the chain to which the atom belongs, as a String
-
GetResidueName
()¶ Returns the name of the residue to which the atom belongs, as a String
-
GetAtomName
()¶ Returns the name of the atom, as a String
-
GetQualifiedAtomName
()¶ Returns the qualified name of the atom (the chain name, followed by a unique residue identifier and the atom name. For example: “A.GLY2.CA”)
-
-
class
ResidueRDMap
¶ Dictionary-like object containing the list of interatomic distances that originate from a single residue to be checked during a run of the Local Distance Difference Test algorithm (key = pair of
UniqueAtomIdentifier
, value = pair of floats)
-
class
GlobalRDMap
¶ Dictionary-like object containing all the
ResidueRDMap
objects related to residues of a single structure (key =ResNum
, value =ResidueRDMap
)
-
PrintResidueRDMap
(residue_distance_list)¶ Prints to standard output all the distances contained in a
ResidueRDMap
object
-
PrintGlobalRDMap
(global_distance_list)¶ Prints to standard output all the distances contained in each of the
ResidueRDMap
objects that make up aGlobalRDMap
object.
qsscoring
– Quaternary Structure (QS) scores¶
Scoring of quaternary structures as in Martino’s 2017 paper.
Note
Requirements for use:
- A default
compound library
must be defined and accessible viaGetDefaultLib()
. This is set by default when executing scripts withost
. Otherwise, you must set this withSetDefaultLib()
. - ClustalW must be installed (unless you provide chain mappings)
- Python modules numpy and scipy must be installed and available
(e.g. use
pip install scipy numpy
)
Authors: Gerardo Tauriello, Martino Bertoni
-
exception
QSscoreError
¶ Exception to be raised for “acceptable” exceptions in QS scoring.
Those are cases we might want to capture for default behavior.
-
class
QSscorer
(ent_1, ent_2)¶ Object to compute QS scores.
Simple usage without any precomputed contacts, symmetries and mappings:
import ost from ost.mol.alg import qsscoring # load two biounits to compare ent_full = ost.io.LoadPDB('3ia3', remote=True) ent_1 = ent_full.Select('cname=A,D') ent_2 = ent_full.Select('cname=B,C') # get score ost.PushVerbosityLevel(3) try: qs_scorer = qsscoring.QSscorer(ent_1, ent_2) ost.LogScript('QSscore:', str(qs_scorer.global_score)) ost.LogScript('Chain mapping used:', str(qs_scorer.chain_mapping)) except qsscoring.QSscoreError as ex: # default handling: report failure and set score to 0 ost.LogError('QSscore failed:', str(ex)) qs_score = 0
For maximal performance when computing QS scores of the same entity with many others, it is advisable to construct and reuse
QSscoreEntity
objects.Any known / precomputed information can be filled into the appropriate attribute here (no checks done!). Otherwise most quantities are computed on first access and cached (lazy evaluation). Setters are provided to set values with extra checks (e.g.
SetSymmetries()
).All necessary seq. alignments are done by global BLOSUM62-based alignment. A multiple sequence alignment is performed with ClustalW unless
chain_mapping
is provided manually. You will need to have an executableclustalw
orclustalw2
in yourPATH
or you must setclustalw_bin
accordingly. Otherwise an exception (ost.settings.FileNotFound
) is thrown.Formulas for QS scores:
- QS_best = weighted_scores / (weight_sum + weight_extra_mapped) - QS_global = weighted_scores / (weight_sum + weight_extra_all) -> weighted_scores = sum(w(min(d1,d2)) * (1 - abs(d1-d2)/12)) for shared -> weight_sum = sum(w(min(d1,d2))) for shared -> weight_extra_mapped = sum(w(d)) for all mapped but non-shared -> weight_extra_all = sum(w(d)) for all non-shared -> w(d) = 1 if d <= 5, exp(-2 * ((d-5.0)/4.28)^2) else
Parameters: - ent_1 (
QSscoreEntity
,EntityHandle
orEntityView
) – First structure to be scored. - ent_2 (
QSscoreEntity
,EntityHandle
orEntityView
) – Second structure to be scored.
Raises: QSscoreError
if input structures are invalid or are monomers or have issues that make it impossible for a QS score to be computed.-
qs_ent_1
¶ QSscoreEntity
object for ent_1 given at construction. If entity names (original_name
) are not unique, we set it to ‘pdb_1’ usingSetName()
.
-
qs_ent_2
¶ QSscoreEntity
object for ent_2 given at construction. If entity names (original_name
) are not unique, we set it to ‘pdb_2’ usingSetName()
.
-
calpha_only
¶ True if any of the two structures is CA-only (after cleanup).
Type: bool
-
max_ca_per_chain_for_cm
¶ Maximal number of CA atoms to use in each chain to determine chain mappings. Setting this to -1 disables the limit. Limiting it speeds up determination of symmetries and chain mappings. By default it is set to 100.
Type: int
-
SetSymmetries
(symm_1, symm_2)¶ Set user-provided symmetry groups.
These groups are restricted to chain names appearing in
ent_to_cm_1
andent_to_cm_2
respectively. They are only valid if they cover all chains and both symm_1 and symm_2 have same lengths of symmetry group tuples. Otherwise trivial symmetry group used (seesymm_1
).Parameters:
-
alignments
¶ List of successful sequence alignments using
chain_mapping
.There will be one alignment for each mapped chain and they are ordered by their chain names in
qs_ent_1
.The sequences of the alignments have views attached into
QSscoreEntity.ent
ofqs_ent_1
andqs_ent_2
.Getter: Computed on first use (cached) Type: list
ofAlignmentHandle
-
best_score
¶ QS-score without penalties.
Like
global_score
, but neglecting additional residues or chains in one of the biounits (i.e. the score is calculated considering only mapped chains and residues).Getter: Computed on first use (cached) Type: float
-
chain_mapping
¶ Mapping from
ent_to_cm_1
toent_to_cm_2
.Properties:
- Mapping is between chains of same chem. group (see
chem_mapping
) - Each chain can appear only once in mapping
- All chains of complex with less chains are mapped
- Symmetry (
symm_1
,symm_2
) is taken into account
Details on algorithms used to find mapping:
- We try all pairs of chem. mapped chains within symmetry group and get superpose-transformation for them
- First option: check for “sufficient overlap” of other chain-pairs
- For each chain-pair defined above: apply superposition to full oligomer and map chains based on structural overlap
- Structural overlap = X% of residues in second oligomer covered within Y Angstrom of a (chem. mapped) chain in first oligomer. We successively try (X,Y) = (80,4), (40,6) and (20,8) to be less and less strict in mapping (warning shown for most permissive one).
- If multiple possible mappings are found, we choose the one which leads to the lowest multi-chain-RMSD given the superposition
- Fallback option: try all mappings to find minimal multi-chain-RMSD
(warning shown)
- For each chain-pair defined above: apply superposition, try all (!) possible chain mappings (within symmetry group) and keep mapping with lowest multi-chain-RMSD
- Repeat procedure above to resolve symmetry. Within the symmetry group we can use the chain mapping computed before and we just need to find which symmetry group in first oligomer maps to which in the second one. We again try all possible combinations...
- Limitations:
- Trying all possible mappings is a combinatorial nightmare (factorial). We throw an exception if too many combinations (e.g. octomer vs octomer with no usable symmetry)
- The mapping is forced: the “best” mapping will be chosen independently of how badly they fit in terms of multi-chain-RMSD
- As a result, such a forced mapping can lead to a large range of
resulting QS scores. An extreme example was observed between 1on3.1
and 3u9r.1, where
global_score
can range from 0.12 to 0.43 for mappings with very similar multi-chain-RMSD.
Getter: Computed on first use (cached) Type: dict
with key / value =str
(chain names, key forent_to_cm_1
, value forent_to_cm_2
)Raises: QSscoreError
if there are too many combinations to check to find a chain mapping.- Mapping is between chains of same chem. group (see
-
chem_mapping
¶ Inter-complex mapping of chemical groups.
Each group (see
QSscoreEntity.chem_groups
) is mapped according to highest sequence identity. Alignment is between longest sequences in groups.Limitations:
- If different numbers of groups, we map only the groups for the complex with less groups (rest considered unmapped and shown as warning)
- The mapping is forced: the “best” mapping will be chosen independently of how low the seq. identity may be
Getter: Computed on first use (cached) Type: dict
with key =tuple
of chain names inqs_ent_1
and value =tuple
of chain names inqs_ent_2
.Raises: QSscoreError
if we end up having less than 2 chains for either entity in the mapping (can happen if chains do not have CA atoms).
-
clustalw_bin
¶ Full path to
clustalw
orclustalw2
executable to use for multiple sequence alignments (unlesschain_mapping
is provided manually).Getter: Located in path on first use (cached) Type: str
-
ent_to_cm_1
¶ Subset of
qs_ent_1
used to compute chain mapping and symmetries.Properties:
- Includes only residues aligned according to
chem_mapping
- Includes only 1 CA atom per residue
- Has at least 5 and at most
max_ca_per_chain_for_cm
atoms per chain - All chains of the same chemical group have the same number of atoms
(also in
ent_to_cm_2
according tochem_mapping
) - All chains appearing in
chem_mapping
appear in this entity (so the two can be safely used together)
This entity might be transformed (i.e. all positions rotated/translated by same transformation matrix) if this can speed up computations. So do not assume fixed global positions (but relative distances will remain fixed).
Getter: Computed on first use (cached) Type: EntityHandle
Raises: QSscoreError
if any chain ends up having less than 5 res.- Includes only residues aligned according to
-
ent_to_cm_2
¶ Subset of
qs_ent_1
used to compute chain mapping and symmetries (seeent_to_cm_1
for details).
-
global_score
¶ QS-score with penalties.
The range of the score is between 0 (i.e. no interface residues are shared between biounits) and 1 (i.e. the interfaces are identical).
The global QS-score is computed applying penalties when interface residues or entire chains are missing (i.e. anything that is not mapped in
mapped_residues
/chain_mapping
) in one of the biounits.Getter: Computed on first use (cached) Type: float
-
lddt_mdl
¶ The model entity used for lDDT scoring (
lddt_score
) and annotated with local scores.Local scores are available as residue properties named ‘lddt’ and on each atom as a B-factor. Only CA atoms are considered if
calpha_only
is True, otherwise this is an all-atom score.Since, the lDDT computation requires a single chain with mapped residue numbering, all chains are appended into a single chain X with unique residue numbers according to the column-index in the alignment. The alignments are in the same order as they appear in
alignments
. Additional residues are appended at the end of the chain with unique residue numbers.Getter: Computed on first use (cached) Type: EntityHandle
-
lddt_ref
¶ The reference entity used for lDDT scoring (
lddt_score
).This is a single chain X with residue numbers matching ones in
lddt_mdl
where aligned and unique numbers for additional residues.Getter: Computed on first use (cached) Type: EntityHandle
-
lddt_score
¶ The multi-chain lDDT score.
Note
lDDT is not considering over-prediction (i.e. extra chains) and hence is not symmetric. Here, we consider
qs_ent_1
as the reference andqs_ent_2
as the model. The alignments fromalignments
are used to map residue numbers and chains.The score is computed with OST’s
LocalDistDiffTest()
function with a single distance threshold of 2 A and an inclusion radius of 8 A. You can uselddt_mdl
andlddt_ref
to get entities on which you can call any other lDDT function with any other set of parameters.Getter: Computed on first use (cached) Type: float
-
mapped_residues
¶ Mapping of shared residues in
alignments
.Getter: Computed on first use (cached) Type: dict
mapped_residues[c1][r1] = r2 with: c1 = Chain name in first entity (= first sequence in aln), r1 = Residue number in first entity, r2 = Residue number in second entity
-
superposition
¶ Superposition result based on shared CA atoms in
alignments
.The superposition can be used to map
QSscoreEntity.ent
ofqs_ent_1
onto the one ofqs_ent_2
. Useost.geom.Invert()
if you need the opposite transformation.Getter: Computed on first use (cached) Type: ost.mol.alg.SuperpositionResult
-
symm_1
¶ Symmetry groups for
qs_ent_1
used to speed up chain mapping.This is a list of chain-lists where each chain-list can be used reconstruct the others via cyclic C or dihedral D symmetry. The first chain-list is used as a representative symmetry group. For heteromers, the group-members must contain all different seqres in oligomer.
Example: symm. groups [(A,B,C), (D,E,F), (G,H,I)] means that there are symmetry transformations to get (D,E,F) and (G,H,I) from (A,B,C).
Properties:
- All symmetry group tuples have the same length (num. of chains)
- All chains in
ent_to_cm_1
appear (w/o duplicates) - For heteros: symmetry group tuples have all different chem. groups
- Trivial symmetry group = one tuple with all chains (used if inconsistent data provided or if no symmetry is found)
- Either compatible to
symm_2
or trivial symmetry groups used. Compatibility requires same lengths of symmetry group tuples and it must be possible to get an overlap (80% of residues covered within 6 A of a (chem. mapped) chain) of all chains in representative symmetry groups by superposing one pair of chains.
Getter: Computed on first use (cached) Type: list
oftuple
ofstr
(chain names)
- ent_1 (
-
class
QSscoreEntity
(ent)¶ Entity with cached entries for QS scoring.
Any known / precomputed information can be filled into the appropriate attribute here as long as they are labelled as read/write. Otherwise the quantities are computed on first access and cached (lazy evaluation). The heaviest load is expected when computing
contacts
andcontacts_ca
.Parameters: ent ( EntityHandle
orEntityView
) – Entity to be used for QS scoring. A copy of it will be processed.-
is_valid
¶ True, if successfully initialized. False, if input structure is monomer or has less than 2 protein chains with >= 20 residues.
Type: bool
-
original_name
¶ Name set for ent when object was created.
Type: str
-
ent
¶ Cleaned version of ent passed at construction. Hydrogens are removed, the entity is processed with a
RuleBasedProcessor
and chains listed inremoved_chains
have been removed. The name of this entity might change during scoring (seeGetName()
). Otherwise, this will be fixed.Type: EntityHandle
-
removed_chains
¶ Chains removed from ent passed at construction. These are ligand and water chains as well as small (< 20 res.) peptides or chains with no amino acids (determined by chem. type, which is set by rule based processor).
Type: list
ofstr
-
calpha_only
¶ Whether entity is CA-only (i.e. it has 0 CB atoms)
Type: bool
-
GetAlignment
(c1, c2)¶ Get sequence alignment of chain c1 with chain c2. Computed on first use based on
ca_chains
(cached).Parameters: - c1 (
str
) – Chain name for first chain to align. - c2 (
str
) – Chain name for second chain to align.
Return type: AlignmentHandle
or None if it failed.- c1 (
-
GetAngles
(c1, c2)¶ Get Euler angles from superposition of chain c1 with chain c2. Computed on first use based on
ca_chains
(cached).Parameters: - c1 (
str
) – Chain name for first chain to superpose. - c2 (
str
) – Chain name for second chain to superpose.
Returns: 3 Euler angles (may contain nan if something fails).
Return type: numpy.array
- c1 (
-
GetAxis
(c1, c2)¶ Get axis of symmetry from superposition of chain c1 with chain c2. Computed on first use based on
ca_chains
(cached).Parameters: - c1 (
str
) – Chain name for first chain to superpose. - c2 (
str
) – Chain name for second chain to superpose.
Returns: Rotational axis (may contain nan if something fails).
Return type: numpy.array
- c1 (
-
GetName
()¶ Wrapper to
GetName()
ofent
. This is used to uniquely identify the entity while scoring. The name may therefore change whileoriginal_name
remains fixed.
-
SetName
(new_name)¶ Wrapper to
SetName()
ofent
. Use this to change unique identifier while scoring (seeGetName()
).
-
ca_chains
¶ Map of chain names in
ent
to sequences with attached view to CA-only chains (intoca_entity
). Useful for alignments and superpositions.Getter: Computed on first use (cached) Type: dict
(key =str
, value =SequenceHandle
)
-
ca_entity
¶ Reduced representation of
ent
with only CA atoms. This guarantees that each included residue has exactly one atom.Getter: Computed on first use (cached) Type: EntityHandle
-
chem_groups
¶ Intra-complex group of chemically identical (seq. id. > 95%) polypeptide chains as extracted from
ca_chains
. First chain in group is the one with the longest sequence.Getter: Computed on first use (cached) Type: list
oflist
ofstr
(chain names)
-
contacts
¶ Connectivity dictionary (read/write). As given by
GetContacts()
with calpha_only = False onent
.Getter: Computed on first use (cached) Setter: Uses FilterContacts()
to ensure that we only keep contacts for chains in the cleaned entity.Type: See return type of GetContacts()
-
contacts_ca
¶ CA-only connectivity dictionary (read/write). Like
contacts
but with calpha_only = True inGetContacts()
.
-
-
FilterContacts
(contacts, chain_names)¶ Filter contacts to contain only contacts for chains in chain_names.
Parameters: - contacts (
dict
) – Connectivity dictionary as produced byGetContacts()
. - chain_names (
list
or (better)set
) – Chain names to keep.
Returns: New connectivity dictionary (format as in
GetContacts()
)Return type: dict
- contacts (
-
GetContacts
(entity, calpha_only, dist_thr=12.0)¶ Get inter-chain contacts of a macromolecular entity.
Contacts are pairs of residues within a given distance belonging to different chains. They are stored once per pair and include the CA/CB-CA/CB distance.
Parameters: - entity (
EntityHandle
orEntityView
) – An entity to check connectivity for. - calpha_only (
bool
) – If True, we only consider CA-CA distances. Else, we use CB unless the residue is a GLY. - dist_thr (
float
) – Maximal CA/CB-CA/CB distance to be considered in contact.
Returns: A connectivity dictionary. A pair of residues with chain names ch_name1 & ch_name2 (ch_name1 < ch_name2), residue numbers res_num1 & res_num2 and distance dist (<= dist_thr) are stored as result[ch_name1][ch_name2][res_num1][res_num2] = dist.
Return type: dict
- entity (
Steric Clashes¶
The following function detects steric clashes in atomic structures. Two atoms are clashing if their euclidian distance is smaller than a threshold value (minus a tolerance offset).
-
FilterClashes
(entity, clashing_distances, always_remove_bb=False)¶ This function filters out residues with non-bonded clashing atoms. If the clashing atom is a backbone atom, the complete residue is removed from the structure, if the atom is part of the sidechain, only the sidechain atoms are removed. This behavior is changed by the always_remove_bb flag: when the flag is set to True the whole residue is removed even if a clash is just detected in the side-chain.
The function returns a view containing all elements (residues, atoms) that have not been removed from the input structure, plus a
ClashingInfo
object containing information about the detected clashes.Two atoms are defined as clashing if their distance is shorter than the reference distance minus a tolerance threshold. The information about the clashing distances and the tolerance thresholds for all possible pairs of atoms is passed to the function as a parameter.
Hydrogen and deuterium atoms are ignored by this function.
Parameters: - entity (
EntityView
orEntityHandle
) – The input entity - clashing_distances (
ClashingDistances
) – information about the clashing distances - always_remove_bb (
bool
) – if set to True, the whole residue is removed even if the clash happens in the side-chain
Returns: A tuple of two elements: The filtered
EntityView
, and aClashingInfo
object- entity (
-
CheckStereoChemistry
(entity, bond_stats, angle_stats, bond_tolerance, angle_tolerance, always_remove_bb=False)¶ This function filters out residues with severe stereo-chemical violations. If the violation involves a backbone atom, the complete residue is removed from the structure, if it involves an atom that is part of the sidechain, only the sidechain is removed. This behavior is changed by the always_remove_bb flag: when the flag is set to True the whole residue is removed even if a violation is just detected in the side-chain.
The function returns a view containing all elements (residues, atoms) that have not been removed from the input structure, plus a
StereoChemistryInfo
object containing information about the detected stereo-chemical violations.A violation is defined as a bond length that lies outside of the range: [mean_length-std_dev*bond_tolerance, mean_length+std_dev*bond_tolerance] or an angle width outside of the range [mean_width-std_dev*angle_tolerance, mean_width+std_dev*angle_tolerance ]. The information about the mean lengths and widths and the corresponding standard deviations is passed to the function using two parameters.
Hydrogen and deuterium atoms are ignored by this function.
Parameters: - entity (
EntityView
orEntityHandle
) – The input entity - bond_stats (
StereoChemicalParams
) – statistics about bond lengths - angle_stats (
StereoChemicalParams
) – statistics about angle widths - bond_tolerance (
float
) – tolerance for bond lengths (in standard deviations) - angle_tolerance (
float
) – tolerance for angle widths (in standard deviations) - always_remove_bb (
bool
) – if set to True, the whole residue is removed even if the clash happens in the side-chain
Returns: A tuple of two elements: The filtered
EntityView
, and aStereoChemistryInfo
object- entity (
-
class
ClashingInfo
¶ This object is returned by the
FilterClashes()
function, and contains information about the clashes detected by the function.-
GetClashCount
()¶ Returns: number of clashes between non-bonded atoms detected in the input structure
-
GetAverageOffset
()¶ Returns: a value in Angstroms representing the average offset by which clashing atoms lie closer than the minimum acceptable distance (which of course differs for each possible pair of elements)
-
GetClashList
()¶ Returns: list of detected inter-atomic clashes Return type: list
ofClashEvent
-
-
class
ClashEvent
¶ This object contains all the information relative to a single clash detected by the
FilterClashes()
function-
GetFirstAtom
()¶ -
GetSecondAtom
()¶ Returns: atoms which clash Return type: UniqueAtomIdentifier
-
GetModelDistance
()¶ Returns: distance (in Angstroms) between the two clashing atoms as observed in the model
-
GetAdjustedReferenceDistance
()¶ Returns: minimum acceptable distance (in Angstroms) between the two atoms involved in the clash, as defined in ClashingDistances
-
-
class
StereoChemistryInfo
¶ This object is returned by the
CheckStereoChemistry()
function, and contains information about bond lengths and planar angle widths in the structure that diverge from the parameters tabulated by Engh and Huber in the International Tables of Crystallography. Only elements that diverge from the tabulated value by a minimumnumber of standard deviations (defined when the CheckStereoChemistry function is called) are reported.-
GetBadBondCount
()¶ Returns: number of bonds where a serious violation was detected
-
GetBondCount
()¶ Returns: total number of bonds in the structure checked by the CheckStereoChemistry function
-
GetAvgZscoreBonds
()¶ Returns: average z-score of all the bond lengths in the structure, computed using Engh and Huber’s mean and standard deviation values
-
GetBadAngleCount
()¶ Returns: number of planar angles where a serious violation was detected
-
GetAngleCount
()¶ Returns: total number of planar angles in the structure checked by the CheckStereoChemistry function
-
GetAvgZscoreAngles
()¶ Returns: average z-score of all the planar angle widths, computed using Engh and Huber’s mean and standard deviation values.
-
GetBondViolationList
()¶ Returns: list of bond length violations detected in the structure Return type: list
ofStereoChemicalBondViolation
-
GetAngleViolationList
()¶ Returns: list of angle width violations detected in the structure Return type: list
ofStereoChemicalAngleViolation
-
-
class
StereoChemicalBondViolation
¶ This object contains all the information relative to a single detected violation of stereo-chemical parameters in a bond length
-
GetFirstAtom
()¶ -
GetSecondAtom
()¶ Returns: first / second atom of the bond Return type: UniqueAtomIdentifier
-
GetBondLength
()¶ Returns: length of the bond (in Angstroms) as observed in the model
-
GetAllowedRange
()¶ Returns: allowed range of bond lengths (in Angstroms), according to Engh and Huber’s tabulated parameters and the tolerance threshold used when the CheckStereoChemistry()
function was calledReturn type: tuple
(minimum and maximum)
-
-
class
StereoChemicalAngleViolation
¶ This object contains all the information relative to a single detected violation of stereo-chemical parameters in a planar angle width
-
GetFirstAtom
()¶ -
GetSecondAtom
()¶ -
GetThirdAtom
()¶ Returns: first / second (vertex) / third atom that defines the planar angle Return type: UniqueAtomIdentifier
-
GetAngleWidth
()¶ Returns: width of the planar angle (in degrees) as observed in the model
-
GetAllowedRange
()¶ Returns: allowed range of angle widths (in degrees), according to Engh and Huber’s tabulated parameters and the tolerance threshold used when the CheckStereoChemistry()
function was calledReturn type: tuple
(minimum and maximum)
-
-
class
ClashingDistances
¶ Object containing information about clashing distances between non-bonded atoms
-
ClashingDistances
()¶ Creates an empty distance list
-
SetClashingDistance
(ele1, ele2, clash_distance, tolerance)¶ Adds or replaces an entry in the list
Parameters: - ele1 – string containing the first element’s name
- ele2 – string containing the second element’s name
- clash_distance – minimum clashing distance (in Angstroms)
- tolerance – tolerance threshold (in Angstroms)
-
GetClashingDistance
(ele1, ele2)¶ Returns: reference distance and a tolerance threshold (both in Angstroms) for two elements
Return type: tuple
(minimum clashing distance, tolerance threshold)Parameters: - ele1 – string containing the first element’s name
- ele2 – string containing the second element’s name
-
GetAdjustedClashingDistance
(ele1, ele2)¶ Returns: reference distance (in Angstroms) for two elements, already adjusted by the tolerance threshold
Parameters: - ele1 – string containing the first element’s name
- ele2 – string containing the second element’s name
-
GetMaxAdjustedDistance
()¶ Returns: longest clashing distance (in Angstroms) in the list, after adjustment with tolerance threshold
-
IsEmpty
()¶ Returns: True if the list is empty (i.e. in an invalid, useless state)
-
PrintAllDistances
()¶ Prints all distances in the list to standard output
-
-
class
StereoChemicalParams
¶ Object containing stereo-chemical information about bonds and angles. For each item (bond or angle in a specific residue), stores the mean and standard deviation
-
StereoChemicalParams
()¶ Creates an empty parameter list
-
SetParam
(item, residue, mean, standard_dev)¶ Adds or replaces an entry in the list
Parameters: - item – string defining a bond (format: X-Y) or an angle (format: X-Y-Z), where X,Y an Z are atom names
- residue – string containing the residue type for this entry
- mean – mean bond length (in Angstroms) or angle width (in degrees)
- standard_dev – standard deviation of the bond length (in Angstroms) or of the angle width (in degrees)
-
IsEmpty
()¶ Returns: True if the list is empty (i.e. in an invalid, useless state)
-
PrintAllParameters
()¶ Prints all entries in the list to standard output
-
-
FillClashingDistances
(file_content)¶ -
FillBondStereoChemicalParams
(file_content)¶ -
FillAngleStereoChemicalParams
(file_content)¶ These three functions fill a list of reference clashing distances, a list of stereo-chemical parameters for bonds and a list of stereo-chemical parameters for angles, respectively, starting from the content of a parameter file.
Parameters: file_content ( list
ofstr
) – list of lines from the parameter fileReturn type: ClashingDistances
orStereoChemicalParams
-
FillClashingDistancesFromFile
(filename)¶ -
FillBondStereoChemicalParamsFromFile
(filename)¶ -
FillAngleStereoChemicalParamsFromFile
(filename)¶ These three functions fill a list of reference clashing distances, a list of stereo-chemical parameters for bonds and a list of stereo-chemical parameters for angles, respectively, starting from a file path.
Parameters: filename ( str
) – path to parameter fileReturn type: ClashingDistances
orStereoChemicalParams
-
DefaultClashingDistances
()¶ -
DefaultBondStereoChemicalParams
()¶ -
DefaultAngleStereoChemicalParams
()¶ These three functions fill a list of reference clashing distances, a list of stereo-chemical parameters for bonds and a list of stereo-chemical parameters for angles, respectively, using the default parameter files distributed with OpenStructure.
Return type: ClashingDistances
orStereoChemicalParams
-
ResidueNamesMatch
(probe, reference)¶ The function requires a reference structure and a probe structure. The function checks that all the residues in the reference structure that appear in the probe structure (i.e., that have the same ResNum) are of the same residue type. Chains are comapred by order, not by chain name (i.e.: the first chain of the reference will be compared with the first chain of the probe structure, etc.)
Parameters: - probe (
EntityView
) – the structure to test - reference (
EntityView
) – the reference structure
Returns: True if the residue names are the same, False otherwise
- probe (
Superposing structures¶
-
Superpose
(ent_a, ent_b, match='number', atoms='all', iterative=False, max_iterations=5, distance_threshold=3.0)¶ Superposes the model entity onto the reference. To do so, two views are created, returned with the result. atoms describes what goes into these views and match the selection method. For superposition,
SuperposeSVD()
orIterativeSuperposeSVD()
are called (depending on iterative). For matching, the following methods are recognised:number
- select residues by residue number, includes atoms, callsMatchResidueByNum()
index
- select residues by index in chain, includes atoms, callsMatchResidueByIdx()
local-aln
- select residues from a Smith/Waterman alignment, includes atoms, callsMatchResidueByLocalAln()
global-aln
- select residues from a Needleman/Wunsch alignment, includes atoms, callsMatchResidueByGlobalAln()
Parameters: - ent_a (
EntityView
orEntityHandle
) – The model entity (superposition transform is applied on full entity handle here) - ent_b (
EntityView
orEntityHandle
) – The reference entity - match (
str
) – Method to gather residues/ atoms - atoms (
str
,list
,set
) – The subset of atoms to be used in the superposition - iterative (
bool
) – Whether or not to use iterative superpositon. - max_iterations (
int
) – Max. number of iterations forIterativeSuperposeSVD()
(only if iterative = True) - distance_threshold (
float
) – Distance threshold forIterativeSuperposeSVD()
(only if iterative = True)
Returns: An instance of
SuperpositionResult
.
-
ParseAtomNames
(atoms)¶ Parses different representations of a list of atom names and returns a
set
, understandable byMatchResidueByNum()
. In essence, this function translates- None to
None
- ‘all’ to
None
- ‘backbone’ to
set(['N', 'CA', 'C', 'O'])
- ‘aname1, aname2’ to
set(['aname1', 'aname2'])
['aname1', 'aname2']
toset(['aname1', 'aname2'])
Parameters: atoms ( str
,list
,set
) – Identifier or list of atomsReturns: A set
of atoms.- None to
-
MatchResidueByNum
(ent_a, ent_b, atoms='all')¶ Returns a tuple of views containing exactly the same number of atoms. Residues are matched by residue number. A subset of atoms to be included in the views can be specified in the atoms argument. Regardless of what the list of atoms says, only those present in two matched residues will be included in the views. Chains are processed in the order they occur in the entities. If ent_a and ent_b contain a different number of chains, processing stops with the lower count.
Parameters: - ent_a (
EntityView
orEntityHandle
) – The first entity - ent_b (
EntityView
orEntityHandle
) – The second entity - atoms (
str
,list
,set
) – The subset of atoms to be included in the two views.
Returns: Two
EntityView
instances with the same amount of residues matched by number. Each residue will have the same number & type of atoms.- ent_a (
-
MatchResidueByIdx
(ent_a, ent_b, atoms='all')¶ Returns a tuple of views containing exactly the same number of atoms. Residues are matched by position in the chains of an entity. A subset of atoms to be included in the views can be specified in the atoms argument. Regardless of what the list of atoms says, only those present in two matched residues will be included in the views. Chains are processed in order of appearance. If ent_a and ent_b contain a different number of chains, processing stops with the lower count. The number of residues per chain is supposed to be the same.
Parameters: - ent_a (
EntityView
orEntityHandle
) – The first entity - ent_b (
EntityView
orEntityHandle
) – The second entity - atoms (
str
,list
,set
) – The subset of atoms to be included in the two views.
Returns: Two
EntityView
instances with the same amount of residues matched by position. Each residue will have the same number & type of atoms.- ent_a (
-
MatchResidueByLocalAln
(ent_a, ent_b, atoms='all')¶ Match residues by local alignment. Takes ent_a and ent_b, extracts the sequences chain-wise and aligns them in Smith/Waterman manner using the BLOSUM62 matrix for scoring. Only residues which are marked as
peptide linking
are considered for alignment. The residues of the entities are then matched based on this alignment. Only atoms present in both residues are included in the views. Chains are processed in order of appearance. If ent_a and ent_b contain a different number of chains, processing stops with the lower count.Parameters: - ent_a (
EntityView
orEntityHandle
) – The first entity - ent_b (
EntityView
orEntityHandle
) – The second entity - atoms (
str
,list
,set
) – The subset of atoms to be included in the two views.
Returns: Two
EntityView
instances with the same number of residues. Each residue will have the same number & type of atoms.- ent_a (
-
MatchResidueByGlobalAln
(ent_a, ent_b, atoms='all')¶ Match residues by global alignment. Same as
MatchResidueByLocalAln()
but performs a global Needleman/Wunsch alignment of the sequences using the BLOSUM62 matrix for scoring.Parameters: - ent_a (
EntityView
orEntityHandle
) – The first entity - ent_b (
EntityView
orEntityHandle
) – The second entity - atoms (
str
,list
,set
) – The subset of atoms to be included in the two views.
Returns: Two
EntityView
instances with the same number of residues. Each residue will have the same number & type of atoms.- ent_a (
-
class
SuperpositionResult
¶ -
rmsd
¶ RMSD of the superposed entities.
-
view1
¶ -
view2
¶ Two
EntityView
used in superposition (not set if methods withVec3List
used).
-
fraction_superposed
¶ -
rmsd_superposed_atoms
¶ -
ncycles
¶ For iterative superposition (
IterativeSuperposeSVD()
): fraction and RMSD of atoms that were superposed with a distance below the given threshold and the number of iteration cycles performed.
-
-
SuperposeSVD
(view1, view2, apply_transform=True)¶ -
SuperposeSVD
(list1, list2) Superposition of two sets of atoms minimizing RMSD using a classic SVD based algorithm.
Note that the atom positions in the view are taken blindly in the order in which the atoms appear.
Parameters: - view1 (
EntityView
) – View on the model entity - view2 (
EntityView
) – View on the reference entity - list1 (
Vec3List
) – List of atom positions for model entity - list2 (
Vec3List
) – List of atom positions for reference entity - apply_transform (
bool
) – If True, the superposition transform is applied to the (full!) entity handle linked to view1.
Returns: An instance of
SuperpositionResult
.- view1 (
-
IterativeSuperposeSVD
(view1, view2, max_iterations=5, distance_threshold=3.0, apply_transform=True)¶ -
IterativeSuperposeSVD
(list1, list2, max_iterations=5, distance_threshold=3.0) Iterative superposition of two sets of atoms. In each iteration cycle, we keep a fraction of atoms with distances below distance_threshold and get the superposition considering only those atoms.
Note that the atom positions in the view are taken blindly in the order in which the atoms appear.
Parameters: - view1 (
EntityView
) – View on the model entity - view2 (
EntityView
) – View on the reference entity - list1 (
Vec3List
) – List of atom positions for model entity - list2 (
Vec3List
) – List of atom positions for reference entity - max_iterations (
int
) – Max. number of iterations to be performed - distance_threshold (
float
) – Distance threshold defining superposed atoms - apply_transform (
bool
) – If True, the superposition transform is applied to the (full!) entity handle linked to view1.
Returns: An instance of
SuperpositionResult
.Raises: Exception if atom counts do not match or if less than 3 atoms.
- view1 (
-
CalculateRMSD
(view1, view2, transformation=geom.Mat4())¶ Returns: RMSD of atom positions (taken blindly in the order in which the atoms appear) in the two given views.
Return type: float
Parameters: - view1 (
EntityView
) – View on the model entity - view2 (
EntityView
) – View on the reference entity - transformation (
Mat4
) – Optional transformation to apply on each atom position of view1.
- view1 (
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 surface - probe_radius (
float
) – Radius of probe to determine accessible surface area - include_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.
- 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
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.
-
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.
- 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.
- frames (
-
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.
- traj (
-
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 computed - stride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- traj (
-
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.
- traj (
-
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.
- traj (
-
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.
- traj (
-
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 computed - sele2 (
EntityView
.) – The selection from which the second center of mass is computed - stride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- traj (
-
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 structure - sele_view (
EntityView
.) – The selection compared to the reference_view - stride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- traj (
-
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 atoms - view2 (
EntityView
.) – The second group of atoms - stride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- traj (
-
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 taken - view_atoms (
EntityView
.) – The second group of atoms - stride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- traj (
-
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 atoms - view_ring2 (
EntityView
.) – Second group of atoms - stride – Size of the increment of the frame’s index between two consecutive frames analyzed.
- traj (
helix_kinks
– Algorithms to calculate Helix Kinks¶
Functions to calculate helix kinks: bend, face shift and wobble angles
Author: Niklaus Johner
-
AnalyzeHelixKink
(t, sele, proline=False)¶ This function calculates the bend, wobble and face-shift angles in an alpha- helix over a trajectory. The determination is more stable if there are at least 4 residues on each side (8 is even better) of the proline around which the helix is kinked. The selection should contain all residues in the correct order and with no gaps and no missing C-alphas.
Parameters: - t (
CoordGroup
) – The trajectory to be analyzed - sele (
EntityView
) – A selection containing the alpha helix to be analyzed - proline (
ost.mol.EntityView
) – A selection containing only the proline (or another residue) around which the helix is kinked. If False, the proline will be searched for automatically
Returns: A tuple (bend_angle, face_shift, wobble_angle).
Return type: (FloatList, FLoatList, FloatList)
- t (
-
CalculateHelixKink
(sele, proline=False)¶ This function calculates the bend, wobble and face-shift angles in an alpha- helix of an EntityView. The determination is more stable if there are at least 4 residues on each side (8 is even better) of the proline around which the helix is kinked. The selection should contain all residues in the correct order and with no gaps and no missing C-alphas.
Parameters: - sele (
EntityView
) – A selection containing the alpha helix to be analyzed - proline (
ost.mol.EntityView
) – A selection containing only the proline (or another residue) around which the helix is kinked. If False, the proline will be searched for automatically
Returns: A tuple (bend_angle, face_shift, wobble_angle).
Return type: (float, float, float)
- sele (
trajectory_analysis
– DRMSD, pairwise distances and more¶
This Module requires numpy
This module contains functions to analyze trajectories, mainly similiraty measures baed on RMSDS and pairwise distances.
Author: Niklaus Johner (niklaus.johner@unibas.ch)
-
AverageDistanceMatrixFromTraj
(t, sele, first=0, last=-1)¶ This function calcultes the distance between each pair of atoms in sele, averaged over the trajectory t.
Parameters: - t (
CoordGroupHandle
) – the trajectory - sele (
EntityView
) – the selection used to determine the atom pairs - first (
int
) – the first frame of t to be used - last (
int
) – the last frame of t to be used
Returns: a numpy NpairsxNpairs matrix, where Npairs is the number of atom pairs in sele.
- t (
-
DistRMSDFromTraj
(t, sele, ref_sele, radius=7.0, average=False, seq_sep=4, first=0, last=-1)¶ This function calculates the distance RMSD from a trajectory. The distances selected for the calculation are all the distances between pair of atoms from residues that are at least seq_sep apart in the sequence and that are smaller than radius in ref_sel. The number and order of atoms in ref_sele and sele should be the same.
Parameters: - t (
CoordGroupHandle
) – the trajectory - sele (
EntityView
) – the selection used to calculate the distance RMSD - ref_sele (
EntityView
) – the reference selection used to determine the atom pairs and reference distances - radius (
float
) – the upper limit of distances in ref_sele considered for the calculation - seq_sep (
int
) – the minimal sequence separation between atom pairs considered for the calculation - average (
bool
) – use the average distance in the trajectory as reference instead of the distance obtained from ref_sele - first (
int
) – the first frame of t to be used - last (
int
) – the last frame of t to be used
Returns: a numpy vecor dist_rmsd(Nframes).
- t (
-
DistanceMatrixFromPairwiseDistances
(distances, p=2)¶ This function calculates an distance matrix M(NframesxNframes) from the pairwise distances matrix D(NpairsxNframes), where Nframes is the number of frames in the trajectory and Npairs the number of atom pairs. M[i,j] is the distance between frame i and frame j calculated as a p-norm of the differences in distances from the two frames (distance-RMSD for p=2).
Parameters: - distances – a pairwise distance matrix as obtained from
PairwiseDistancesFromTraj()
- p – exponent used for the p-norm.
Returns: a numpy NframesxNframes matrix, where Nframes is the number of frames.
- distances – a pairwise distance matrix as obtained from
-
PairwiseDistancesFromTraj
(t, sele, first=0, last=-1, seq_sep=1)¶ This function calculates the distances between any pair of atoms in sele with sequence separation larger than seq_sep from a trajectory t. It return a matrix containing one line for each atom pair and Nframes columns, where Nframes is the number of frames in the trajectory.
Parameters: - t (
CoordGroupHandle
) – the trajectory - sele (
EntityView
) – the selection used to determine the atom pairs - first (
int
) – the first frame of t to be used - last (
int
) – the last frame of t to be used - seq_sep (
int
) – The minimal sequence separation between atom pairs
Returns: a numpy NpairsxNframes matrix.
- t (
-
RMSD_Matrix_From_Traj
(t, sele, first=0, last=-1, align=True, align_sele=None)¶ This function calculates a matrix M such that M[i,j] is the RMSD (calculated on sele) between frames i and j of the trajectory t aligned on sele.
Parameters: - t (
CoordGroupHandle
) – the trajectory - sele (
EntityView
) – the selection used for alignment and RMSD calculation - first (
int
) – the first frame of t to be used - last (
int
) – the last frame of t to be used
Returns: Returns a numpy NframesxNframes matrix, where Nframes is the number of frames.
- t (
structure_analysis
– Functions to analyze structures¶
Some functions for analyzing structures
Author: Niklaus Johner (Niklaus.Johner@unibas.ch)
-
CalculateBestFitLine
(sele1)¶ This function calculates the best fit line to the atoms in sele1.
Parameters: sele1 ( EntityView
) –Returns: Line3
-
CalculateBestFitPlane
(sele1)¶ This function calculates the best fit plane to the atoms in sele1.
Parameters: sele1 ( EntityView
) –Returns: Plane
-
CalculateDistanceDifferenceMatrix
(sele1, sele2)¶ This function calculates the pairwise distance differences between two selections (
EntityView
). The two selections should have the same number of atoms It returns an NxN DistanceDifferenceMatrix M (where N is the number of atoms in sele1) where M[i,j]=||(sele2.atoms[i].pos-sele2.atoms[j].pos)||-||(sele1.atoms[i].pos-sele1.atoms[j].pos)||Parameters: - sele1 (
EntityView
) – - sele2 (
EntityView
) –
Returns: NxN numpy matrix
- sele1 (
-
CalculateHelixAxis
(sele1)¶ This function calculates the best fit cylinder to the CA atoms in sele1, and returns its axis. Residues should be ordered correctly in sele1.
Parameters: sele1 ( EntityView
) –Returns: Line3
-
GetAlphaHelixContent
(sele1)¶ This function calculates the content of alpha helix in a view. All residues in the view have to ordered and adjacent (no gaps allowed)
Parameters: sele1 ( EntityView
) –Returns: float
-
GetDistanceBetwCenterOfMass
(sele1, sele2)¶ This function calculates the distance between the centers of mass of sele1 and sele2, two selections from the same Entity.
Parameters: - sele1 (
EntityView
) – - sele2 (
EntityView
) –
Returns: float
- sele1 (
-
GetFrameFromEntity
(eh)¶ This function returns a CoordFrame from an EntityHandle
Parameters: eh ( EntityHandle
) –Returns: ost.mol.CoordFrame
-
GetMinDistBetwCenterOfMassAndView
(sele1, sele2)¶ This function calculates the minimal distance between sele2 and the center of mass of sele1, two selections from the same Entity.
Parameters: - sele1 (
EntityView
) – The selection from which the center of mass is taken - sele2 (
EntityView
) –
Returns: distance (
float
)- sele1 (
-
GetMinDistanceBetweenViews
(sele1, sele2)¶ This function calculates the minimal distance between sele1 and sele2, two selections from the same Entity.
Parameters: - sele1 (
EntityView
) – - sele2 (
EntityView
) –
Returns: float
- sele1 (