2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
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"
61 c = toupper(fgetc(stdin));
63 while ((c != 'Y') && (c != 'N'));
68 t_specbond *get_specbonds(int *nspecbond)
70 const char *sbfile = "specbond.dat";
72 t_specbond *sb = nullptr;
73 char r1buf[32], r2buf[32], a1buf[32], a2buf[32], nr1buf[32], nr2buf[32];
79 nlines = get_lines(sbfile, &lines);
86 for (i = 0; (i < nlines); i++)
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)
91 fprintf(stderr, "Invalid line '%s' in %s\n", lines[i], sbfile);
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);
103 sb[n].length = length;
112 fprintf(stderr, "%d out of %d lines of %s converted successfully\n",
120 void done_specbonds(int nsb, t_specbond sb[])
124 for (i = 0; (i < nsb); i++)
130 sfree(sb[i].newres1);
131 sfree(sb[i].newres2);
135 static bool is_special(int nsb, t_specbond sb[], char *res, char *atom)
139 for (i = 0; (i < nsb); i++)
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)))
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)
156 char *at1, *at2, *res1, *res2;
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;
163 for (i = 0; (i < nsb); 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)))
172 if ((0.9*sb[i].length < d) && (1.1*sb[i].length > d))
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)))
183 if ((0.9*sb[i].length < d) && (1.1*sb[i].length > d))
192 static void rename_1res(t_atoms *pdba, int resind, char *newres, bool bVerbose)
196 printf("Using rtp entry %s for %s %d\n",
198 *pdba->resinfo[resind].name,
199 pdba->resinfo[resind].nr);
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);
206 int mk_specbonds(t_atoms *pdba, rvec x[], bool bInteractive,
207 t_ssbond **specbonds, bool bVerbose)
209 t_specbond *sb = nullptr;
210 t_ssbond *bonds = nullptr;
216 int ai, aj, index_sb;
221 sb = get_specbonds(&nsb);
225 snew(specp, pdba->nr);
229 for (i = 0; (i < pdba->nr); i++)
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.
234 if (is_special(nsb, sb, *pdba->resinfo[pdba->atom[i].resind].name,
235 *pdba->atomname[i]) &&
237 pdba->atom[sgp[nspec-1]].resind == pdba->atom[i].resind &&
238 gmx_strcasecmp(*pdba->atomname[sgp[nspec-1]],
239 *pdba->atomname[i]) == 0))
241 specp[nspec] = pdba->atom[i].resind;
246 /* distance matrix d[nspec][nspec] */
248 for (i = 0; (i < nspec); i++)
253 for (i = 0; (i < nspec); i++)
256 for (j = 0; (j < nspec); j++)
259 d[i][j] = std::sqrt(distance2(x[ai], x[aj]));
265 fprintf(stderr, "Special Atom Distance matrix:\n");
266 for (b = 0; (b < nspec); b += MAXCOL)
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++)
273 sprintf(buf, "%s%d", *pdba->resinfo[pdba->atom[sgp[i]].resind].name,
274 pdba->resinfo[specp[i]].nr);
275 fprintf(stderr, "%8s", buf);
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++)
283 sprintf(buf, "%s%d", *pdba->atomname[sgp[i]], sgp[i]+1);
284 fprintf(stderr, "%8s", buf);
286 fprintf(stderr, "\n");
288 e = std::min(b+MAXCOL, nspec);
289 for (i = b+1; (i < nspec); i++)
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 sprintf(buf, "%s%d", *pdba->atomname[sgp[i]], sgp[i]+1);
295 fprintf(stderr, "%8s", buf);
297 for (j = b; (j < e2); j++)
299 fprintf(stderr, " %7.3f", d[i][j]);
301 fprintf(stderr, "\n");
308 for (i = 0; (i < nspec); i++)
311 for (j = i+1; (j < nspec); 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))
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;
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 */
338 rename_1res(pdba, specp[i], sb[index_sb].newres2, bVerbose);
339 rename_1res(pdba, specp[j], sb[index_sb].newres1, bVerbose);
343 rename_1res(pdba, specp[i], sb[index_sb].newres1, bVerbose);
344 rename_1res(pdba, specp[j], sb[index_sb].newres2, bVerbose);
352 for (i = 0; (i < nspec); i++)
360 done_specbonds(nsb, sb);