PBxplore API cookbook — PB assignation

Note

This page is initialy a jupyter notebook. You can see a notebook HTML render of it or download the notebook itself.

We hereby demonstrate how to use the PBxplore API to assign PB sequences.

from __future__ import print_function, division
from pprint import pprint
%cd ../../../
/Users/jon/dev/PBxplore
import pbxplore as pbx

Use the built-in structure parser

Assign PB for a single structure

The pbxplore.chains_from_files() function is the prefered way to read PDB and PDBx/mmCIF files using PBxplore. This function takes a list of file path as argument, and yield each chain it can read from these files. It provides a single interface to read PDB and PDBx/mmCIF files, to read single model and multimodel files, and to read a single file of a collection of files.

Here we want to read a single file with a single model and a single chain. Therefore, we need the first and only record that is yield by pbxplore.chains_from_files(). This record contains a name for the chain, and the chain itself as a pbxplore.structure.Chain object. Note that, even if we want to read a single file, we need to provide it as a list to pbxplore.chains_from_files().

structure_reader = pbx.chains_from_files(['demo1_assignation/1BTA.pdb'])
chain_name, chain = next(structure_reader)
print(chain_name)
print(chain)
Read 1 chain(s) in demo1_assignation/1BTA.pdb
demo1_assignation/1BTA.pdb | chain A
Chain A / model : 1434 atoms

Protein Blocks are assigned based on the dihedral angles of the backbone. So we need to calculate them. The pbxplore.structure.Chain.get_phi_psi_angles() methods calculate these angles and return them in a form that can be directly provided to the assignement function.

The dihedral angles are returned as a dictionnary. Each key of this dictionary is a residue number, and each value is a dictionary with the phi and psi angles.

dihedrals = chain.get_phi_psi_angles()
pprint(dihedrals)
{1: {'phi': None, 'psi': -171.6556313444855},
 2: {'phi': -133.80467711845586, 'psi': 153.74322760775027},
 3: {'phi': -134.66175688926953, 'psi': 157.30476083095584},
 4: {'phi': -144.49159910635186, 'psi': 118.59706956501037},
 5: {'phi': -100.12866913978127, 'psi': 92.98634825528089},
 6: {'phi': -83.48980457968895, 'psi': 104.23730726195485},
 7: {'phi': -64.77163869310709, 'psi': -43.25159835828049},
 8: {'phi': -44.47885842536948, 'psi': -25.89184262616925},
 9: {'phi': -94.90790101955957, 'psi': -47.182577907117775},
 10: {'phi': -41.31267169232996, 'psi': 133.73743399231304},
 11: {'phi': -119.15122785547305, 'psi': -11.82789586402356},
 12: {'phi': -174.21196552933984, 'psi': 175.87239770676175},
 13: {'phi': -56.61341695443224, 'psi': -45.74767617535588},
 14: {'phi': -50.78226415072095, 'psi': -45.3742585970337},
 15: {'phi': -57.93584481869442, 'psi': -43.329444361460844},
 16: {'phi': -55.209603541130434, 'psi': -56.47559202715399},
 17: {'phi': -64.51979885245254, 'psi': -18.577118068149446},
 18: {'phi': -70.24273354141468, 'psi': -55.153744337676926},
 19: {'phi': -65.20648546633561, 'psi': -41.28370221159946},
 20: {'phi': -58.98821952110768, 'psi': -35.78957701447905},
 21: {'phi': -66.8659714296852, 'psi': -42.14634696303375},
 22: {'phi': -67.34201665142825, 'psi': -57.40438549689628},
 23: {'phi': -52.29793609141382, 'psi': -66.09120830346023},
 24: {'phi': -61.19010445362886, 'psi': -14.807316930892569},
 25: {'phi': 54.951586944206355, 'psi': 47.59528477656777},
 26: {'phi': -69.51531755580697, 'psi': 161.10806531443862},
 27: {'phi': -57.36300935545188, 'psi': -179.66365615297644},
 28: {'phi': -79.91369005407893, 'psi': -18.494472196394668},
 29: {'phi': -93.51717329199727, 'psi': 13.80253054655975},
 30: {'phi': -38.40653214238887, 'psi': 105.85297788366393},
 31: {'phi': -64.01559307951965, 'psi': -5.507357757886837},
 32: {'phi': 50.06519606710964, 'psi': 3.604730286754302},
 33: {'phi': -84.83560576923662, 'psi': -176.04877012701309},
 34: {'phi': -76.65985981150652, 'psi': -36.89428882663367},
 35: {'phi': -66.20745817863622, 'psi': -36.19018119951471},
 36: {'phi': -80.76844188891471, 'psi': -55.88509876949212},
 37: {'phi': -45.0995601497454, 'psi': -50.82304368319501},
 38: {'phi': -58.512419169182465, 'psi': -56.4318511704347},
 39: {'phi': -44.00775783983471, 'psi': -26.06209153795419},
 40: {'phi': -79.6799641005731, 'psi': -51.3827703817916},
 41: {'phi': -58.80532943671335, 'psi': -49.46425322450557},
 42: {'phi': -75.73059711071141, 'psi': 3.9162670655634235},
 43: {'phi': -177.14613562249534, 'psi': 60.46495675947551},
 44: {'phi': 177.12658169328853, 'psi': -66.62887199130637},
 45: {'phi': -58.436100708193806, 'psi': 149.59997847317612},
 46: {'phi': -102.66050573267097, 'psi': 132.43212727859543},
 47: {'phi': -114.52132755246623, 'psi': 169.33012343233455},
 48: {'phi': -61.39150617820462, 'psi': 136.7035538929314},
 49: {'phi': -113.17589693608565, 'psi': 156.54195412530404},
 50: {'phi': -117.26440335376822, 'psi': 138.51305036902693},
 51: {'phi': -120.03410170277817, 'psi': 81.75707989178757},
 52: {'phi': -77.60981590398819, 'psi': 83.18451037698443},
 53: {'phi': -79.65858964180552, 'psi': 111.40143302647459},
 54: {'phi': -100.37011629225776, 'psi': 150.03395825502497},
 55: {'phi': 49.87330458406237, 'psi': 68.74199803405018},
 56: {'phi': -73.87938409722335, 'psi': -66.7355521840301},
 57: {'phi': -56.20534388077749, 'psi': -35.207843043514686},
 58: {'phi': -66.38284564180043, 'psi': -32.21866387324769},
 59: {'phi': -94.6778115344365, 'psi': 17.686140221665553},
 60: {'phi': -111.48538994784963, 'psi': -38.09776457861392},
 61: {'phi': -70.64502750557983, 'psi': -62.8582975880629},
 62: {'phi': -33.50588994671665, 'psi': -32.02546270762559},
 63: {'phi': -128.57384077349852, 'psi': 62.57927537310066},
 64: {'phi': -12.365761900396365, 'psi': 106.99327496259977},
 65: {'phi': 73.68588813063124, 'psi': 32.131558860201714},
 66: {'phi': -89.05260862755028, 'psi': -69.16778908477181},
 67: {'phi': -77.83088301001709, 'psi': -21.564910924673597},
 68: {'phi': -71.32122280651765, 'psi': -21.859413182600065},
 69: {'phi': -81.4118653034867, 'psi': -55.2935117883826},
 70: {'phi': -52.047970110313145, 'psi': -43.22593946145588},
 71: {'phi': -59.215594114973726, 'psi': -45.283196644537554},
 72: {'phi': -52.67186926130671, 'psi': -38.127901315075064},
 73: {'phi': -75.00963018964649, 'psi': -30.83999517691734},
 74: {'phi': -74.69878930178584, 'psi': -35.042954979175136},
 75: {'phi': -80.22740138668189, 'psi': -37.2721834868002},
 76: {'phi': -63.3253002341084, 'psi': -46.736848174955014},
 77: {'phi': -62.577975558265166, 'psi': -38.836376804396195},
 78: {'phi': -58.4371262613883, 'psi': -30.932534133630554},
 79: {'phi': -77.25603045197096, 'psi': -28.810984281581455},
 80: {'phi': -65.77402807318447, 'psi': -6.587861693755428},
 81: {'phi': 113.27162201541087, 'psi': -14.067924223417435},
 82: {'phi': -63.856071155072016, 'psi': 160.46313493362334},
 83: {'phi': -109.29442965951228, 'psi': 65.33016925110071},
 84: {'phi': -94.29902268445335, 'psi': 87.93029438989075},
 85: {'phi': -52.91938395571083, 'psi': 98.897475962567},
 86: {'phi': -73.44769372512917, 'psi': 114.6488125441093},
 87: {'phi': -114.16119204550668, 'psi': 101.24805765454327},
 88: {'phi': -96.78933556699712, 'psi': 106.74340425527281},
 89: {'phi': -109.02775603395975, 'psi': None}}

The dihedral angles can be provided to the pbxplore.assign() function that assigns a Protein Block to each residue, and that returns the PB sequence as a string. Note that the first and last two residues are assigned to the Z jocker block as some dihedral angles cannot be calculated.

pb_seq = pbx.assign(dihedrals)
print(pb_seq)
ZZdddfklonbfklmmmmmmmmnopafklnoiaklmmmmmnoopacddddddehkllmmmmngoilmmmmmmmmmmmmnopacdcddZZ

Assign PB for several models of a single file

A single PDB file can contain several models. Then, we do not want to read only the first chain. Instead, we want to iterate over all the chains.

for chain_name, chain in pbx.chains_from_files(['demo1_assignation/2LFU.pdb']):
    dihedrals = chain.get_phi_psi_angles()
    pb_seq = pbx.assign(dihedrals)
    print('* {}'.format(chain_name))
    print('  {}'.format(pb_seq))
Read 10 chain(s) in demo1_assignation/2LFU.pdb
* demo1_assignation/2LFU.pdb | model 1 | chain A
  ZZbghiacfkbccdddddehiadddddddddddfklggcdddddddddddddehifbdcddddddddddfklopadddddfhpamlnopcddddddehjadddddehjacbddddddddfklmaccddddddfbgniaghiapaddddddfklnoambZZ
* demo1_assignation/2LFU.pdb | model 2 | chain A
  ZZpcfblcffbccdddddeehjacdddddddddfklggcddddddddddddddfblghiadddddddddfklopadddddehpmmmnopcddddddeehiacdddfblopadcddddddfklpaccdddddfklmlmgcdehiaddddddfklmmgopZZ
* demo1_assignation/2LFU.pdb | model 3 | chain A
  ZZmgghiafbbccdddddehjbdcdddddddddfklggcddddddddddddddfbfghpacddddddddfklopadddddehiaklmmmgcdddddeehiaddddfkbgciacdddddefklpaccddddddfkgojbdfehpaddddddfkbccfbgZZ
* demo1_assignation/2LFU.pdb | model 4 | chain A
  ZZcghiacfkbacdddddfbhpacdddddddddfklmcfdddddddddddddehiacddddddddddddfknopadddddfkpamlnopaddddddehjaccdddfklnopacddddddfklmpccdddddddehiabghehiaddddddfklpccfkZZ
* demo1_assignation/2LFU.pdb | model 5 | chain A
  ZZpaehiehkaccdddddehjbccdddddddddfklggcddddddddddddddfbhpadddddddddddfklopadddddehiamlmmpccdddddeehiadddddfbacddcddddddfklmaccddddddfbgghiafehiadddddddfklpacfZZ
* demo1_assignation/2LFU.pdb | model 6 | chain A
  ZZmghbacfkbccdddddeehpacdddddddddfklggcdddddddddddddehiacadddddddddddfklopadddddehiaklnopcddddddeehiadddehjlnopacddddddfklmaccddddehiaehbgcdehiadddddddfehjlpcZZ
* demo1_assignation/2LFU.pdb | model 7 | chain A
  ZZcchbacfkbccdddddfehpacdddddddddfklggcdddddddddddddddehjapadddddddddfknopadddddfklmmmnopcddddddehjiddddddfknopacddddddfklpaccdddddfklmaacdfehpadddddehjblckknZZ
* demo1_assignation/2LFU.pdb | model 8 | chain A
  ZZcehjdeehiacdjdddedjbdcdddddddddfklggcdddddddddddddddbfblbacddddddddfklopacddddehiamlnopaddddddehjacddddfehpaaccdddddefklpaccdddddfklmbfbehehiaddddddffkgoiehZZ
* demo1_assignation/2LFU.pdb | model 9 | chain A
  ZZpccdjdfkbccdddddehhpacdddddddddfklggcdddddddddddddehiacbdcdddddddddfklopadddddehiammnopcddddddeejiadddehjlgobacddddddfklmpccddddehiacbcbdfehpadddddehjklmklmZZ
* demo1_assignation/2LFU.pdb | model 10 | chain A
  ZZccfklcfkbccdddddehjbdcdddddddddfklggcdddddddddddddehiapaccdddddddddfklopadddddehjamlnopaddddddehjddcdddfbfghpacddddddfklpaccddddddfbcfbacfehpadddddddekpghiaZZ

Assign PB for a set of structures

The pbxplore.chains_from_files() function can also handle several chains from several files.

import glob
files = glob.glob('demo1_assignation/*.pdb')
for chain_name, chain in pbx.chains_from_files(files):
    dihedrals = chain.get_phi_psi_angles()
    pb_seq = pbx.assign(dihedrals)
    print('* {}'.format(chain_name))
    print('  {}'.format(pb_seq))
Read 2 chain(s) in demo1_assignation/1AY7.pdb
* demo1_assignation/1AY7.pdb | chain A
  ZZbjadfklmcfklmmmmmmmmnnpaafbfkgopacehlnomaccddehjaccdddddehklpnbjadcdddfbehiacddfegolaccdddfkZZ
* demo1_assignation/1AY7.pdb | chain B
  ZZcddfklpcbfklmmmmmmmmnopafklgoiaklmmmmmmmmpacddddddehkllmmmmnnommmmmmmmmmmmmmnopacddddZZ
Read 1 chain(s) in demo1_assignation/1BTA.pdb
* demo1_assignation/1BTA.pdb | chain A
  ZZdddfklonbfklmmmmmmmmnopafklnoiaklmmmmmnoopacddddddehkllmmmmngoilmmmmmmmmmmmmnopacdcddZZ
Read 10 chain(s) in demo1_assignation/2LFU.pdb
* demo1_assignation/2LFU.pdb | model 1 | chain A
  ZZbghiacfkbccdddddehiadddddddddddfklggcdddddddddddddehifbdcddddddddddfklopadddddfhpamlnopcddddddehjadddddehjacbddddddddfklmaccddddddfbgniaghiapaddddddfklnoambZZ
* demo1_assignation/2LFU.pdb | model 2 | chain A
  ZZpcfblcffbccdddddeehjacdddddddddfklggcddddddddddddddfblghiadddddddddfklopadddddehpmmmnopcddddddeehiacdddfblopadcddddddfklpaccdddddfklmlmgcdehiaddddddfklmmgopZZ
* demo1_assignation/2LFU.pdb | model 3 | chain A
  ZZmgghiafbbccdddddehjbdcdddddddddfklggcddddddddddddddfbfghpacddddddddfklopadddddehiaklmmmgcdddddeehiaddddfkbgciacdddddefklpaccddddddfkgojbdfehpaddddddfkbccfbgZZ
* demo1_assignation/2LFU.pdb | model 4 | chain A
  ZZcghiacfkbacdddddfbhpacdddddddddfklmcfdddddddddddddehiacddddddddddddfknopadddddfkpamlnopaddddddehjaccdddfklnopacddddddfklmpccdddddddehiabghehiaddddddfklpccfkZZ
* demo1_assignation/2LFU.pdb | model 5 | chain A
  ZZpaehiehkaccdddddehjbccdddddddddfklggcddddddddddddddfbhpadddddddddddfklopadddddehiamlmmpccdddddeehiadddddfbacddcddddddfklmaccddddddfbgghiafehiadddddddfklpacfZZ
* demo1_assignation/2LFU.pdb | model 6 | chain A
  ZZmghbacfkbccdddddeehpacdddddddddfklggcdddddddddddddehiacadddddddddddfklopadddddehiaklnopcddddddeehiadddehjlnopacddddddfklmaccddddehiaehbgcdehiadddddddfehjlpcZZ
* demo1_assignation/2LFU.pdb | model 7 | chain A
  ZZcchbacfkbccdddddfehpacdddddddddfklggcdddddddddddddddehjapadddddddddfknopadddddfklmmmnopcddddddehjiddddddfknopacddddddfklpaccdddddfklmaacdfehpadddddehjblckknZZ
* demo1_assignation/2LFU.pdb | model 8 | chain A
  ZZcehjdeehiacdjdddedjbdcdddddddddfklggcdddddddddddddddbfblbacddddddddfklopacddddehiamlnopaddddddehjacddddfehpaaccdddddefklpaccdddddfklmbfbehehiaddddddffkgoiehZZ
* demo1_assignation/2LFU.pdb | model 9 | chain A
  ZZpccdjdfkbccdddddehhpacdddddddddfklggcdddddddddddddehiacbdcdddddddddfklopadddddehiammnopcddddddeejiadddehjlgobacddddddfklmpccddddehiacbcbdfehpadddddehjklmklmZZ
* demo1_assignation/2LFU.pdb | model 10 | chain A
  ZZccfklcfkbccdddddehjbdcdddddddddfklggcdddddddddddddehiapaccdddddddddfklopadddddehjamlnopaddddddehjddcdddfbfghpacddddddfklpaccddddddfbcfbacfehpadddddddekpghiaZZ
Read 1 chain(s) in demo1_assignation/3ICH.pdb
* demo1_assignation/3ICH.pdb | chain A
  ZZccdfbdcdddddehjbdebjcdddddfklmmmlmmmmmmmmnopnopajeopacfbdcehibacehiamnonopgocdfkbjbdcdfblmbccfbghiacdddebehiafkbccddfbdcfklgokaccfbdcfbhklmmmmmmmpccdfkopafbacddfbgcddddfbacddddZZ

Assign PB for frames in a trajectory

PB sequences can be assigned from a trajectory. To do so, we use the pbxplore.chains_from_trajectory() function that takes the path to a trajectory and the path to the corresponding topology as argument. Any file formats readable by MDAnalysis can be used. Except for its arguments, pbxplore.chains_from_trajectory() works the same as pbxplore.chains_from_files().

** Note that MDAnalysis is required to use this feature. **

trajectory = 'demo2/barstar_md_traj.xtc'
topology = 'demo2/barstar_md_traj.gro'
for chain_name, chain in pbx.chains_from_trajectory(trajectory, topology):
    dihedrals = chain.get_phi_psi_angles()
    pb_seq = pbx.assign(dihedrals)
    print('* {}'.format(chain_name))
    print('  {}'.format(pb_seq))
* demo2/barstar_md_traj.xtc | frame 0
  ZZdddfklpmbfklmmmmmmmmnopafklgoiaklmmmmmmmmpacddddddehklmmmmmoghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 1
  ZZdddfklpcbfklmmmmmmmmnopafkbghiaklmmmmmmmmpccddddddehklmmmmmcehilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 2
  ZZdddfklpcbfklmmmmmmmmnopafklgoiaklmmmmmmmmpccddddddehklmmmmnpghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 3
  ZZdddfklccbfklmmmmmmmmnopafkbghiaklmmmmmnopaccddddddehkllmmmmpghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 4
  ZZdddfklpmbfklmmmmmmmmnopafklgoiaklmmmmmnopaccddddddehklmmmmmpghjllmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 5
  ZZdddfklpmbfklmmmmmmmmnopafkbgoiaklmmmmmnopbacddddddehklmmmmmpghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 6
  ZZdddfklpmblmlmmmmmmmmnopafkbghiaklmmmmmnopbacddddddehjllmmmnoghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 7
  ZZdddfklcfblmlmmmmmmmmnopafkbgoiaklmmmmmmmmcacddddddehklmmmmnpghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 8
  ZZdddfklpgbfklmmmmmmmmnopafklgoiaklmmmmmnojaccddddddehklmmmmmpghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 9
  ZZdddfklpcbfklmmmmmmmmnopafklgoiaklmmmmmmombacddddddehkllmmmnbghilkmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 10
  ZZdddfklmmblklmmmmmmmmnopafkbgoiaklmmmmmmmmppcddddddehkllmmmmbghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 11
  ZZddefklpcbfklmmmmmmmmnopafklghiaklmmmmmnopbacddddddehklmmmmmpghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 12
  ZZdddfklcfbfklmmmmmmmmnopafklgoiaklmmmmmmmmpacddddddehklmmmmmpghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 13
  ZZdddfklpmbfklmmmmmmmmnopafklgoiaklmmmmmnopbacddddddehklmmmmmpghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 14
  ZZdddfklpmbfklmmmmmmmmnopafklghiaklmmmmmmmmpccddddddehklmmmmmpghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 15
  ZZdddfklpcbfklmmmmmmmmnopafklghiaklmmmmmmmmpccddddddehklmmmmmbghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 16
  ZZdddfklpmbfklmmmmmmmmnopafkbghiaklmmmmmnoobacddddddehkllmmmmpghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 17
  ZZdddfklombfklmmmmmmmmnopafkbghiaklmmmmmnopbacddddddehklmmmmmpghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 18
  ZZdddfklpmbfklmmmmmmmmnopafklghiaklmmmmmmmmpacddddddehklmmmmmpghilmmmmmmmmmmmmnopcdddddZZ
* demo2/barstar_md_traj.xtc | frame 19
  ZZdddfklpmblmlmmmmmmmmnopafkbghiaklmmmmmnopbacddddddehklmmmmnbghijmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 20
  ZZdddfklcfbfklmmmmmmmmnopafklgoiaklmmmmmmmmpacddddddehkllmmmmmghijklmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 21
  ZZdddfklpcbfklmmmmmmmmnopafklghiaklmmmmmmmmpccddddddehklmmmmmpghijklmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 22
  ZZdddfklpmbfklmmmmmmmmnopafkbghiaklmmmmmmoopacddddddehklmmmmmmghijmlmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 23
  ZZdddfklpfblmlmmmmmmmmnopafkbgoiaklmmmmmmombacddddddehklmmmmmmghijmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 24
  ZZdddfklpmbfklmmmmmmmmnopafkbccdfklmmmmmnopaccddddddehklmmmmmbghijklmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 25
  ZZdddfklpgbfklmmmmmmmmnopafkbccdfklmmmmmmmmpacddddddehklmmmmmpghklmmmmmmmmmmmmnopadddddZZ
* demo2/barstar_md_traj.xtc | frame 26
  ZZdddfklpmbfklmmmmmmmmnopafkbccdfklmmmmmmnmpacddddddehklmmmmmbghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 27
  ZZdddfklpcbfklmmmmmmmmnopafkbccbfklmmmmmmmmpccddddddehklmmmmmbghillmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 28
  ZZdddfklmmbfmlmmmmmmmmnopafkbccdfklmmmmmmmmpacddddddehkllmmmmpghklmmmmmmmmmmmmnopadddddZZ
* demo2/barstar_md_traj.xtc | frame 29
  ZZdddfkbpcbfklmmmmmmmmnopafkbccbfklmmmmmnombacddddddehkllmmmmoghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 30
  ZZdddfklpmbfklmmmmmmmmnopafkbccbfklmmmmmnombccddddddehjlmmmmmpghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 31
  ZZdddfklpmbfklmmmmmmmmnopafklccdfklmmmmmnopaacddddddehklmmmmmpghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 32
  ZZdddfklcfbfklmmmmmmmmnopafkbcbdfklmmmmmmmmpacddddddehklmmmmmpghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 33
  ZZdddfklpcbfmlmmmmmmmmnopafkbckbfklmmmmmmmmmccddddddehjlmmmmmoghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 34
  ZZdddfklcfbfklmmmmmmmmnopafklccdfklmmmmmmmmpccddddddehklmmmmmpghilmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 35
  ZZdddfklcfbfklmmmmmmmmnopafkbccdfklmmmmmmmmppcddddddehklmmmmmpghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 36
  ZZdddfklpgbjklmmmmmmmmnopafkbccbfklmmmmmmmmpacddddddehklmmmmmoghjlmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 37
  ZZdddfklpfbfklmmmmmmmmnopafklccbfklmmmmmnopbacddddddehklmmmmmoghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 38
  ZZdddfklccbfklmmmmmmmmnopafkbckbfklmmmmmmompacddddddehklmmmmmoghjlmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 39
  ZZdddfklpmblmlmmmmmmmmnopafkbccdfklmmmmmmmmcacddddddehklmmmmmoghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 40
  ZZdddfklcfbfmlmmmmmmmmnopafkbccbfklmmmmmmmmmccddddddehklmmmmmpghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 41
  ZZdddfklpmbfklmmmmmmmmnopafkbccdfklmmmmmmnopacddddddehklmmmmmpghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 42
  ZZdddfklpmblmlmmmmmmmmnopafkbccdfklmmmmmmgoiacddddddehklmmmmmmgoklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 43
  ZZdddfklpcbfklmmmmmmmmnopafklccdfklmmmmmmmmpacddddddehklmmmmmmghklmmmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 44
  ZZdddfklcfbfklmmmmmmmmnopafkbccdfklmmmmmmmmgccddddddehjllmmmmmghiaklmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 45
  ZZdddfklcfbfmlmmmmmmmmnopafkbccbfklmmmmmmmmcccddddddehklmmmmmpghiamlmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 46
  ZZdddfkbcfbfklmmmmmmmmnopafkbccdfklmmmmmmmmpacddddddehkllmmmmpghijmlmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 47
  ZZdddfkbcfbfklmmmmmmmmnopafkbccdfklmmmmmmmmpccddddddehklmmmmmpghijklmmmmmmmmmmnopadddddZZ
* demo2/barstar_md_traj.xtc | frame 48
  ZZdddfkbcfbfklmmmmmmmmnopafkbccdfklmmmmmmomcacddddddehkllmmmmpghilmlmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 49
  ZZdddffbcfbfklmmmmmmmmnopafkbccbfklmmmmmmnmpacddddddehklmmmmmpghijklmmmmmmmmmmnopacddddZZ
* demo2/barstar_md_traj.xtc | frame 50
  ZZdddfklgfbfklmmmmmmmmnopafklccdfklmmmmmmmmmccddddddehklmmmmmpghijmmmmmmmmmmmmnopacddddZZ

Use a different structure parser

Providing the dihedral angles can be formated as expected by pbxplore.assign(), the source of these angles does not matter. For instance, other PDB parser can be used with PBxplore.

BioPython

import Bio.PDB
import math

for model in Bio.PDB.PDBParser().get_structure("demo1_assignation/2LFU",
                                               "demo1_assignation/2LFU.pdb"):
    for chain in model:
        polypeptides = Bio.PDB.PPBuilder().build_peptides(chain)
        for poly_index, poly in enumerate(polypeptides):
            dihedral_list = poly.get_phi_psi_list()
            dihedrals = {}
            for resid, (phi, psi) in enumerate(dihedral_list, start=1):
                if not phi is None:
                    phi = 180 * phi / math.pi
                if not psi is None:
                    psi = 180 * psi / math.pi
                dihedrals[resid] = {'phi': phi, 'psi': psi}
        print(model, chain)
        pb_seq = pbx.assign(dihedrals)
        print(pb_seq)
<Model id=0> <Chain id=A>
ZZbghiacfkbccdddddehiadddddddddddfklggcdddddddddddddehifbdcddddddddddfklopadddddfhpamlnopcddddddehjadddddehjacbddddddddfklmaccddddddfbgniaghiapaddddddfklnoambZZ
<Model id=1> <Chain id=A>
ZZpcfblcffbccdddddeehjacdddddddddfklggcddddddddddddddfblghiadddddddddfklopadddddehpmmmnopcddddddeehiacdddfblopadcddddddfklpaccdddddfklmlmgcdehiaddddddfklmmgopZZ
<Model id=2> <Chain id=A>
ZZmgghiafbbccdddddehjbdcdddddddddfklggcddddddddddddddfbfghpacddddddddfklopadddddehiaklmmmgcdddddeehiaddddfkbgciacdddddefklpaccddddddfkgojbdfehpaddddddfkbccfbgZZ
<Model id=3> <Chain id=A>
ZZcghiacfkbacdddddfbhpacdddddddddfklmcfdddddddddddddehiacddddddddddddfknopadddddfkpamlnopaddddddehjaccdddfklnopacddddddfklmpccdddddddehiabghehiaddddddfklpccfkZZ
<Model id=4> <Chain id=A>
ZZpaehiehkaccdddddehjbccdddddddddfklggcddddddddddddddfbhpadddddddddddfklopadddddehiamlmmpccdddddeehiadddddfbacddcddddddfklmaccddddddfbgghiafehiadddddddfklpacfZZ
<Model id=5> <Chain id=A>
ZZmghbacfkbccdddddeehpacdddddddddfklggcdddddddddddddehiacadddddddddddfklopadddddehiaklnopcddddddeehiadddehjlnopacddddddfklmaccddddehiaehbgcdehiadddddddfehjlpcZZ
<Model id=6> <Chain id=A>
ZZcchbacfkbccdddddfehpacdddddddddfklggcdddddddddddddddehjapadddddddddfknopadddddfklmmmnopcddddddehjiddddddfknopacddddddfklpaccdddddfklmaacdfehpadddddehjblckknZZ
<Model id=7> <Chain id=A>
ZZcehjdeehiacdjdddedjbdcdddddddddfklggcdddddddddddddddbfblbacddddddddfklopacddddehiamlnopaddddddehjacddddfehpaaccdddddefklpaccdddddfklmbfbehehiaddddddffkgoiehZZ
<Model id=8> <Chain id=A>
ZZpccdjdfkbccdddddehhpacdddddddddfklggcdddddddddddddehiacbdcdddddddddfklopadddddehiammnopcddddddeejiadddehjlgobacddddddfklmpccddddehiacbcbdfehpadddddehjklmklmZZ
<Model id=9> <Chain id=A>
ZZccfklcfkbccdddddehjbdcdddddddddfklggcdddddddddddddehiapaccdddddddddfklopadddddehjamlnopaddddddehjddcdddfbfghpacddddddfklpaccddddddfbcfbacfehpadddddddekpghiaZZ