1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
48 #include "gromacs/fileio/futil.h"
49 #include "fflibutil.h"
50 #include "hackblock.h"
51 #include "gmx_fatal.h"
61 static void get_xlatoms(const char *fn, FILE *fp,
62 int *nptr, t_xlate_atom **xlptr)
64 char filebase[STRLEN];
66 char abuf[1024], rbuf[1024], repbuf[1024], dumbuf[1024];
71 fflib_filename_base(fn, filebase, STRLEN);
76 while (get_a_line(fp, line, STRLEN))
78 na = sscanf(line, "%s%s%s%s", rbuf, abuf, repbuf, dumbuf);
79 /* Check if we are reading an old format file with the number of items
82 if (na == 1 && n == *nptr && sscanf(rbuf, "%d", &idum) == 1)
88 gmx_fatal(FARGS, "Expected a residue name and two atom names in file '%s', not '%s'", fn, line);
92 xl[n].filebase = strdup(filebase);
94 /* Use wildcards... */
95 if (strcmp(rbuf, "*") != 0)
97 xl[n].res = strdup(rbuf);
104 /* Replace underscores in the string by spaces */
105 while ((_ptr = strchr(abuf, '_')) != 0)
110 xl[n].atom = strdup(abuf);
111 xl[n].replace = strdup(repbuf);
119 static void done_xlatom(int nxlate, t_xlate_atom *xlatom)
123 for (i = 0; (i < nxlate); i++)
125 sfree(xlatom[i].filebase);
126 if (xlatom[i].res != NULL)
128 sfree(xlatom[i].res);
130 sfree(xlatom[i].atom);
131 sfree(xlatom[i].replace);
136 void rename_atoms(const char *xlfile, const char *ffdir,
137 t_atoms *atoms, t_symtab *symtab, const t_restp *restp,
138 gmx_bool bResname, gmx_residuetype_t rt, gmx_bool bReorderNum,
142 int nxlate, a, i, resind;
143 t_xlate_atom *xlatom;
146 char c, *rnm, atombuf[32], *ptr0, *ptr1;
147 gmx_bool bReorderedNum, bRenamed, bMatch;
153 fp = libopen(xlfile);
154 get_xlatoms(xlfile, fp, &nxlate, &xlatom);
159 nf = fflib_search_file_end(ffdir, ".arn", FALSE, &f);
160 for (i = 0; i < nf; i++)
162 fp = fflib_open(f[i]);
163 get_xlatoms(f[i], fp, &nxlate, &xlatom);
170 for (a = 0; (a < atoms->nr); a++)
172 resind = atoms->atom[a].resind;
175 rnm = *(atoms->resinfo[resind].name);
179 rnm = *(atoms->resinfo[resind].rtp);
182 strcpy(atombuf, *(atoms->atomname[a]));
183 bReorderedNum = FALSE;
186 if (isdigit(atombuf[0]))
189 for (i = 0; ((size_t)i < strlen(atombuf)-1); i++)
191 atombuf[i] = atombuf[i+1];
194 bReorderedNum = TRUE;
198 for (i = 0; (i < nxlate) && !bRenamed; i++)
200 /* Check if the base file name of the rtp and arn entry match */
202 gmx_strcasecmp(restp[resind].filebase, xlatom[i].filebase) == 0)
204 /* Match the residue name */
205 bMatch = (xlatom[i].res == NULL ||
206 (gmx_strcasecmp("protein", xlatom[i].res) == 0 &&
207 gmx_residuetype_is_protein(rt, rnm)) ||
208 (gmx_strcasecmp("DNA", xlatom[i].res) == 0 &&
209 gmx_residuetype_is_dna(rt, rnm)) ||
210 (gmx_strcasecmp("RNA", xlatom[i].res) == 0 &&
211 gmx_residuetype_is_rna(rt, rnm)));
215 ptr1 = xlatom[i].res;
216 while (ptr0[0] != '\0' && ptr1[0] != '\0' &&
217 (ptr0[0] == ptr1[0] || ptr1[0] == '?'))
222 bMatch = (ptr0[0] == '\0' && ptr1[0] == '\0');
224 if (bMatch && strcmp(atombuf, xlatom[i].atom) == 0)
226 /* We have a match. */
227 /* Don't free the old atomname,
228 * since it might be in the symtab.
230 ptr0 = strdup(xlatom[i].replace);
233 printf("Renaming atom '%s' in residue %d %s to '%s'\n",
235 atoms->resinfo[resind].nr,
236 *atoms->resinfo[resind].name,
239 atoms->atomname[a] = put_symtab(symtab, ptr0);
244 if (bReorderedNum && !bRenamed)
246 atoms->atomname[a] = put_symtab(symtab, atombuf);
250 done_xlatom(nxlate, xlatom);