-
Notifications
You must be signed in to change notification settings - Fork 6
/
system.py
431 lines (380 loc) · 15.5 KB
/
system.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
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
import logging
import os
import re
from collections import defaultdict
from typing import Tuple
import networkx as nx
import parmed as pm
from rdkit import Chem
from transformato.utils import get_toppar_dir
logger = logging.getLogger(__name__)
class SystemStructure(object):
def __init__(self, configuration: dict, structure: str):
"""
A class that contains all informations for a single small ligand in different environments.
Everything is constructed automatically from the charmm-gui folder.
Parameters
----------
configuration: dict
the configuration dictionary obtained with utils.load_config_yaml
structure: str
either 'structure1' or 'structure2'
"""
self.structure: str = structure
self.name: str = configuration["system"][structure]["name"]
self.tlc: str = configuration["system"][structure]["tlc"]
self.charmm_gui_base: str = configuration["system"][structure]["charmm_gui_dir"]
self.psfs: defaultdict = defaultdict(pm.charmm.CharmmPsfFile)
self.offset: defaultdict = defaultdict(int)
self.parameter = self._read_parameters("waterbox")
self.cgenff_version: float
self.envs = set()
# running a binding-free energy calculation?
if configuration["simulation"]["free-energy-type"] == "rbfe":
self.envs = set(["complex", "waterbox"])
for env in self.envs:
parameter = self._read_parameters(env)
# set up system
self.psfs[env] = self._initialize_system(configuration, env)
# load parameters
self.psfs[env].load_parameters(parameter)
# get offset
self.offset[
env
] = self._determine_offset_and_set_possible_dummy_properties(
self.psfs[env]
)
# generate rdkit mol object of small molecule
self.mol: Chem.Mol = self._generate_rdkit_mol(
"complex", self.psfs["complex"][f":{self.tlc}"]
)
self.graph: nx.Graph = self.mol_to_nx(self.mol)
elif (
configuration["simulation"]["free-energy-type"] == "rsfe"
or configuration["simulation"]["free-energy-type"] == "asfe"
):
self.envs = set(["waterbox", "vacuum"])
for env in self.envs:
parameter = self._read_parameters(env)
# set up system
self.psfs[env] = self._initialize_system(configuration, env)
# load parameters
self.psfs[env].load_parameters(parameter)
# get offset
self.offset[
env
] = self._determine_offset_and_set_possible_dummy_properties(
self.psfs[env]
)
# generate rdkit mol object of small molecule
self.mol: Chem.Mol = self._generate_rdkit_mol(
"waterbox", self.psfs["waterbox"][f":{self.tlc}"]
)
self.graph: nx.Graph = self.mol_to_nx(self.mol)
else:
raise NotImplementedError(
"only binding and solvation free energy implemented."
)
def _set_hmr(self, configuration: dict, env: str):
# check if HMR is set
if configuration["simulation"].get("HMR", False):
logger.info("Using HMR.")
for bond in self.psfs[env].bonds:
# Only take the ones with at least one hydrogen
atom1, atom2 = bond.atom1, bond.atom2
# print(dir(atom1))
if atom1.atomic_number != 1 and atom2.atomic_number != 1:
continue
if "TIP3" in [atom1.residue.name, atom2.residue.name]:
continue
if atom2.atomic_number == 1:
atom1, atom2 = (
atom2,
atom1,
) # now atom1 is hydrogen for sure
if atom2.atomic_number != 1:
mass1 = atom1.mass
mass2 = atom2.mass
new_mass1 = mass1 * 3.0
new_mass2 = mass2 - new_mass1 + mass1
atom1.mass = new_mass1
atom2.mass = new_mass2
@staticmethod
def mol_to_nx(mol: Chem.Mol):
try:
from tf_routes import preprocessing
return preprocessing._mol_to_nx_full_weight(mol)
except ModuleNotFoundError:
G = nx.Graph()
for atom in mol.GetAtoms():
G.add_node(
atom.GetIdx(),
atomic_num=atom.GetAtomicNum(),
formal_charge=atom.GetFormalCharge(),
chiral_tag=atom.GetChiralTag(),
hybridization=atom.GetHybridization(),
num_explicit_hs=atom.GetNumExplicitHs(),
is_aromatic=atom.GetIsAromatic(),
)
for bond in mol.GetBonds():
G.add_edge(
bond.GetBeginAtomIdx(),
bond.GetEndAtomIdx(),
bond_type=bond.GetBondType(),
)
return G
def _read_parameters(self, env: str) -> pm.charmm.CharmmParameterSet:
"""
Reads in topparameters from a toppar dir and ligand specific parameters.
Parameters
----------
env: str
waterbox,complex or vacuum
Returns
----------
parameter : pm.charmm.CharmmParameterSet
parameters obtained from the CHARMM-GUI output dir.
"""
# the parameters for the vacuum system is parsed from the waterbox charmm-gui directory
if env == "vacuum":
env = "waterbox"
charmm_gui_env = self.charmm_gui_base + env
tlc = self.tlc
tlc_lower = str(tlc).lower()
toppar_dir = f"{charmm_gui_env}/toppar"
if os.path.isdir(toppar_dir):
pass
else:
toppar_dir = get_toppar_dir()
# if custom parameter are added they are located in l1,l2
parameter_files = tuple()
l1 = f"{charmm_gui_env}/{tlc_lower}/{tlc_lower}.rtf"
l2 = f"{charmm_gui_env}/{tlc_lower}/{tlc_lower}.prm"
l3 = f"{charmm_gui_env}/{tlc_lower}/{tlc_lower}.str"
for file_path in [l1, l2, l3]:
if os.path.isfile(file_path):
parameter_files += (file_path,)
else:
logger.critical(
f"Custom ligand parameters are not present in {file_path}"
)
# check cgenff versions
if parameter_files:
with open(parameter_files[0]) as f:
_ = f.readline()
cgenff = f.readline().rstrip()
logger.info(f"CGenFF version: {cgenff}")
cgenff_version = re.findall("\d+\.\d+", cgenff)[0]
self.cgenff_version = float(cgenff_version)
parameter_files += (f"{toppar_dir}/top_all36_prot.rtf",)
parameter_files += (f"{toppar_dir}/par_all36m_prot.prm",)
parameter_files += (f"{toppar_dir}/par_all36_na.prm",)
parameter_files += (f"{toppar_dir}/top_all36_na.rtf",)
parameter_files += (f"{toppar_dir}/top_all36_cgenff.rtf",)
parameter_files += (f"{toppar_dir}/par_all36_cgenff.prm",)
parameter_files += (f"{toppar_dir}/par_all36_lipid.prm",)
parameter_files += (f"{toppar_dir}/top_all36_lipid.rtf",)
parameter_files += (f"{toppar_dir}/toppar_water_ions.str",)
parameter_files += (
f"{toppar_dir}/toppar_all36_prot_na_combined.str",
) # if modified aminoacids are needed
if os.path.isfile(f"{toppar_dir}/toppar_all36_moreions.str"):
parameter_files += (f"{toppar_dir}/toppar_all36_moreions.str",)
# set up parameter objec
parameter = pm.charmm.CharmmParameterSet(*parameter_files)
return parameter
def _initialize_system(
self, configuration: dict, env: str
) -> pm.charmm.CharmmPsfFile:
"""
Generates the psf file and sets the coordinates from the CHARMM-GUI files.
Parameters
----------
configuration: dict
the configuration dictionary obtained with utils.load_config_yaml
env: str
waterbox,complex or vacuum
Returns
----------
psf : pm.charmm.CharmmPsfFile
"""
if env == "vacuum":
# take the structures from the waterbox system and extract only the ligand
taken_from = "waterbox"
psf_file_name = configuration["system"][self.structure][taken_from][
"psf_file_name"
]
crd_file_name = configuration["system"][self.structure][taken_from][
"crd_file_name"
]
psf_file_path = (
f"{self.charmm_gui_base}/{taken_from}/openmm/{psf_file_name}.psf"
)
crd_file_path = (
f"{self.charmm_gui_base}/{taken_from}/openmm/{crd_file_name}.crd"
)
# load psf
psf = pm.charmm.CharmmPsfFile(psf_file_path)
coord = pm.charmm.CharmmCrdFile(crd_file_path)
psf.coordinates = coord.coordinates
lp = False
for atom in psf.atoms:
if hasattr(atom, "frame_type"):
lp = True
# this is used for creating the vacuum structure,
# unfortunatly parmed forgets afterwards about the frame_type
# which is necessary for the check_for_lp function
if lp:
g = psf.groups
frame_idx = []
frame_frame = []
for atom in psf.atoms:
if hasattr(atom, "frame_type"):
frame_idx.append(atom.idx)
frame_frame.append(atom.frame_type)
psf = psf[f":{self.tlc}"] # the important part
psf.groups = g
for atom in psf.atoms:
if atom.idx in frame_idx:
atom.frame_type = frame_frame[frame_idx.index(atom.idx)]
else:
psf = psf[f":{self.tlc}"]
else:
psf_file_name = configuration["system"][self.structure][env][
"psf_file_name"
]
crd_file_name = configuration["system"][self.structure][env][
"crd_file_name"
]
psf_file_path = f"{self.charmm_gui_base}/{env}/openmm/{psf_file_name}.psf"
crd_file_path = f"{self.charmm_gui_base}/{env}/openmm/{crd_file_name}.crd"
psf = pm.charmm.CharmmPsfFile(psf_file_path)
coord = pm.charmm.CharmmCrdFile(crd_file_path)
psf.coordinates = coord.coordinates
return psf
def _determine_offset_and_set_possible_dummy_properties(
self, psf: pm.charmm.CharmmPsfFile
) -> int:
"""
Determines the offset and sets possible properties on the psf.
Parameters
----------
psf : pm.charmm.CharmmPsfFile
env: str
waterbox,complex or vacuum
Returns
----------
"""
assert type(psf) == pm.charmm.CharmmPsfFile
if len(psf.view[f":{self.tlc}"].atoms) < 1:
raise RuntimeError(f"No ligand selected for tlc: {self.tlc}")
psf.number_of_dummys = 0
psf.mutations_to_default = 0
idx_list = []
for atom in psf.view[f":{self.tlc}"].atoms:
idx_list.append(int(atom.idx))
# charge, epsilon and rmin are directly modiefied
atom.initial_charge = atom.charge
atom.initial_epsilon = atom.epsilon
atom.initial_rmin = atom.rmin
return min(idx_list)
def _return_small_molecule(self, env: str) -> Chem.rdchem.Mol:
import glob
charmm_gui_env = self.charmm_gui_base + env
possible_files = []
for ending in ["sdf", "mol", "mol2"]:
possible_files.extend(glob.glob(f"{charmm_gui_env}/*/*{ending}"))
# looking for small molecule
# start with sdf
mol2_detected = False
mol_detected = False
found_small_molecule = False
for f in possible_files:
if f.endswith(".sdf"):
suppl = Chem.SDMolSupplier(f, removeHs=False)
mol = next(suppl)
logger.info(f"SDF file loaded: {f}")
found_small_molecule = True
return mol
elif f.endswith(".mol2"):
mol2_detected = True
elif f.endswith(".mol"):
mol_detected = True
if mol2_detected == True or mol_detected == True:
raise RuntimeError(
"Please convert mol2 or mol file to sdf using obabel: {possible_files}."
)
if not found_small_molecule:
raise FileNotFoundError(
"Could not load small molecule sdf file in {charmm_gui_env}. Aborting."
)
def _generate_rdkit_mol(
self, env: str, psf: pm.charmm.CharmmPsfFile
) -> Chem.rdchem.Mol:
"""
Generates the rdkit mol object.
Parameters
----------
env: str
waterbox,complex or vacuum
psf: pm.charmm.CharmmPsfFile
Returns
----------
mol: rdkit.Chem.rdchem.Mol
"""
assert type(psf) == pm.charmm.CharmmPsfFile
mol = self._return_small_molecule(env)
(
atom_idx_to_atom_name,
_,
atom_name_to_atom_type,
atom_idx_to_atom_partial_charge,
) = self.generate_atom_tables_from_psf(psf)
for atom in mol.GetAtoms():
atom.SetProp("atom_name", atom_idx_to_atom_name[atom.GetIdx()])
atom.SetProp(
"atom_type",
atom_name_to_atom_type[atom_idx_to_atom_name[atom.GetIdx()]],
)
atom.SetProp("atom_index", str(atom.GetIdx()))
atom.SetProp(
"atom_charge", str(atom_idx_to_atom_partial_charge[atom.GetIdx()])
)
# check if psf and sdf have same indeces
for a in mol.GetAtoms():
if str(psf[a.GetIdx()].element_name) == str(a.GetSymbol()):
pass
else:
raise RuntimeError("PSF to mol conversion did not work! Aborting.")
return mol
def generate_atom_tables_from_psf(
self, psf: pm.charmm.CharmmPsfFile
) -> Tuple[dict, dict, dict, dict]:
"""
Generate mapping dictionaries for a molecule in a psf.
Parameters
----------
psf: pm.charmm.CharmmPsfFile
Returns
----------
dict's
"""
atom_idx_to_atom_name = dict()
atom_name_to_atom_idx = dict()
atom_name_to_atom_type = dict()
atom_idx_to_atom_partial_charge = dict()
for atom in psf.view[f":{self.tlc}"].atoms:
atom_name = atom.name
atom_index = atom.idx
atom_type = atom.type
atom_charge = atom.charge
atom_idx_to_atom_name[atom_index] = atom_name
atom_name_to_atom_idx[atom_name] = atom_index
atom_name_to_atom_type[atom_name] = atom_type
atom_idx_to_atom_partial_charge[atom_index] = atom_charge
return (
atom_idx_to_atom_name,
atom_name_to_atom_idx,
atom_name_to_atom_type,
atom_idx_to_atom_partial_charge,
)