Bio.PDB.internal_coords module¶
Classes to support internal coordinates for protein structures.
Internal coordinates comprise Psi, Phi and Omega dihedral angles along the protein backbone, Chi angles along the sidechains, and all 3-atom angles and bond lengths comprising a protein chain. These routines can compute internal coordinates from atom XYZ coordinates, and compute atom XYZ coordinates from internal coordinates.
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.
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.
Internal coordinates may also be exported as OpenSCAD data arrays for generating 3D printed protein models. OpenSCAD software is provided as proof-of-concept for generating such models.
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; methods generally apply IC_Residue methods along chain. 
- IC_Residue: Extends for Biopython Residue on .internal_coord attribute.
- Most control and methods of interest are in this class, see API. 
- Dihedron: four joined atoms forming a dihedral angle.
- Dihedral angle, homogeneous atom coordinates in local coordinate space, references to relevant Hedra and IC_Residue. Methods to compute 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: biopython Chain object reference
- The Chain object this extends 
- initNCaC: AtomKey indexed dictionary of N, Ca, C atom coordinates.
- 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. 
- MaxPeptideBond: **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 fully disordered (missing) residues. 
- ordered_aa_ic_list: list of IC_Residue objects
- IC_Residue objects ic algorithms can process (e.g. no waters) 
- hedra: dict indexed by 3-tuples of AtomKeys
- Hedra forming residues in this chain 
- hedraLen: int length of hedra dict
- hedraNdx: dict mapping hedra AtomKeys to numpy array data
- dihedra: dict indexed by 4-tuples of AtomKeys
- Dihedra forming (overlapping) this residue 
- dihedraLen: int length of dihedra dict
- dihedraNdx: dict mapping dihedra AtomKeys to numpy array data
- atomArray: numpy array of homogeneous atom coords for chain
- atomArrayIndex: dict mapping AtomKeys to atomArray indexes
- numpy arrays for vector processing of chain di/hedra:
- hedraIC: length-angle-length entries for each hedron
- hAtoms: homogeneous atom coordinates (3x4) of hedra, central atom at origin
- hAtomsR: hAtoms in reverse order
- hAtoms_needs_update: booleans indicating whether hAtoms represent hedraIC
- dihedraIC: dihedral angles for each dihedron
- dAtoms: homogeneous atom coordinates (4x4) of dihedra, second atom at origin
- dAtoms_needs_update: booleans indicating whether dAtoms represent dihedraIC
 
 - Methods - internal_to_atom_coordinates(verbose, start, fin) - Process ic data to Residue/Atom coordinates; calls assemble_residues() followed by coords_to_structure() - assemble_residues(verbose, start, fin) - Generate IC_Residue atom coords from internal coordinates - coords_to_structure() - update Biopython Residue.Atom coords from IC_Residue coords for all Residues with IC_Residue attributes - atom_to_internal_coordinates(verbose) - Calculate dihedrals, angles, bond lengths (internal coordinates) for Atom data - link_residues() - Call link_dihedra() on each IC_Residue (needs rprev, rnext set) - set_residues() - Add .internal_coord attribute for all Residues in parent Chain, populate ordered_aa_ic_list, set IC_Residue rprev, rnext or initNCaC coordinates - write_SCAD() - Write OpenSCAD matrices for internal coordinate data comprising chain - 
MaxPeptideBond= 1.4¶
 - 
__init__(self, parent: ‘Chain’, verbose: bool = False) → None¶
- Initialize IC_Chain object, with or without residue/Atom data. - Parameters
- parent – Biopython Chain object Chain object this extends 
 
 - 
clear_ic(self)¶
- Clear residue internal_coord settings for this chain. 
 - 
set_residues(self, verbose: bool = False) → None¶
- Initialize internal_coord data for loaded Residues. - Add IC_Residue as .internal_coord attribute for each Residue in parent Chain; populate ordered_aa_ic_list with IC_Residue references for residues which can be built (amino acids and some hetatms); set rprev and rnext on each sequential IC_Residue, populate initNCaC at start and after chain breaks. 
 - 
link_residues(self) → None¶
- link_dihedra() for each IC_Residue; needs rprev, rnext set. - Called by PICIO:read_PIC() after finished reading chain 
 - 
assemble_residues(self, verbose: bool = False, start: Union[int, NoneType] = None, fin: Union[int, NoneType] = None) → None¶
- Generate IC_Residue atom coords from internal coordinates. - Filter positions between start and fin if set, find appropriate start coordinates for each residue and pass to IC_Residue.assemble() - Parameters
- bool (verbose) – default False describe runtime problems 
- Param
- start, fin lists sequence position, insert code for begin, end of subregion to process 
 
 - 
coords_to_structure(self) → None¶
- Promote all ic atom_coords to Biopython Residue/Atom coords. - IC atom_coords are homogeneous [4], Biopython atom coords are XYZ [3]. 
 - 
init_edra(self) → None¶
- Create chain level di/hedra arrays. - If called by read_PIC, self.di/hedra = {} and object tree has IC data. -> build chain arrays from IC data - If called at start of atom_to_internal_coords, self.di/hedra fully populated. -> create empty chain numpy arrays - In both cases, fix di/hedra object attributes to be views on chain-level array data 
 - 
init_atom_coords(self) → None¶
- Set chain level di/hedra initial coord arrays from IC_Residue data. 
 - 
internal_to_atom_coordinates(self, verbose: bool = False, start: Union[int, NoneType] = None, fin: Union[int, NoneType] = None, promote: Union[bool, NoneType] = True) → None¶
- Process, IC data to Residue/Atom coords. - Not yet vectorized. - Parameters
- bool (promote) – default False describe runtime problems 
- bool – default True If True (the default) copy result atom XYZ coordinates to Biopython Atom objects for access by other Biopython methods; otherwise, updated atom coordinates must be accessed through IC_Residue and hedron objects. 
 
- Param
- start, fin lists sequence position, insert code for begin, end of subregion to process 
 
 - 
atom_to_internal_coordinates(self, verbose: bool = False) → None¶
- Calculate dihedrals, angles, bond lengths for Atom data. - Parameters
- bool (verbose) – default False describe runtime problems 
 
 - 
write_SCAD(self, fp: <class 'TextIO'>, backboneOnly: bool) → None¶
- Write self to file fp as OpenSCAD data matrices. - Works with write_SCAD() and embedded OpenSCAD routines in SCADIO.py. The OpenSCAD code explicitly creates spheres and cylinders to represent atoms and bonds in a 3D model. Options are available to support rotatable bonds and magnetic hydrogen bonds. - Matrices are written to link, enumerate and describe residues, dihedra, hedra, and chains, mirroring contents of the relevant IC_* data structures. - The OpenSCAD matrix of hedra has additional information as follows: - the atom and bond state (single, double, resonance) are logged so that covalent radii may be used for atom spheres in the 3D models 
- bonds and atoms are tracked so that each is only created once 
- bond options for rotation and magnet holders for hydrogen bonds may be specified 
 - Note the application of IC_Chain attribute MaxPeptideBond: missing residues may be linked (joining chain segments with arbitrarily long bonds) by setting this to a large value. - All ALTLOC (disordered) residues and atoms are written to the output model. 
 
- 
class Bio.PDB.internal_coords.IC_Residue(parent: Residue, NO_ALTLOC: bool = False)¶
- Bases: - object- Class to extend Biopython Residue with internal coordinate data. - Attributes
- 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
- Local copy of atom homogeneous coordinates [4] for work distinct from Bopython Residue/Atom values 
- 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. Set by link_dihedra() 
- is20AA: bool
- True if residue is one of 20 standard amino acids, based on Residue resname 
- accept_atoms: tuple
- list of PDB atom names to use when generatiing 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 backbone 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 - There is currently no option to output internal coordinates with D instead of H 
- accept_resnames: tuple
- 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.py and support in atom_to_internal_coordinates() 
- gly_Cbeta: bool default False
- override class variable to True to generate internal coordinates for glycine CB atoms in atom_to_internal_coordinates(). - IC_Residue.gly_Cbeta = True 
- allBonds: bool default False
- whereas a PDB file just specifies atoms, OpenSCAD output for 3D printing needs all bonds specified explicitly - otherwise, e.g. PHE rings will not be closed. This variable is managed by the Write_SCAD() code and enables this. 
- cic: IC_Chain default None
- parent chain IC_Chain object 
- scale: optional float
- used for OpenSCAD output to generate gly_Cbeta bond length 
 
 - Methods - applyMtx() - multiply all IC_Residue atom_cords by passed matrix - assemble(atomCoordsIn, resetLocation, verbose) - Compute atom coordinates for this residue from internal coordinates - atm241(coord) - Convert 1x3 cartesian coords to 4x1 homogeneous coords - coords_to_residue() - Convert homogeneous atom_coords to Biopython cartesian Atom coords - atom_to_internal_coordinates(verbose) - Create hedra and dihedra for atom coordinates - get_angle() - Return angle for passed key - get_length() - Return bond length for specified pair - link_dihedra() - Link dihedra to this residue, form id3_dh_index - load_PIC(edron) - Process parsed (di-/h-)edron data from PIC file - pick_angle() - Find Hedron or Dihedron for passed key - pick_length() - Find hedra for passed AtomKey pair - rak(atom info) - Residue AtomKey - per residue AtomKey result cache - set_angle() - Set angle for passed key (no position updates) - set_length() - Set bond length in all relevant hedra for specified pair - write_PIC(pdbid, chainId, s) - Generate PIC format strings for this residue - 
accept_resnames= ('CYG', 'YCM', 'UNK')¶
 - 
AllBonds: bool = False¶
 - 
__init__(self, parent: ‘Residue’, NO_ALTLOC: bool = False) → None¶
- Initialize IC_Residue with parent Biopython Residue. - Parameters
- parent – Biopython Residue object The Biopython Residue this object extends 
- NO_ALTLOC – bool default False Option to disable processing altloc disordered atoms, use selected. 
 
 
 - 
rak(self, atm: Union[str, Bio.PDB.Atom.Atom]) → ’AtomKey’¶
- Cache calls to AtomKey for this residue. 
 - 
build_rak_cache(self) → None¶
- Create explicit entries for for atoms so don’t miss altlocs. 
 - 
accept_mainchain= ('N', 'CA', 'C', 'O', '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', 'OXT')¶
 - 
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', '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', 'OXT', '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')¶
 - 
gly_Cbeta= False¶
 - 
static atm241(coord: <built-in function array>) → <built-in function array>¶
- Convert 1x3 cartesian coordinates to 4x1 homogeneous coordinates. 
 - 
__repr__(self) → str¶
- Print string is parent Residue ID. 
 - 
pretty_str(self) → str¶
- Nice string for residue ID. 
 - 
load_PIC(self, edron: Dict[str, str]) → None¶
- Process parsed (di-/h-)edron data from PIC file. - Parameters
- edron – parse dictionary from Edron.edron_re 
 
 - 
link_dihedra(self, verbose: bool = False) → None¶
- Housekeeping after loading all residues and dihedra. - Link dihedra to this residue 
- form id3_dh_index 
- form ak_set 
- set NCaCKey to be available AtomKeys 
 
 - 
set_flexible(self) → None¶
- For OpenSCAD, mark N-CA and CA-C bonds to be flexible joints. 
 - 
set_hbond(self) → None¶
- For OpenSCAD, mark H-N and C-O bonds to be hbonds (magnets). 
 - 
default_startpos(self) → Dict[ForwardRef(‘AtomKey’), <built-in function array>]¶
- Generate default N-Ca-C coordinates to build this residue from. 
 - 
get_startpos(self) → Dict[ForwardRef(‘AtomKey’), <built-in function array>]¶
- Find N-Ca-C coordinates to build this residue from. 
 - 
clear_transforms(self)¶
- Set cst and rcst attributes to none before assemble(). 
 - 
assemble(self, resetLocation: bool = False, verbose: bool = False) → Union[Dict[ForwardRef(‘AtomKey’), <built-in function array>], Dict[Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)], <built-in function array>], NoneType]¶
- Compute atom coordinates for this residue from internal coordinates. - Join dihedrons 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 coord_space()) - Not vectorized (yet). - 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 N-Ca-C hedron at origin. 
- Returns
- Dict of AtomKey -> homogeneous atom coords for residue in protein space relative to previous residue 
 
 - 
atom_to_internal_coordinates(self, verbose: bool = False) → None¶
- Create hedra and dihedra for atom coordinates. - Parameters
- verbose – bool default False warn about missing N, Ca, C backbone atoms. 
 
 - 
build_glyCB(self, gCBd: ‘Dihedron’)¶
- Populate values for Gly C-beta, rest of chain complete. - Data averaged from Sep 2019 Dunbrack cullpdb_pc20_res2.2_R1.0 restricted to structures with amide protons. - 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 avg_rslt | sd_rslt | count |-122.682194862932 | 5.04403040513919 | 14098 |
 - 
write_PIC(self, pdbid: str = '0PDB', chainid: str = 'A', s: str = '') → str¶
- Write PIC format lines for this residue. - Parameters
- pdbid (str) – PDB idcode string; default 0PDB 
- chainid (str) – PDB Chain ID character; default A 
- s (str) – result string to add to 
 
 
 - 
coords_to_residue(self, rnext: bool = False) → None¶
- Convert self.atom_coords to biopython Residue Atom coords. - Copy homogeneous IC_Residue atom_coords to self.residue cartesian Biopython Atom coords. - Parameters
- rnext – bool default False next IC_Residue has no atom coords due to missing atoms, so try to populate with any available coords calculated for this residue di/hedra extending into next 
 
 - 
relative_atom_re= re.compile('^(-?[10])([A-Z]+)$')¶
 - 
pick_angle(self, angle_key: Union[Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)], Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)], str]) → Union[ForwardRef(‘Hedron’), ForwardRef(‘Dihedron’), NoneType]¶
- 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 
- except for tuples of AtomKeys, no option to access 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 is handled here, but may be an issue if accessing directly. - The following 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 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(self, angle_key: Union[Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)], Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)], str]) → Union[float, NoneType]¶
- Get dihedron or hedron angle for specified key. - See pick_angle() for key specifications. 
 - 
set_angle(self, angle_key: Union[Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)], Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)], str], v: float)¶
- Set dihedron or hedron angle for specified key. - See pick_angle() for key specifications. 
 - 
pick_length(self, ak_spec: Union[str, Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)]]) → Tuple[Union[List[ForwardRef(‘Hedron’)], NoneType], Union[Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)], NoneType]]¶
- 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. 
 
 - 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, tuple of atom keys 
 
 - 
get_length(self, ak_spec: Union[str, Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)]]) → Union[float, NoneType]¶
- Get bond length for specified atom pair. - See pick_length() for ak_spec. 
 - 
set_length(self, ak_spec: Union[str, Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)]], val: float) → None¶
- Set bond length for specified atom pair. - See pick_length() for ak_spec. 
 - 
applyMtx(self, mtx: <built-in function array>) → None¶
- Apply matrix to atom_coords for this IC_Residue. 
 - 
__annotations__= {'AllBonds': <class 'bool'>}¶
 
- 
class Bio.PDB.internal_coords.Edron(*args: Union[List[AtomKey], Tuple[AtomKey, AtomKey, AtomKey], Tuple[AtomKey, AtomKey, AtomKey, 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) AtomKeys 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 
- 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 
 
 - Methods - gen_key([AtomKey, …] or AtomKey, …) (Static Method) - generate a ‘:’-joined string of AtomKey Ids - gen_acs(atom_coords) - generate tuple of atom coords for keys in self.aks - 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<)¶
 - 
static gen_key(lst: Union[List[str], List[ForwardRef(‘AtomKey’)]]) → str¶
- Generate string of ‘:’-joined AtomKey strings from input. - Parameters
- lst – list of AtomKey objects or id strings 
 
 - 
__init__(self, *args: Union[List[ForwardRef(‘AtomKey’)], Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)], Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘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’ … 
 - 
gen_acs(self, atom_coords: Dict[ForwardRef('AtomKey'), <built-in function array>]) → Tuple[<built-in function array>, …]¶
- Generate tuple of atom coord arrays for keys in self.aks. - Parameters
- atom_coords – AtomKey dict of atom coords for residue 
- Raises
- MissingAtomError any atoms in self.aks missing coordinates 
 
 - 
is_backbone(self) → bool¶
- Report True for contains only N, C, CA, O, H atoms. 
 - 
__repr__(self) → str¶
- Tuple of AtomKeys is default repr string. 
 - 
__hash__(self) → int¶
- Hash calculated at init from aks tuple. 
 - 
__eq__(self, other: object) → bool¶
- Test for equality. 
 - 
__ne__(self, other: object) → bool¶
- Test for inequality. 
 - 
__gt__(self, other: object) → bool¶
- Test greater than. 
 - 
__ge__(self, other: object) → bool¶
- Test greater or equal. 
 - 
__lt__(self, other: object) → bool¶
- Test less than. 
 - 
__le__(self, other: object) → bool¶
- Test less or equal. 
 
- 
class Bio.PDB.internal_coords.Hedron(*args: Union[List[AtomKey], Tuple[AtomKey, AtomKey, 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. Stored in two orientations, with the 3rd (forward) or first (reversed) atom on the +Z axis. - Attributes
- lal: numpy array of len12, angle, len23
- len12 = distance between 1st and 2nd atom angle = angle (degrees) formed by 3 atoms len23 = distance between 2nd and 3rd atoms 
- atoms: 3x4 numpy arrray (view on chain array)
- 3 homogeneous atoms comprising hedron, 1st on XZ, 2nd at origin, 3rd on +Z 
- atomsR: 3x4 numpy array (view on chain array)
- atoms reversed, 1st on +Z, 2nd at origin, 3rd on XZ plane 
 
 - Methods - get_length() - get bond length for specified atom pair - set_length() - set bond length for specified atom pair - angle(), len12(), len23() - gettters and setters for relevant attributes (angle in degrees) - 
__init__(self, *args: Union[List[ForwardRef(‘AtomKey’)], Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)]], **kwargs: str) → None¶
- Initialize Hedron with sequence of AtomKeys, kwargs. - Acceptable input:
- As for Edron, plus optional ‘len12’, ‘angle’, ‘len23’ keyworded values. 
 
 - 
__repr__(self) → str¶
- Print string for Hedron object. 
 - 
property angle¶
- Get this hedron angle. 
 - 
property len12¶
- Get first length for Hedron. 
 - 
property len23¶
- Get second length for Hedron. 
 - 
get_length(self, ak_tpl: Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)]) → Union[float, NoneType]¶
- Get bond length for specified atom pair. - Parameters
- ak_tpl – tuple of AtomKeys pair of atoms in this Hedron 
 
 - 
set_length(self, ak_tpl: Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)], newLength: float)¶
- Set bond length for specified atom pair; sets needs_update. - Parameters
- ak_tpl – tuple of AtomKeys pair of atoms in this Hedron 
 
 
- 
class Bio.PDB.internal_coords.Dihedron(*args: Union[List[AtomKey], Tuple[AtomKey, AtomKey, AtomKey, 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 
- 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 
- initial_coords: tuple[4] of numpy arrays [4]
- Local atom coords for 4 atoms, [0] on XZ plane, [1] at origin, [2] on +Z, [3] rotated by dihedral 
- a4_pre_rotation: numpy array [4]
- 4th atom of dihedral aligned to XZ plane (angle not applied) 
- ic_residue: 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 (configured by set_hedra()) 
- cst, rcst: numpy array [4][4]
- transforms to and from coordinate space defined by first hedron. set by IC_Residue.assemble(). defined by id3 order NOT h1key order (atoms may be reversed between these two) 
 
 - Methods - set_hedra() - work out hedra keys and orientation for this dihedron - angle() - getter/setter for dihdral angle in degrees - 
__init__(self, *args: Union[List[ForwardRef(‘AtomKey’)], Tuple[ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’), ForwardRef(‘AtomKey’)]], **kwargs: str) → None¶
- Initialize Dihedron with sequence of AtomKeys and optional dihedral angle. - Acceptable input:
- As for Edron, plus optional ‘dihedral’ keyworded angle value. 
 
 - 
__repr__(self) → str¶
- Print string for Dihedron object. 
 - 
set_hedra(self) → Tuple[bool, Bio.PDB.internal_coords.Hedron, Bio.PDB.internal_coords.Hedron]¶
- Work out hedra keys and set rev flag. 
 - 
property angle¶
- Get dihedral angle. 
 
- 
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, insertion code, 1 or 3 char residue name, atom name, altloc, and occupancy 
 - 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 
- endnum_re: compiled regex (Class Attribute)
- A compiled regular expresion capturing digits at end of a string 
- d2h: bool (Class Attribute)
- Convert D atoms to H on input; 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 (not set here) 
 
 - Methods - altloc_match(other) - Returns True if this AtomKey matches other AtomKey excluding altloc and occupancy fields - 
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?))?$')¶
 - 
endnum_re= re.compile('\\D+(\\d+)$')¶
 - 
fieldNames= ('respos', 'icode', 'resname', 'atm', 'altloc', 'occ')¶
 - 
class fieldsDef(respos, icode, resname, atm, altloc, occ)¶
- Bases: - tuple- 
__getnewargs__(self)¶
- Return self as a plain tuple. Used by copy and pickle. 
 - 
static __new__(_cls, respos, icode, resname, atm, altloc, occ)¶
- Create new instance of fieldsDef(respos, icode, resname, atm, altloc, occ) 
 - 
__repr__(self)¶
- Return a nicely formatted representation string 
 - 
__slots__= ()¶
 - 
property altloc¶
- Alias for field number 4 
 - 
property atm¶
- Alias for field number 3 
 - 
property icode¶
- Alias for field number 1 
 - 
property occ¶
- Alias for field number 5 
 - 
property resname¶
- Alias for field number 2 
 - 
property respos¶
- Alias for field number 0 
 
- 
 - 
fields= fieldsDef(respos=0, icode=1, resname=2, atm=3, altloc=4, occ=5)¶
 - 
d2h= False¶
 - 
__init__(self, *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 
 
 - 
__repr__(self) → str¶
- Repr string from id. 
 - 
__hash__(self) → int¶
- Hash calculated at init from akl tuple. 
 - 
altloc_match(self, other: ‘AtomKey’) → bool¶
- Test AtomKey match other discounting occupancy and altloc. 
 - 
__ne__(self, other: object) → bool¶
- Test for inequality. 
 - 
__eq__(self, other: object) → bool¶
- Test for equality. 
 - 
__gt__(self, other: object) → bool¶
- Test greater than. 
 - 
__ge__(self, other: object) → bool¶
- Test greater or equal. 
 - 
__lt__(self, other: object) → bool¶
- Test less than. 
 - 
__le__(self, 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. :param float num: 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.