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