Add AWH biasing module + tests
[gromacs.git] / src / gromacs / gmxpreprocess / xlate.cpp
blob760fd0c1d73ecbcc42f018a386483f02b5dfff46
1 /*
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,2017, 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.
37 #include "gmxpre.h"
39 #include "xlate.h"
41 #include <ctype.h>
42 #include <string.h>
44 #include "gromacs/gmxpreprocess/fflibutil.h"
45 #include "gromacs/gmxpreprocess/hackblock.h"
46 #include "gromacs/topology/residuetypes.h"
47 #include "gromacs/topology/symtab.h"
48 #include "gromacs/utility/cstringutil.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/futil.h"
51 #include "gromacs/utility/smalloc.h"
52 #include "gromacs/utility/strdb.h"
54 typedef struct {
55 char *filebase;
56 char *res;
57 char *atom;
58 char *replace;
59 } t_xlate_atom;
61 static void get_xlatoms(const char *fn, FILE *fp,
62 int *nptr, t_xlate_atom **xlptr)
64 char filebase[STRLEN];
65 char line[STRLEN];
66 char abuf[1024], rbuf[1024], repbuf[1024], dumbuf[1024];
67 char *_ptr;
68 int n, na, idum;
69 t_xlate_atom *xl;
71 fflib_filename_base(fn, filebase, STRLEN);
73 n = *nptr;
74 xl = *xlptr;
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
80 * on the first line.
82 if (na == 1 && n == *nptr && sscanf(rbuf, "%d", &idum) == 1)
84 continue;
86 if (na != 3)
88 gmx_fatal(FARGS, "Expected a residue name and two atom names in file '%s', not '%s'", fn, line);
91 srenew(xl, n+1);
92 xl[n].filebase = gmx_strdup(filebase);
94 /* Use wildcards... */
95 if (strcmp(rbuf, "*") != 0)
97 xl[n].res = gmx_strdup(rbuf);
99 else
101 xl[n].res = nullptr;
104 /* Replace underscores in the string by spaces */
105 while ((_ptr = strchr(abuf, '_')) != nullptr)
107 *_ptr = ' ';
110 xl[n].atom = gmx_strdup(abuf);
111 xl[n].replace = gmx_strdup(repbuf);
112 n++;
115 *nptr = n;
116 *xlptr = xl;
119 static void done_xlatom(int nxlate, t_xlate_atom *xlatom)
121 int i;
123 for (i = 0; (i < nxlate); i++)
125 sfree(xlatom[i].filebase);
126 if (xlatom[i].res != nullptr)
128 sfree(xlatom[i].res);
130 sfree(xlatom[i].atom);
131 sfree(xlatom[i].replace);
133 sfree(xlatom);
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,
139 gmx_bool bVerbose)
141 FILE *fp;
142 int nxlate, a, i, resind;
143 t_xlate_atom *xlatom;
144 int nf;
145 char **f;
146 char c, *rnm, atombuf[32], *ptr0, *ptr1;
147 gmx_bool bReorderedNum, bRenamed, bMatch;
148 gmx_bool bStartTerm, bEndTerm;
150 nxlate = 0;
151 xlatom = nullptr;
152 if (xlfile != nullptr)
154 fp = libopen(xlfile);
155 get_xlatoms(xlfile, fp, &nxlate, &xlatom);
156 fclose(fp);
158 else
160 nf = fflib_search_file_end(ffdir, ".arn", FALSE, &f);
161 for (i = 0; i < nf; i++)
163 fp = fflib_open(f[i]);
164 get_xlatoms(f[i], fp, &nxlate, &xlatom);
165 gmx_ffclose(fp);
166 sfree(f[i]);
168 sfree(f);
171 for (a = 0; (a < atoms->nr); a++)
173 resind = atoms->atom[a].resind;
175 bStartTerm = (resind == 0) || atoms->resinfo[resind].chainnum != atoms->resinfo[resind-1].chainnum;
176 bEndTerm = (resind >= atoms->nres-1) || atoms->resinfo[resind].chainnum != atoms->resinfo[resind+1].chainnum;
178 if (bResname)
180 rnm = *(atoms->resinfo[resind].name);
182 else
184 rnm = *(atoms->resinfo[resind].rtp);
187 strcpy(atombuf, *(atoms->atomname[a]));
188 bReorderedNum = FALSE;
189 if (bReorderNum)
191 if (isdigit(atombuf[0]))
193 c = atombuf[0];
194 for (i = 0; ((size_t)i < strlen(atombuf)-1); i++)
196 atombuf[i] = atombuf[i+1];
198 atombuf[i] = c;
199 bReorderedNum = TRUE;
202 bRenamed = FALSE;
203 for (i = 0; (i < nxlate) && !bRenamed; i++)
205 /* Check if the base file name of the rtp and arn entry match */
206 if (restp == nullptr ||
207 gmx_strcasecmp(restp[resind].filebase, xlatom[i].filebase) == 0)
209 /* Match the residue name */
210 bMatch = (xlatom[i].res == nullptr ||
211 (gmx_strcasecmp("protein-nterm", xlatom[i].res) == 0 &&
212 gmx_residuetype_is_protein(rt, rnm) && bStartTerm) ||
213 (gmx_strcasecmp("protein-cterm", xlatom[i].res) == 0 &&
214 gmx_residuetype_is_protein(rt, rnm) && bEndTerm) ||
215 (gmx_strcasecmp("protein", xlatom[i].res) == 0 &&
216 gmx_residuetype_is_protein(rt, rnm)) ||
217 (gmx_strcasecmp("DNA", xlatom[i].res) == 0 &&
218 gmx_residuetype_is_dna(rt, rnm)) ||
219 (gmx_strcasecmp("RNA", xlatom[i].res) == 0 &&
220 gmx_residuetype_is_rna(rt, rnm)));
221 if (!bMatch)
223 ptr0 = rnm;
224 ptr1 = xlatom[i].res;
225 while (ptr0[0] != '\0' && ptr1[0] != '\0' &&
226 (ptr0[0] == ptr1[0] || ptr1[0] == '?'))
228 ptr0++;
229 ptr1++;
231 bMatch = (ptr0[0] == '\0' && ptr1[0] == '\0');
233 if (bMatch && strcmp(atombuf, xlatom[i].atom) == 0)
235 /* We have a match. */
236 /* Don't free the old atomname,
237 * since it might be in the symtab.
239 ptr0 = gmx_strdup(xlatom[i].replace);
240 if (bVerbose)
242 printf("Renaming atom '%s' in residue %d %s to '%s'\n",
243 *atoms->atomname[a],
244 atoms->resinfo[resind].nr,
245 *atoms->resinfo[resind].name,
246 ptr0);
248 atoms->atomname[a] = put_symtab(symtab, ptr0);
249 bRenamed = TRUE;
253 if (bReorderedNum && !bRenamed)
255 atoms->atomname[a] = put_symtab(symtab, atombuf);
259 done_xlatom(nxlate, xlatom);