torsion module

Protein backbone torsion angle parameterization.

Module author: Brynmor K. Chapman.

Module author: Michael S. Chapman <chapmami@ohsu.edu> (secondary).

Authors:

Michael S. Chapman <chapmami@ohsu.edu>,

Oregon Health & Science University

Version:

0.5, March 23, 2016

Changed in version 08/19/2011.

Changed in version 0.3.7: 11/1/12 (MSC) Consolidating from TAProtein.py & TAProtein_copy.py which will be deprecated. Some parts to be moved to transform.py. Improving documentation, consistency with other parts of RSRef.

Changed in version 0.5.0: 05/04/15 ReStructuredText docs

Changed in version 0.5.2: 08/26/15 prior_main, next_main replaced by bonded_main Should now be direction-neutral, and therefore robust to different ordering of atoms in the PDB file, or to arbitrary directions in the edge determination in the bond graph.

Changed in version 0.5.3: 09/27/15 get_all_psi() & other locations now tolerate missing OXT, issuing warning.

Changed in version 0.5.5: 10/13/15 Anisotropic Us can ride with atom rotations.

Todo

Brynmor needs to edit entire documentation.

Todo

Connect HETATMs in same residue, noting that would have to change code near backbone_must_be_polymer.

Warning

Requires atom records in order, N before CA... Else throws an AssertionError exception like:

  • Setting N-side neighbor of atom 7 to 5, but previously set to 0.
  • Setting C-side neighbor of atom 4 to 7, but previously set to 8.

The reset value will correct, but the assertion is checking that it has not been assigned previously. (Line ~401).

Error

This error may have been resolved 10/31/12, not completely sure. The last angle: (a) prior to C-terminal carboxylate; (b) in a selected sub-part of a chain; is not refining properly. The torsion angle changes are tiny such that in (a) the carboxylate is not rotated back into density; and in (b) so that there is not a clean break with the adjoining structure. Brynmor suspects that the masking of moving atoms might be the issue. Could also be that near-last atoms within the refining selection are being fixed (no derivatives). 10/31/12 - Bug may have been resolved. (perhaps only if -0.2 < psi < -1.0?)

Warning

If a chain is not terminated with OXT or O1/O2, then neither phi nor psi of the final residue is refined.

Warning

The following problem was likely fixed by recoding the bonded neighbor code, Summer 2015, but here is the original warning: At/after carboxy terminus, garbage connections are sometimes being added to atoms that follow, even if outside selection. It may be when HETATMS happen to contain backbone-like atom names. dihedral_name() returns “ERROR”, but otherwise program is allowed to continue, because the torsion angles appear not to be set to refine, so perhaps benign, but should try to figure it out.

Todo

Pair() requires an exact match, including of alternate location labels, so alternate conformer may be ignored and neighboring torsion angles excluded. Need more general Pair() that allows: (1) user to specify equivalence of selections, and (2) allows lower stringency to make sensible pairings.

class torsion.Graph(edges=None, vertices=None)

Bases: object

Container for directed/undirected graphs.

Optionally directed graph, if edges provided, or empty if size not None.

Creates a directed graph with the specified edges (if provided) or an empty undirected graph with the specified number of vertices

Parameters:
  • edges ({int:{int}} OR {(int, int)} OR [(int, int)] OR [[int, int]] OR [{int, int}]) – optional map: int -> sets of vertices (edge from u to v represented by v in edges[u]) OR set of pairs of vertices (edge from u to v represented by (u, v) in edges)
  • vertices ({int} OR [int] OR int (interpreted as range(vertices))) – the set of vertices in the graph (optional)
addDirectedEdge(u, v)
addEdge(u, v)
addVertex(v)
dfs(vertexSet='all', forceDistinctComponents=False, componentList=None)

Depth-first search, augmented to track order of edge traversal.

Parameters:
  • vertexSet ({int} OR [int]) – the set of vertices to consider when picking a new root (optional, defaults to ‘all’).
  • forceDistinctComponents (bool) – forces search to stop upon reaching a vertex in the given vertexSet
  • componentList ({int: int}) – if provided, creates a list of trees in DFS forest (optional, starts empty, ‘destroyed’ by method)
Returns:

visit_times, finish_times, edge_list

Return type:

([int], [int], [(int, int)])

dfsVisit(v, edge, edges, visitTimes, finishTimes, stack, nstack, time)
find_articulation_edges()

Finds the articulation edges in an undirected graph.

Returns:articulation_edges, a list of edges which subdivide a connected component when removed.
Return type:[(int, int)]
find_bridgeless_components(startPoints='all')

Finds the bridgeless components.

These are the maximal connected subgraphs (without articulation edges) of an undirected graph.

Returns:bridgeless_components, a list of components, each represented as a list of vertices.
Return type:[[int]]
get_reverse()

Non-destructive reversal of edge directions in a directed graph.

Returns:reverse the new graph with the same vertices and edges, but with edge orientation reversed.
Return type:Graph
removeDirectedEdge(u, v)
removeEdge(u, v)
scc(startPoints='all', forceDistinctComponents=False)

Strongly connected component algorithm.

A connected undirected graph is one in which any pair of vertices is connected by a path. Finds all sub-graphs which are maximal with the connected property in the sense that no additional vertices keep it connected.

Parameters:startPoints (list(int)) – restrict search to given vertices, restricting returned components to those containing these vertices. Expected to work only with undirected graphs. QQQ
Returns:list of the strongly connected components of the graph (each represented as a list of vertices)
Return type:[[int]]
version_ID = None

version number to check whether any changes have been made

class torsion.ImpactRestraint

Bases: object

Refinement restraint to limit change of lowest impact dihedrals.

class torsion.Limit(limit, weight, small=0.0001, ramp=1.01, reference=None, torsions=None, verbose=True)

Bases: object

Refinement restraint on total change in torsion angles.

Changed in version 04/23/14,: restraint can be from optional reference structure instead of current structure.

Tested by regression. Regression on restarted refinements inexact, because cannot control number of evaluations on final cycle of first batch.

Initialize a lasso or harmonic restraint on dihedral changes.

This should be called at the start of a refinement when the objective function is initialized. Without a reference, the restraint will be on changes from the start of a batch of refinement. The reference can be used to restrain always to the input coordinates.

Parameters:
  • limit (float) – sum(abs(torsion changes)) where harmonic restraint starts (degrees).
  • weight (float) – empirical (user-adjusted) weight for the restraint.
  • small (float) – dampen the contributions of torsion changes < small (degrees) with a linear weight that falls to zero at change=0, where there is a discontinuity in the derivative. (Rotations < small will be less restrained. Rotations of 0.1 deg correspond to to atomic shifts O(0.03Å) in a typical protein.)
  • ramp (float) – for L1Norm only, the restraint ramps up from zero at limit to the full restraint at ramp * limit.
  • reference (Atoms or NoneType) – if torsion angle changes are to be restrained an optional reference state of the moving structure (same atoms, order etc.) can be provided (eg. input structure); else any restraint will apply from the start of each refinement batch.
  • torsions (TorsionGp) – torsion angle parameterization needed if reference.
L1Norm = True

Restraint is a modified L1-Norm (def.), else custom L1sq function. With appropriately different weights, L1-Norm & L1sq can give ~identical structure. L1-Norm is more conventional, less sensitive to exact choice of weight, near-zeros the changes on more parameters and therefore the default. (To change interactively, try EXEC torsion.Limit.L1Norm=False or EXEC torsion.Limit.L1Norm=False).

change(limit=None, weight=None, small=None, continuation=True)

Change the parameters of the torsional L1 restraint.

Parameters:
  • limit (float) – sum(abs(torsion changes)) where harmonic restraint starts (degrees).
  • weight (float) – empirical (user-adjusted) weight for the restraint.
  • small (float) – dampen the contributions of torsion changes < small (degrees) with a linear weight that falls to zero at change=0, where there is a discontinuity in the derivative. (Rotations < small will be less restrained. Rotations of 0.1 deg correspond to to atomic shifts O(0.03Å) in a typical protein.)
  • continuation (bool) – if within an ongoing batch of refinement, a warning will be issued.
Returns:

whether the parameters are actually changed.

Return type:

bool

dfL1Norm(changes, verbose=False)

Partial derivatives for a soft limit on total torsional change.

Let f = W(Σii - α0,i|.di - L).D

f = W(Σii - α0,i|.di).D - WLD

where di = 0 when Σ < L, and increases linearly to 1.0 at Σ = ramp * L

Differentiation is much simplified with the assumption that ∂D/∂αn = 0. This is exactly true in most cases (where Σ < L or Σ > ramp * L), and approximately true for L < Σ < ramp * L, because D depends on hundreds of α. So, proceeding with this approximation...

∂D/∂αn = WD(|αn - α0,n|.∂dn)/∂αn + sign(αn - α0,n).dn

Parameters:changes (ndarray(dtype=float)) – Changes in torsion angles (current - start), degrees. These are passed as an argument (rather than calculated from the structures) on the assumption that they can be more efficiently extracted directly from the refinement parameter array.
Returns:∂f()/∂αn, units = /deg.
Return type:ndarray(dtype=float, size=#torsion_angles)
dfL1sq(changes, verbose=False)

Partial derivatives for a soft limit on total torsional change.

Let Y = (Σii - α0,i|.Di - L)

∂Y/∂αn = 2Y∂Y/∂αn

∂Y/∂αn = |αn - α0,n|.∂Dn)/∂αn + sign(αn - α0,n).Dn

if |αn - α0,n| > S (S being a small rotation):

∂Y/∂αn = sign(αn - α0,n).Dn

else:

∂Y/∂αn = |αn - α0,n| sign(αn - α0,n)/S + Dnsign(αn - α0,n)

∂Y/∂αn = 2Dnsign(αn - α0,n)

Warning

dfL1Norm is very similar, more conventional and easier to adjust weight.

Parameters:changes (ndarray(dtype=float)) – Changes in torsion angles (current - start), degrees. These are passed as an argument (rather than calculated from the structures) on the assumption that they can be more efficiently extracted directly from the refinement parameter array.
Returns:(0 if self.f() <= 0 else ∂f()/∂αn, units = /deg.
Return type:ndarray(dtype=float, size=#torsion_angles)
fL1Norm(changes)

Objective function for a soft limit on total torsional change.

Based on a modified L1-Norm, i.e. max(0, Σii - α0,i|.di - L).DW where di is zero for |αi - α0,i| = 0, avoiding the non-differentiability of the L1-norm at this point. Thus, di avoids the need for more complicated algorithms [Schmidt-2005], at a cost of incompletely restraining small changes, so that they are likely to survive as non-zero. The approach is a modification of the LASSO method, which is similar to Basis Persuit. Limit (target) L is an additional modification of the basic LASSO as is the attenuator D that is required for it to be smoothly differentiable.

[Schmidt-2005]Schmidt M. (2005) Least Squares Optimization with L1-Norm Regularization CS542B Project Report, Dec. 2005.
Parameters:changes (ndarray(dtype=float)) – Changes in torsion angles (current - start), degrees. These are passed as an argument (rather than calculated from the structures) on the assumption that they can be more efficiently extracted directly from the refinement parameter array.
Returns:Σii - α0,i|.di - L).DW where D, W are a damping constant and weight and runs over all i torsion angles, where... Di = max(1, |αi - α0,i|/S) where S is the limit below which a change is considered small.
Return type:float

Changed in version 4/30/14: - introduction of L, D

fL1sq(changes)

Objective function for a soft limit on total torsional change.

Based on square of a modified L1-Norm.

Warning

fL1Norm is very similar, more conventional and easier to adjust weight.

Parameters:changes (ndarray(dtype=float)) – Changes in torsion angles (current - start), degrees. These are passed as an argument (rather than calculated from the structures) on the assumption that they can be more efficiently extracted directly from the refinement parameter array.
Returns:W.max(0, Σii - α0,i|.Di - L)2, where D, L, W are a damping constant, soft limit and weight and runs over all i torsion angles, where... Di = max(1, |αi - α0,i|/S) where S is the limit below which a change is considered small.
Return type:float
static mod360(a)
stats(changes, verbose=False)

Calculate total torsion angle change (current - start).

i - α0,i|.Di)

Returns:total, total (undamped, D=0), number torsion angles, number damped, total change from reference
Return type:float, float, int, int, float
class torsion.TAProtein(atoms, moving_selection=None, link=None, unlink=None)

Bases: object

Protein structure parameterized by backbone dihedrals: φ, ψ (& ω).

Assuming that side chain torsions are invariant, simplifies greatly with a graph that connects all side-chain atoms to CA for rigid rotations. (Would require extensive re-write to extend to all dihedral angles.)

There is currently no special treatment of Prolines. Are we going to restrain possible values, or do we make invariant by selecting out? Currently fully variable.

Parameterize protein as connected atoms with variable φ, ψ.

Parameters:
  • moving_selection (atoms.Atoms.Selection) – Atoms that can be moved within the structure.
  • link ((Group of) 2-atom atoms.Selection) – Additional pair(s) of atoms to be linked for joint rotations about (other) torsion bonds.
  • unlink ((Group of) 2-atom atoms.Selection) – Disconnect pair(s) of atoms from those automatically linked so that not rotated together about torsion bonds.
all_dependents

Atoms dependent on any variable torsion angle.

Requires generate_rotation_masks to have been run with current configuration or new calculation triggered with self._all_dependents set to False.

Returns:True for atoms that depend on torsion(s).
Return type:boolean-like array (atoms.Selection)
apply_torsion(parameters, atoms, refineGroup, backboneOnly=True, update_anisotropic=True, deriv=False)

Applies a set of torsion angle changes.

Aligns the transformed molecule on the original so that there is not an overall rotation. However, this will not achieve a rigid group refinement, so would likely have to combine torsional refinement with overall or group rigid refinement.

Parameters:
  • parameters (array-like, size = #variable torsion angles) – torsion angle changes (from initial structure).
  • atoms (atoms.Atoms) – Initial structure (where parameters = 0).
  • refineGroup (boolean nd.array) – rotatable bonds, as returned in torsion_bonds of phi_psi_selection, get_all_bonds_in_selection etc. Should be consistent with the atom selection used to initialize TAProtein. (documented by MSC).
  • backboneOnly (bool) – True only if refineGroup includes only backbone bonds, then exploit an efficiency.
  • update_anisotropic (bool) – rotate anisotropic Us, else leave (out of date) for a later calculation (if intermediate steps do not depend on them).
  • deriv (bool) – calculate partial derivatives of atomic parameters with respect to torsion angles.
Returns:

updated coordinates, u-tensors, derivatives of coordinates with respect to the torsion angles]

Return type:

ndarray(shape=(3,n_atom)), Anisotropics or NoneType, ndarray(shape=(n_bonds,3,n_atoms))

Changed in version 0.3.7: dXbydP.shape[3] is now moving_atoms.size, not len(atoms). This reduces the size of matrix multiplications at a cost of more complicated indexing. Further savings are possible by similarly redimensioning dXfbydP. This would require changes to calling code and likely lesser savings, because only used in one multiplication.

Warning

refineGroup could be a source of fragility if calculated with a different atom selection than with which TAProtein initialized.

apply_torsion_deprecated(parameters, atoms, refineGroup, backboneOnly=True, deriv=False)

Applies a set of torsion angle changes.

Aligns the transformed molecule on the original so that there is not an overall rotation. However, this will not achieve a rigid group refinement, so would likely have to combine torsional refinement with overall or group rigid refinement.

Deprecated since version 0.3.7.

Parameters:
  • parameters (array-like, size = #variable torsion angles) – torsion angle changes (from initial structure).
  • atoms (atoms.Atoms) – Initial structure (where parameters = 0).
  • refineGroup (boolean nd.array) – rotatable bonds, as returned in torsion_bonds of phi_psi_selection, get_all_bonds_in_selection etc. Should be consistent with the atom selection used to initialize TAProtein. (documented by MSC).
  • backboneOnly (bool) – True only if refineGroup includes only backbone bonds, then exploit an efficiency.
Returns:

updated coordinates [, derivatives of coordinates with respect to the torsion angles]

Return type:

ndarray(shape=(3,n_atom))[, ndarray(shape=(n_bonds,3,n_atoms))]

Warning

refineGroup could be a source of fragility if calculated with a different atom selection than with which TAProtein initialized.

bond_dependents

Atoms dependent on each (variable) torsion angle.

Requires generate_rotation_masks to have been run with current configuration or new calculation triggered with self._all_dependents set to False.

Returns:True for atoms that depend on torsion.
Return type:boolean-like array (atoms.Selection)
compare(target, map2target, plot=True, color=False, pseudo=False, normalize=True, hinge=False, threshold=3.0, gap=1, start=[], end=[], verbose=True)

Compare φ, ψ angles in self to those in target atoms.

Todo

map2target suggests that no longer need target to be ordered like atoms - check that this is true. (2/13/13).

Parameters:
  • target (Atoms) – coordinates within which to compare.
  • map2target (dict) – keys are the indices of self.atoms previously paired and values are the indices of corresponding atoms in target.
  • plot (bool) – graph of phi, psi, combined vs. residue number
  • color (bool) – return copy of self.atoms, B-factors representing differences in nearest φ, ψ (suitable for coloring by PyMol).
  • pseudo (bool) – sum φi with ψi-1 when coloring atoms.
  • normalize (bool) – color and plot magnitude of differences, and normalize the color scale by relative magnitude.
  • hinge (bool) – analyze possible hinges using threshold and gap arguments.
  • threshold (float) – rotation above which a dihedral can be considered part of a hinge. If above 1.0, units are degrees, else it is the fraction of the sum of dihedral differences.
  • gap (int) – residues between hinge dihedrals that can be bridged when consolidating.
  • start ([‘C-Ni’, ...] where C is char chain, N is int residue number and i is an optional insertion designator (char).) – starting residues for manually designated hinges.
  • end (see start) – ending residues to match those of start.
Returns:

(atoms, mean_change, std_change, each mean_combined, std_combined, combined, combined_ids) - atoms for 3D display (if color; differences stored as B-factors); mean_change and std_change - statistics for individual φ and ψ; each - changes in each variable dihedral; doc - identifying string for each angle; mean_combined, std_combined - statistics for pseudo dihedrals; combined - changes in pseudo dihedrals (φi + ψi-1); combined_ids - tuple of “each” indices for (ψi-1, φi) of each pseudo-dihedral; doc_combined - identifying string for each pseudo-dihedral.

Return type:

Atoms or NoneType, float, float, ndarray, list(str), float, float, ndarray, list(tuple), list(str)

copy(pairing, significance=-0.001, pseudo=False, confirm=False)

Set selected dihedrals to those of the target.

This is a very simple algorithm that applies to all dihdreals with differences greater than significance. There is not yet any way to restrict which regions of the molecule will be changed.

Parameters:
  • pairing (Pair) – paired atomic structures, the torsion angles of “a” to acquire the values of “b”.
  • significance (float) – reset dihedrals whose changes exceed this number of standard deviations.
  • pseudo (bool) – Selection on combined φi, ψi-1 pseudo-dihedral, else individual dihedrals.
  • confirm (bool) – Each change must be confirmed interactively.
debug = False

Edit which atoms move as connected groups about each torsion angle.

Modifies self.bond_graph in-place.

Parameters:
  • moving_selection (atoms.Atoms) – Atoms that can be moved within the structure.
  • link ((Group of) 2-atom atoms.Selection) – Additional pair(s) of atoms to be linked for joint rotations about (other) torsion bonds.
  • unlink ((Group of) 2-atom atoms.Selection) – Disconnect pair(s) of atoms from those automatically linked so that not rotated together about torsion bonds.
  • check (bool) – Harsh check that all atoms to be linked / unlinked are within moving_selection.
Returns:

None

generate_rotation_masks(fixOmega=True)

Identify rotatable bonds and atom positions that depend on them.

Determined from the graph of atom connectivity. Updates instance variables bond_atoms and rotation_mask.

Parameters:fixOmega (bool) – skip omega (peptide) dihedrals to accelerate this method. Note that fixOmega=True would block refinement of omega if requested in other program options.
Returns:None.

Todo

could fix Pro phi here (or elsewhere) if desired.

Todo

make properties for bond_atoms, rotation masks and components to minimize recalculation here.

generate_structure(atoms, moving_selection)

Structure atom connectivity.

Calculates instance variables bond_graph, next_residue, prior_main, next_main.

Parameters:atoms (Atoms) – Coordinates.

Todo

needs documenting (Brynmor).

get_all_backbone_atoms()

Identify all backbone atoms.

Returns:atom_selection
Return type:atoms.Selection
get_all_backbone_bonds()

Identify all backbone bonds that are rotatable in principle.

Computes a boolean mask representing all backbone bonds in the structure.

Returns:bond_mask
Return type:ndarray(shape=n_bonds, dtype=bool)
get_all_bonds_crossing_cut(atom_selection)

Mask all bonds with exactly one atom within selection.

Computes a boolean mask representing all bonds exactly one endpoint of which is contained in the given selection

Parameters:atom_selection (atoms.Selection OR atoms.Group) – selection identifying the cut
Returns:bond_mask
Return type:ndarray(shape=n_atoms, dtype=bool)
get_all_bonds_in_selection(atom_selection)

Mask all bonds with both atoms within an atom selection.

Computes a boolean mask representing all bonds, both endpoints of which are contained in the given selection.

Parameters:atom_selection (atoms.Selection OR atoms.Group) – selection in which to find bonds
Returns:bond_mask
Return type:ndarray(shape=n_bonds, dtype=bool)
get_all_omega()

Identify all ω.

Computes boolean mask representing all ω torsion angles.

Returns:bond_mask
Return type:ndarray(shape=n_bonds, dtype=bool)
get_all_phi()

Identify rotatable φ.

Computes a boolean mask representing all φ bonds in the structure (which can be rotated).

Returns:bond_mask
Return type:ndarray(shape=n_bonds, dtype=bool)
get_all_phi_psi()

Identify all rotatable φ and ψ.

Computes boolean mask representing all φ & ψ bonds in the structure (which can be rotated)

Returns:bond_mask
Return type:ndarray(shape=n_bonds, dtype=bool)
get_all_psi()

Identify all rotatable ψ.

Computes a boolean mask representing all ψ bonds in the structure (which can be rotated)

Returns:bond_mask
Return type:ndarray(shape=n_bonds, dtype=bool)
get_edge_set(bond_mask)

Computes the set of graph edges corresponding to the boolean mask.

Parameters:bond_mask (ndarray(shape=n_bonds, dtype=bool)) – boolean mask (e.g. as returned by self.get_all_bonds_in_selection)
Returns:edge_set
Return type:set((int, int))
get_rotation_mask(bondID)

Identify atoms whose positions are dependent on each torsion angle.

Mask specifying atoms rotating as a function of a given torsion angle.

Parameters:bondID (int) – index of bond to rotate (index of the atom in the bond which is closer to the N-terminus).
Returns:rotation_mask
Return type:N.ndarray(shape = (n, ), dtype = bool)
get_selected_phi(selection)

Identify φ bonds within atom selection.

Parameters:selection (Atoms.Selection) – boolean array of atoms.
Returns:bond_mask
Return type:ndarray(shape=n_bonds, dtype=bool)
get_selected_psi(selection)

Identify ψ bonds within atom selection.

Parameters:selection (Atoms.Selection) – boolean array of atoms.
Returns:bond_mask
Return type:ndarray(shape=n_bonds, dtype=bool)
join_all(bond_group)
join_atoms(selection)

Add a connection between 2 atoms to the graph.

Modifies self.bond_graph in-place.

When one atom is the N-terminal N of a chain, TAProtein fails to exclude its phi. Perhaps not critical, because it is listed as a non-variable angle. A kluge in dihedral_name is needed which should be removed when the problem can be corrected.

Parameters:selection (atoms.Selection) – boolean True for exactly 2 atoms to be connected.
phi_psi_selection(selection_criteria='all', verbose=False)

Identify which torsion angles are variable.

Among all backbone dihedrals within the structure (or the sub-structure defined by moving_structure in TAProtein.__init__), identify the φ, ψ wholy within selection_criteria.

Returns:Of all backbone dihedrals (listed N to C), True for φ, ψ within atom selection as True, others False.
Return type:ndarray(dtype=bool, size=3*sub-structure residues)

Note

Excludes first φ of (sub)structure. In a structure without hydrogens, alignment of the transformed stucture on the original means that the first φ has no effect on the structure.

Todo

Michael to test what happens to different chains etc. and change documentation accordingly.

print_bonds(atom_number)
print_connectivity()
quartet(pair)

Add immediate backbone neighbors to a backbone pair.

Parameters:pair (list) – atomic indices of a pair of bonded backbone atoms.
Return quartet:indices of atoms that define the torsion angle about pair.
Rtype quartet:(int*4)
save_selection(name, mask, verbose=True)

Save bond selection.

Parameters:
  • name (str) – key for storage dictionary.
  • mask (list or ndarray) – boolean mask
  • verbose (bool) – warnings when over-writing previous save.
saved_selection(name, verbose=True)

Recall previously saved bond selection.

Parameters:
  • name (str) – key for storage dictionary.
  • verbose (bool) – warning message when name not found.
Returns:

logical array

Return type:

list or ndarray

Raise:

KeyError (in addition to warning) if name not found.

separate_all(bond_group)
separate_atoms(selection)

Remove a connection between 2 atoms to the graph.

Modifies self.bond_graph in-place.

Parameters:selection (atoms.Selection) – boolean True for exactly 2 atoms to be connected.
vs_old_apply = False
warnings = False
torsion.dihedral_atoms(atoms, indices)

Indentify atoms of a dihedral angle.

Parameters:
  • atoms (atoms.Atoms) – coordinate set
  • indices ((int, )*4) – numerical indices of 4 atoms defining dihedral angle.
Returns:

atom identifiers for 4 atoms defining dihedral or None if missing.

Return type:

tuple of str/NoneType

torsion.dihedral_name(atoms, indices)

Name of an indexed dihedral angle.

Parameters:
  • atoms (atoms.Atoms) – coordinate set
  • indices ((int, )*4) – numerical indices of 4 atoms defining dihedral angle.
Returns:

name of backbone dihedral, “phi”, “psi”, “omega” or “unknown”

Return type:

string

Note

todo Delete kluge for N-phi when understand why a link to chain N-terminus is causing 1st phi to be re-instated as a fixed angle.

Warning

At/after carboxy terminus, sometimes called with garbage atoms. Returning “ERROR”, as likely not critical, but should follow up.

torsion.dihedral_value(atoms, indices, degrees=False)

Value of an indexed dihedral.

See also

get_dihedral() for similar functionality.

Module author: Michael S. Chapman

Parameters:
  • atoms (atoms.Atoms) – coordinate set
  • indices ((int, )*4) – numerical indices of 4 atoms defining dihedral angle.
  • degrees (bool) – report units of degrees, else radians
Returns:

float dihedral angle.

torsion.getResIDKeys(atoms, selection=None)

Creates list of tuple residue IDs sorted on chain & residue number.

Sorts the keys first by chain ID number, then by residue number, then by residue insert char

Parameters:
  • atoms (atoms.Atoms) – relevant data from PDB file
  • selection (atoms.Selection or NoneType) – use only selected atoms (None for all).
Returns:

ResIDKeys list of distinct (chain, resnum, insert) tuples in ascending lexicographic order as above.

Return type:

[(str, int, str)]

torsion.get_dihedral(xyz, indices)

Calculate the value of a dihedral (radians).

Module author: Brynmor K. Chapman.

See also

dihedral_value() for similar method

torsion.hinges(residue_label, residue_number, delta, dihedral=None, doc='current structure & target', threshold=3.0, gap=1, atoms=None, start=[], end=[])

Find & document hinges where dihedral differences > threshold.

Parameters:
  • residue_label (ndarray(dtype=str)) – unique identifier w/ chain & residue to label residue.
  • residue_number (int or float) – numeric
  • delta (ndarray(dtype=float)) – difference in each paired dihedral angle (deg.)
  • dihedral (ndarray(dtype=str)) – name of dihedral angle, or None if combined per residue.
  • doc (str) – annotation of what delta is (for printing).
  • threshold (float) – absolute value of diheral (deg) above which is hinge.
  • gap (int) – number of residues that can be bridged to consolidate hinge. Use 0 to consolidate only consecutive residues; Use < 0 to avoid consolidation (even over an individual dihedral).
  • atoms (Atoms or NoneType) – atomic coordinates needed to calculate compound rotations.
  • start ([‘C-Ni’, ...] where C is char chain, N is int residue number and i is an optional insertion designator (char).) – starting residues for manually designated hinges.
  • end (see start) – ending residues to match those of start.
Returns:

hinge identifier for each dihedral (-1 if None).

Return type:

ndarray(dtype=int)

torsion.warn(msg)