3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
42 #include "gromacs/fileio/pdbio.h"
57 c = toupper(fgetc(stdin));
59 while ((c != 'Y') && (c != 'N'));
64 t_specbond *get_specbonds(int *nspecbond)
66 const char *sbfile = "specbond.dat";
68 t_specbond *sb = NULL;
69 char r1buf[32], r2buf[32], a1buf[32], a2buf[32], nr1buf[32], nr2buf[32];
75 nlines = get_lines(sbfile, &lines);
82 for (i = 0; (i < nlines); i++)
84 if (sscanf(lines[i], "%s%s%d%s%s%d%lf%s%s",
85 r1buf, a1buf, &nb1, r2buf, a2buf, &nb2, &length, nr1buf, nr2buf) != 9)
87 fprintf(stderr, "Invalid line '%s' in %s\n", lines[i], sbfile);
91 sb[n].res1 = strdup(r1buf);
92 sb[n].res2 = strdup(r2buf);
93 sb[n].newres1 = strdup(nr1buf);
94 sb[n].newres2 = strdup(nr2buf);
95 sb[n].atom1 = strdup(a1buf);
96 sb[n].atom2 = strdup(a2buf);
99 sb[n].length = length;
108 fprintf(stderr, "%d out of %d lines of %s converted successfully\n",
116 void done_specbonds(int nsb, t_specbond sb[])
120 for (i = 0; (i < nsb); i++)
126 sfree(sb[i].newres1);
127 sfree(sb[i].newres2);
131 static gmx_bool is_special(int nsb, t_specbond sb[], char *res, char *atom)
135 for (i = 0; (i < nsb); i++)
137 if (((strncmp(sb[i].res1, res, 3) == 0) &&
138 (gmx_strcasecmp(sb[i].atom1, atom) == 0)) ||
139 ((strncmp(sb[i].res2, res, 3) == 0) &&
140 (gmx_strcasecmp(sb[i].atom2, atom) == 0)))
148 static gmx_bool is_bond(int nsb, t_specbond sb[], t_atoms *pdba, int a1, int a2,
149 real d, int *index_sb, gmx_bool *bSwap)
152 char *at1, *at2, *res1, *res2;
154 at1 = *pdba->atomname[a1];
155 at2 = *pdba->atomname[a2];
156 res1 = *pdba->resinfo[pdba->atom[a1].resind].name;
157 res2 = *pdba->resinfo[pdba->atom[a2].resind].name;
161 fprintf(stderr, "Checking %s-%d %s-%d and %s-%d %s-%d: %g ",
162 res1, pdba->resinfo[pdba->atom[a1].resind].nr, at1, a1+1,
163 res2, pdba->resinfo[pdba->atom[a2].resind].nr, at2, a2+1, d);
166 for (i = 0; (i < nsb); i++)
169 if (((strncmp(sb[i].res1, res1, 3) == 0) &&
170 (gmx_strcasecmp(sb[i].atom1, at1) == 0) &&
171 (strncmp(sb[i].res2, res2, 3) == 0) &&
172 (gmx_strcasecmp(sb[i].atom2, at2) == 0)))
175 if ((0.9*sb[i].length < d) && (1.1*sb[i].length > d))
179 fprintf(stderr, "%g\n", sb[i].length);
184 if (((strncmp(sb[i].res1, res2, 3) == 0) &&
185 (gmx_strcasecmp(sb[i].atom1, at2) == 0) &&
186 (strncmp(sb[i].res2, res1, 3) == 0) &&
187 (gmx_strcasecmp(sb[i].atom2, at1) == 0)))
190 if ((0.9*sb[i].length < d) && (1.1*sb[i].length > d))
194 fprintf(stderr, "%g\n", sb[i].length);
202 fprintf(stderr, "\n");
207 static void rename_1res(t_atoms *pdba, int resind, char *newres, gmx_bool bVerbose)
211 printf("Using rtp entry %s for %s %d\n",
213 *pdba->resinfo[resind].name,
214 pdba->resinfo[resind].nr);
216 /* this used to free *resname, which fucks up the symtab! */
217 snew(pdba->resinfo[resind].rtp, 1);
218 *pdba->resinfo[resind].rtp = strdup(newres);
221 int mk_specbonds(t_atoms *pdba, rvec x[], gmx_bool bInteractive,
222 t_ssbond **specbonds, gmx_bool bVerbose)
224 t_specbond *sb = NULL;
225 t_ssbond *bonds = NULL;
229 gmx_bool bDoit, bSwap;
231 int ai, aj, index_sb;
236 sb = get_specbonds(&nsb);
240 snew(specp, pdba->nr);
244 for (i = 0; (i < pdba->nr); i++)
246 /* Check if this atom is special and if it is not a double atom
247 * in the input that still needs to be removed.
249 if (is_special(nsb, sb, *pdba->resinfo[pdba->atom[i].resind].name,
250 *pdba->atomname[i]) &&
252 pdba->atom[sgp[nspec-1]].resind == pdba->atom[i].resind &&
253 gmx_strcasecmp(*pdba->atomname[sgp[nspec-1]],
254 *pdba->atomname[i]) == 0))
256 specp[nspec] = pdba->atom[i].resind;
261 /* distance matrix d[nspec][nspec] */
263 for (i = 0; (i < nspec); i++)
268 for (i = 0; (i < nspec); i++)
271 for (j = 0; (j < nspec); j++)
274 d[i][j] = sqrt(distance2(x[ai], x[aj]));
280 fprintf(stderr, "Special Atom Distance matrix:\n");
281 for (b = 0; (b < nspec); b += MAXCOL)
283 /* print resname/number column headings */
284 fprintf(stderr, "%8s%8s", "", "");
285 e = min(b+MAXCOL, nspec-1);
286 for (i = b; (i < e); i++)
288 sprintf(buf, "%s%d", *pdba->resinfo[pdba->atom[sgp[i]].resind].name,
289 pdba->resinfo[specp[i]].nr);
290 fprintf(stderr, "%8s", buf);
292 fprintf(stderr, "\n");
293 /* print atomname/number column headings */
294 fprintf(stderr, "%8s%8s", "", "");
295 e = min(b+MAXCOL, nspec-1);
296 for (i = b; (i < e); i++)
298 sprintf(buf, "%s%d", *pdba->atomname[sgp[i]], sgp[i]+1);
299 fprintf(stderr, "%8s", buf);
301 fprintf(stderr, "\n");
303 e = min(b+MAXCOL, nspec);
304 for (i = b+1; (i < nspec); i++)
306 sprintf(buf, "%s%d", *pdba->resinfo[pdba->atom[sgp[i]].resind].name,
307 pdba->resinfo[specp[i]].nr);
308 fprintf(stderr, "%8s", buf);
309 sprintf(buf, "%s%d", *pdba->atomname[sgp[i]], sgp[i]+1);
310 fprintf(stderr, "%8s", buf);
312 for (j = b; (j < e2); j++)
314 fprintf(stderr, " %7.3f", d[i][j]);
316 fprintf(stderr, "\n");
323 for (i = 0; (i < nspec); i++)
326 for (j = i+1; (j < nspec); j++)
329 if (is_bond(nsb, sb, pdba, ai, aj, d[i][j], &index_sb, &bSwap))
331 fprintf(stderr, "%s %s-%d %s-%d and %s-%d %s-%d%s",
332 bInteractive ? "Link" : "Linking",
333 *pdba->resinfo[pdba->atom[ai].resind].name,
334 pdba->resinfo[specp[i]].nr,
335 *pdba->atomname[ai], ai+1,
336 *pdba->resinfo[pdba->atom[aj].resind].name,
337 pdba->resinfo[specp[j]].nr,
338 *pdba->atomname[aj], aj+1,
339 bInteractive ? " (y/n) ?" : "...\n");
340 bDoit = bInteractive ? yesno() : TRUE;
344 /* Store the residue numbers in the bonds array */
345 bonds[nbonds].res1 = specp[i];
346 bonds[nbonds].res2 = specp[j];
347 bonds[nbonds].a1 = strdup(*pdba->atomname[ai]);
348 bonds[nbonds].a2 = strdup(*pdba->atomname[aj]);
349 /* rename residues */
352 rename_1res(pdba, specp[i], sb[index_sb].newres2, bVerbose);
353 rename_1res(pdba, specp[j], sb[index_sb].newres1, bVerbose);
357 rename_1res(pdba, specp[i], sb[index_sb].newres1, bVerbose);
358 rename_1res(pdba, specp[j], sb[index_sb].newres2, bVerbose);
366 for (i = 0; (i < nspec); i++)
374 done_specbonds(nsb, sb);