-
Notifications
You must be signed in to change notification settings - Fork 8
/
make_septop.py
358 lines (312 loc) · 13.9 KB
/
make_septop.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
"""
Functions to build the separated topology of two ligands.
"""
import parmed as pmd
def combine_ligands_top(top_A, top_B, septop, ligand='LIG'):
"""Create a combined topology file of the complex A and ligand B. Insert ligand B into the same molecule section as ligand A.
Parameters
----------
top_A : str
Topology file of complex A
top_B : str
Topology file of complex B
septop : str
Topology file with ligand B inserted in complex A topology
"""
complex_A = pmd.load_file(top_A)
complex_B = pmd.load_file(top_B)
# Grab ligand B from complex B
lig2 = complex_B[ligand, :]
#Find residue number ligand A
for res_nr, res in enumerate(complex_A.residues):
if res.name == ligand:
ligand_res_nr = res_nr
#Split complex into individial molecules ( [ moleculetype ] sections GROMACS)
split_A = complex_A.split()
#Find index of the [ moleculetype ] section of ligand A
count = 0
for ind,i in enumerate(split_A):
count += len(i[1])
if i[0].residues[0].name == ligand:
lig1_pos = ind
break
#Create new topology where ligand B is inserted directly after ligand A
sep_top = complex_A[:ligand_res_nr+1, :] + lig2 + complex_A[ligand_res_nr+1:, :]
# combine ligA and ligB into a single molecule entry based on index [ moleculetype ] secion ligand
# sep_top.write(septop, combine = [[lig1_pos, lig1_pos+1]])
sep_top.write(septop, combine = [[count-1, count]])
# sep_top.write(septop, combine = [[ligand_res_nr, ligand_res_nr+1]])
return
def make_section_dictionaries(text):
# Create dictionary of different section of the topology file
diclist = []
dic = {}
start_section = False
index=0
while index < len(text):
line = text[index]
# '[ ' marks the start of a section
if '[ ' in line:
if not start_section: #and index == 0:
start_section = True
key = line
dic[key] = []
index+=1
# reset for next section dictionary
else:
# break if already next section
diclist.append(dic)
dic = {}
start_section=False
elif start_section:
# append all the lines that belong to that section
dic[key].append(line)
index+=1
else:
index+=1
# append the last section to the diclist
diclist.append(dic)
return diclist
def create_A_and_B_state_ligand(line, A_B_state='vdwq_q', lig = 1):
"""Create A and B state topology for a ligand.
Parameters
----------
line : str
'Atom line': with atomtype, mass, charge,...
A_B_state : str
Interactions in the A state and in the B state.
vdwq_vdwq: ligand fully interacting in A and B state
vdwq_vdw: vdw interactions and electrostatics in the A_state, only vdw in the B_state
vdw_vdwq: charge
vdw_dummy
dummy_vdw
vdwq_dummy
vdwq_scaled-vdw
scaled-vdw_dummy
dummy_scaled-vdwq
scaled-vdwq_vdwq
lig: int
if first or second ligand
Returns
-------
text : str
Atoms line for topology file with A and B state parameters
"""
atom_number = line.split()[0]
atom_type = line.split()[1]
if lig == 1:
atom_type = "".join(('LIG1_', atom_type))
if lig == 2:
atom_type = "".join(('LIG2_', atom_type))
scaled_atomtype = "".join(('scaled_', atom_type))
residue_nr = line.split()[2]
residue_name = line.split()[3]
atom_name = line.split()[4]
cgnr = line.split()[5]
charge = line.split()[6]
mass = line.split()[7]
# A and B state are the same
if A_B_state == 'vdwq_vdwq':
text = ' ' + atom_number + ' ' + atom_type + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + charge + ' ' + mass + ' ' + \
atom_type + ' ' + charge + ' ' + mass + '\n'
# Turn on vdw
elif A_B_state == 'dummy_vdw':
charge = str(0.0)
text = ' ' + atom_number + ' d%s ' % atom_type + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + charge + ' ' + mass + ' ' + \
atom_type + ' ' + charge + ' ' + mass + '\n'
# Turn vdw off
elif A_B_state == 'vdw_dummy':
charge = str(0.0)
text = ' ' + atom_number + ' ' + atom_type + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + charge + ' ' + mass + \
' d%s ' % atom_type + ' ' + charge + ' ' + mass + '\n'
# Turn vdw and electrostatics off
elif A_B_state == 'vdwq_dummy':
text = ' ' + atom_number + ' ' + atom_type + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + charge + ' ' + mass + \
' d%s ' % atom_type + ' ' + ' 0.0 ' + ' ' + mass + '\n'
# uncharge
elif A_B_state == 'vdwq_vdw':
text = ' ' + atom_number + ' ' + atom_type + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + charge + ' ' + mass + ' ' + \
atom_type + ' ' + ' 0.0 ' + ' ' + mass + '\n'
# charge
elif A_B_state == 'vdw_vdwq':
text = ' ' + atom_number + ' ' + atom_type + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + str(0.0) + ' ' + \
mass + ' ' + atom_type + ' ' + charge + ' ' + mass + '\n'
# Posre off
elif A_B_state == 'dummy':
charge = str(0.0)
text = ' ' + atom_number + ' d%s ' % atom_type + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + charge + ' ' + mass + ' ' + '\n'
# Turn vdw and electrostatics off
elif A_B_state == 'vdwq':
text = ' ' + atom_number + ' ' + atom_type + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + charge + ' ' + mass + ' ' + '\n'
# vdw till gamma and charges
elif A_B_state == 'vdwq_scaled-vdw':
text = ' ' + atom_number + ' ' + atom_type + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + charge + ' ' + mass + ' ' + \
scaled_atomtype + ' ' + ' 0.0 ' + ' ' + mass + '\n'
# Turn vdw off from scaled ligand
elif A_B_state == 'scaled-vdw_dummy':
charge = str(0.0)
text = ' ' + atom_number + ' ' + scaled_atomtype + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + charge + ' ' + mass + \
' d%s ' % atom_type + ' ' + charge + ' ' + mass + '\n'
# Turn on vdw
elif A_B_state == 'dummy_scaled-vdwq':
text = ' ' + atom_number + ' d%s ' % atom_type + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + str(0.0) + ' ' + mass + ' ' + \
scaled_atomtype + ' ' + charge + ' ' + mass + '\n'
# turn on rest of vdw
elif A_B_state == 'scaled-vdwq_vdwq':
text = ' ' + atom_number + ' ' + scaled_atomtype + ' ' + residue_nr + ' ' + \
residue_name + ' ' + atom_name + ' ' + cgnr + ' ' + charge + ' ' + mass + ' ' + \
atom_type + ' ' + charge + ' ' + mass + '\n'
else:
print('Transformation not implemented yet')
return text
def atom_types_ligand(in_top, ligand='LIG'):
"""Store atom types ligand in a list. Atom types have to come from individual ligand files, because parmed just
combines different atom types to a single one if they have the same name.
Parameters
----------
in_top : str
topology file of complex and both ligands (generated from combine_ligands_top)
ligand = str
three letter code for the ligand residue, default = 'LIG'
"""
file = open(in_top, 'r')
text = file.readlines()
file.close()
end_text = len(text)
atomtype = []
for line in text:
if ligand in line and not line.startswith(';') and line.split()[0].isdigit():
at = line.split()[1]
for l in text:
if l.startswith(';') or len(l.split()) == 0:
continue
if l.split()[0] == at and l not in atomtype:
atomtype.append(l)
return atomtype
def create_top(in_top, out_top, gamma, A_B_state_ligA, A_B_state_ligB, in_top_A, in_top_B, ligand='LIG'):
"""Create separated topology
Parameters
----------
in_top : str
topology file of complex and both ligands (generated from combine_ligands_top)
out_top: str
name for output topology file
gamma: int
scaling parameter for REST scaling
A_B_state_ligA : str
Interactions in the A state and in the B state for ligand A.
vdwq_vdwq: ligand fully interacting in A and B state
vdwq_vdw: vdw interactions and electrostatics in the A_state, only vdw in the B_state
vdw_vdwq: charge
vdw_dummy
dummy_vdw
vdwq_dummy
A_B_state_ligB: str
Interactions in the A state and in the B state for ligand B
ligand = str
three letter code for the ligand residue, default = 'LIG'
Examples
--------
Turn on vdw ligand B: create_top(top, turnon_vdw_B, 'vdwq_vdwq', 'dummy_vdw')
Charge ligand B while uncharging ligand A: create_top(top, charge_uncharge, 'vdwq_vdw', 'vdw_vdwq')
Turn off vdw ligand A: create_top(top, turnoff_vdw_A, 'vdw_dummy', 'vdwq_vdwq')
"""
atomtype_i = atom_types_ligand(in_top_A, ligand=ligand)
atomtype_j = atom_types_ligand(in_top_B, ligand=ligand)
file = open(in_top, 'r')
text = file.readlines()
file.close()
file = open(out_top, 'w')
end_text = len(text)
count = 0
outtext = []
section = 0
# print text does print out the proper file
#print(text)
diclist = make_section_dictionaries(text)
for dic in diclist:
# Create dictionary of different sections
#dic = make_section_dictionary(text[count:])
#print(dic)
count += 1
# Loop through sections
for key, value in dic.items():
# For every atomtype add a dummy-atomtype with no vdW interactions
if 'atomtypes' in key:
outtext.append(key)
for v in value:
if v.startswith(';') or v.startswith('\n'):
continue
if v.split()[0] in [l.split()[0] for l in atomtype_i]:
continue
if v.split()[0] in [l.split()[0] for l in atomtype_j]:
continue
else:
outtext.append(v)
for inx,atomtype in enumerate([atomtype_i, atomtype_j]):
for l in atomtype:
if inx == 0:
new_at = "".join(('LIG1_', l))
elif inx == 1:
new_at = "".join(('LIG2_', l))
at_line = new_at.split()
#Scale the epsilon of the vdw terms by a constant gamma
scaled_at = 'scaled_' + at_line[0] + ' ' + at_line[1] + ' ' + at_line[
2] + ' 0.0 ' + ' A ' + at_line[5] + ' ' + str(float(at_line[6])*gamma) + '\n'
###To do: decide if I only allow those 4 endstates or make it more flexible!!!
dummy = 'd' + at_line[0] + ' ' + at_line[1] + ' ' + at_line[
2] + ' 0.0 ' + ' A ' + ' 0.0 ' + ' 0.0\n'
outtext.append(new_at)
outtext.append(scaled_at)
outtext.append(dummy)
outtext.append('\n\n')
outtext.append('[ nonbond_params ]\n')
for i in atomtype_i:
i = "".join(('LIG1_', i.split()[0]))
scaled_i = "".join(('scaled_', i.split()[0]))
for j in atomtype_j:
j = "".join(('LIG2_', j.split()[0]))
scaled_j = "".join(('scaled_', j.split()[0]))
nb = i + ' ' + j + ' 1 0 0\n'
scaled_nb = scaled_i + ' ' + scaled_j + ' 1 0 0\n'
outtext.append(nb)
outtext.append(scaled_nb)
outtext.append('\n\n')
# Modify the ligand system
elif 'atoms' in key and ligand in value[2].split():
outtext.append(key)
atomindex_i, atomindex_j = [], []
for v in value:
if v.startswith(';') or v.startswith('\n'):
outtext.append(v)
# For ligand A create an A and a B state according to A_B_state_ligA
if ligand in v and int(v.split()[2]) == 1 and not v.startswith(';'):
atomindex_i.append(v.split()[0])
line = create_A_and_B_state_ligand(v, A_B_state_ligA, lig = 1)
outtext.append(line)
# For ligand B create an A and a B state according to A_B_state_ligB
if ligand in v and int(v.split()[2]) == 2 and not v.startswith(';'):
atomindex_j.append(v.split()[0])
line = create_A_and_B_state_ligand(v, A_B_state_ligB, lig = 2)
outtext.append(line)
outtext.append('\n\n')
else:
outtext.append(key)
outtext.append(value)
section += 1
for sec in outtext:
for line in sec:
file.write(line)
file.close()
return out_top