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