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, 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.
44 #include "gromacs/fileio/pdbio.h"
45 #include "gromacs/fileio/strdb.h"
46 #include "gromacs/utility/cstringutil.h"
47 #include "gromacs/utility/smalloc.h"
50 #include "gromacs/math/vec.h"
59 c = toupper(fgetc(stdin));
61 while ((c != 'Y') && (c != 'N'));
66 t_specbond *get_specbonds(int *nspecbond)
68 const char *sbfile = "specbond.dat";
70 t_specbond *sb = NULL;
71 char r1buf[32], r2buf[32], a1buf[32], a2buf[32], nr1buf[32], nr2buf[32];
77 nlines = get_lines(sbfile, &lines);
84 for (i = 0; (i < nlines); i++)
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)
89 fprintf(stderr, "Invalid line '%s' in %s\n", lines[i], sbfile);
93 sb[n].res1 = strdup(r1buf);
94 sb[n].res2 = strdup(r2buf);
95 sb[n].newres1 = strdup(nr1buf);
96 sb[n].newres2 = strdup(nr2buf);
97 sb[n].atom1 = strdup(a1buf);
98 sb[n].atom2 = strdup(a2buf);
101 sb[n].length = length;
110 fprintf(stderr, "%d out of %d lines of %s converted successfully\n",
118 void done_specbonds(int nsb, t_specbond sb[])
122 for (i = 0; (i < nsb); i++)
128 sfree(sb[i].newres1);
129 sfree(sb[i].newres2);
133 static gmx_bool is_special(int nsb, t_specbond sb[], char *res, char *atom)
137 for (i = 0; (i < nsb); i++)
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)))
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)
154 char *at1, *at2, *res1, *res2;
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;
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);
168 for (i = 0; (i < nsb); 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)))
177 if ((0.9*sb[i].length < d) && (1.1*sb[i].length > d))
181 fprintf(stderr, "%g\n", sb[i].length);
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)))
192 if ((0.9*sb[i].length < d) && (1.1*sb[i].length > d))
196 fprintf(stderr, "%g\n", sb[i].length);
204 fprintf(stderr, "\n");
209 static void rename_1res(t_atoms *pdba, int resind, char *newres, gmx_bool bVerbose)
213 printf("Using rtp entry %s for %s %d\n",
215 *pdba->resinfo[resind].name,
216 pdba->resinfo[resind].nr);
218 /* this used to free *resname, which messes up the symtab! */
219 snew(pdba->resinfo[resind].rtp, 1);
220 *pdba->resinfo[resind].rtp = strdup(newres);
223 int mk_specbonds(t_atoms *pdba, rvec x[], gmx_bool bInteractive,
224 t_ssbond **specbonds, gmx_bool bVerbose)
226 t_specbond *sb = NULL;
227 t_ssbond *bonds = NULL;
231 gmx_bool bDoit, bSwap;
233 int ai, aj, index_sb;
238 sb = get_specbonds(&nsb);
242 snew(specp, pdba->nr);
246 for (i = 0; (i < pdba->nr); i++)
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.
251 if (is_special(nsb, sb, *pdba->resinfo[pdba->atom[i].resind].name,
252 *pdba->atomname[i]) &&
254 pdba->atom[sgp[nspec-1]].resind == pdba->atom[i].resind &&
255 gmx_strcasecmp(*pdba->atomname[sgp[nspec-1]],
256 *pdba->atomname[i]) == 0))
258 specp[nspec] = pdba->atom[i].resind;
263 /* distance matrix d[nspec][nspec] */
265 for (i = 0; (i < nspec); i++)
270 for (i = 0; (i < nspec); i++)
273 for (j = 0; (j < nspec); j++)
276 d[i][j] = sqrt(distance2(x[ai], x[aj]));
282 fprintf(stderr, "Special Atom Distance matrix:\n");
283 for (b = 0; (b < nspec); b += MAXCOL)
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++)
290 sprintf(buf, "%s%d", *pdba->resinfo[pdba->atom[sgp[i]].resind].name,
291 pdba->resinfo[specp[i]].nr);
292 fprintf(stderr, "%8s", buf);
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++)
300 sprintf(buf, "%s%d", *pdba->atomname[sgp[i]], sgp[i]+1);
301 fprintf(stderr, "%8s", buf);
303 fprintf(stderr, "\n");
305 e = min(b+MAXCOL, nspec);
306 for (i = b+1; (i < nspec); i++)
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);
314 for (j = b; (j < e2); j++)
316 fprintf(stderr, " %7.3f", d[i][j]);
318 fprintf(stderr, "\n");
325 for (i = 0; (i < nspec); i++)
328 for (j = i+1; (j < nspec); j++)
331 if (is_bond(nsb, sb, pdba, ai, aj, d[i][j], &index_sb, &bSwap))
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;
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 = strdup(*pdba->atomname[ai]);
350 bonds[nbonds].a2 = strdup(*pdba->atomname[aj]);
351 /* rename residues */
354 rename_1res(pdba, specp[i], sb[index_sb].newres2, bVerbose);
355 rename_1res(pdba, specp[j], sb[index_sb].newres1, bVerbose);
359 rename_1res(pdba, specp[i], sb[index_sb].newres1, bVerbose);
360 rename_1res(pdba, specp[j], sb[index_sb].newres2, bVerbose);
368 for (i = 0; (i < nspec); i++)
376 done_specbonds(nsb, sb);