/
docking_utils.py
308 lines (258 loc) · 9.49 KB
/
docking_utils.py
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
"""
This file contains utilities for molecular docking.
"""
from typing import List, Optional, Tuple
import os
import numpy as np
from deepchem.utils.typing import RDKitMol
from deepchem.utils.pdbqt_utils import pdbqt_to_pdb
def write_vina_conf(protein_filename: str,
ligand_filename: str,
centroid: np.ndarray,
box_dims: np.ndarray,
conf_filename: str,
num_modes: int = 9,
exhaustiveness: int = None) -> None:
"""Writes Vina configuration file to disk.
Autodock Vina accepts a configuration file which provides options
under which Vina is invoked. This utility function writes a vina
configuration file which directs Autodock vina to perform docking
under the provided options.
Parameters
----------
protein_filename: str
Filename for protein
ligand_filename: str
Filename for the ligand
centroid: np.ndarray
A numpy array with shape `(3,)` holding centroid of system
box_dims: np.ndarray
A numpy array of shape `(3,)` holding the size of the box to dock
conf_filename: str
Filename to write Autodock Vina configuration to.
num_modes: int, optional (default 9)
The number of binding modes Autodock Vina should find
exhaustiveness: int, optional
The exhaustiveness of the search to be performed by Vina
"""
with open(conf_filename, "w") as f:
f.write("receptor = %s\n" % protein_filename)
f.write("ligand = %s\n\n" % ligand_filename)
f.write("center_x = %f\n" % centroid[0])
f.write("center_y = %f\n" % centroid[1])
f.write("center_z = %f\n\n" % centroid[2])
f.write("size_x = %f\n" % box_dims[0])
f.write("size_y = %f\n" % box_dims[1])
f.write("size_z = %f\n\n" % box_dims[2])
f.write("num_modes = %d\n\n" % num_modes)
if exhaustiveness is not None:
f.write("exhaustiveness = %d\n" % exhaustiveness)
def write_gnina_conf(protein_filename: str,
ligand_filename: str,
conf_filename: str,
num_modes: int = 9,
exhaustiveness: int = None,
**kwargs) -> None:
"""Writes GNINA configuration file to disk.
GNINA accepts a configuration file which provides options
under which GNINA is invoked. This utility function writes a
configuration file which directs GNINA to perform docking
under the provided options.
Parameters
----------
protein_filename: str
Filename for protein
ligand_filename: str
Filename for the ligand
conf_filename: str
Filename to write Autodock Vina configuration to.
num_modes: int, optional (default 9)
The number of binding modes GNINA should find
exhaustiveness: int, optional
The exhaustiveness of the search to be performed by GNINA
kwargs:
Args supported by GNINA documented here
https://github.com/gnina/gnina#usage
"""
with open(conf_filename, "w") as f:
f.write("receptor = %s\n" % protein_filename)
f.write("ligand = %s\n\n" % ligand_filename)
f.write("autobox_ligand = %s\n\n" % protein_filename)
if exhaustiveness is not None:
f.write("exhaustiveness = %d\n" % exhaustiveness)
f.write("num_modes = %d\n\n" % num_modes)
for k, v in kwargs.items():
f.write("%s = %s\n" % (str(k), str(v)))
def read_gnina_log(log_file: str) -> np.array:
"""Read GNINA logfile and get docking scores.
GNINA writes computed binding affinities to a logfile.
Parameters
----------
log_file: str
Filename of logfile generated by GNINA.
Returns
-------
scores: np.array, dimension (num_modes, 3)
Array of binding affinity (kcal/mol), CNN pose score,
and CNN affinity for each binding mode.
"""
scores = []
lines = open(log_file).readlines()
mode_start = np.inf
for idx, line in enumerate(lines):
if line[:6] == '-----+':
mode_start = idx
if idx > mode_start:
mode = line.split()
score = [float(x) for x in mode[1:]]
scores.append(score)
scores = np.array(scores)
return scores
def load_docked_ligands(
pdbqt_output: str) -> Tuple[List[RDKitMol], List[float]]:
"""This function loads ligands docked by autodock vina.
Autodock vina writes outputs to disk in a PDBQT file format. This
PDBQT file can contain multiple docked "poses". Recall that a pose
is an energetically favorable 3D conformation of a molecule. This
utility function reads and loads the structures for multiple poses
from vina's output file.
Parameters
----------
pdbqt_output: str
Should be the filename of a file generated by autodock vina's
docking software.
Returns
-------
Tuple[List[rdkit.Chem.rdchem.Mol], List[float]]
Tuple of `molecules, scores`. `molecules` is a list of rdkit
molecules with 3D information. `scores` is the associated vina
score.
Notes
-----
This function requires RDKit to be installed.
"""
try:
from rdkit import Chem
except ModuleNotFoundError:
raise ImportError("This function requires RDKit to be installed.")
lines = open(pdbqt_output).readlines()
molecule_pdbqts = []
scores = []
current_pdbqt: Optional[List[str]] = None
for line in lines:
if line[:5] == "MODEL":
current_pdbqt = []
elif line[:19] == "REMARK VINA RESULT:":
words = line.split()
# the line has format
# REMARK VINA RESULT: score ...
# There is only 1 such line per model so we can append it
scores.append(float(words[3]))
elif line[:6] == "ENDMDL":
molecule_pdbqts.append(current_pdbqt)
current_pdbqt = None
else:
# FIXME: Item "None" of "Optional[List[str]]" has no attribute "append"
current_pdbqt.append(line) # type: ignore
molecules = []
for pdbqt_data in molecule_pdbqts:
pdb_block = pdbqt_to_pdb(pdbqt_data=pdbqt_data)
mol = Chem.MolFromPDBBlock(str(pdb_block), sanitize=False, removeHs=False)
molecules.append(mol)
return molecules, scores
def prepare_inputs(protein: str,
ligand: str,
replace_nonstandard_residues: bool = True,
remove_heterogens: bool = True,
remove_water: bool = True,
add_hydrogens: bool = True,
pH: float = 7.0,
optimize_ligand: bool = True,
pdb_name: Optional[str] = None) -> Tuple[RDKitMol, RDKitMol]:
"""This prepares protein-ligand complexes for docking.
Autodock Vina requires PDB files for proteins and ligands with
sensible inputs. This function uses PDBFixer and RDKit to ensure
that inputs are reasonable and ready for docking. Default values
are given for convenience, but fixing PDB files is complicated and
human judgement is required to produce protein structures suitable
for docking. Always inspect the results carefully before trying to
perform docking.
Parameters
----------
protein: str
Filename for protein PDB file or a PDBID.
ligand: str
Either a filename for a ligand PDB file or a SMILES string.
replace_nonstandard_residues: bool (default True)
Replace nonstandard residues with standard residues.
remove_heterogens: bool (default True)
Removes residues that are not standard amino acids or nucleotides.
remove_water: bool (default True)
Remove water molecules.
add_hydrogens: bool (default True)
Add missing hydrogens at the protonation state given by `pH`.
pH: float (default 7.0)
Most common form of each residue at given `pH` value is used.
optimize_ligand: bool (default True)
If True, optimize ligand with RDKit. Required for SMILES inputs.
pdb_name: Optional[str]
If given, write sanitized protein and ligand to files called
"pdb_name.pdb" and "ligand_pdb_name.pdb"
Returns
-------
Tuple[RDKitMol, RDKitMol]
Tuple of `protein_molecule, ligand_molecule` with 3D information.
Note
----
This function requires RDKit and OpenMM to be installed.
Read more about PDBFixer here: https://github.com/openmm/pdbfixer.
Examples
--------
>>> p, m = prepare_inputs('3cyx', 'CCC')
>>> p.GetNumAtoms()
1415
>>> m.GetNumAtoms()
11
>>> p, m = prepare_inputs('3cyx', 'CCC', remove_heterogens=False)
>>> p.GetNumAtoms()
1720
"""
try:
from rdkit import Chem
from rdkit.Chem import AllChem
from pdbfixer import PDBFixer
from simtk.openmm.app import PDBFile
except ModuleNotFoundError:
raise ImportError(
"This function requires RDKit and OpenMM to be installed.")
if protein.endswith('.pdb'):
fixer = PDBFixer(protein)
else:
fixer = PDBFixer(url='https://files.rcsb.org/download/%s.pdb' % (protein))
if ligand.endswith('.pdb'):
m = Chem.MolFromPDBFile(ligand)
else:
m = Chem.MolFromSmiles(ligand, sanitize=True)
# Apply common fixes to PDB files
if replace_nonstandard_residues:
fixer.findMissingResidues()
fixer.findNonstandardResidues()
fixer.replaceNonstandardResidues()
if remove_heterogens and not remove_water:
fixer.removeHeterogens(True)
if remove_heterogens and remove_water:
fixer.removeHeterogens(False)
if add_hydrogens:
fixer.addMissingHydrogens(pH)
PDBFile.writeFile(fixer.topology, fixer.positions, open('tmp.pdb', 'w'))
p = Chem.MolFromPDBFile('tmp.pdb', sanitize=True)
os.remove('tmp.pdb')
# Optimize ligand
if optimize_ligand:
m = Chem.AddHs(m) # need hydrogens for optimization
AllChem.EmbedMolecule(m)
AllChem.MMFFOptimizeMolecule(m)
if pdb_name:
Chem.rdmolfiles.MolToPDBFile(p, '%s.pdb' % (pdb_name))
Chem.rdmolfiles.MolToPDBFile(m, 'ligand_%s.pdb' % (pdb_name))
return (p, m)