Code beautification with uncrustify
[alexxy/gromacs.git] / src / programs / pdb2gmx / specbond.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 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <ctype.h>
40 #include <math.h>
41 #include "typedefs.h"
42 #include "pdbio.h"
43 #include "strdb.h"
44 #include "string2.h"
45 #include "smalloc.h"
46 #include "specbond.h"
47 #include "pdb2top.h"
48 #include "vec.h"
49 #include "macros.h"
50
51 gmx_bool yesno(void)
52 {
53     char c;
54
55     do
56     {
57         c = toupper(fgetc(stdin));
58     }
59     while ((c != 'Y') && (c != 'N'));
60
61     return (c == 'Y');
62 }
63
64 t_specbond *get_specbonds(int *nspecbond)
65 {
66     const char  *sbfile = "specbond.dat";
67
68     t_specbond  *sb = NULL;
69     char         r1buf[32], r2buf[32], a1buf[32], a2buf[32], nr1buf[32], nr2buf[32];
70     double       length;
71     int          nb1, nb2;
72     char       **lines;
73     int          nlines, i, n;
74
75     nlines = get_lines(sbfile, &lines);
76     if (nlines > 0)
77     {
78         snew(sb, nlines);
79     }
80
81     n = 0;
82     for (i = 0; (i < nlines); i++)
83     {
84         if (sscanf(lines[i], "%s%s%d%s%s%d%lf%s%s",
85                    r1buf, a1buf, &nb1, r2buf, a2buf, &nb2, &length, nr1buf, nr2buf) != 9)
86         {
87             fprintf(stderr, "Invalid line '%s' in %s\n", lines[i], sbfile);
88         }
89         else
90         {
91             sb[n].res1    = strdup(r1buf);
92             sb[n].res2    = strdup(r2buf);
93             sb[n].newres1 = strdup(nr1buf);
94             sb[n].newres2 = strdup(nr2buf);
95             sb[n].atom1   = strdup(a1buf);
96             sb[n].atom2   = strdup(a2buf);
97             sb[n].nbond1  = nb1;
98             sb[n].nbond2  = nb2;
99             sb[n].length  = length;
100             n++;
101         }
102         sfree(lines[i]);
103     }
104     if (nlines > 0)
105     {
106         sfree(lines);
107     }
108     fprintf(stderr, "%d out of %d lines of %s converted successfully\n",
109             n, nlines, sbfile);
110
111     *nspecbond = n;
112
113     return sb;
114 }
115
116 void done_specbonds(int nsb, t_specbond sb[])
117 {
118     int i;
119
120     for (i = 0; (i < nsb); i++)
121     {
122         sfree(sb[i].res1);
123         sfree(sb[i].res2);
124         sfree(sb[i].atom1);
125         sfree(sb[i].atom2);
126         sfree(sb[i].newres1);
127         sfree(sb[i].newres2);
128     }
129 }
130
131 static gmx_bool is_special(int nsb, t_specbond sb[], char *res, char *atom)
132 {
133     int i;
134
135     for (i = 0; (i < nsb); i++)
136     {
137         if (((strncmp(sb[i].res1, res, 3) == 0) &&
138              (gmx_strcasecmp(sb[i].atom1, atom) == 0)) ||
139             ((strncmp(sb[i].res2, res, 3) == 0) &&
140              (gmx_strcasecmp(sb[i].atom2, atom) == 0)))
141         {
142             return TRUE;
143         }
144     }
145     return FALSE;
146 }
147
148 static gmx_bool is_bond(int nsb, t_specbond sb[], t_atoms *pdba, int a1, int a2,
149                         real d, int *index_sb, gmx_bool *bSwap)
150 {
151     int   i;
152     char *at1, *at2, *res1, *res2;
153
154     at1  = *pdba->atomname[a1];
155     at2  = *pdba->atomname[a2];
156     res1 = *pdba->resinfo[pdba->atom[a1].resind].name;
157     res2 = *pdba->resinfo[pdba->atom[a2].resind].name;
158
159     if (debug)
160     {
161         fprintf(stderr, "Checking %s-%d %s-%d and %s-%d %s-%d: %g ",
162                 res1, pdba->resinfo[pdba->atom[a1].resind].nr, at1, a1+1,
163                 res2, pdba->resinfo[pdba->atom[a2].resind].nr, at2, a2+1, d);
164     }
165
166     for (i = 0; (i < nsb); i++)
167     {
168         *index_sb = i;
169         if (((strncmp(sb[i].res1, res1, 3) == 0)  &&
170              (gmx_strcasecmp(sb[i].atom1, at1) == 0) &&
171              (strncmp(sb[i].res2, res2, 3) == 0)  &&
172              (gmx_strcasecmp(sb[i].atom2, at2) == 0)))
173         {
174             *bSwap = FALSE;
175             if ((0.9*sb[i].length < d) && (1.1*sb[i].length > d))
176             {
177                 if (debug)
178                 {
179                     fprintf(stderr, "%g\n", sb[i].length);
180                 }
181                 return TRUE;
182             }
183         }
184         if (((strncmp(sb[i].res1, res2, 3) == 0)  &&
185              (gmx_strcasecmp(sb[i].atom1, at2) == 0) &&
186              (strncmp(sb[i].res2, res1, 3) == 0)  &&
187              (gmx_strcasecmp(sb[i].atom2, at1) == 0)))
188         {
189             *bSwap = TRUE;
190             if ((0.9*sb[i].length < d) && (1.1*sb[i].length > d))
191             {
192                 if (debug)
193                 {
194                     fprintf(stderr, "%g\n", sb[i].length);
195                 }
196                 return TRUE;
197             }
198         }
199     }
200     if (debug)
201     {
202         fprintf(stderr, "\n");
203     }
204     return FALSE;
205 }
206
207 static void rename_1res(t_atoms *pdba, int resind, char *newres, gmx_bool bVerbose)
208 {
209     if (bVerbose)
210     {
211         printf("Using rtp entry %s for %s %d\n",
212                newres,
213                *pdba->resinfo[resind].name,
214                pdba->resinfo[resind].nr);
215     }
216     /* this used to free *resname, which fucks up the symtab! */
217     snew(pdba->resinfo[resind].rtp, 1);
218     *pdba->resinfo[resind].rtp = strdup(newres);
219 }
220
221 int mk_specbonds(t_atoms *pdba, rvec x[], gmx_bool bInteractive,
222                  t_ssbond **specbonds, gmx_bool bVerbose)
223 {
224     t_specbond *sb    = NULL;
225     t_ssbond   *bonds = NULL;
226     int         nsb;
227     int         nspec, nbonds;
228     int        *specp, *sgp;
229     gmx_bool    bDoit, bSwap;
230     int         i, j, b, e, e2;
231     int         ai, aj, index_sb;
232     real      **d;
233     char        buf[10];
234
235     nbonds = 0;
236     sb     = get_specbonds(&nsb);
237
238     if (nsb > 0)
239     {
240         snew(specp, pdba->nr);
241         snew(sgp, pdba->nr);
242
243         nspec = 0;
244         for (i = 0; (i < pdba->nr); i++)
245         {
246             /* Check if this atom is special and if it is not a double atom
247              * in the input that still needs to be removed.
248              */
249             if (is_special(nsb, sb, *pdba->resinfo[pdba->atom[i].resind].name,
250                            *pdba->atomname[i]) &&
251                 !(nspec > 0 &&
252                   pdba->atom[sgp[nspec-1]].resind == pdba->atom[i].resind &&
253                   gmx_strcasecmp(*pdba->atomname[sgp[nspec-1]],
254                                  *pdba->atomname[i]) == 0))
255             {
256                 specp[nspec] = pdba->atom[i].resind;
257                 sgp[nspec]   = i;
258                 nspec++;
259             }
260         }
261         /* distance matrix d[nspec][nspec] */
262         snew(d, nspec);
263         for (i = 0; (i < nspec); i++)
264         {
265             snew(d[i], nspec);
266         }
267
268         for (i = 0; (i < nspec); i++)
269         {
270             ai = sgp[i];
271             for (j = 0; (j < nspec); j++)
272             {
273                 aj      = sgp[j];
274                 d[i][j] = sqrt(distance2(x[ai], x[aj]));
275             }
276         }
277         if (nspec > 1)
278         {
279 #define MAXCOL 7
280             fprintf(stderr, "Special Atom Distance matrix:\n");
281             for (b = 0; (b < nspec); b += MAXCOL)
282             {
283                 /* print resname/number column headings */
284                 fprintf(stderr, "%8s%8s", "", "");
285                 e = min(b+MAXCOL, nspec-1);
286                 for (i = b; (i < e); i++)
287                 {
288                     sprintf(buf, "%s%d", *pdba->resinfo[pdba->atom[sgp[i]].resind].name,
289                             pdba->resinfo[specp[i]].nr);
290                     fprintf(stderr, "%8s", buf);
291                 }
292                 fprintf(stderr, "\n");
293                 /* print atomname/number column headings */
294                 fprintf(stderr, "%8s%8s", "", "");
295                 e = min(b+MAXCOL, nspec-1);
296                 for (i = b; (i < e); i++)
297                 {
298                     sprintf(buf, "%s%d", *pdba->atomname[sgp[i]], sgp[i]+1);
299                     fprintf(stderr, "%8s", buf);
300                 }
301                 fprintf(stderr, "\n");
302                 /* print matrix */
303                 e = min(b+MAXCOL, nspec);
304                 for (i = b+1; (i < nspec); i++)
305                 {
306                     sprintf(buf, "%s%d", *pdba->resinfo[pdba->atom[sgp[i]].resind].name,
307                             pdba->resinfo[specp[i]].nr);
308                     fprintf(stderr, "%8s", buf);
309                     sprintf(buf, "%s%d", *pdba->atomname[sgp[i]], sgp[i]+1);
310                     fprintf(stderr, "%8s", buf);
311                     e2 = min(i, e);
312                     for (j = b; (j < e2); j++)
313                     {
314                         fprintf(stderr, " %7.3f", d[i][j]);
315                     }
316                     fprintf(stderr, "\n");
317                 }
318             }
319         }
320
321         snew(bonds, nspec);
322
323         for (i = 0; (i < nspec); i++)
324         {
325             ai = sgp[i];
326             for (j = i+1; (j < nspec); j++)
327             {
328                 aj = sgp[j];
329                 if (is_bond(nsb, sb, pdba, ai, aj, d[i][j], &index_sb, &bSwap))
330                 {
331                     fprintf(stderr, "%s %s-%d %s-%d and %s-%d %s-%d%s",
332                             bInteractive ? "Link" : "Linking",
333                             *pdba->resinfo[pdba->atom[ai].resind].name,
334                             pdba->resinfo[specp[i]].nr,
335                             *pdba->atomname[ai], ai+1,
336                             *pdba->resinfo[pdba->atom[aj].resind].name,
337                             pdba->resinfo[specp[j]].nr,
338                             *pdba->atomname[aj], aj+1,
339                             bInteractive ? " (y/n) ?" : "...\n");
340                     bDoit = bInteractive ? yesno() : TRUE;
341
342                     if (bDoit)
343                     {
344                         /* Store the residue numbers in the bonds array */
345                         bonds[nbonds].res1 = specp[i];
346                         bonds[nbonds].res2 = specp[j];
347                         bonds[nbonds].a1   = strdup(*pdba->atomname[ai]);
348                         bonds[nbonds].a2   = strdup(*pdba->atomname[aj]);
349                         /* rename residues */
350                         if (bSwap)
351                         {
352                             rename_1res(pdba, specp[i], sb[index_sb].newres2, bVerbose);
353                             rename_1res(pdba, specp[j], sb[index_sb].newres1, bVerbose);
354                         }
355                         else
356                         {
357                             rename_1res(pdba, specp[i], sb[index_sb].newres1, bVerbose);
358                             rename_1res(pdba, specp[j], sb[index_sb].newres2, bVerbose);
359                         }
360                         nbonds++;
361                     }
362                 }
363             }
364         }
365
366         for (i = 0; (i < nspec); i++)
367         {
368             sfree(d[i]);
369         }
370         sfree(d);
371         sfree(sgp);
372         sfree(specp);
373
374         done_specbonds(nsb, sb);
375         sfree(sb);
376     }
377
378     *specbonds = bonds;
379
380     return nbonds;
381 }