1 /*
2 * This file is part of the GROMACS molecular simulation package.
3 *
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2017,2018 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
11 *
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
16 *
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
21 *
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 *
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
34 *
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
37 */
38 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
40
41 #include "hizzie.h"
42
43 #include <cmath>
44 #include <cstdio>
45 #include <cstring>
46
47 #include "gromacs/fileio/pdbio.h"
48 #include "gromacs/gmxpreprocess/pdb2top.h"
49 #include "gromacs/gmxpreprocess/toputil.h"
50 #include "gromacs/math/functions.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/topology/block.h"
54 #include "gromacs/topology/symtab.h"
55 #include "gromacs/utility/arraysize.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/smalloc.h"
59
in_strings(char * key,int nstr,const char ** str)60 static int in_strings(char* key, int nstr, const char** str)
61 {
62 int j;
63
64 for (j = 0; (j < nstr); j++)
65 {
66 if (strcmp(str[j], key) == 0)
67 {
68 return j;
69 }
70 }
71
72 return -1;
73 }
74
hbond(rvec x[],int i,int j,real distance)75 static bool hbond(rvec x[], int i, int j, real distance)
76 {
77 real tol = distance * distance;
78 rvec tmp;
79
80 rvec_sub(x[i], x[j], tmp);
81
82 return (iprod(tmp, tmp) < tol);
83 }
84
chk_allhb(t_atoms * pdba,rvec x[],t_blocka * hb,const bool donor[],const bool accept[],real dist)85 static void chk_allhb(t_atoms* pdba, rvec x[], t_blocka* hb, const bool donor[], const bool accept[], real dist)
86 {
87 int i, j, k, ii, natom;
88
89 natom = pdba->nr;
90 snew(hb->index, natom + 1);
91 snew(hb->a, 6 * natom);
92 hb->nr = natom;
93 hb->nra = 6 * natom;
94
95 k = ii = 0;
96 hb->index[ii++] = 0;
97 for (i = 0; (i < natom); i++)
98 {
99 if (donor[i])
100 {
101 for (j = i + 1; (j < natom); j++)
102 {
103 if ((accept[j]) && (hbond(x, i, j, dist)))
104 {
105 hb->a[k++] = j;
106 }
107 }
108 }
109 else if (accept[i])
110 {
111 for (j = i + 1; (j < natom); j++)
112 {
113 if ((donor[j]) && (hbond(x, i, j, dist)))
114 {
115 hb->a[k++] = j;
116 }
117 }
118 }
119 hb->index[ii++] = k;
120 }
121 hb->nra = k;
122 }
123
chk_hbonds(int i,t_atoms * pdba,rvec x[],const bool ad[],bool hbond[],rvec xh,real angle,real dist)124 static bool chk_hbonds(int i, t_atoms* pdba, rvec x[], const bool ad[], bool hbond[], rvec xh, real angle, real dist)
125 {
126 bool bHB;
127 int j, aj, ri, natom;
128 real d2, dist2, a;
129 rvec nh, oh;
130
131 natom = pdba->nr;
132 bHB = FALSE;
133 ri = pdba->atom[i].resind;
134 dist2 = gmx::square(dist);
135 for (j = 0; (j < natom); j++)
136 {
137 /* Check whether the other atom is a donor/acceptor and not i */
138 if ((ad[j]) && (j != i))
139 {
140 /* Check whether the other atom is on the same ring as well */
141 if ((pdba->atom[j].resind != ri)
142 || ((strcmp(*pdba->atomname[j], "ND1") != 0) && (strcmp(*pdba->atomname[j], "NE2") != 0)))
143 {
144 aj = j;
145 d2 = distance2(x[i], x[j]);
146 rvec_sub(x[i], xh, nh);
147 rvec_sub(x[aj], xh, oh);
148 a = RAD2DEG * acos(cos_angle(nh, oh));
149 if ((d2 < dist2) && (a > angle))
150 {
151 hbond[i] = TRUE;
152 bHB = TRUE;
153 }
154 }
155 }
156 }
157 return bHB;
158 }
159
calc_ringh(rvec xattach,rvec xb,rvec xc,rvec xh)160 static void calc_ringh(rvec xattach, rvec xb, rvec xc, rvec xh)
161 {
162 rvec tab, tac;
163 real n;
164
165 /* Add a proton on a ring to atom attach at distance 0.1 nm */
166 rvec_sub(xattach, xb, tab);
167 rvec_sub(xattach, xc, tac);
168 rvec_add(tab, tac, xh);
169 n = 0.1 / norm(xh);
170 svmul(n, xh, xh);
171 rvec_inc(xh, xattach);
172 }
173
set_histp(t_atoms * pdba,rvec * x,t_symtab * symtab,real angle,real dist)174 void set_histp(t_atoms* pdba, rvec* x, t_symtab* symtab, real angle, real dist)
175 {
176 static const char* prot_acc[] = { "O", "OD1", "OD2", "OE1", "OE2", "OG", "OG1", "OH", "OW" };
177 #define NPA asize(prot_acc)
178 static const char* prot_don[] = { "N", "NH1", "NH2", "NE", "ND1", "ND2", "NE2",
179 "NZ", "OG", "OG1", "OH", "NE1", "OW" };
180 #define NPD asize(prot_don)
181
182 bool * donor, *acceptor;
183 bool* hbond;
184 bool bHDd, bHEd;
185 rvec xh1, xh2;
186 int natom;
187 int i, j, nd, na, hisind, type = -1;
188 int nd1, ne2, cg, cd2, ce1;
189 t_blocka* hb;
190 char* atomnm;
191
192 natom = pdba->nr;
193
194 i = 0;
195 while (i < natom && gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
196 {
197 i++;
198 }
199 if (natom == i)
200 {
201 return;
202 }
203
204 /* A histidine residue exists that requires automated assignment, so
205 * doing the analysis of donors and acceptors is worthwhile. */
206 fprintf(stderr,
207 "Analysing hydrogen-bonding network for automated assignment of histidine\n"
208 " protonation.");
209
210 snew(donor, natom);
211 snew(acceptor, natom);
212 snew(hbond, natom);
213 snew(hb, 1);
214
215 nd = na = 0;
216 for (j = 0; (j < natom); j++)
217 {
218 if (in_strings(*pdba->atomname[j], NPA, prot_acc) != -1)
219 {
220 acceptor[j] = TRUE;
221 na++;
222 }
223 if (in_strings(*pdba->atomname[j], NPD, prot_don) != -1)
224 {
225 donor[j] = TRUE;
226 nd++;
227 }
228 }
229 fprintf(stderr, " %d donors and %d acceptors were found.\n", nd, na);
230 chk_allhb(pdba, x, hb, donor, acceptor, dist);
231 fprintf(stderr, "There are %d hydrogen bonds\n", hb->nra);
232
233 /* Now do the HIS stuff */
234 hisind = -1;
235 while (i < natom)
236 {
237 if (gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
238 {
239 i++;
240 }
241 else
242 {
243 if (pdba->atom[i].resind != hisind)
244 {
245 hisind = pdba->atom[i].resind;
246
247 /* Find the atoms in the ring */
248 nd1 = ne2 = cg = cd2 = ce1 = -1;
249 while (i < natom && pdba->atom[i].resind == hisind)
250 {
251 atomnm = *pdba->atomname[i];
252 if (strcmp(atomnm, "CD2") == 0)
253 {
254 cd2 = i;
255 }
256 else if (strcmp(atomnm, "CG") == 0)
257 {
258 cg = i;
259 }
260 else if (strcmp(atomnm, "CE1") == 0)
261 {
262 ce1 = i;
263 }
264 else if (strcmp(atomnm, "ND1") == 0)
265 {
266 nd1 = i;
267 }
268 else if (strcmp(atomnm, "NE2") == 0)
269 {
270 ne2 = i;
271 }
272
273 i++;
274 }
275
276 if (!((cg == -1) || (cd2 == -1) || (ce1 == -1) || (nd1 == -1) || (ne2 == -1)))
277 {
278 calc_ringh(x[nd1], x[cg], x[ce1], xh1);
279 calc_ringh(x[ne2], x[ce1], x[cd2], xh2);
280
281 bHDd = chk_hbonds(nd1, pdba, x, acceptor, hbond, xh1, angle, dist);
282 chk_hbonds(nd1, pdba, x, donor, hbond, xh1, angle, dist);
283 bHEd = chk_hbonds(ne2, pdba, x, acceptor, hbond, xh2, angle, dist);
284 chk_hbonds(ne2, pdba, x, donor, hbond, xh2, angle, dist);
285
286 if (bHDd)
287 {
288 if (bHEd)
289 {
290 type = ehisH;
291 }
292 else
293 {
294 type = ehisA;
295 }
296 }
297 else
298 {
299 type = ehisB;
300 }
301 fprintf(stderr, "Will use %s for residue %d\n", hh[type], pdba->resinfo[hisind].nr);
302 }
303 else
304 {
305 gmx_fatal(FARGS, "Incomplete ring in HIS%d", pdba->resinfo[hisind].nr);
306 }
307
308 pdba->resinfo[hisind].rtp = put_symtab(symtab, hh[type]);
309 }
310 }
311 }
312 done_blocka(hb);
313 sfree(hb);
314 sfree(donor);
315 sfree(acceptor);
316 sfree(hbond);
317 }
318