Bio.PDB.internal_coords module

Classes to support internal coordinates for protein structures.

Internal coordinates comprise Psi, Omega and Phi dihedral angles along the protein backbone, Chi angles along the sidechains, and all 3-atom angles and bond lengths defining a protein chain. These routines can compute internal coordinates from atom XYZ coordinates, and compute atom XYZ coordinates from internal coordinates.

Secondary benefits include the ability to align and compare residue environments in 3D structures, support for 2D atom distance plots, converting a distance plot plus chirality information to a structure, generating an OpenSCAD description of a structure for 3D printing, and reading/writing structures as internal coordinate data files.

Terms and key data structures: Internal coordinates are defined on sequences of atoms which span residues or follow accepted nomenclature along sidechains. To manage these sequences and support Biopython’s disorder mechanisms, AtomKey specifiers are implemented to capture residue, atom and variant identification in a single object. A Hedron object is specified as three sequential AtomKeys, comprising two bond lengths and the bond angle between them. A Dihedron consists of four sequential AtomKeys, linking two Hedra with a dihedral angle between them.

Algorithmic overview: The Internal Coordinates module combines a specification of connected atoms as hedra and dihedra in the ic_data file with routines here to transform XYZ coordinates of these atom sets between a local coordinate system and the world coordinates supplied in e.g. a PDB or mmCif data file. The local coordinate system places the center atom of a hedron at the origin (0,0,0), one leg on the +Z axis, and the other leg on the XZ plane (see Hedron). Measurement and creation or manipulation of hedra and dihedra in the local coordinate space is straightforward, and the calculated transformation matrices enable assembling these subunits into a protein chain starting from supplied (PDB) coordinates for the initial N-Ca-C atoms.

Psi and Phi angles are defined on atoms from adjacent residues in a protein chain, see e.g. pick_angle() and ic_data for the relevant mapping between residues and backbone dihedral angles.

Transforms to and from the dihedron local coordinate space described above are accessible via IC_Chain.dCoordSpace and Dihedron attributes .cst and .rcst, and may be applied in the alignment and comparison of residues and their environments with code along the lines of:

chi1 = ric0.pick_angle("chi1") # chi1 space defined with CA at origin
cst = np.transpose(chi1.cst) # transform TO chi1 local space
newAtomCoords = oldAtomCoords.dot(cst)

The core algorithms were developed independently during 1993-4 for ‘’Development and Application of a Three-dimensional Description of Amino Acid Environments in Protein,’’ Miller, Douthart, and Dunker, Advances in Molecular Bioinformatics, IOS Press, 1994, ISBN 90 5199 172 x, pp. 9-30.

A Protein Internal Coordinate (.pic) file format is defined to capture sufficient detail to reproduce a PDB file from chain starting coordinates (first residue N, Ca, C XYZ coordinates) and remaining internal coordinates. These files are used internally to verify that a given structure can be regenerated from its internal coordinates. See PICIO for reading and writing .pic files and structure_rebuild_test() to determine if a specific PDB or mmCif datafile has sufficient information to interconvert between cartesian and internal coordinates.

Internal coordinates may also be exported as OpenSCAD data arrays for generating 3D printed protein models. OpenSCAD software is provided as a starting point and proof-of-concept for generating such models. See SCADIO and this Thingiverse project for a more advanced example.

Refer to distance_plot() and distance_to_internal_coordinates() for converting structure data to/from 2D distance plots.

The following classes comprise the core functionality for processing internal coordinates and are sufficiently related and coupled to place them together in this module:

IC_Chain: Extends Biopython Chain on .internal_coord attribute.

Manages connected sequence of residues and chain breaks; holds numpy arrays for all atom coordinates and bond geometries. For ‘parallel’ processing IC_Chain methods operate on these arrays with single numpy commands.

IC_Residue: Extends Biopython Residue on .internal_coord attribute.

Access for per residue views on internal coordinates and methods for serial (residue by residue) assembly.

Dihedron: four joined atoms forming a dihedral angle.

Dihedral angle, homogeneous atom coordinates in local coordinate space, references to relevant Hedra and IC_Residue. Getter methods for residue dihedral angles, bond angles and bond lengths.

Hedron: three joined atoms forming a plane.

Contains homogeneous atom coordinates in local coordinate space as well as bond lengths and angle between them.

Edron: base class for Hedron and Dihedron classes.

Tuple of AtomKeys comprising child, string ID, mainchain membership boolean and other routines common for both Hedra and Dihedra. Implements rich comparison.

AtomKey: keys (dictionary and string) for referencing atom sequences.

Capture residue and disorder/occupancy information, provides a no-whitespace key for .pic files, and implements rich comparison.

Custom exception classes: HedronMatchError and MissingAtomError

class Bio.PDB.internal_coords.IC_Chain(parent: Chain, verbose: bool = False)

Bases: object

Class to extend Biopython Chain with internal coordinate data.

Attributes
chain: object reference

The Biopython Bio.PDB.Chain object this extends

MaxPeptideBond: float

Class attribute to detect chain breaks. Override for fully contiguous chains with some very long bonds - e.g. for 3D printing (OpenSCAD output) a structure with missing residues. MaxPeptideBond

ParallelAssembleResidues: bool

Class attribute affecting internal_to_atom_coords. Short (50 residue and less) chains are faster to assemble without the overhead of creating numpy arrays, and the algorithm is easier to understand and trace processing a single residue at a time. Clearing (set to False) this flag will switch to the serial algorithm

ordered_aa_ic_list: list

IC_Residue objects internal_coords algorithms can process (e.g. no waters)

initNCaC: List of N, Ca, C AtomKey tuples (NCaCKeys).

NCaCKeys start chain segments (first residue or after chain break). These 3 atoms define the coordinate space for a contiguous chain segment, as initially specified by PDB or mmCIF file.

AAsiz = int

AtomArray size, number of atoms in this chain

atomArray: numpy array

homogeneous atom coords ([x,, y, z, 1.0]) for every atom in chain

atomArrayIndex: dict

maps AtomKeys to atomArray indexes

hedra: dict

Hedra forming residues in this chain; indexed by 3-tuples of AtomKeys.

hedraLen: int

length of hedra dict

hedraNdx: dict

maps hedra AtomKeys to numeric index into hedra data arrays e.g. hedraL12 below

a2ha_map: [hedraLen x 3]

atom indexes in hedraNdx order

dihedra: dict

Dihedra forming residues in this chain; indexed by 4-tuples of AtomKeys.

dihedraLen: int

length of dihedra dict

dihedraNdx: dict

maps dihedra AtomKeys to dihedra data arrays e.g. dihedraAngle

a2da_map[dihedraLen x 4]

AtomNdx’s in dihedraNdx order

d2a_map[dihedraLen x [4]]

AtomNdx’s for each dihedron (reshaped a2da_map)

Numpy arrays for vector processing of chain di/hedra:
hedraL12: numpy array

bond length between hedron 1st and 2nd atom

hedraAngle: numpy array

bond angle for each hedron, in degrees

hedraL23: numpy array

bond length between hedron 2nd and 3rd atom

id3_dh_index: dict

maps hedron key to list of dihedra starting with hedron, used by assemble and bond_rotate to find dihedra with h1 key

id32_dh_index: dict

like id3_dh_index, find dihedra from h2 key

hAtoms: numpy array

homogeneous atom coordinates (3x4) of hedra, central atom at origin

hAtomsR: numpy array

hAtoms in reverse orientation

hAtoms_needs_update: numpy array of bool

indicates whether hAtoms represent hedraL12/A/L23

dihedraAngle: numpy array

dihedral angles for each dihedron

dAtoms: numpy array

homogeneous atom coordinates (4x4) of dihedra, second atom at origin

dAtoms_needs_update: numpy array of bool

indicates whether dAtoms represent dihedraAngle

dCoordSpace: numpy array

forward and reverse transform matrices standardising positions of first hedron. See dCoordSpace.

dcs_valid: bool

indicates dCoordSpace up to date

See also attributes generated bymeth:build_edraArrays for indexing
di/hedra data elements.

Methods

internal_to_atom_coordinates:

Process ic data to Residue/Atom coordinates; calls assemble_residues()

assemble_residues:

Generate IC_Chain atom coords from internal coordinates (parallel)

assemble_residues_ser:

Generate IC_Residue atom coords from internal coordinates (serial)

atom_to_internal_coordinates:

Calculate dihedrals, angles, bond lengths (internal coordinates) for Atom data

write_SCAD:

Write OpenSCAD matrices for internal coordinate data comprising chain; this is a support routine, see SCADIO.write_SCAD() to generate OpenSCAD description of a protein chain.

distance_plot:

Generate 2D plot of interatomic distances with optional filter

distance_to_internal_coordinates:

Compute internal coordinates from distance plot and array of dihedral angle signs.

MaxPeptideBond = 1.4

Larger C-N distance than this will be chain break

ParallelAssembleResidues = True

Enable parallel internal_to_atom algorithm, is slower for short chains

AAsiz = 0

Number of atoms in this chain (size of atomArray)

atomArray: numpy.array = None

AAsiz x [4] of float np.float64 homogeneous atom coordinates, all atoms in chain.

dCoordSpace = None

[2][dihedraLen][4][4] : 2 arrays of 4x4 coordinate space transforms for each dihedron. The first [0] converts TO standard space with first atom on the XZ plane, the second atom at the origin, the third on the +Z axis, and the fourth placed according to the dihedral angle. The second [1] transform returns FROM the standard space to world coordinates (PDB file input or whatever is current). Also accessible as .cst (forward transform) and .rcst (reverse transform) in Dihedron.

dcsValid = None

True if dCoordSpace is up to date. Use update_dCoordSpace() if needed.

__init__(parent: Chain, verbose: bool = False) None

Initialize IC_Chain object, with or without residue/Atom data.

Parameters

parent (Bio.PDB.Chain) – Biopython Chain object Chain object this extends

__deepcopy__(memo) Bio.PDB.internal_coords.IC_Chain

Implement deepcopy for IC_Chain.

clear_ic()

Clear residue internal_coord settings for this chain.

build_atomArray() None

Build IC_Chain numpy coordinate array from biopython atoms.

See also init_edra() for more complete initialization of IC_Chain.

Inputs:
self.aksetset

AtomKey s in this chain

Generates:
self.AAsizint

number of atoms in chain (len(akset))

self.aktupleAAsiz x AtomKeys

sorted akset AtomKeys

self.atomArrayIndex[AAsiz] of int

numerical index for each AtomKey in aktuple

self.atomArrayValidAAsiz x bool

atomArray coordinates current with internal coordinates if True

self.atomArrayAAsiz x np.float64[4]

homogeneous atom coordinates; Biopython Atom coordinates are view into this array after execution

rak_cachedict

lookup cache for AtomKeys for each residue

build_edraArrays() None

Build chain level hedra and dihedra arrays.

Used by init_edra() and _hedraDict2chain(). Should be private method but exposed for documentation.

Inputs:
self.dihedraLenint

number of dihedra needed

self.hedraLenint

number of hedra needed

self.AAsizint

length of atomArray

self.hedraNdxdict

maps hedron keys to range(hedraLen)

self.dihedraNdxdict

maps dihedron keys to range(dihedraLen)

self.hedradict

maps Hedra keys to Hedra for chain

self.atomArrayAAsiz x np.float64[4]

homogeneous atom coordinates for chain

self.atomArrayIndexdict

maps AtomKeys to atomArray

self.atomArrayValidAAsiz x bool

indicates coord is up-to-date

Generates:
self.dCoordSpace[2][dihedraLen][4][4]

transforms to/from dihedron coordinate space

self.dcsValiddihedraLen x bool

indicates dCoordSpace is current

self.hAtomshedraLen x 3 x np.float64[4]

atom coordinates in hCoordSpace

self.hAtomsRhedraLen x 3 x np.float64[4]

hAtoms in reverse order (trading space for time)

self.hAtoms_needs_updatehedraLen x bool

indicates hAtoms, hAtoms current

self.a2h_mapAAsiz x [int …]

maps atomArrayIndex to hedraNdx’s with that atom

self.a2ha_map[hedraLen x 3]

AtomNdx’s in hedraNdx order

self.h2aahedraLen x [int …]

maps hedraNdx to atomNdx’s in hedron (reshaped later)

Hedron.ndxint

self.hedraNdx value stored inside Hedron object

self.dRevdihedraLen x bool

dihedron reversed if true

self.dH1ndx, dH2ndx[dihedraLen]

hedraNdx’s for 1st and 2nd hedra

self.h1d_maphedraLen x []

hedraNdx -> [dihedra using hedron]

Dihedron.h1key, h2key[AtomKey …]

hedron keys for dihedron, reversed as needed

Dihedron.hedron1, hedron2Hedron

references inside dihedron to hedra

Dihedron.ndxint

self.dihedraNdx info inside Dihedron object

Dihedron.cst, rcstnp.float64p4][4]

dCoordSpace references inside Dihedron

self.a2da_map[dihedraLen x 4]

AtomNdx’s in dihedraNdx order

self.d2a_map[dihedraLen x [4]]

AtomNdx’s for each dihedron (reshaped a2da_map)

self.dFwdbool

dihedron is not Reversed if True

self.a2d_mapAAsiz x [[dihedraNdx]

[atom ndx 0-3 of atom in dihedron]], maps atom indexes to dihedra and atoms in them

self.dAtoms_needs_updatedihedraLen x bool

atoms in h1, h2 are current if False

assemble_residues(verbose: bool = False) None

Generate atom coords from internal coords (vectorised).

This is the ‘Numpy parallel’ version of assemble_residues_ser().

Starting with dihedra already formed by init_atom_coords(), transform each from dihedron local coordinate space into protein chain coordinate space. Iterate until all dependencies satisfied.

Does not update dCoordSpace as assemble_residues_ser() does. Call update_dCoordSpace() if needed. Faster to do in single operation once all atom coordinates finished.

Parameters

verbose (bool) – default False. Report number of iterations to compute changed dihedra

generates:
self.dSet: AAsiz x dihedraLen x 4

maps atoms in dihedra to atomArray

self.dSetValid[dihedraLen][4] of bool

map of valid atoms into dihedra to detect 3 or 4 atoms valid

Output coordinates written to atomArray. Biopython Bio.PDB.Atom coordinates are a view on this data.

assemble_residues_ser(verbose: bool = False, start: Optional[int] = None, fin: Optional[int] = None) None

Generate IC_Residue atom coords from internal coordinates (serial).

See assemble_residues() for ‘numpy parallel’ version.

Filter positions between start and fin if set, find appropriate start coordinates for each residue and pass to assemble()

Parameters
  • verbose (bool) – default False. Describe runtime problems

  • start,fin (int) – default None. Sequence position for begin, end of subregion to generate coords for.

init_edra(verbose: bool = False) None

Create chain and residue di/hedra structures, arrays, atomArray.

Inputs:

self.ordered_aa_ic_list : list of IC_Residue

Generates:
  • edra objects, self.di/hedra (executes _create_edra())

  • atomArray and support (executes build_atomArray())

  • self.hedraLen : number of hedra in structure

  • hedraL12 : numpy arrays for lengths, angles (empty)

  • hedraAngle ..

  • hedraL23 ..

  • self.hedraNdx : dict mapping hedrakeys to hedraL12 etc

  • self.dihedraLen : number of dihedra in structure

  • dihedraAngle ..

  • dihedraAngleRads : np arrays for angles (empty)

  • self.dihedraNdx : dict mapping dihedrakeys to dihedraAngle

init_atom_coords() None

Set chain level di/hedra initial coords from angles and distances.

Initializes atom coordinates in local coordinate space for hedra and dihedra, will be transformed appropriately later by dCoordSpace matrices for assembly.

update_dCoordSpace(workSelector: Optional[numpy.ndarray] = None) None

Compute/update coordinate space transforms for chain dihedra.

Requires all atoms updated so calls assemble_residues() (returns immediately if all atoms already assembled).

Parameters

workSelector ([bool]) – Optional mask to select dihedra for update

propagate_changes() None

Track through di/hedra to invalidate dependent atoms.

internal_to_atom_coordinates(verbose: bool = False, start: Optional[int] = None, fin: Optional[int] = None) None

Process IC data to Residue/Atom coords.

Parameters
  • verbose (bool) – default False. Describe runtime problems

  • start,fin (int) – Optional sequence positions for begin, end of subregion to process.

Note

Setting start or fin activates serial assemble_residues_ser() instead of (Numpy parallel) assemble_residues(). Start C-alpha will be at origin.

atom_to_internal_coordinates(verbose: bool = False) None

Calculate dihedrals, angles, bond lengths for Atom data.

Generates atomArray (through init_edra), value arrays for hedra and dihedra, and coordinate space transforms for dihedra.

Generates Gly C-beta if specified, see IC_Residue.gly_Cbeta

Parameters

verbose (bool) – default False. describe runtime problems

distance_plot(filter: Optional[numpy.ndarray] = None) numpy.ndarray

Generate 2D distance plot from atomArray.

Default is to calculate distances for all atoms. To generate the classic C-alpha distance plot, pass a boolean mask array like:

atmNameNdx = internal_coords.AtomKey.fields.atm
CaSelect = [
    atomArrayIndex.get(k)
    for k in atomArrayIndex.keys()
    if k.akl[atmNameNdx] == "CA"
]
plot = distance_plot(CaSelect)

Alternatively, this will select all backbone atoms:

backboneSelect = [
    atomArrayIndex.get(k)
    for k in atomArrayIndex.keys()
    if k.is_backbone()
]
Parameters

filter ([bool]) – restrict atoms for calculation

See also

distance_to_internal_coordinates(), which requires the default all atom distance plot.

dihedral_signs() numpy.ndarray

Get sign array (+1/-1) for each element of chain dihedraAngle array.

Required for distance_to_internal_coordinates()

distplot_to_dh_arrays(distplot: numpy.ndarray) None

Load di/hedra distance arays from distplot.

Fill IC_Chain arrays hedraL12, L23, L13 and dihedraL14 distance value arrays from input distplot. Distplot and di/hedra distance arrays must index according to AtomKey mappings in IC_Chain .hedraNdx and .dihedraNdx (created in IC_Chain.init_edra())

Call atom_to_internal_coordinates() (or at least init_edra()) to generate a2ha_map and d2a_map before running this.

Explcitly removed from distance_to_internal_coordinates() so user may populate these chain di/hedra distance arrays by other methods.

distance_to_internal_coordinates(dihedra_signs: numpy.ndarray) None

Compute chain di/hedra from from distance and chirality data.

Distance properties on hedra L12, L23, L13 and dihedra L14 configured by distplot_to_dh_arrays() or alternative loader.

dihedraAngles result multiplied by dihedra_signs at final step to recover chirality information lost in distance plot (mirror image of structure has same distances but opposite sign dihedral angles).

Note that chain breaks will cause errors in rebuilt structure, use copy_initNCaCs() to avoid this

Based on Blue’s answer to The dihedral angles of a tetrahedron in terms of its edge lengths on math.stackexchange.com.

Other values from that analysis included here as comments for completeness:

  • oa = hedron1 L12 if reverse else hedron1 L23

  • ob = hedron1 L23 if reverse else hedron1 L12

  • ac = hedron2 L12 if reverse else hedron2 L23

  • ab = hedron1 L13 = law of cosines on OA, OB (hedron1 L12, L23)

  • oc = hedron2 L13 = law of cosines on OA, AC (hedron2 L12, L23)

  • bc = dihedron L14

target is OA, the dihedral angle along edge oa

copy_initNCaCs(other: Bio.PDB.internal_coords.IC_Chain) None

Copy atom coordinates for initNCaC atoms from other IC_Chain.

Needed for distance_to_internal_coordinates() if target has chain breaks (otherwise each fragment will start at origin)

__annotations__ = {'atomArray': <built-in function array>, 'atomArrayIndex': "Dict['AtomKey', int]", 'bpAtomArray': "List['Atom']", 'ordered_aa_ic_list': 'List[IC_Residue]'}
class Bio.PDB.internal_coords.IC_Residue(parent: Residue)

Bases: object

Class to extend Biopython Residue with internal coordinate data.

Parameters
parent: biopython Residue object this class extends
Attributes
no_altloc: bool default False

Class variable, disable processing of ALTLOC atoms if True, use only selected atoms.

accept_atoms: tuple

Class variable accept_atoms, list of PDB atom names to use when generating internal coordinates. Default is:

accept_atoms = accept_mainchain + accept_hydrogens

to exclude hydrogens in internal coordinates and generated PDB files, override as:

IC_Residue.accept_atoms = IC_Residue.accept_mainchain

to get only mainchain atoms plus amide proton, use:

IC_Residue.accept_atoms = IC_Residue.accept_mainchain + ('H',)

to convert D atoms to H, set AtomKey.d2h = True and use:

IC_Residue.accept_atoms = (
    accept_mainchain + accept_hydrogens + accept_deuteriums
)

Note that accept_mainchain = accept_backbone + accept_sidechain. Thus to generate sequence-agnostic conformational data for e.g. structure alignment in dihedral angle space, use:

IC_Residue.accept_atoms = accept_backbone

or set gly_Cbeta = True and use:

IC_Residue.accept_atoms = accept_backbone + ('CB',)

There is currently no option to output internal coordinates with D instead of H.

accept_resnames: tuple

Class variable accept_resnames, list of 3-letter residue names for HETATMs to accept when generating internal coordinates from atoms. HETATM sidechain will be ignored, but normal backbone atoms (N, CA, C, O, CB) will be included. Currently only CYG, YCM and UNK; override at your own risk. To generate sidechain, add appropriate entries to ic_data_sidechains in ic_data and support in IC_Chain.atom_to_internal_coordinates().

gly_Cbeta: bool default False

Class variable !data:gly_Cbeta, override to True to generate internal coordinates for glycine CB atoms in IC_Chain.atom_to_internal_coordinates()

IC_Residue.gly_Cbeta = True
residue: Biopython Residue object reference

The Residue object this extends

hedra: dict indexed by 3-tuples of AtomKeys

Hedra forming this residue

dihedra: dict indexed by 4-tuples of AtomKeys

Dihedra forming (overlapping) this residue

rprev, rnext: lists of IC_Residue objects

References to adjacent (bonded, not missing, possibly disordered) residues in chain

atom_coords: AtomKey indexed dict of numpy [4] arrays

View into IC_Chain’s atomArray of homogeneous coordinates

alt_ids: list of char

AltLoc IDs from PDB file

bfactors: dict

AtomKey indexed B-factors as read from PDB file

NCaCKey: List of tuples of AtomKeys

List of tuples of N, Ca, C backbone atom AtomKeys; usually only 1 but more if backbone altlocs.

is20AA: bool

True if residue is one of 20 standard amino acids, based on Residue resname

isAccept: bool

True if is20AA or in accept_resnames below

rbase: tuple

residue position, insert code or none, resname (1 letter if standard amino acid)

cic: IC_Chain default None

parent chain IC_Chain object

scale: optional float

used for OpenSCAD output to generate gly_Cbeta bond length

Methods

assemble(atomCoordsIn, resetLocation, verbose)

Compute atom coordinates for this residue from internal coordinates

get_angle()

Return angle for passed key

get_length()

Return bond length for specified pair

pick_angle()

Find Hedron or Dihedron for passed key

pick_length()

Find hedra for passed AtomKey pair

set_angle()

Set angle for passed key (no position updates)

set_length()

Set bond length in all relevant hedra for specified pair

bond_rotate(delta)

adjusts related dihedra angles by delta, e.g. rotating psi (N-Ca-C-N) will adjust the adjacent N-Ca-C-O by the same amount to avoid clashes

bond_set(angle)

uses bond_rotate to set specified dihedral to angle and adjust related dihedra accordingly

rak(atom info)

cached AtomKeys for this residue

accept_resnames = ('CYG', 'YCM', 'UNK')

Add 3-letter residue name here for non-standard residues with normal backbone. CYG included for test case 4LGY (1305 residue contiguous chain). Safe to add more names for N-CA-C-O backbones, any more complexity will need additions to accept_atoms, ic_data_sidechains in ic_data and support in IC_Chain.atom_to_internal_coordinates()

no_altloc: bool = False

Set True to filter altloc atoms on input and only work with Biopython default Atoms

gly_Cbeta: bool = False

Create beta carbons on all Gly residues.

Setting this to True will generate internal coordinates for Gly C-beta carbons in atom_to_internal_coordinates().

Data averaged from Sep 2019 Dunbrack cullpdb_pc20_res2.2_R1.0 restricted to structures with amide protons. Please see

PISCES: A Protein Sequence Culling Server

‘G. Wang and R. L. Dunbrack, Jr. PISCES: a protein sequence culling server. Bioinformatics, 19:1589-1591, 2003.’

Ala avg rotation of OCCACB from NCACO query:

select avg(g.rslt) as avg_rslt, stddev(g.rslt) as sd_rslt, count(*)
from
(select f.d1d, f.d2d,
(case when f.rslt > 0 then f.rslt-360.0 else f.rslt end) as rslt
from (select d1.angle as d1d, d2.angle as d2d,
(d2.angle - d1.angle) as rslt from dihedron d1,
dihedron d2 where d1.rdh_class='AOACACAACB' and
d2.rdh_class='ANACAACAO' and d1.pdb=d2.pdb and d1.chn=d2.chn
and d1.res=d2.res) as f) as g

results:

| avg_rslt          | sd_rslt          | count   |
| -122.682194862932 | 5.04403040513919 | 14098   |
accept_backbone = ('N', 'CA', 'C', 'O', 'OXT')
accept_sidechain = ('CB', 'CG', 'CG1', 'OG1', 'OG', 'SG', 'CG2', 'CD', 'CD1', 'SD', 'OD1', 'ND1', 'CD2', 'ND2', 'CE', 'CE1', 'NE', 'OE1', 'NE1', 'CE2', 'OE2', 'NE2', 'CE3', 'CZ', 'NZ', 'CZ2', 'CZ3', 'OD2', 'OH', 'CH2', 'NH1', 'NH2')
accept_mainchain = ('N', 'CA', 'C', 'O', 'OXT', 'CB', 'CG', 'CG1', 'OG1', 'OG', 'SG', 'CG2', 'CD', 'CD1', 'SD', 'OD1', 'ND1', 'CD2', 'ND2', 'CE', 'CE1', 'NE', 'OE1', 'NE1', 'CE2', 'OE2', 'NE2', 'CE3', 'CZ', 'NZ', 'CZ2', 'CZ3', 'OD2', 'OH', 'CH2', 'NH1', 'NH2')
accept_hydrogens = ('H', 'H1', 'H2', 'H3', 'HA', 'HA2', 'HA3', 'HB', 'HB1', 'HB2', 'HB3', 'HG2', 'HG3', 'HD2', 'HD3', 'HE2', 'HE3', 'HZ1', 'HZ2', 'HZ3', 'HG11', 'HG12', 'HG13', 'HG21', 'HG22', 'HG23', 'HZ', 'HD1', 'HE1', 'HD11', 'HD12', 'HD13', 'HG', 'HG1', 'HD21', 'HD22', 'HD23', 'NH1', 'NH2', 'HE', 'HH11', 'HH12', 'HH21', 'HH22', 'HE21', 'HE22', 'HE2', 'HH', 'HH2')
accept_deuteriums = ('D', 'D1', 'D2', 'D3', 'DA', 'DA2', 'DA3', 'DB', 'DB1', 'DB2', 'DB3', 'DG2', 'DG3', 'DD2', 'DD3', 'DE2', 'DE3', 'DZ1', 'DZ2', 'DZ3', 'DG11', 'DG12', 'DG13', 'DG21', 'DG22', 'DG23', 'DZ', 'DD1', 'DE1', 'DD11', 'DD12', 'DD13', 'DG', 'DG1', 'DD21', 'DD22', 'DD23', 'ND1', 'ND2', 'DE', 'DH11', 'DH12', 'DH21', 'DH22', 'DE21', 'DE22', 'DE2', 'DH', 'DH2')
accept_atoms = ('N', 'CA', 'C', 'O', 'OXT', 'CB', 'CG', 'CG1', 'OG1', 'OG', 'SG', 'CG2', 'CD', 'CD1', 'SD', 'OD1', 'ND1', 'CD2', 'ND2', 'CE', 'CE1', 'NE', 'OE1', 'NE1', 'CE2', 'OE2', 'NE2', 'CE3', 'CZ', 'NZ', 'CZ2', 'CZ3', 'OD2', 'OH', 'CH2', 'NH1', 'NH2', 'H', 'H1', 'H2', 'H3', 'HA', 'HA2', 'HA3', 'HB', 'HB1', 'HB2', 'HB3', 'HG2', 'HG3', 'HD2', 'HD3', 'HE2', 'HE3', 'HZ1', 'HZ2', 'HZ3', 'HG11', 'HG12', 'HG13', 'HG21', 'HG22', 'HG23', 'HZ', 'HD1', 'HE1', 'HD11', 'HD12', 'HD13', 'HG', 'HG1', 'HD21', 'HD22', 'HD23', 'NH1', 'NH2', 'HE', 'HH11', 'HH12', 'HH21', 'HH22', 'HE21', 'HE22', 'HE2', 'HH', 'HH2')

Change accept_atoms to restrict atoms processed. See IC_Residue for usage.

__init__(parent: Residue) None

Initialize IC_Residue with parent Biopython Residue.

Parameters

parent (Residue) – Biopython Residue object. The Biopython Residue this object extends

__deepcopy__(memo)

Deep copy implementation for IC_Residue.

rak(atm: Union[str, Bio.PDB.Atom.Atom]) Bio.PDB.internal_coords.AtomKey

Cache calls to AtomKey for this residue.

__repr__() str

Print string is parent Residue ID.

pretty_str() str

Nice string for residue ID.

set_flexible() None

For OpenSCAD, mark N-CA and CA-C bonds to be flexible joints.

See SCADIO.write_SCAD()

set_hbond() None

For OpenSCAD, mark H-N and C-O bonds to be hbonds (magnets).

See SCADIO.write_SCAD()

clear_transforms()

Invalidate dihedra coordinate space attributes before assemble().

Coordinate space attributes are Dihedron.cst and .rcst, and IC_Chain.dCoordSpace

assemble(resetLocation: bool = False, verbose: bool = False) Optional[Union[Dict[Bio.PDB.internal_coords.AtomKey, numpy.array], Dict[Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], numpy.array]]]

Compute atom coordinates for this residue from internal coordinates.

This is the IC_Residue part of the assemble_residues_ser() serial version, see assemble_residues() for numpy vectorized approach which works at the IC_Chain level.

Join prepared dihedra starting from N-CA-C and N-CA-CB hedrons, computing protein space coordinates for backbone and sidechain atoms

Sets forward and reverse transforms on each Dihedron to convert from protein coordinates to dihedron space coordinates for first three atoms (see IC_Chain.dCoordSpace)

Call init_atom_coords() to update any modified di/hedra before coming here, this only assembles dihedra into protein coordinate space.

Algorithm

Form double-ended queue, start with c-ca-n, o-c-ca, n-ca-cb, n-ca-c.

if resetLocation=True, use initial coords from generating dihedron for n-ca-c initial positions (result in dihedron coordinate space)

while queue not empty

get 3-atom hedron key

for each dihedron starting with hedron key (1st hedron of dihedron)

if have coordinates for all 4 atoms already

add 2nd hedron key to back of queue

else if have coordinates for 1st 3 atoms

compute forward and reverse transforms to take 1st 3 atoms to/from dihedron initial coordinate space

use reverse transform to get position of 4th atom in current coordinates from dihedron initial coordinates

add 2nd hedron key to back of queue

else

ordering failed, put hedron key at back of queue and hope next time we have 1st 3 atom positions (should not happen)

loop terminates (queue drains) as hedron keys which do not start any dihedra are removed without action

Parameters

resetLocation (bool) – default False. - Option to ignore start location and orient so initial N-Ca-C hedron at origin.

Returns

Dict of AtomKey -> homogeneous atom coords for residue in protein space relative to previous residue

Also directly updates IC_Chain.atomArray as assemble_residues() does.

split_akl(lst: Union[Tuple[Bio.PDB.internal_coords.AtomKey, ...], List[Bio.PDB.internal_coords.AtomKey]], missingOK: bool = False) List[Tuple[Bio.PDB.internal_coords.AtomKey, ...]]

Get AtomKeys for this residue (ak_set) for generic list of AtomKeys.

Changes and/or expands a list of ‘generic’ AtomKeys (e.g. ‘N, C, C’) to be specific to this Residue’s altlocs etc., e.g. ‘(N-Ca_A_0.3-C, N-Ca_B_0.7-C)’

Given a list of AtomKeys (aks) for a Hedron or Dihedron,
return:

list of matching aks that have id3_dh in this residue (ak may change if occupancy != 1.00)

or

multiple lists of matching aks expanded for all atom altlocs

or

empty list if any of atom_coord(ak) missing and not missingOK

Parameters
  • lst (list) – list[3] or [4] of AtomKeys. Non-altloc AtomKeys to match to specific AtomKeys for this residue

  • missingOK (bool) – default False, see above.

atom_sernum = None
atom_chain = None
pdb_residue_string() str

Generate PDB ATOM records for this residue as string.

Convenience method for functionality not exposed in PDBIO.py. Increments IC_Residue.atom_sernum if not None

Parameters
  • IC_Residue.atom_sernum – Class variable default None. Override and increment atom serial number if not None

  • IC_Residue.atom_chain – Class variable. Override atom chain id if not None

Todo

move to PDBIO

pic_flags = _pfDef(psi=1, omg=2, phi=4, tau=8, chi1=16, chi2=32, chi3=64, chi4=128, chi5=256, pomg=512, chi=496, classic_b=525, classic=1021, hedra=1024, primary=2048, secondary=4096, all=7168, initAtoms=8192, bFactors=16384)

Used by PICIO.write_PIC() to control classes of values to be defaulted.

picFlagsDefault = 31744

Default is all dihedra + initial tau atoms + bFactors.

picFlagsDict = {'all': 7168, 'bFactors': 16384, 'chi': 496, 'chi1': 16, 'chi2': 32, 'chi3': 64, 'chi4': 128, 'chi5': 256, 'classic': 1021, 'classic_b': 525, 'hedra': 1024, 'initAtoms': 8192, 'omg': 2, 'phi': 4, 'pomg': 512, 'primary': 2048, 'psi': 1, 'secondary': 4096, 'tau': 8}

Dictionary of pic_flags values to use as needed.

pick_angle(angle_key: Union[Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], str]) Optional[Union[Bio.PDB.internal_coords.Hedron, Bio.PDB.internal_coords.Dihedron]]

Get Hedron or Dihedron for angle_key.

Parameters

angle_key

  • tuple of 3 or 4 AtomKeys

  • string of atom names (‘CA’) separated by :’s

  • string of [-1, 0, 1]<atom name> separated by ‘:’s. -1 is previous residue, 0 is this residue, 1 is next residue

  • psi, phi, omg, omega, chi1, chi2, chi3, chi4, chi5

  • tau (N-CA-C angle) see Richardson1981

  • tuples of AtomKeys is only access for alternate disordered atoms

Observe that a residue’s phi and omega dihedrals, as well as the hedra comprising them (including the N:Ca:C tau hedron), are stored in the n-1 di/hedra sets; this overlap is handled here, but may be an issue if accessing directly.

The following print commands are equivalent (except for sidechains with non-carbon atoms for chi2):

ric = r.internal_coord
print(
    r,
    ric.get_angle("psi"),
    ric.get_angle("phi"),
    ric.get_angle("omg"),
    ric.get_angle("tau"),
    ric.get_angle("chi2"),
)
print(
    r,
    ric.get_angle("N:CA:C:1N"),
    ric.get_angle("-1C:N:CA:C"),
    ric.get_angle("-1CA:-1C:N:CA"),
    ric.get_angle("N:CA:C"),
    ric.get_angle("CA:CB:CG:CD"),
)

See ic_data.py for detail of atoms in the enumerated sidechain angles and the backbone angles which do not span the peptide bond. Using ‘s’ for current residue (‘self’) and ‘n’ for next residue, the spanning (overlapping) angles are:

(sN, sCA, sC, nN)   # psi
(sCA, sC, nN, nCA)  # omega i+1
(sC, nN, nCA, nC)   # phi i+1
(sCA, sC, nN)
(sC, nN, nCA)
(nN, nCA, nC)       # tau i+1
Returns

Matching Hedron, Dihedron, or None.

get_angle(angle_key: Union[Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], str]) Optional[float]

Get dihedron or hedron angle for specified key.

See pick_angle() for key specifications.

set_angle(angle_key: Union[Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], str], v: float)

Set dihedron or hedron angle for specified key.

See pick_angle() for key specifications.

bond_rotate(angle_key: Union[Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], str], delta: float)

Rotate set of overlapping dihedrals by delta degrees.

See pick_angle() for key specifications.

bond_set(angle_key: Union[Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], str], val: float)

Set dihedron to val, update overlapping dihedra by same amount.

See pick_angle() for key specifications.

pick_length(ak_spec: Union[str, Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey]]) Tuple[Optional[List[Bio.PDB.internal_coords.Hedron]], Optional[Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey]]]

Get list of hedra containing specified atom pair.

Parameters

ak_spec

  • tuple of two AtomKeys

  • string: two atom names separated by ‘:’, e.g. ‘N:CA’ with optional position specifier relative to self, e.g. ‘-1C:N’ for preceding peptide bond. Position specifiers are -1, 0, 1.

The following are equivalent:

ric = r.internal_coord
print(
    r,
    ric.get_length("0C:1N"),
)
print(
    r,
    None
    if not ric.rnext
    else ric.get_length((ric.rak("C"), ric.rnext[0].rak("N"))),
)
Returns

list of hedra containing specified atom pair as tuples of AtomKeys

get_length(ak_spec: Union[str, Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey]]) Optional[float]

Get bond length for specified atom pair.

See pick_length() for ak_spec.

set_length(ak_spec: Union[str, Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey]], val: float) None

Set bond length for specified atom pair.

See pick_length() for ak_spec.

applyMtx(mtx: numpy.array) None

Apply matrix to atom_coords for this IC_Residue.

__annotations__ = {'_AllBonds': <class 'bool'>, 'ak_set': 'Set[AtomKey]', 'akc': 'Dict[Union[str, Atom], AtomKey]', 'alt_ids': 'Union[List[str], None]', 'bfactors': 'Dict[str, float]', 'cic': 'IC_Chain', 'dihedra': 'Dict[DKT, Dihedron]', 'gly_Cbeta': <class 'bool'>, 'hedra': 'Dict[HKT, Hedron]', 'no_altloc': <class 'bool'>, 'rnext': 'List[IC_Residue]', 'rprev': 'List[IC_Residue]'}
class Bio.PDB.internal_coords.Edron(*args: Union[List[Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey]], **kwargs: str)

Bases: object

Base class for Hedron and Dihedron classes.

Supports rich comparison based on lists of AtomKeys.

Attributes
aks: tuple

3 (hedron) or 4 (dihedron) AtomKey s defining this di/hedron

id: str

‘:’-joined string of AtomKeys for this di/hedron

needs_update: bool

indicates di/hedron local atom_coords do NOT reflect current di/hedron angle and length values in hedron local coordinate space

dh_class: str

sequence of atoms (no position or residue) comprising di/hedron for statistics

rdh_class: str

sequence of residue, atoms comprising di/hedron for statistics

crdh_class: tuple

tuple of covalent radii classses comprising di/hedron for statistics

edron_re: compiled regex (Class Attribute)

A compiled regular expression matching string IDs for Hedron and Dihedron objects

cic: IC_Chain reference

Chain internal coords object containing this hedron

ndx: int

index into IC_Chain level numpy data arrays for di/hedra. Set in IC_Chain.init_edra()

rc: int

number of residues involved in this edron

Methods

gen_key([AtomKey, …] or AtomKey, …) (Static Method)

generate a ‘:’-joined string of AtomKey Ids

is_backbone()

Return True if all aks atoms are N, Ca, C or O

edron_re = re.compile('^(?P<pdbid>\\w+)?\\s(?P<chn>[\\w|\\s])?\\s(?P<a1>[\\w\\-\\.]+):(?P<a2>[\\w\\-\\.]+):(?P<a3>[\\w\\-\\.]+)(:(?P<a4>[\\w\\-\\.]+))?\\s+(((?P<len12>\\S+)\\s+(?P<angle>\\S+)\\s+(?P<len23>\\S+)\\s*$)|((?P<)

A compiled regular expression matching string IDs for Hedron and Dihedron objects

static gen_key(lst: List[Bio.PDB.internal_coords.AtomKey]) str

Generate string of ‘:’-joined AtomKey strings from input.

Parameters

lst (list) – list of AtomKey objects

__init__(*args: Union[List[Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey]], **kwargs: str) None

Initialize Edron with sequence of AtomKeys.

Acceptable input:

[ AtomKey, … ] : list of AtomKeys AtomKey, … : sequence of AtomKeys as args {‘a1’: str, ‘a2’: str, … } : dict of AtomKeys as ‘a1’, ‘a2’ …

__deepcopy__(memo)

Deep copy implementation for Edron.

is_backbone() bool

Report True for contains only N, C, CA, O, H atoms.

__repr__() str

Tuple of AtomKeys is default repr string.

__hash__() int

Hash calculated at init from aks tuple.

__eq__(other: object) bool

Test for equality.

__ne__(other: object) bool

Test for inequality.

__gt__(other: object) bool

Test greater than.

__ge__(other: object) bool

Test greater or equal.

__lt__(other: object) bool

Test less than.

__le__(other: object) bool

Test less or equal.

class Bio.PDB.internal_coords.Hedron(*args: Union[List[Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey]], **kwargs: str)

Bases: Bio.PDB.internal_coords.Edron

Class to represent three joined atoms forming a plane.

Contains atom coordinates in local coordinate space: central atom at origin, one terminal atom on XZ plane, and the other on the +Z axis. Stored in two orientations, with the 3rd (forward) or first (reversed) atom on the +Z axis. See Dihedron for use of forward and reverse orientations.

Attributes
len12: float

distance between first and second atoms

len23: float

distance between second and third atoms

angle: float

angle (degrees) formed by three atoms in hedron

xrh_class: string

only for hedron spanning 2 residues, will have ‘X’ for residue contributing only one atom

Methods

get_length()

get bond length for specified atom pair

set_length()

set bond length for specified atom pair

angle(), len12(), len23()

setters for relevant attributes (angle in degrees)

__init__(*args: Union[List[Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey]], **kwargs: str) None

Initialize Hedron with sequence of AtomKeys, kwargs.

Acceptable input:

As for Edron, plus optional ‘len12’, ‘angle’, ‘len23’ keyworded values.

__repr__() str

Print string for Hedron object.

property angle: float

Get this hedron angle.

property len12

Get first length for Hedron.

property len23: float

Get second length for Hedron.

get_length(ak_tpl: Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey]) Optional[float]

Get bond length for specified atom pair.

Parameters

ak_tpl (tuple) – tuple of AtomKeys. Pair of atoms in this Hedron

set_length(ak_tpl: Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey], newLength: float)

Set bond length for specified atom pair; sets needs_update.

Parameters

ak_tpl (tuple) – tuple of AtomKeys Pair of atoms in this Hedron

cic: Bio.PDB.internal_coords.IC_Chain
class Bio.PDB.internal_coords.Dihedron(*args: Union[List[Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey]], **kwargs: str)

Bases: Bio.PDB.internal_coords.Edron

Class to represent four joined atoms forming a dihedral angle.

Attributes
angle: float

Measurement or specification of dihedral angle in degrees; prefer IC_Residue.bond_set() to set

hedron1, hedron2: Hedron object references

The two hedra which form the dihedral angle

h1key, h2key: tuples of AtomKeys

Hash keys for hedron1 and hedron2

id3,id32: tuples of AtomKeys

First 3 and second 3 atoms comprising dihedron; hxkey orders may differ

ric: IC_Residue object reference

IC_Residue object containing this dihedral

reverse: bool

Indicates order of atoms in dihedron is reversed from order of atoms in hedra

primary: bool

True if this is psi, phi, omega or a sidechain chi angle

pclass: string (primary angle class)

rdh_class with X for adjacent residue according to nomenclature (psi, omega, phi)

cst, rcst: numpy [4][4] arrays

transformations to (cst) and from (rcst) Dihedron coordinate space defined with atom 2 (Hedron 1 center atom) at the origin. Views on IC_Chain.dCoordSpace.

Methods

angle()

getter/setter for dihdral angle in degrees; prefer IC_Residue.bond_set()

bits()

return IC_Residue.pic_flags bitmask for dihedron psi, omega, etc

__init__(*args: Union[List[Bio.PDB.internal_coords.AtomKey], Tuple[Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey, Bio.PDB.internal_coords.AtomKey]], **kwargs: str) None

Init Dihedron with sequence of AtomKeys and optional dihedral angle.

Acceptable input:

As for Edron, plus optional ‘dihedral’ keyworded angle value.

__repr__() str

Print string for Dihedron object.

property angle: float

Get dihedral angle.

static angle_dif(a1: Union[float, numpy.ndarray], a2: Union[float, numpy.ndarray])

Get angle difference between two +/- 180 angles.

https://stackoverflow.com/a/36001014/2783487

static angle_avg(alst: List, in_rads: bool = False, out_rads: bool = False)

Get average of list of +/-180 angles.

Parameters
  • alst (List) – list of angles to average

  • in_rads (bool) – input values are in radians

  • out_rads (bool) – report result in radians

static angle_pop_sd(alst: List, avg: float)

Get population standard deviation for list of +/-180 angles.

should be sample std dev but avoid len(alst)=1 -> div by 0

difference(other: Bio.PDB.internal_coords.Dihedron) float

Get angle difference between this and other +/- 180 angles.

bits() int

Get IC_Residue.pic_flags bitmasks for self is psi, omg, phi, pomg, chiX.

cic: Bio.PDB.internal_coords.IC_Chain
class Bio.PDB.internal_coords.AtomKey(*args: Union[Bio.PDB.internal_coords.IC_Residue, Bio.PDB.Atom.Atom, List, Dict, str], **kwargs: str)

Bases: object

Class for dict keys to reference atom coordinates.

AtomKeys capture residue and disorder information together, and provide a no-whitespace string key for .pic files.

Supports rich comparison and multiple ways to instantiate.

AtomKeys contain:

residue position (respos), insertion code (icode), 1 or 3 char residue name (resname), atom name (atm), altloc (altloc), and occupancy (occ)

Use AtomKey.fields to get the index to the component of interest by name:

Get C-alpha atoms from IC_Chain atomArray and atomArrayIndex with AtomKeys:

atmNameNdx = internal_coords.AtomKey.fields.atm
CaSelection = [
    atomArrayIndex.get(k)
    for k in atomArrayIndex.keys()
    if k.akl[atmNameNdx] == "CA"
]
AtomArrayCa = atomArray[CaSelection]

Get all phenylalanine atoms in a chain:

resNameNdx = internal_coords.AtomKey.fields.resname
PheSelection = [
    atomArrayIndex.get(k)
    for k in atomArrayIndex.keys()
    if k.akl[resNameNdx] == "F"
]
AtomArrayPhe = atomArray[PheSelection]

‘resname’ will be the uppercase 1-letter amino acid code if one of the 20 standard residues, otherwise the supplied 3-letter code. Supplied as input or read from .rbase attribute of IC_Residue.

Attributes
akl: tuple

All six fields of AtomKey

fieldNames: tuple (Class Attribute)

Mapping of key index positions to names

fields: namedtuple (Class Attribute)

Mapping of field names to index positions.

id: str

‘_’-joined AtomKey fields, excluding ‘None’ fields

atom_re: compiled regex (Class Attribute)

A compiled regular expression matching the string form of the key

d2h: bool (Class Attribute) default False

Convert D atoms to H on input if True; must also modify IC_Residue.accept_atoms

missing: bool default False

AtomKey __init__’d from string is probably missing, set this flag to note the issue. Set by IC_Residue.rak()

ric: IC_Residue default None

If initialised with IC_Residue, this references the IC_residue

Methods

altloc_match(other)

Returns True if this AtomKey matches other AtomKey excluding altloc and occupancy fields

is_backbone()

Returns True if atom is N, CA, C, O or H

atm()

Returns atom name, e.g. N, CA, CB, etc.

atom_re = re.compile('^(?P<respos>-?\\d+)(?P<icode>[A-Za-z])?_(?P<resname>[a-zA-Z]+)_(?P<atm>[A-Za-z0-9]+)(?:_(?P<altloc>\\w))?(?:_(?P<occ>-?\\d\\.\\d?\\d?))?$')

Pre-compiled regular expression to match an AtomKey string.

fieldNames = ('respos', 'icode', 'resname', 'atm', 'altloc', 'occ')
fields = _fieldsDef(respos=0, icode=1, resname=2, atm=3, altloc=4, occ=5)

Use this namedtuple to access AtomKey fields. See AtomKey

d2h = False

Set True to convert D Deuterium to H Hydrogen on input.

__init__(*args: Union[Bio.PDB.internal_coords.IC_Residue, Bio.PDB.Atom.Atom, List, Dict, str], **kwargs: str) None

Initialize AtomKey with residue and atom data.

Examples of acceptable input:

(<IC_Residue>, 'CA', ...)    : IC_Residue with atom info
(<IC_Residue>, <Atom>)       : IC_Residue with Biopython Atom
([52, None, 'G', 'CA', ...])  : list of ordered data fields
(52, None, 'G', 'CA', ...)    : multiple ordered arguments
({respos: 52, icode: None, atm: 'CA', ...}) : dict with fieldNames
(respos: 52, icode: None, atm: 'CA', ...) : kwargs with fieldNames
52_G_CA, 52B_G_CA, 52_G_CA_0.33, 52_G_CA_B_0.33  : id strings
__deepcopy__(memo)

Deep copy implementation for AtomKey.

__repr__() str

Repr string from id.

__hash__() int

Hash calculated at init from akl tuple.

altloc_match(other: Bio.PDB.internal_coords.AtomKey) bool

Test AtomKey match to other discounting occupancy and altloc.

is_backbone() bool

Return True if is N, C, CA, O, or H.

atm() str

Return atom name : N, CA, CB, O etc.

__ne__(other: object) bool

Test for inequality.

__eq__(other: object) bool

Test for equality.

__gt__(other: object) bool

Test greater than.

__ge__(other: object) bool

Test greater or equal.

__lt__(other: object) bool

Test less than.

__le__(other: object) bool

Test less or equal.

Bio.PDB.internal_coords.set_accuracy_95(num: float) float

Reduce floating point accuracy to 9.5 (xxxx.xxxxx).

Used by Hedron and Dihedron classes writing PIC and SCAD files.

Parameters

num (float) – input number

Returns

float with specified accuracy

exception Bio.PDB.internal_coords.HedronMatchError

Bases: Exception

Cannot find hedron in residue for given key.

exception Bio.PDB.internal_coords.MissingAtomError

Bases: Exception

Missing atom coordinates for hedron or dihedron.