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.
48 #include "gromacs/fileio/pdbio.h"
49 #include "gromacs/gmxpreprocess/pdb2top.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/smalloc.h"
54 #include "gromacs/utility/strdb.h"
62 c = toupper(fgetc(stdin));
64 while ((c != 'Y') && (c != 'N'));
69 t_specbond *get_specbonds(int *nspecbond)
71 const char *sbfile = "specbond.dat";
73 t_specbond *sb = nullptr;
74 char r1buf[32], r2buf[32], a1buf[32], a2buf[32], nr1buf[32], nr2buf[32];
80 nlines = get_lines(sbfile, &lines);
87 for (i = 0; (i < nlines); i++)
89 if (sscanf(lines[i], "%s%s%d%s%s%d%lf%s%s",
90 r1buf, a1buf, &nb1, r2buf, a2buf, &nb2, &length, nr1buf, nr2buf) != 9)
92 fprintf(stderr, "Invalid line '%s' in %s\n", lines[i], sbfile);
96 sb[n].res1 = gmx_strdup(r1buf);
97 sb[n].res2 = gmx_strdup(r2buf);
98 sb[n].newres1 = gmx_strdup(nr1buf);
99 sb[n].newres2 = gmx_strdup(nr2buf);
100 sb[n].atom1 = gmx_strdup(a1buf);
101 sb[n].atom2 = gmx_strdup(a2buf);
104 sb[n].length = length;
113 fprintf(stderr, "%d out of %d lines of %s converted successfully\n",
121 void done_specbonds(int nsb, t_specbond sb[])
125 for (i = 0; (i < nsb); i++)
131 sfree(sb[i].newres1);
132 sfree(sb[i].newres2);
136 static bool is_special(int nsb, t_specbond sb[], char *res, char *atom)
140 for (i = 0; (i < nsb); i++)
142 if (((strncmp(sb[i].res1, res, 3) == 0) &&
143 (gmx_strcasecmp(sb[i].atom1, atom) == 0)) ||
144 ((strncmp(sb[i].res2, res, 3) == 0) &&
145 (gmx_strcasecmp(sb[i].atom2, atom) == 0)))
153 static bool is_bond(int nsb, t_specbond sb[], t_atoms *pdba, int a1, int a2,
154 real d, int *index_sb, bool *bSwap)
157 char *at1, *at2, *res1, *res2;
159 at1 = *pdba->atomname[a1];
160 at2 = *pdba->atomname[a2];
161 res1 = *pdba->resinfo[pdba->atom[a1].resind].name;
162 res2 = *pdba->resinfo[pdba->atom[a2].resind].name;
164 for (i = 0; (i < nsb); i++)
167 if (((strncmp(sb[i].res1, res1, 3) == 0) &&
168 (gmx_strcasecmp(sb[i].atom1, at1) == 0) &&
169 (strncmp(sb[i].res2, res2, 3) == 0) &&
170 (gmx_strcasecmp(sb[i].atom2, at2) == 0)))
173 if ((0.9*sb[i].length < d) && (1.1*sb[i].length > d))
178 if (((strncmp(sb[i].res1, res2, 3) == 0) &&
179 (gmx_strcasecmp(sb[i].atom1, at2) == 0) &&
180 (strncmp(sb[i].res2, res1, 3) == 0) &&
181 (gmx_strcasecmp(sb[i].atom2, at1) == 0)))
184 if ((0.9*sb[i].length < d) && (1.1*sb[i].length > d))
193 static void rename_1res(t_atoms *pdba, int resind, char *newres, bool bVerbose)
197 printf("Using rtp entry %s for %s %d\n",
199 *pdba->resinfo[resind].name,
200 pdba->resinfo[resind].nr);
202 /* this used to free *resname, which messes up the symtab! */
203 snew(pdba->resinfo[resind].rtp, 1);
204 *pdba->resinfo[resind].rtp = gmx_strdup(newres);
207 int mk_specbonds(t_atoms *pdba, rvec x[], bool bInteractive,
208 t_ssbond **specbonds, bool bVerbose)
210 t_specbond *sb = nullptr;
211 t_ssbond *bonds = nullptr;
217 int ai, aj, index_sb;
222 sb = get_specbonds(&nsb);
226 snew(specp, pdba->nr);
230 for (i = 0; (i < pdba->nr); i++)
232 /* Check if this atom is special and if it is not a double atom
233 * in the input that still needs to be removed.
235 if (is_special(nsb, sb, *pdba->resinfo[pdba->atom[i].resind].name,
236 *pdba->atomname[i]) &&
238 pdba->atom[sgp[nspec-1]].resind == pdba->atom[i].resind &&
239 gmx_strcasecmp(*pdba->atomname[sgp[nspec-1]],
240 *pdba->atomname[i]) == 0))
242 specp[nspec] = pdba->atom[i].resind;
247 /* distance matrix d[nspec][nspec] */
249 for (i = 0; (i < nspec); i++)
254 for (i = 0; (i < nspec); i++)
257 for (j = 0; (j < nspec); j++)
260 d[i][j] = std::sqrt(distance2(x[ai], x[aj]));
266 fprintf(stderr, "Special Atom Distance matrix:\n");
267 for (b = 0; (b < nspec); b += MAXCOL)
269 /* print resname/number column headings */
270 fprintf(stderr, "%8s%8s", "", "");
271 e = std::min(b+MAXCOL, nspec-1);
272 for (i = b; (i < e); i++)
274 sprintf(buf, "%s%d", *pdba->resinfo[pdba->atom[sgp[i]].resind].name,
275 pdba->resinfo[specp[i]].nr);
276 fprintf(stderr, "%8s", buf);
278 fprintf(stderr, "\n");
279 /* print atomname/number column headings */
280 fprintf(stderr, "%8s%8s", "", "");
281 e = std::min(b+MAXCOL, nspec-1);
282 for (i = b; (i < e); i++)
284 sprintf(buf, "%s%d", *pdba->atomname[sgp[i]], sgp[i]+1);
285 fprintf(stderr, "%8s", buf);
287 fprintf(stderr, "\n");
289 e = std::min(b+MAXCOL, nspec);
290 for (i = b+1; (i < nspec); i++)
292 sprintf(buf, "%s%d", *pdba->resinfo[pdba->atom[sgp[i]].resind].name,
293 pdba->resinfo[specp[i]].nr);
294 fprintf(stderr, "%8s", buf);
295 sprintf(buf, "%s%d", *pdba->atomname[sgp[i]], sgp[i]+1);
296 fprintf(stderr, "%8s", buf);
298 for (j = b; (j < e2); j++)
300 fprintf(stderr, " %7.3f", d[i][j]);
302 fprintf(stderr, "\n");
309 for (i = 0; (i < nspec); i++)
312 for (j = i+1; (j < nspec); j++)
315 /* Ensure creation of at most nspec special bonds to avoid overflowing bonds[] */
316 if (nbonds < nspec && is_bond(nsb, sb, pdba, ai, aj, d[i][j], &index_sb, &bSwap))
318 fprintf(stderr, "%s %s-%d %s-%d and %s-%d %s-%d%s",
319 bInteractive ? "Link" : "Linking",
320 *pdba->resinfo[pdba->atom[ai].resind].name,
321 pdba->resinfo[specp[i]].nr,
322 *pdba->atomname[ai], ai+1,
323 *pdba->resinfo[pdba->atom[aj].resind].name,
324 pdba->resinfo[specp[j]].nr,
325 *pdba->atomname[aj], aj+1,
326 bInteractive ? " (y/n) ?" : "...\n");
327 bDoit = bInteractive ? yesno() : TRUE;
331 /* Store the residue numbers in the bonds array */
332 bonds[nbonds].res1 = specp[i];
333 bonds[nbonds].res2 = specp[j];
334 bonds[nbonds].a1 = gmx_strdup(*pdba->atomname[ai]);
335 bonds[nbonds].a2 = gmx_strdup(*pdba->atomname[aj]);
336 /* rename residues */
339 rename_1res(pdba, specp[i], sb[index_sb].newres2, bVerbose);
340 rename_1res(pdba, specp[j], sb[index_sb].newres1, bVerbose);
344 rename_1res(pdba, specp[i], sb[index_sb].newres1, bVerbose);
345 rename_1res(pdba, specp[j], sb[index_sb].newres2, bVerbose);
353 for (i = 0; (i < nspec); i++)
361 done_specbonds(nsb, sb);