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