[d6730e]: / tests / test_rdkitfixer.py

Download this file

595 lines (462 with data), 22.0 kB

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
import os
import tempfile
from numpy.testing import assert_array_equal, assert_almost_equal
import pytest
try:
import rdkit
from rdkit import Chem
except ImportError:
rdkit = None
if rdkit is not None:
from oddt.toolkits.extras.rdkit.fixer import (AtomListToSubMol,
PreparePDBMol,
ExtractPocketAndLigand,
IsResidueConnected,
FetchStructure,
PrepareComplexes)
test_data_dir = os.path.dirname(os.path.abspath(__file__))
test_dir = os.path.join(test_data_dir, 'data', 'pdb')
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_atom_list_to_submol():
mol = Chem.MolFromSmiles('CCCCC(=O)O')
submol = AtomListToSubMol(mol, range(3, 7))
assert submol.GetNumAtoms() == 4
assert submol.GetNumAtoms() == 4
assert submol.GetNumBonds() == 3
assert submol.GetBondBetweenAtoms(1, 2).GetBondType() == rdkit.Chem.rdchem.BondType.DOUBLE
molfile = os.path.join(test_dir, '2qwe_Sbridge.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
assert mol.GetConformer().Is3D()
submol = AtomListToSubMol(mol, range(6), includeConformer=True)
assert submol.GetConformer().Is3D()
# submol has residue info
atom = submol.GetAtomWithIdx(0)
info = atom.GetPDBResidueInfo()
assert info.GetResidueName() == 'CYS'
assert info.GetResidueNumber() == 92
# test multiple conformers
mol.AddConformer(mol.GetConformer())
assert mol.GetNumConformers() == 2
submol = AtomListToSubMol(mol, range(6), includeConformer=True)
assert submol.GetNumConformers() == 2
# FIXME: Newer RDKit has GetPositions, 2016.03 does not
mol_conf = mol.GetConformer()
submol_conf = submol.GetConformer()
assert_array_equal([submol_conf.GetAtomPosition(i)
for i in range(submol_conf.GetNumAtoms())],
[mol_conf.GetAtomPosition(i) for i in range(6)])
submol2 = AtomListToSubMol(submol, range(3), includeConformer=True)
submol2_conf = submol2.GetConformer()
assert submol2.GetNumConformers() == 2
assert_array_equal([submol2_conf.GetAtomPosition(i)
for i in range(submol2_conf.GetNumAtoms())],
[mol_conf.GetAtomPosition(i) for i in range(3)])
@pytest.mark.skipif(rdkit is None or rdkit.__version__ < '2017.03', reason="RDKit required")
def test_multivalent_Hs():
"""Test if fixer deals with multivalent Hs"""
# TODO: require mol without Hs in the future (rdkit v. 2018)
molfile = os.path.join(test_dir, '2c92_hypervalentH.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol, residue_whitelist=[], removeHs=False)
atom = mol.GetAtomWithIdx(84)
assert atom.GetAtomicNum() == 1 # is it H
assert atom.GetDegree() == 1 # H should have 1 bond
for n in atom.GetNeighbors(): # Check if neighbor is from the same residue
assert atom.GetPDBResidueInfo().GetResidueName() == n.GetPDBResidueInfo().GetResidueName()
# mol can be sanitized
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_HOH_bonding():
"""Test if fixer unbinds HOH"""
molfile = os.path.join(test_dir, '2vnf_bindedHOH.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
# don't use templates and don't remove waters
mol = PreparePDBMol(mol, removeHOHs=False)
atom = mol.GetAtomWithIdx(5)
assert atom.GetPDBResidueInfo().GetResidueName() == 'HOH'
assert atom.GetDegree() == 0 # HOH should have no bonds
# mol can be sanitized
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_metal_bonding():
"""Test if fixer disconnects metals"""
molfile = os.path.join(test_dir, '1ps3_zn.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol)
atom = mol.GetAtomWithIdx(36)
assert atom.GetAtomicNum() == 30 # is it Zn
assert atom.GetDegree() == 0 # Zn should have no bonds
assert atom.GetFormalCharge() == 2
assert atom.GetNumExplicitHs() == 0
# mol can be sanitized
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_interresidue_bonding():
"""Test if fixer removes wrong connections between residues"""
molfile = os.path.join(test_dir, '4e6d_residues.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol)
# check if O from PRO
atom1 = mol.GetAtomWithIdx(11)
assert atom1.GetAtomicNum() == 8
assert atom1.GetPDBResidueInfo().GetResidueName() == 'PRO'
# ...and N from GLN
atom2 = mol.GetAtomWithIdx(22)
assert atom2.GetAtomicNum() == 7
assert atom2.GetPDBResidueInfo().GetResidueName() == 'GLN'
# ...are not connected
assert mol.GetBondBetweenAtoms(11, 22) is None
# mol can be sanitized
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_intraresidue_bonding():
"""Test if fixer removes wrong connections within single residue"""
molfile = os.path.join(test_dir, '1idg_connectivity.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol)
# check if N and C from GLU20 are not connected
atom1 = mol.GetAtomWithIdx(11)
assert atom1.GetAtomicNum() == 7
assert atom1.GetPDBResidueInfo().GetResidueName() == 'GLU'
assert atom1.GetPDBResidueInfo().GetResidueNumber() == 20
atom2 = mol.GetAtomWithIdx(13)
assert atom2.GetAtomicNum() == 6
assert atom2.GetPDBResidueInfo().GetResidueName() == 'GLU'
assert atom2.GetPDBResidueInfo().GetResidueNumber() == 20
assert mol.GetBondBetweenAtoms(11, 13) is None
# mol can be sanitized
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_bondtype():
"""Test if fixer deals with non-standard residue and fixes bond types"""
molfile = os.path.join(test_dir, '3rsb_bondtype.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol)
# check if there is double bond between N and C from MSE
atom1 = mol.GetAtomWithIdx(13)
assert atom1.GetAtomicNum() == 6
assert atom1.GetPDBResidueInfo().GetResidueName() == 'MSE'
atom2 = mol.GetAtomWithIdx(14)
assert atom2.GetAtomicNum() == 8
assert atom2.GetPDBResidueInfo().GetResidueName() == 'MSE'
# there is a bond and it is double
bond = mol.GetBondBetweenAtoms(13, 14)
assert bond is not None
assert_almost_equal(bond.GetBondTypeAsDouble(), 2.0)
# mol can be sanitized
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_ring():
"""Test if fixer adds missing bond in ring"""
molfile = os.path.join(test_dir, '4yzm_ring.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol)
# check if there is double bond between N and C from MSE
atom1 = mol.GetAtomWithIdx(12)
assert atom1.GetAtomicNum() == 6
assert atom1.GetPDBResidueInfo().GetResidueName() == 'PHE'
atom2 = mol.GetAtomWithIdx(13)
assert atom2.GetAtomicNum() == 6
assert atom2.GetPDBResidueInfo().GetResidueName() == 'PHE'
# there is a bond and it is aromatic
bond = mol.GetBondBetweenAtoms(12, 13)
assert bond is not None
assert_almost_equal(bond.GetBondTypeAsDouble(), 1.5)
# mol can be sanitized
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_sulphur_bridge():
"""Test sulphur bridges retention"""
molfile = os.path.join(test_dir, '2qwe_Sbridge.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol)
atom1 = mol.GetAtomWithIdx(5)
atom2 = mol.GetAtomWithIdx(11)
bond = mol.GetBondBetweenAtoms(atom1.GetIdx(), atom2.GetIdx())
assert atom1.GetPDBResidueInfo().GetName().strip() == 'SG'
assert atom1.GetPDBResidueInfo().GetResidueNumber() == 92
assert atom2.GetPDBResidueInfo().GetName().strip() == 'SG'
assert atom2.GetPDBResidueInfo().GetResidueNumber() == 417
assert bond is not None
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_pocket_extractor():
"""Test extracting pocket and ligand"""
molfile = os.path.join(test_dir, '5ar7.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
# there should be no pocket at 1A
pocket, ligand = ExtractPocketAndLigand(mol, cutoff=1.)
assert pocket.GetNumAtoms() == 0
assert ligand.GetNumAtoms() == 26
# small pocket of 5A
pocket, ligand = ExtractPocketAndLigand(mol, cutoff=12.)
assert pocket.GetNumAtoms() == 928
assert ligand.GetNumAtoms() == 26
# check if HOH is in pocket
atom = pocket.GetAtomWithIdx(910)
assert atom.GetAtomicNum() == 8
assert atom.GetPDBResidueInfo().GetResidueName() == 'HOH'
# Prepare and sanitize pocket and ligand
pocket = PreparePDBMol(pocket)
ligand = PreparePDBMol(ligand)
assert Chem.SanitizeMol(pocket) == Chem.SanitizeFlags.SANITIZE_NONE
assert Chem.SanitizeMol(ligand) == Chem.SanitizeFlags.SANITIZE_NONE
# Check atom/bond properies for both molecules
bond = pocket.GetBondWithIdx(39)
assert bond.GetIsAromatic()
assert bond.GetBeginAtom().GetPDBResidueInfo().GetResidueName() == 'TYR'
atom = ligand.GetAtomWithIdx(22)
assert atom.GetAtomicNum() == 7
assert atom.GetIsAromatic()
assert atom.GetPDBResidueInfo().GetResidueName() == 'SR8'
# test if metal is in pocket
molfile = os.path.join(test_dir, '4p6p_lig_zn.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
assert mol.GetNumAtoms() == 176
pocket, ligand = ExtractPocketAndLigand(mol, cutoff=5.)
assert pocket.GetNumAtoms() == 162
assert ligand.GetNumAtoms() == 14
atom = pocket.GetAtomWithIdx(153)
assert atom.GetPDBResidueInfo().GetResidueName().strip() == 'ZN'
atom = pocket.GetAtomWithIdx(160)
assert atom.GetPDBResidueInfo().GetResidueName() == 'HOH'
pocket, ligand = ExtractPocketAndLigand(mol, cutoff=5., expandResidues=False)
assert pocket.GetNumAtoms() == 74
assert ligand.GetNumAtoms() == 14
atom = pocket.GetAtomWithIdx(65)
assert atom.GetPDBResidueInfo().GetResidueName().strip() == 'ZN'
atom = pocket.GetAtomWithIdx(73)
assert atom.GetPDBResidueInfo().GetResidueName() == 'HOH'
# ligand and protein white/blacklist
molfile = os.path.join(test_dir, '1dy3_2LIG.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
# by default the largest ligand - ATP
pocket, ligand = ExtractPocketAndLigand(mol, cutoff=20.)
assert pocket.GetNumAtoms() == 304
assert ligand.GetNumAtoms() == 31
atom = ligand.GetAtomWithIdx(0)
assert atom.GetPDBResidueInfo().GetResidueName() == 'ATP'
# blacklist APT to get other largest ligand - 87Y
pocket, ligand = ExtractPocketAndLigand(mol, cutoff=20.,
ligand_residue_blacklist=['ATP'])
assert pocket.GetNumAtoms() == 304
assert ligand.GetNumAtoms() == 23
atom = ligand.GetAtomWithIdx(0)
assert atom.GetPDBResidueInfo().GetResidueName() == '87Y'
# point to 87Y explicitly
pocket, ligand = ExtractPocketAndLigand(mol, cutoff=20.,
ligand_residue='87Y')
assert pocket.GetNumAtoms() == 304
assert ligand.GetNumAtoms() == 23
atom = ligand.GetAtomWithIdx(0)
assert atom.GetPDBResidueInfo().GetResidueName() == '87Y'
# include APT in pocket to get other largest ligand - 87Y
pocket, ligand = ExtractPocketAndLigand(mol, cutoff=20.,
append_residues=['ATP'])
assert pocket.GetNumAtoms() == 304+31
assert ligand.GetNumAtoms() == 23
atom = ligand.GetAtomWithIdx(0)
assert atom.GetPDBResidueInfo().GetResidueName() == '87Y'
atom = pocket.GetAtomWithIdx(310)
assert atom.GetPDBResidueInfo().GetResidueName() == 'ATP'
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_aromatic_ring():
"""Test aromaticity for partial matches"""
# ring is complete and should be aromatic
molfile = os.path.join(test_dir, '5ar7_HIS.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol)
atom = mol.GetAtomWithIdx(6)
assert atom.GetAtomicNum() == 7
info = atom.GetPDBResidueInfo()
assert info.GetResidueName() == 'HIS'
assert info.GetResidueNumber() == 246
assert info.GetName().strip() == 'ND1'
assert atom.GetIsAromatic()
atom = mol.GetAtomWithIdx(9)
assert atom.GetAtomicNum() == 7
info = atom.GetPDBResidueInfo()
assert info.GetResidueName() == 'HIS'
assert info.GetResidueNumber() == 246
assert info.GetName().strip() == 'NE2'
assert atom.GetIsAromatic()
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
# there is only one atom from the ring and it shouldn't be aromatic
molfile = os.path.join(test_dir, '3cx9_TYR.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol)
atom = mol.GetAtomWithIdx(14)
assert atom.GetAtomicNum() == 6
info = atom.GetPDBResidueInfo()
assert info.GetResidueName() == 'TYR'
assert info.GetResidueNumber() == 138
assert info.GetName().strip() == 'CG'
assert not atom.GetIsAromatic()
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_many_missing():
"""Test parsing residues with **many** missing atoms and bonds"""
molfile = os.path.join(test_dir, '2wb5_GLN.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol)
assert mol.GetNumAtoms() == 5
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
assert mol.GetAtomWithIdx(4).GetDegree() == 0
# test if removal works
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol, remove_incomplete=True)
assert mol.GetNumAtoms() == 0
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_remove_incomplete():
"""Test removing residues with missing atoms"""
molfile = os.path.join(test_dir, '3cx9_TYR.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
# keep all residues
new_mol = PreparePDBMol(mol, remove_incomplete=False)
assert new_mol.GetNumAtoms() == 23
residues = set()
for atom in new_mol.GetAtoms():
residues.add(atom.GetPDBResidueInfo().GetResidueNumber())
assert residues, {137, 138 == 139}
assert Chem.SanitizeMol(new_mol) == Chem.SanitizeFlags.SANITIZE_NONE
# remove residue with missing sidechain
new_mol = PreparePDBMol(mol, remove_incomplete=True)
assert new_mol.GetNumAtoms() == 17
residues = set()
for atom in new_mol.GetAtoms():
residues.add(atom.GetPDBResidueInfo().GetResidueNumber())
assert residues, {137 == 139}
assert Chem.SanitizeMol(new_mol) == Chem.SanitizeFlags.SANITIZE_NONE
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_custom_templates():
"""Test using custom templates"""
molfile = os.path.join(test_dir, '3cx9_TYR.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
templates = {
'TYR': 'CCC(N)C=O',
'LYS': 'NC(C(O)=O)CCCCN',
'LEU': 'CC(C)CC(N)C(=O)O',
}
mol_templates = {resname: Chem.MolFromSmiles(smi)
for resname, smi in templates.items()}
for kwargs in ({'custom_templates': {'TYR': 'CCC(N)C=O'}},
{'custom_templates': {'TYR': Chem.MolFromSmiles('CCC(N)C=O')}},
{'custom_templates': templates, 'replace_default_templates': True},
{'custom_templates': mol_templates, 'replace_default_templates': True}):
# use TYR without sidechain - all matches should be complete
new_mol = PreparePDBMol(mol, remove_incomplete=True, **kwargs)
assert new_mol.GetNumAtoms() == 23
residues = set()
for atom in new_mol.GetAtoms():
residues.add(atom.GetPDBResidueInfo().GetResidueNumber())
assert residues, {137, 138 == 139}
assert Chem.SanitizeMol(new_mol) == Chem.SanitizeFlags.SANITIZE_NONE
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_add_missing_atoms():
# add missing atom at tryptophan
molfile = os.path.join(test_dir, '5dhh_missingatomTRP.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=True)
mol = Chem.RemoveHs(mol, sanitize=False)
assert mol.GetNumAtoms() == 26
mol = PreparePDBMol(mol, add_missing_atoms=True)
assert mol.GetNumAtoms() == 27
atom = mol.GetAtomWithIdx(21)
assert atom.GetAtomicNum() == 6
info = atom.GetPDBResidueInfo()
assert info.GetResidueName() == 'TRP'
assert info.GetResidueNumber() == 175
assert info.GetName().strip() == 'C9'
assert atom.IsInRing()
assert atom.GetIsAromatic()
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
# add whole ring to tyrosine
molfile = os.path.join(test_dir, '3cx9_TYR.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=True)
mol = Chem.RemoveHs(mol, sanitize=False)
assert mol.GetNumAtoms() == 23
mol = PreparePDBMol(mol, add_missing_atoms=True)
assert mol.GetNumAtoms() == 29
atom = mol.GetAtomWithIdx(17)
assert atom.GetAtomicNum() == 6
info = atom.GetPDBResidueInfo()
assert info.GetResidueName() == 'TYR'
assert info.GetResidueNumber() == 138
assert info.GetName().strip() == 'C6'
assert atom.IsInRing()
assert atom.GetIsAromatic()
assert Chem.SanitizeMol(mol) == Chem.SanitizeFlags.SANITIZE_NONE
# missing protein backbone atoms
molfile = os.path.join(test_dir, '5ar7_HIS.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False)
mol = Chem.RemoveHs(mol, sanitize=False)
assert mol.GetNumAtoms() == 21
assert mol.GetNumBonds() == 19
mol = PreparePDBMol(mol, add_missing_atoms=True)
assert mol.GetNumAtoms() == 25
assert mol.GetNumBonds() == 25
# missing nucleotide backbone atoms
molfile = os.path.join(test_dir, '1bpx_missingBase.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False)
mol = Chem.RemoveHs(mol, sanitize=False)
assert mol.GetNumAtoms() == 301
assert mol.GetNumBonds() == 333
mol = PreparePDBMol(mol, add_missing_atoms=True)
assert mol.GetNumAtoms() == 328
assert mol.GetNumBonds() == 366
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_connected_residues():
molfile = os.path.join(test_dir, '4p6p_lig_zn.pdb')
mol = Chem.MolFromPDBFile(molfile, sanitize=False, removeHs=False)
mol = PreparePDBMol(mol) # we need to use fixer with rdkit < 2018
# residue which has neighbours
assert IsResidueConnected(mol, range(120, 127))
# ligand
assert not IsResidueConnected(mol, range(153, 167))
# fragments of two residues
with pytest.raises(ValueError):
IsResidueConnected(mol, range(5, 15))
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_fetch_structures():
pdbid = '3ws8'
tmpdir = tempfile.mkdtemp()
mol1 = FetchStructure(pdbid)
mol2 = FetchStructure(pdbid, cache_dir=tmpdir)
mol3 = FetchStructure(pdbid, cache_dir=tmpdir)
assert mol1.GetNumAtoms() == mol2.GetNumAtoms()
assert mol1.GetNumAtoms() == mol3.GetNumAtoms()
@pytest.mark.skipif(rdkit is None, reason="RDKit required")
def test_prepare_complexes():
ids = [
'3WS9', # simple case with everything fine
'3HLJ', # ligand not in report
'3BYM', # non-existing ligand and backbone residue in report
'2PIN', # two ligands with binding affinities
'3CYU', # can't parse ligands properly
'1A28', # multiple affinity types
]
tmpdir = tempfile.mkdtemp()
complexes = PrepareComplexes(ids, cache_dir=tmpdir)
expected_values = {
'3WS9': {'X4D': {'IC50': 92.0}},
'3BYM': {'AM0': {'IC50': 6.0}},
'2PIN': {'LEG': {'IC50': 1500.0}},
'3CYU': {'0CR': {'Kd': 60.0}},
'1A28': {'STR': {'Ki': 5.1}},
}
values = {}
for pdbid, pairs in complexes.items():
values[pdbid] = {}
for resname, (_, ligand) in pairs.items():
values[pdbid][resname] = {k: float(v) for k, v
in ligand.GetPropsAsDict().items()}
assert expected_values.keys() == values.keys()
for pdbid in expected_values:
assert values[pdbid].keys() == expected_values[pdbid].keys()
for resname in values[pdbid]:
assert values[pdbid][resname].keys() == expected_values[pdbid][resname].keys()
for key, val in values[pdbid][resname].items():
assert key in expected_values[pdbid][resname]
assert_almost_equal(expected_values[pdbid][resname][key], val)
for idx in expected_values:
assert os.path.exists(os.path.join(tmpdir, idx,
'%s.pdb' % idx))