Remove unused thole polarization rfac parameter
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / hizzie.cpp
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,2021, 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/utilities.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/topology/block.h"
55 #include "gromacs/topology/symtab.h"
56 #include "gromacs/utility/arraysize.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/smalloc.h"
60
61 static int in_strings(char* key, int nstr, const char** str)
62 {
63     int j;
64
65     for (j = 0; (j < nstr); j++)
66     {
67         if (strcmp(str[j], key) == 0)
68         {
69             return j;
70         }
71     }
72
73     return -1;
74 }
75
76 static bool hbond(rvec x[], int i, int j, real distance)
77 {
78     real tol = distance * distance;
79     rvec tmp;
80
81     rvec_sub(x[i], x[j], tmp);
82
83     return (iprod(tmp, tmp) < tol);
84 }
85
86 static void chk_allhb(t_atoms* pdba, rvec x[], t_blocka* hb, const bool donor[], const bool accept[], real dist)
87 {
88     int i, j, k, ii, natom;
89
90     natom = pdba->nr;
91     snew(hb->index, natom + 1);
92     snew(hb->a, 6 * natom);
93     hb->nr  = natom;
94     hb->nra = 6 * natom;
95
96     k = ii          = 0;
97     hb->index[ii++] = 0;
98     for (i = 0; (i < natom); i++)
99     {
100         if (donor[i])
101         {
102             for (j = i + 1; (j < natom); j++)
103             {
104                 if ((accept[j]) && (hbond(x, i, j, dist)))
105                 {
106                     hb->a[k++] = j;
107                 }
108             }
109         }
110         else if (accept[i])
111         {
112             for (j = i + 1; (j < natom); j++)
113             {
114                 if ((donor[j]) && (hbond(x, i, j, dist)))
115                 {
116                     hb->a[k++] = j;
117                 }
118             }
119         }
120         hb->index[ii++] = k;
121     }
122     hb->nra = k;
123 }
124
125 static bool chk_hbonds(int i, t_atoms* pdba, rvec x[], const bool ad[], bool hbond[], rvec xh, real angle, real dist)
126 {
127     bool bHB;
128     int  j, aj, ri, natom;
129     real d2, dist2, a;
130     rvec nh, oh;
131
132     natom = pdba->nr;
133     bHB   = FALSE;
134     ri    = pdba->atom[i].resind;
135     dist2 = gmx::square(dist);
136     for (j = 0; (j < natom); j++)
137     {
138         /* Check whether the other atom is a donor/acceptor and not i */
139         if ((ad[j]) && (j != i))
140         {
141             /* Check whether the other atom is on the same ring as well */
142             if ((pdba->atom[j].resind != ri)
143                 || ((strcmp(*pdba->atomname[j], "ND1") != 0) && (strcmp(*pdba->atomname[j], "NE2") != 0)))
144             {
145                 aj = j;
146                 d2 = distance2(x[i], x[j]);
147                 rvec_sub(x[i], xh, nh);
148                 rvec_sub(x[aj], xh, oh);
149                 a = gmx::c_rad2Deg * acos(cos_angle(nh, oh));
150                 if ((d2 < dist2) && (a > angle))
151                 {
152                     hbond[i] = TRUE;
153                     bHB      = TRUE;
154                 }
155             }
156         }
157     }
158     return bHB;
159 }
160
161 static void calc_ringh(rvec xattach, rvec xb, rvec xc, rvec xh)
162 {
163     rvec tab, tac;
164     real n;
165
166     /* Add a proton on a ring to atom attach at distance 0.1 nm */
167     rvec_sub(xattach, xb, tab);
168     rvec_sub(xattach, xc, tac);
169     rvec_add(tab, tac, xh);
170     n = 0.1 / norm(xh);
171     svmul(n, xh, xh);
172     rvec_inc(xh, xattach);
173 }
174
175 void set_histp(t_atoms* pdba, rvec* x, t_symtab* symtab, real angle, real dist)
176 {
177     static const char* prot_acc[] = { "O", "OD1", "OD2", "OE1", "OE2", "OG", "OG1", "OH", "OW" };
178 #define NPA asize(prot_acc)
179     static const char* prot_don[] = { "N",  "NH1", "NH2", "NE", "ND1", "ND2", "NE2",
180                                       "NZ", "OG",  "OG1", "OH", "NE1", "OW" };
181 #define NPD asize(prot_don)
182
183     bool *          donor, *acceptor;
184     bool*           hbond;
185     bool            bHDd, bHEd;
186     rvec            xh1, xh2;
187     int             natom;
188     int             i, j, nd, na, hisind;
189     HistidineStates type = HistidineStates::Count;
190     int             nd1, ne2, cg, cd2, ce1;
191     t_blocka*       hb;
192     char*           atomnm;
193
194     natom = pdba->nr;
195
196     i = 0;
197     while (i < natom && gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
198     {
199         i++;
200     }
201     if (natom == i)
202     {
203         return;
204     }
205
206     /* A histidine residue exists that requires automated assignment, so
207      * doing the analysis of donors and acceptors is worthwhile. */
208     fprintf(stderr,
209             "Analysing hydrogen-bonding network for automated assignment of histidine\n"
210             " protonation.");
211
212     snew(donor, natom);
213     snew(acceptor, natom);
214     snew(hbond, natom);
215     snew(hb, 1);
216
217     nd = na = 0;
218     for (j = 0; (j < natom); j++)
219     {
220         if (in_strings(*pdba->atomname[j], NPA, prot_acc) != -1)
221         {
222             acceptor[j] = TRUE;
223             na++;
224         }
225         if (in_strings(*pdba->atomname[j], NPD, prot_don) != -1)
226         {
227             donor[j] = TRUE;
228             nd++;
229         }
230     }
231     fprintf(stderr, " %d donors and %d acceptors were found.\n", nd, na);
232     chk_allhb(pdba, x, hb, donor, acceptor, dist);
233     fprintf(stderr, "There are %d hydrogen bonds\n", hb->nra);
234
235     /* Now do the HIS stuff */
236     hisind = -1;
237     while (i < natom)
238     {
239         if (gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
240         {
241             i++;
242         }
243         else
244         {
245             if (pdba->atom[i].resind != hisind)
246             {
247                 hisind = pdba->atom[i].resind;
248
249                 /* Find the  atoms in the ring */
250                 nd1 = ne2 = cg = cd2 = ce1 = -1;
251                 while (i < natom && pdba->atom[i].resind == hisind)
252                 {
253                     atomnm = *pdba->atomname[i];
254                     if (strcmp(atomnm, "CD2") == 0)
255                     {
256                         cd2 = i;
257                     }
258                     else if (strcmp(atomnm, "CG") == 0)
259                     {
260                         cg = i;
261                     }
262                     else if (strcmp(atomnm, "CE1") == 0)
263                     {
264                         ce1 = i;
265                     }
266                     else if (strcmp(atomnm, "ND1") == 0)
267                     {
268                         nd1 = i;
269                     }
270                     else if (strcmp(atomnm, "NE2") == 0)
271                     {
272                         ne2 = i;
273                     }
274
275                     i++;
276                 }
277
278                 if (!((cg == -1) || (cd2 == -1) || (ce1 == -1) || (nd1 == -1) || (ne2 == -1)))
279                 {
280                     calc_ringh(x[nd1], x[cg], x[ce1], xh1);
281                     calc_ringh(x[ne2], x[ce1], x[cd2], xh2);
282
283                     bHDd = chk_hbonds(nd1, pdba, x, acceptor, hbond, xh1, angle, dist);
284                     chk_hbonds(nd1, pdba, x, donor, hbond, xh1, angle, dist);
285                     bHEd = chk_hbonds(ne2, pdba, x, acceptor, hbond, xh2, angle, dist);
286                     chk_hbonds(ne2, pdba, x, donor, hbond, xh2, angle, dist);
287
288                     if (bHDd)
289                     {
290                         if (bHEd)
291                         {
292                             type = HistidineStates::H;
293                         }
294                         else
295                         {
296                             type = HistidineStates::A;
297                         }
298                     }
299                     else
300                     {
301                         type = HistidineStates::B;
302                     }
303                     fprintf(stderr,
304                             "Will use %s for residue %d\n",
305                             enumValueToString(type),
306                             pdba->resinfo[hisind].nr);
307                 }
308                 else
309                 {
310                     gmx_fatal(FARGS, "Incomplete ring in HIS%d", pdba->resinfo[hisind].nr);
311                 }
312
313                 pdba->resinfo[hisind].rtp = put_symtab(symtab, enumValueToString(type));
314             }
315         }
316     }
317     done_blocka(hb);
318     sfree(hb);
319     sfree(donor);
320     sfree(acceptor);
321     sfree(hbond);
322 }