allpy
changeset 701:8ca932b9e951
structure/SequenceMixin: add resi_begin and resi_end parameters to PDB loader
Optional arguments resi_begin and resi_end where added to methods
set_pdb_chain() and auto_pdb() of SequenceMixin.
These arguments can be used to restrict range of pdb residues,
that can be connected with monomers from sequence.
These arguments are tuples like (hetero_flag, sequence_id, insertion_code).
See biopdb documantation for more information.
closes #61
author | boris <bnagaev@gmail.com> |
---|---|
date | Wed, 06 Jul 2011 20:49:31 +0200 |
parents | bf97a2bc9d08 |
children | eb2a902b1e0a |
files | allpy/structure.py |
diffstat | 1 files changed, 19 insertions(+), 3 deletions(-) [+] |
line diff
1.1 --- a/allpy/structure.py Wed Jul 06 18:05:51 2011 +0400 1.2 +++ b/allpy/structure.py Wed Jul 06 20:49:31 2011 +0200 1.3 @@ -68,7 +68,8 @@ 1.4 xyz_only=False, 1.5 min_continuous_match=config.min_continuous_match, 1.6 min_match=config.min_match, 1.7 - max_waste_in_pdb=config.max_waste_in_pdb): 1.8 + max_waste_in_pdb=config.max_waste_in_pdb, 1.9 + resi_begin=None, resi_end=None): 1.10 """ Read Pdb chain from file 1.11 1.12 and align each Monomer with PDB.Residue 1.13 @@ -77,6 +78,11 @@ 1.14 * min_continuous_match 1.15 * min_match 1.16 * max_waste_in_pdb 1.17 + * resi_begin and resi_end -- identifiers of first and last pdb 1.18 + residues of residue range to be used. If only one is provided, 1.19 + another parameter is set to first or last residue of chain. 1.20 + Id of residue is (hetero_flag, sequence_id, insertion_code), 1.21 + see biopdb documantation for more information. 1.22 1.23 The sequence is aligned with pdb_sequence using muscle. 1.24 The alignment is splitted to continuous_blocks 1.25 @@ -93,6 +99,12 @@ 1.26 if not xyz_only: 1.27 self.pdb_chain = chain 1.28 pdb_sequence = self.__class__.from_pdb_chain(chain) 1.29 + if resi_begin: 1.30 + while pdb_sequence[0].pdb_residue.get_id() != resi_begin: 1.31 + pdb_sequence.pop(0) 1.32 + if resi_end: 1.33 + while pdb_sequence[-1].pdb_residue.get_id() != resi_end: 1.34 + pdb_sequence.pop() 1.35 Alignment = self.types.Alignment 1.36 a = Alignment() 1.37 a.append_sequence(self) 1.38 @@ -158,9 +170,12 @@ 1.39 sequence.append(monomer) 1.40 return sequence 1.41 1.42 - def auto_pdb(self, conformity=None, pdb_getter=download_pdb, xyz_only=False): 1.43 + def auto_pdb(self, conformity=None, pdb_getter=download_pdb, xyz_only=False, 1.44 + resi_begin=None, resi_end=None): 1.45 """ Add pdb information to each monomer 1.46 1.47 + * resi_begin and resi_end -- see set_pdb_chain 1.48 + 1.49 raise AssertionError if info has not been successfully added 1.50 1.51 TODO: conformity_file 1.52 @@ -171,7 +186,8 @@ 1.53 chain = match['chain'] 1.54 model = match['model'] 1.55 pdb_file = pdb_getter(code) 1.56 - self.set_pdb_chain(pdb_file, code, chain, model, xyz_only=xyz_only) 1.57 + self.set_pdb_chain(pdb_file, code, chain, model, xyz_only=xyz_only, 1.58 + resi_begin=resi_begin, resi_end=resi_end) 1.59 1.60 def save_pdb(self, out_filename, new_chain=None, new_model=None): 1.61 """ Save pdb_chain to out_file """