Bio.PDB.FragmentMapper module

Classify protein backbone structure with Kolodny et al’s fragment libraries.

It can be regarded as a form of objective secondary structure classification. Only fragments of length 5 or 7 are supported (ie. there is a ‘central’ residue).

Full reference:

Kolodny R, Koehl P, Guibas L, Levitt M. Small libraries of protein fragments model native protein structures accurately. J Mol Biol. 2002 323(2):297-307.

The definition files of the fragments can be obtained from:

http://github.com/csblab/fragments/

You need these files to use this module.

The following example uses the library with 10 fragments of length 5. The library files can be found in directory ‘fragment_data’.

>>> from Bio.PDB.PDBParser import PDBParser
>>> from Bio.PDB.FragmentMapper import FragmentMapper
>>> parser = PDBParser()
>>> structure = parser.get_structure("1a8o", "PDB/1A8O.pdb")
>>> model = structure[0]
>>> fm = FragmentMapper(model, lsize=10, flength=5, fdir="PDB")
>>> chain = model['A']
>>> res152 = chain[152]
>>> res157 = chain[157]
>>> res152 in fm # is res152 mapped? (fragment of a C-alpha polypeptide)
False
>>> res157 in fm # is res157 mapped? (fragment of a C-alpha polypeptide)
True
class Bio.PDB.FragmentMapper.Fragment(length, fid)

Bases: object

Represent a polypeptide C-alpha fragment.

__init__(self, length, fid)

Initialize fragment object.

Parameters
  • length (int) – length of the fragment

  • fid (int) – id for the fragment

get_resname_list(self)

Get residue list.

Returns

the residue names

Return type

[string, string,..]

get_id(self)

Get identifier for the fragment.

Returns

id for the fragment

Return type

int

get_coords(self)

Get the CA coordinates in the fragment.

Returns

the CA coords in the fragment

Return type

Numeric (Nx3) array

add_residue(self, resname, ca_coord)

Add a residue.

Parameters
  • resname (string) – residue name (eg. GLY).

  • ca_coord (Numeric array with length 3) – the c-alpha coorinates of the residues

__len__(self)

Return lengt of the fragment.

__sub__(self, other)

Return rmsd between two fragments.

Returns

rmsd between fragments

Return type

float

Examples

This is an incomplete but illustrative example:

rmsd = fragment1 - fragment2
__repr__(self)

Represent the fragment object as a string.

Returns <Fragment length=L id=ID> where L=length of fragment and ID the identifier (rank in the library).

class Bio.PDB.FragmentMapper.FragmentMapper(model, lsize=20, flength=5, fdir='.')

Bases: object

Map polypeptides in a model to lists of representative fragments.

__init__(self, model, lsize=20, flength=5, fdir='.')

Create instance of FragmentMapper.

Parameters
  • model (L{Model}) – the model that will be mapped

  • lsize (int) – number of fragments in the library

  • flength (int) – length of fragments in the library

  • fdir (string) – directory where the definition files are found (default=”.”)

__contains__(self, res)

Check if the given residue is in any of the mapped fragments.

__getitem__(self, res)

Get an entry.

Returns

fragment classification

Return type

L{Fragment}