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,2018,2019, 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.
47 #include "gromacs/gmxpreprocess/fflibutil.h"
48 #include "gromacs/gmxpreprocess/grompp-impl.h"
49 #include "gromacs/gmxpreprocess/hackblock.h"
50 #include "gromacs/topology/residuetypes.h"
51 #include "gromacs/topology/symtab.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
56 #include "gromacs/utility/strdb.h"
65 static void get_xlatoms(const std::string &filename, FILE *fp,
66 int *nptr, t_xlate_atom **xlptr)
68 char filebase[STRLEN];
70 char abuf[1024], rbuf[1024], repbuf[1024], dumbuf[1024];
75 fflib_filename_base(filename.c_str(), filebase, STRLEN);
80 while (get_a_line(fp, line, STRLEN))
82 na = sscanf(line, "%s%s%s%s", rbuf, abuf, repbuf, dumbuf);
83 /* Check if we are reading an old format file with the number of items
86 if (na == 1 && n == *nptr && sscanf(rbuf, "%d", &idum) == 1)
92 gmx_fatal(FARGS, "Expected a residue name and two atom names in file '%s', not '%s'", filename.c_str(), line);
96 xl[n].filebase = gmx_strdup(filebase);
98 /* Use wildcards... */
99 if (strcmp(rbuf, "*") != 0)
101 xl[n].res = gmx_strdup(rbuf);
108 /* Replace underscores in the string by spaces */
109 while ((_ptr = strchr(abuf, '_')) != nullptr)
114 xl[n].atom = gmx_strdup(abuf);
115 xl[n].replace = gmx_strdup(repbuf);
123 static void done_xlatom(int nxlate, t_xlate_atom *xlatom)
127 for (i = 0; (i < nxlate); i++)
129 sfree(xlatom[i].filebase);
130 if (xlatom[i].res != nullptr)
132 sfree(xlatom[i].res);
134 sfree(xlatom[i].atom);
135 sfree(xlatom[i].replace);
140 void rename_atoms(const char* xlfile, const char *ffdir,
141 t_atoms *atoms, t_symtab *symtab, const t_restp *restp,
142 bool bResname, ResidueType *rt, bool bReorderNum,
145 int nxlate, a, i, resind;
146 t_xlate_atom *xlatom;
147 char c, *rnm, atombuf[32], *ptr0, *ptr1;
148 bool bReorderedNum, bRenamed, bMatch;
149 bool bStartTerm, bEndTerm;
153 if (xlfile != nullptr)
155 gmx::FilePtr fp = gmx::openLibraryFile(xlfile);
156 get_xlatoms(xlfile, fp.get(), &nxlate, &xlatom);
160 std::vector<std::string> fns = fflib_search_file_end(ffdir, ".arn", FALSE);
161 for (const auto &filename : fns)
163 FILE * fp = fflib_open(filename);
164 get_xlatoms(filename, fp, &nxlate, &xlatom);
169 for (a = 0; (a < atoms->nr); a++)
171 resind = atoms->atom[a].resind;
173 bStartTerm = (resind == 0) || atoms->resinfo[resind].chainnum != atoms->resinfo[resind-1].chainnum;
174 bEndTerm = (resind >= atoms->nres-1) || atoms->resinfo[resind].chainnum != atoms->resinfo[resind+1].chainnum;
178 rnm = *(atoms->resinfo[resind].name);
182 rnm = *(atoms->resinfo[resind].rtp);
185 strcpy(atombuf, *(atoms->atomname[a]));
186 bReorderedNum = FALSE;
189 if (isdigit(atombuf[0]))
192 for (i = 0; (static_cast<size_t>(i) < strlen(atombuf)-1); i++)
194 atombuf[i] = atombuf[i+1];
197 bReorderedNum = TRUE;
201 for (i = 0; (i < nxlate) && !bRenamed; i++)
203 /* Check if the base file name of the rtp and arn entry match */
204 if (restp == nullptr ||
205 gmx_strcasecmp(restp[resind].filebase, xlatom[i].filebase) == 0)
207 /* Match the residue name */
208 bMatch = (xlatom[i].res == nullptr ||
209 (gmx_strcasecmp("protein-nterm", xlatom[i].res) == 0 &&
210 rt->namedResidueHasType(rnm, "Protein") && bStartTerm) ||
211 (gmx_strcasecmp("protein-cterm", xlatom[i].res) == 0 &&
212 rt->namedResidueHasType(rnm, "Protein") && bEndTerm) ||
213 (gmx_strcasecmp("protein", xlatom[i].res) == 0 &&
214 rt->namedResidueHasType(rnm, "Protein")) ||
215 (gmx_strcasecmp("DNA", xlatom[i].res) == 0 &&
216 rt->namedResidueHasType(rnm, "DNA")) ||
217 (gmx_strcasecmp("RNA", xlatom[i].res) == 0 &&
218 rt->namedResidueHasType(rnm, "RNA")));
222 ptr1 = xlatom[i].res;
223 while (ptr0[0] != '\0' && ptr1[0] != '\0' &&
224 (ptr0[0] == ptr1[0] || ptr1[0] == '?'))
229 bMatch = (ptr0[0] == '\0' && ptr1[0] == '\0');
231 if (bMatch && strcmp(atombuf, xlatom[i].atom) == 0)
233 /* We have a match. */
234 /* Don't free the old atomname,
235 * since it might be in the symtab.
237 ptr0 = gmx_strdup(xlatom[i].replace);
240 printf("Renaming atom '%s' in residue %d %s to '%s'\n",
242 atoms->resinfo[resind].nr,
243 *atoms->resinfo[resind].name,
246 atoms->atomname[a] = put_symtab(symtab, ptr0);
251 if (bReorderedNum && !bRenamed)
253 atoms->atomname[a] = put_symtab(symtab, atombuf);
257 done_xlatom(nxlate, xlatom);