Skip to content

Commit

Permalink
respect no_altloc in more places
Browse files Browse the repository at this point in the history
  • Loading branch information
rob-miller committed Jul 2, 2022
1 parent 4812bb5 commit bbb7e0e
Show file tree
Hide file tree
Showing 2 changed files with 35 additions and 8 deletions.
33 changes: 25 additions & 8 deletions Bio/PDB/internal_coords.py
Original file line number Diff line number Diff line change
Expand Up @@ -555,6 +555,7 @@ def setAtomVw(res, atm):

def setResAtmVws(res):
for atm in res.get_atoms():
# copy not filter so ignore no_altloc
if atm.is_disordered():
for altAtom in atm.child_dict.values():
setAtomVw(res, altAtom)
Expand Down Expand Up @@ -641,7 +642,15 @@ def _peptide_check(self, prev: "Residue", curr: "Residue") -> Optional[str]:
if pNatom is None or pCAatom is None:
return "previous residue missing N or Ca"

if not Natom.is_disordered() and not pCatom.is_disordered():
if IC_Residue.no_altloc:
if Natom.is_disordered():
Natom = Natom.selected_child
if pCatom.is_disordered():
pCatom = pCatom.selected_child

if IC_Residue.no_altloc or (
not Natom.is_disordered() and not pCatom.is_disordered()
):
dc = self._atm_dist_chk(
Natom, pCatom, IC_Chain.MaxPeptideBond, self.sqMaxPeptideBond
)
Expand Down Expand Up @@ -747,7 +756,7 @@ def _set_residues(self, verbose: bool = False) -> None:
# select only not hetero or accepted hetero
if res.id[0] == " " or res.id[0] in IC_Residue.accept_resnames:
this_res: List["IC_Residue"] = []
if 2 == res.is_disordered():
if 2 == res.is_disordered() and not IC_Residue.no_altloc:
# print('disordered res:', res.is_disordered(), res)
for r in res.child_dict.values():
if self._add_residue(r, last_res, last_ord_res, verbose):
Expand Down Expand Up @@ -807,8 +816,11 @@ def setAtom(res, atm):
def setResAtms(res):
for atm in res.get_atoms():
if atm.is_disordered():
for altAtom in atm.child_dict.values():
setAtom(res, altAtom)
if IC_Residue.no_altloc:
setAtom(res, atm.selected_child)
else:
for altAtom in atm.child_dict.values():
setAtom(res, altAtom)
else:
setAtom(res, atm)

Expand Down Expand Up @@ -1018,9 +1030,12 @@ def _hedraDict2chain(
initNCaC = []
for atm in ric.residue.get_atoms(): # n.b. only few PIC spec atoms
if 2 == atm.is_disordered():
for altAtom in atm.child_dict.values():
if altAtom.coord is not None:
initNCaC.append(AtomKey(ric, altAtom))
if IC_Residue.no_altloc:
initNCaC.append(AtomKey(ric, atm.selected_child))
else:
for altAtom in atm.child_dict.values():
if altAtom.coord is not None:
initNCaC.append(AtomKey(ric, altAtom))
elif atm.coord is not None:
initNCaC.append(AtomKey(ric, atm))
if initNCaC != []:
Expand Down Expand Up @@ -3263,6 +3278,8 @@ def _pdb_atom_string(atm: Atom, cif_extend: bool = False) -> str:
override atom chain id if not None
"""
if 2 == atm.is_disordered():
if IC_Residue.no_altloc:
return IC_Residue._pdb_atom_string(atm.selected_child, cif_extend)
s = ""
for a in atm.child_dict.values():
s += IC_Residue._pdb_atom_string(a, cif_extend)
Expand Down Expand Up @@ -3511,7 +3528,7 @@ def _write_PIC(
col = 0
for a in sorted(self.residue.get_atoms()):
if 2 == a.is_disordered():
if self.alt_ids is None:
if IC_Residue.no_altloc or self.alt_ids is None:
s, col = self._write_pic_bfac(a.selected_child, s, col)
else:
for atm in a.child_dict.values():
Expand Down
10 changes: 10 additions & 0 deletions Tests/test_PDB_internal_coords.py
Original file line number Diff line number Diff line change
Expand Up @@ -53,6 +53,7 @@ class Rebuild(unittest.TestCase):
pdb_2XHE2 = PDB_parser.get_structure("2XHE", "PDB/2XHE.pdb")
pdb_1A8O = PDB_parser.get_structure("1A8O", "PDB/1A8O.pdb")
cif_3JQH = CIF_parser.get_structure("3JQH", "PDB/3JQH.cif")
cif_3JQH2 = CIF_parser.get_structure("3JQH", "PDB/3JQH.cif")
cif_4CUP = CIF_parser.get_structure("4CUP", "PDB/4CUP.cif")
cif_4CUP2 = CIF_parser.get_structure("4CUP", "PDB/4CUP.cif")
cif_4ZHL = CIF_parser.get_structure("4ZHL", "PDB/4ZHL.cif")
Expand Down Expand Up @@ -105,6 +106,7 @@ def test_rebuild_disordered_atoms_residues(self):
"""Convert disordered protein to internal coordinates and back."""
# 3jqh has both disordered residues
# and disordered atoms in ordered residues

with warnings.catch_warnings(record=True):
warnings.simplefilter("always", PDBConstructionWarning)
r = structure_rebuild_test(self.cif_3JQH, False)
Expand All @@ -118,6 +120,14 @@ def test_rebuild_disordered_atoms_residues(self):
self.assertEqual(len(r["chains"]), 1)
self.assertTrue(r["pass"])

IC_Residue.no_altloc = True
with warnings.catch_warnings(record=True):
warnings.simplefilter("always", PDBConstructionWarning)
r = structure_rebuild_test(self.cif_3JQH2, verbose=False, quick=True)
self.assertEqual(r["aCoordMatchCount"], 167, msg="no_altloc fail")
self.assertTrue(r["pass"], msg="no_altloc fail")
IC_Residue.no_altloc = False

def test_no_crosstalk(self):
"""Deep copy, change few internal coords, test nothing else changes."""
# IC_Chain.ParallelAssembleResidues = False
Expand Down

0 comments on commit bbb7e0e

Please sign in to comment.