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