Code beautification with uncrustify
[alexxy/gromacs.git] / src / programs / pdb2gmx / hizzie.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <stdio.h>
41 #include <string.h>
42 #include "typedefs.h"
43 #include "pdbio.h"
44 #include "smalloc.h"
45 #include "vec.h"
46 #include "physics.h"
47 #include "toputil.h"
48 #include "pdb2top.h"
49 #include "string2.h"
50 #include "macros.h"
51
52 static int in_strings(char *key, int nstr, const char **str)
53 {
54     int j;
55
56     for (j = 0; (j < nstr); j++)
57     {
58         if (strcmp(str[j], key) == 0)
59         {
60             return j;
61         }
62     }
63
64     return -1;
65 }
66
67 static gmx_bool hbond(rvec x[], int i, int j, real distance)
68 {
69     real   tol = distance*distance;
70     rvec   tmp;
71
72     rvec_sub(x[i], x[j], tmp);
73
74     return (iprod(tmp, tmp) < tol);
75 }
76
77 static void chk_allhb(t_atoms *pdba, rvec x[], t_blocka *hb,
78                       gmx_bool donor[], gmx_bool accept[], real dist)
79 {
80     int i, j, k, ii, natom;
81
82     natom = pdba->nr;
83     snew(hb->index, natom+1);
84     snew(hb->a, 6*natom);
85     hb->nr  = natom;
86     hb->nra = 6*natom;
87
88     k               = ii = 0;
89     hb->index[ii++] = 0;
90     for (i = 0; (i < natom); i++)
91     {
92         if (donor[i])
93         {
94             for (j = i+1; (j < natom); j++)
95             {
96                 if ((accept[j]) && (hbond(x, i, j, dist)))
97                 {
98                     hb->a[k++] = j;
99                 }
100             }
101         }
102         else if (accept[i])
103         {
104             for (j = i+1; (j < natom); j++)
105             {
106                 if ((donor[j]) && (hbond(x, i, j, dist)))
107                 {
108                     hb->a[k++] = j;
109                 }
110             }
111         }
112         hb->index[ii++] = k;
113     }
114     hb->nra = k;
115 }
116
117 static void pr_hbonds(FILE *fp, t_blocka *hb, t_atoms *pdba)
118 {
119     int i, j, k, j0, j1;
120
121     fprintf(fp, "Dumping all hydrogen bonds!\n");
122     for (i = 0; (i < hb->nr); i++)
123     {
124         j0 = hb->index[i];
125         j1 = hb->index[i+1];
126         for (j = j0; (j < j1); j++)
127         {
128             k = hb->a[j];
129             fprintf(fp, "%5s%4d%5s - %5s%4d%5s\n",
130                     *pdba->resinfo[pdba->atom[i].resind].name,
131                     pdba->resinfo[pdba->atom[i].resind].nr, *pdba->atomname[i],
132                     *pdba->resinfo[pdba->atom[k].resind].name,
133                     pdba->resinfo[pdba->atom[k].resind].nr, *pdba->atomname[k]);
134         }
135     }
136 }
137
138 static gmx_bool chk_hbonds(int i, t_atoms *pdba, rvec x[],
139                            gmx_bool ad[], gmx_bool hbond[], rvec xh,
140                            real angle, real dist)
141 {
142     gmx_bool bHB;
143     int      j, aj, ri, natom;
144     real     d2, dist2, a;
145     rvec     nh, oh;
146
147     natom = pdba->nr;
148     bHB   = FALSE;
149     ri    = pdba->atom[i].resind;
150     dist2 = sqr(dist);
151     for (j = 0; (j < natom); j++)
152     {
153         /* Check whether the other atom is a donor/acceptor and not i */
154         if ((ad[j]) && (j != i))
155         {
156             /* Check whether the other atom is on the same ring as well */
157             if ((pdba->atom[j].resind != ri) ||
158                 ((strcmp(*pdba->atomname[j], "ND1") != 0) &&
159                  (strcmp(*pdba->atomname[j], "NE2") != 0)))
160             {
161                 aj  = j;
162                 d2  = distance2(x[i], x[j]);
163                 rvec_sub(x[i], xh, nh);
164                 rvec_sub(x[aj], xh, oh);
165                 a  = RAD2DEG * acos(cos_angle(nh, oh));
166                 if ((d2 < dist2) && (a > angle))
167                 {
168                     if (debug)
169                     {
170                         fprintf(debug,
171                                 "HBOND between %s%d-%s and %s%d-%s is %g nm, %g deg\n",
172                                 *pdba->resinfo[pdba->atom[i].resind].name,
173                                 pdba->resinfo[pdba->atom[i].resind].nr, *pdba->atomname[i],
174                                 *pdba->resinfo[pdba->atom[aj].resind].name,
175                                 pdba->resinfo[pdba->atom[aj].resind].nr, *pdba->atomname[aj],
176                                 sqrt(d2), a);
177                     }
178                     hbond[i] = TRUE;
179                     bHB      = TRUE;
180                 }
181             }
182         }
183     }
184     return bHB;
185 }
186
187 static void calc_ringh(rvec xattach, rvec xb, rvec xc, rvec xh)
188 {
189     rvec tab, tac;
190     real n;
191
192     /* Add a proton on a ring to atom attach at distance 0.1 nm */
193     rvec_sub(xattach, xb, tab);
194     rvec_sub(xattach, xc, tac);
195     rvec_add(tab, tac, xh);
196     n = 0.1/norm(xh);
197     svmul(n, xh, xh);
198     rvec_inc(xh, xattach);
199 }
200
201 void set_histp(t_atoms *pdba, rvec *x, real angle, real dist)
202 {
203     static const char *prot_acc[] = {
204         "O", "OD1", "OD2", "OE1", "OE2", "OG", "OG1", "OH", "OW"
205     };
206 #define NPA asize(prot_acc)
207     static const char *prot_don[] = {
208         "N", "NH1", "NH2", "NE", "ND1", "ND2", "NE2", "NZ", "OG", "OG1", "OH", "NE1", "OW"
209     };
210 #define NPD asize(prot_don)
211
212     gmx_bool *donor, *acceptor;
213     gmx_bool *hbond, bHaveH = FALSE;
214     gmx_bool  bHDd, bHEd;
215     rvec      xh1, xh2;
216     int       natom;
217     int       i, j, nd, na, aj, hisind, his0, type = -1;
218     int       nd1, ne2, cg, cd2, ce1;
219     t_blocka *hb;
220     real      d;
221     char     *atomnm;
222
223     natom = pdba->nr;
224
225     i = 0;
226     while (i < natom &&
227            gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
228     {
229         i++;
230     }
231     if (natom == i)
232     {
233         return;
234     }
235
236     /* A histidine residue exists that requires automated assignment, so
237      * doing the analysis of donors and acceptors is worthwhile. */
238     fprintf(stderr,
239             "Analysing hydrogen-bonding network for automated assignment of histidine\n"
240             " protonation.");
241
242     snew(donor, natom);
243     snew(acceptor, natom);
244     snew(hbond, natom);
245     snew(hb, 1);
246
247     nd = na = 0;
248     for (j = 0; (j < natom); j++)
249     {
250         if (in_strings(*pdba->atomname[j], NPA, prot_acc) != -1)
251         {
252             acceptor[j] = TRUE;
253             na++;
254         }
255         if (in_strings(*pdba->atomname[j], NPD, prot_don) != -1)
256         {
257             donor[j] = TRUE;
258             nd++;
259         }
260     }
261     fprintf(stderr, " %d donors and %d acceptors were found.\n", nd, na);
262     chk_allhb(pdba, x, hb, donor, acceptor, dist);
263     if (debug)
264     {
265         pr_hbonds(debug, hb, pdba);
266     }
267     fprintf(stderr, "There are %d hydrogen bonds\n", hb->nra);
268
269     /* Now do the HIS stuff */
270     hisind = -1;
271     while (i < natom)
272     {
273         if (gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
274         {
275             i++;
276         }
277         else
278         {
279             if (pdba->atom[i].resind != hisind)
280             {
281                 hisind = pdba->atom[i].resind;
282
283                 /* Find the  atoms in the ring */
284                 nd1 = ne2 = cg = cd2 = ce1 = -1;
285                 while (i < natom && pdba->atom[i].resind == hisind)
286                 {
287                     atomnm = *pdba->atomname[i];
288                     if (strcmp(atomnm, "CD2") == 0)
289                     {
290                         cd2 = i;
291                     }
292                     else if (strcmp(atomnm, "CG") == 0)
293                     {
294                         cg  = i;
295                     }
296                     else if (strcmp(atomnm, "CE1") == 0)
297                     {
298                         ce1 = i;
299                     }
300                     else if (strcmp(atomnm, "ND1") == 0)
301                     {
302                         nd1 = i;
303                     }
304                     else if (strcmp(atomnm, "NE2") == 0)
305                     {
306                         ne2 = i;
307                     }
308
309                     i++;
310                 }
311
312                 if (!((cg == -1 ) || (cd2 == -1) || (ce1 == -1) ||
313                       (nd1 == -1) || (ne2 == -1)))
314                 {
315                     calc_ringh(x[nd1], x[cg], x[ce1], xh1);
316                     calc_ringh(x[ne2], x[ce1], x[cd2], xh2);
317
318                     bHDd = chk_hbonds(nd1, pdba, x, acceptor, hbond, xh1, angle, dist);
319                     chk_hbonds(nd1, pdba, x, donor, hbond, xh1, angle, dist);
320                     bHEd = chk_hbonds(ne2, pdba, x, acceptor, hbond, xh2, angle, dist);
321                     chk_hbonds(ne2, pdba, x, donor, hbond, xh2, angle, dist);
322
323                     if (bHDd)
324                     {
325                         if (bHEd)
326                         {
327                             type = ehisH;
328                         }
329                         else
330                         {
331                             type = ehisA;
332                         }
333                     }
334                     else
335                     {
336                         type = ehisB;
337                     }
338                     fprintf(stderr, "Will use %s for residue %d\n",
339                             hh[type], pdba->resinfo[hisind].nr);
340                 }
341                 else
342                 {
343                     gmx_fatal(FARGS, "Incomplete ring in HIS%d",
344                               pdba->resinfo[hisind].nr);
345                 }
346
347                 snew(pdba->resinfo[hisind].rtp, 1);
348                 *pdba->resinfo[hisind].rtp = strdup(hh[type]);
349             }
350         }
351     }
352     done_blocka(hb);
353     sfree(hb);
354     sfree(donor);
355     sfree(acceptor);
356     sfree(hbond);
357 }