TEMPy.protein.prot_rep_biopy
- class TEMPy.protein.prot_rep_biopy.Atom(atom_name, coords, record_name='ATOM', serial=None, temp_factor=None, alt_loc=None, icode='', charge='', element='', occ=None, res=None, model=None, chain=None, res_no=None)[source]
A generic atomic class. Subclass should ensure desired fields are set.
- change_init_position()[source]
Change initial position co-ordinates of atom to current position.
- Returns:
new initial position co-ordinates of atom.
- get_pos_mass()[source]
- Returns:
An array containing Vector instances containing 3D coordinates of the atom and and its corresponding mass.
- gridpoints_around_atom(densmap, radius=0.0, dx=0, dy=0, dz=0)[source]
Get the x,y,z box boundaries of a set radius around an atom If radius is set, dx, dy, and dz will be calculated using radius, if not dx, dy, and dz given will be used. :param *densmap*: Map instance containing the grid :param *radius*: Cutoff radius around atom (Angstroms) :param *dx*: gridpoint distance along x (in integers) :param *dy*: gridpoint distance along y (in integers) :param *dz*: gridpoint distance along z (in integers)
- Returns:
Minimum and maximum grid index along each direction of x, y, z e.g. [minx, maxx, miny, maxy, minz, maxz]
- map_grid_position(densMap, index=True)[source]
- Parameters:
*densMap* – EM map object consisting the 3D grid of density values.
*index* – Boolean to return co-ordinates of grid as integer index. If false co-ordinates are returned as floats. Default: True
- Returns:
The co-ordinates and density value of the grid point in a density map closest to this atom. Return 0 if atom is outside of map.
- matrix_transform(rot_mat)[source]
Transform atom using a 3x3 matrix
- Parameters:
*rot_mat* – a 3x3 matrix instance.
- Returns:
Transformed position of atom object
- reset_position()[source]
Translate atom back in its initial position.
- Returns:
initial position co-ordinates of atom.
- set_x(mod)[source]
Change the x co-ordinate of an atom based on the argument.
- Parameters:
*mod* – float value
- Returns:
new x co-ordinate
- set_y(mod)[source]
Change the y co-ordinate of an atom based on the argument.
- Parameters:
*mod* – float value
- Returns:
new y co-ordinate
- set_z(mod)[source]
Change the z co-ordinate of an atom based on the argument.
- Parameters:
*mod* – float value
- Returns:
new x co-ordinate
- class TEMPy.protein.prot_rep_biopy.BioPy_Structure(atomList, filename='Unknown', header='', footer='')[source]
A class representing a Structure object.
- RMSD_from_init_position(CA=False)[source]
Return RMSD of structure from initial position after translation.
- Parameters:
*CA* – True will consider only CA atoms. False will consider all atoms.
- Returns:
RMSD in angstrom
- RMSD_from_same_structure(otherStruct, CA=False, write=True)[source]
Return the RMSD between two structure instances.
- Parameters:
*otherStruct* – Structure instance to compare, containing the same number of atoms as the target instance.
*CA* – True will consider only CA atoms. False will consider all atoms.
- Returns:
RMSD in angstrom
- break_into_segments(rigid_list, chain='')[source]
Return a list of Structure instance based on the rigid body list.
- Parameters:
list* (*rigid) –
- list of rigid body defined as:
[[riA,rfA],..,[riB,rfB]]
- where:
riA is the starting residues number of segment A. rfA is the final residues number of segment A.
- Returns:
List of TEMPy Structure instance
- combine_SSE_structures(structList)[source]
Combine a list of Structure instance into one and return a new Structure instance.
- Parameters:
*structList* – list of Structure instance
- Returns:
New Structure Instance
- combine_structures(structList)[source]
Add a list of Structure instance to the existing structure.
- Parameters:
*structList* – list of Structure instance
- Returns:
New Structure Instance
- find_same_atom(atom_index, otherStruct)[source]
Find if an atom exists in the compared structure, based on atom index.
- Parameters:
*atom_index* – atom number
*otherStruct* – a structure object to compare
- Returns:
If a match is found, it returns the atom object; else it returns a string reporting the mismatch.
- get_atom(index)[source]
Return specific atom in Structure instance.
- Parameters:
*index* – Index of the atom
- Returns:
Returns an Atom instance.
- get_atom_list()[source]
- Returns:
[(RES 1 A: x,y,z), … ,(RES2 1 A: x1,y1,z1)].
- Return type:
An array containing Atom instances of positions of all atoms as
- get_extreme_values()[source]
- Returns:
A 6-tuple containing the minimum and maximum of x, y and z co-ordinates of the structure. Given in order (min_x, max_x, min_y, max_y, min_z, max_z).
- get_next_chain_id(chain_id, chain_list=False)[source]
Renames chains based on mmCIF convention i.e. A -> B, Z -> 0 or g -> h See the ALPHABET list for standard sequence, however any string can be used for chain labelling.
For structures with more than 62 chains, the labels are increased in length, i.e. z -> AA, zz -> AAA, AAA -> BBB
Input: chain_id -> old chain id Output: new_chain_id
- get_pos_mass_list()[source]
- Returns:
Array containing Vector instances of positions of all atoms and mass.
- get_prot_mass_from_atoms()[source]
Calculates Mass (kDa) of the Structure instance, from average mass. Atoms based use get_prot_mass_from_res is more accurate.
- get_prot_mass_from_res(Termini=False)[source]
# ADD by IF 22-4-2013 #from Harpal code calculation of mass from seq. Calculates Mass (kDa) of the Structure instance, from average mass #based on http://web.expasy.org/findmod/findmod_masses.html
- get_residue(resNo)[source]
Get the residue corresponding to the residue number.
- Parameters:
*resNo* – Residues number
- Returns:
Returns a Residues instance.
- get_residue_indexes()[source]
Get the atom indexes for each residues from atomList in the structure
- Returns:
a dictionary containing the index positions for all atoms from each residue in the structure. Items have format: [list, of, atom, indexes, as, ints] Keys have format: (str(chain_label), int(residue_number))
- Return type:
residue_indexes
- get_selection(startRes, finishRes, chain='')[source]
Get a new Structure instance for the selected residues range without considering residues chain.
- Parameters:
*startRes* – Start residue number
*finishRes* – End residue number
- Returns:
New Structure instance
- get_selection_less_than(endRes)[source]
Get a Structure instance comprising all residues with their residue numbers less than endRes.
- Parameters:
*endRes* – a residue number
- Returns:
A Structure instance
- get_selection_more_than(startRes)[source]
Get a Structure instance comprising all residues with their residue numbers greater than startRes.
- Parameters:
*startRes* – a residue number
- Returns:
A Structure instance
- iter_asym_id(asym_id)[source]
Function to increase the chain “numbering” by one e.g. from D -> E or from HA -> IA. Uses the chr() and ord() functions to convert between int and ASCII str One should run rename_chains before to ensure all chain_ids are letters
Input: capital letter from A -> Z Output: Next letter in alphabet e.g. from A -> B. For input Z, ouput is
AA, for input ZZ output is AAA etc..
- matrix_transform(matrix)[source]
Transform Structure using a 3x3 transformation matrix
- Parameters:
*rot_mat* – a 3x3 matrix instance.
- Returns:
Transformed position of Structure object
- randomise_position(max_trans, max_rot, v_grain=30, rad=False, verbose=False)[source]
Randomise the position of the Structure instance using random rotations and translations.
- Parameters:
*max_trans* – Maximum translation permitted
*max_rot* – Maximum rotation permitted (in degree if rad=False)
*v_grain* – Graning Level for the generation of random vetors (default=30)
- Returns:
Transformed position of Structure object
- renumber_atoms(start_num=1)[source]
Renumber the atoms in the structure. After renumbering the starting atom number will be 1 unless start_num
- renumber_residues(startRes=1, missingRes=[])[source]
Renumber the structure starting from startRes. Missing number list to add.
- Parameters:
*startRes* – Starting residue number for renumbering
*missingRes* – A list of missing residue numbers to add
- reorder_residues()[source]
Order residues in atom list by residue number. (NOTE: Does not check for chain information - split by chain first).
- rotate_by_axis_angle(x, y, z, turn, x_trans=0, y_trans=0, z_trans=0, rad=False, com=False)[source]
Rotate the Structure instance around its centre.
- Parameters:
*turn* – angle (in radians if rad == True, else in degrees) to rotate map.
*x – axis to rotate about, ie. x,y,z = 0,0,1 rotates the structure round the xy-plane.
y – axis to rotate about, ie. x,y,z = 0,0,1 rotates the structure round the xy-plane.
z* – axis to rotate about, ie. x,y,z = 0,0,1 rotates the structure round the xy-plane.
*x_trans – extra translational movement if required.
y_trans – extra translational movement if required.
z_trans* – extra translational movement if required.
*com* – centre of mass around which to rotate the structure. If False, rotates around centre of mass of structure.
- rotate_by_euler(x_turn, y_turn, z_turn, x_trans=0, y_trans=0, z_trans=0, rad=False, com=False)[source]
Rotate this Structure instance around its centre.
- Parameters:
*x_turn – Euler angles (in radians if rad == True, else in degrees) used to rotate structure, in order XYZ.
y_turn – Euler angles (in radians if rad == True, else in degrees) used to rotate structure, in order XYZ.
z_turn* – Euler angles (in radians if rad == True, else in degrees) used to rotate structure, in order XYZ.
*x_trans – extra translational movement if required.
y_trans – extra translational movement if required.
z_trans* – extra translational movement if required.
*com* – centre of mass around which to rotate the structure. If False, rotates around centre of mass of structure.
- split_into_chains()[source]
Split the structure into separate chains and returns the list of Structure instance for each chain.
- translate(x, y, z)[source]
Translate the structure.
- Parameters:
*x – distance in Angstroms in respective directions to move structure.
y – distance in Angstroms in respective directions to move structure.
z* – distance in Angstroms in respective directions to move structure.
- Returns:
Translated Structure instance
- class TEMPy.protein.prot_rep_biopy.gemmiAtom(chain, entity, residue, atom)[source]
Class representing an atom as read from a PDB/mmCIF file using gemmi
gemmi docs here: https://gemmi.readthedocs.io/
- class TEMPy.protein.prot_rep_biopy.gemmi_Structure(atomList, gemmi_structure, filename='Unknown', header='', tempy_scores=[])[source]
Class for representing a protein structure.
Extends BioPy_Structure to deal with additional mmCif labels for dividing atoms into entities, asym units and chains.
- add_structure(filename, integrate=False, hetatm=False, water=False)[source]
Add structure from file
- add_structure_instance(new_structure, integrate=False)[source]
Add a seperate structure instance to self
- Arguements:
- new_structure
the new structure instance which will be incorpoated into self
- Optional Arguments:
- integrate
If True, atoms from the new_structure will be incorporated into chains with the same name in self.
- combine_structures(structList)[source]
Combine a list of structure instances with self and returns the new structure as a seperate instance.
Note: this function does not mutate self.
- Arguments:
- structList
List of Structure instances
- Returns:
- combined_struct
New gemmi_Structure instance
- rename_asym_labels(label_list=False, starting_val=False)[source]
Renames asym_id for each atom in alphabetical order
For structures with more than 26 models the ids switch to doulbe lettering (i.e. from Z -> AA and from ZZZ -> AAAA)
- rename_chains(chain_list=False, starting_val=False)[source]
Rename chain labels based on the list of new chain names
- Optional Arguments:
- chain_list
List of chain names If False rename in alphabetical order.
- rename_entity_labels(starting_val=1)[source]
Renames asym_id for each atom in alphabetical order - for structures with more than 26 models the ids switch to doulbe lettering (i.e. from Z -> AA and from ZZZ -> AAAA)
- set_atom_labels()[source]
Sets chain, asym_id and entity_id labels for the atoms in atomList
Important to run this before file writing
- set_global_score(score_name, score)[source]
Function to incorporate global scores into the structure instance so these values can be written into the cif file as label-value pairs
- set_local_score(score_name, scores, verbose=False)[source]
Function to incorporate local scores into structure instance, so that these values can be written into new cif files
This would be efficiently handled using a pandas library