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