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