Multi Coor class

class pdb_manip_py.pdb_manip.Multi_Coor(pdb_in=None, pqr_format=False)

Topologie base on coordinates like pdb or gro.

The coor object containt a a list of dictionnary of atoms indexed on the atom num and the crystal packing info.

Parameters
  • coor_list – list of dictionnary of atom

  • crystal_pack (str) – crystal packing

Note

Files necessary for testing : * ../test/input/1y0m.pdb, ../test/input/1rxz.pdb and ../test/input/4n1m.pdb. To do the unitary test, execute pdb_mani.py (-v for verbose mode)

compute_rmsd_to(atom_sel_2, selec_dict={'name': ['CA']})

Compute RMSD between two atom_dict Then return the RMSD value.

Parameters
  • atom_sel_1 (dict) – atom dictionnary

  • atom_sel_2 (dict) – atom dictionnary

Returns

distance

Return type

float

Example

>>> TEST_OUT = str(getfixture('tmpdir'))
>>> VIP_coor = Multi_Coor(os.path.join(TEST_PATH, '2rri.pdb'))        
Read 20 Model(s)
Succeed to read file ...2rri.pdb, 479 atoms found
>>> aa_seq = VIP_coor.coor_list[0].get_aa_seq()['A']
>>> print(aa_seq)
HSDAVFTDNYTRLRKQMAVKKYLNSILNG
>>> linear_pep = Coor()
>>> pep_out = os.path.join(TEST_OUT, 'tmp_pep.pdb')
>>> pdb_pep = linear_pep.make_peptide(aa_seq, pep_out)        
-Make peptide: HSDAVFTDNYTRLRKQMAVKKYLNSILNG
residue name:X
residue name:H
residue name:S
residue name:D
residue name:A
residue name:V
residue name:F
residue name:T
residue name:D
residue name:N
residue name:Y
residue name:T
residue name:R
residue name:L
residue name:R
residue name:K
residue name:Q
residue name:M
residue name:A
residue name:V
residue name:K
residue name:K
residue name:Y
residue name:L
residue name:N
residue name:S
residue name:I
residue name:L
residue name:N
residue name:G
Succeed to save file ...tmp_pep.pdb
>>> linear_pep.read_pdb(pdb_pep) 
Succeed to read file ...tmp_pep.pdb ,  240 atoms found
>>> rmsd_list = VIP_coor.compute_rmsd_to(linear_pep)
>>> rmsd_str = ['{:.2f}'.format(i) for i in rmsd_list]
>>> rmsd_str
['58.57', '58.40', '58.74', '58.35', '58.60', '58.53', '58.49', '58.40', '58.45', '58.27', '58.52', '58.34', '58.57', '58.33', '58.34', '58.63', '58.61', '58.40', '58.55', '58.32']
read_pdb(pdb_in, pqr_format=False)

Read a pdb file and return atom informations as a dictionnary indexed on the atom num. The fonction can also read pqr files if specified with pqr_format = True, it will only change the column format of beta and occ factors.

Parameters
  • pdb_in (str) – path of the pdb file to read

  • pqr_format (bool, default=False) – Flag for .pqr file format reading.

Example

>>> VIP_coor = Multi_Coor(os.path.join(TEST_PATH, '2rri.pdb'))        
Read 20 Model(s)
Succeed to read file ...2rri.pdb, 479 atoms found
write_pdb(pdb_out, check_file_out=True)

Write a pdb file.

Parameters
  • pdb_out (str) – path of the pdb file to write

  • check_file_out – flag to check if output file already exists

Example

>>> TEST_OUT = str(getfixture('tmpdir'))
>>> VIP_coor = Multi_Coor(os.path.join(TEST_PATH, '2rri.pdb'))        
Read 20 Model(s)
Succeed to read file ...2rri.pdb, 479 atoms found
>>> VIP_coor.write_pdb(os.path.join(TEST_OUT, 'tmp.pdb'))        
Succeed to save file ...tmp.pdb