Fix GPU X buffer ops with empty domain
[gromacs.git] / src / gromacs / gmxpreprocess / xlate.cpp
blob43f232facdfa1b96728257030b1af6b083b0a9c0
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,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.
37 #include "gmxpre.h"
39 #include "xlate.h"
41 #include <cctype>
42 #include <cstring>
44 #include <string>
45 #include <vector>
47 #include "gromacs/gmxpreprocess/fflibutil.h"
48 #include "gromacs/gmxpreprocess/grompp_impl.h"
49 #include "gromacs/topology/residuetypes.h"
50 #include "gromacs/topology/symtab.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/smalloc.h"
55 #include "gromacs/utility/strdb.h"
57 #include "hackblock.h"
59 typedef struct {
60 char *filebase;
61 char *res;
62 char *atom;
63 char *replace;
64 } t_xlate_atom;
66 static void get_xlatoms(const std::string &filename, FILE *fp,
67 int *nptr, t_xlate_atom **xlptr)
69 char filebase[STRLEN];
70 char line[STRLEN];
71 char abuf[1024], rbuf[1024], repbuf[1024], dumbuf[1024];
72 char *_ptr;
73 int n, na, idum;
74 t_xlate_atom *xl;
76 fflib_filename_base(filename.c_str(), filebase, STRLEN);
78 n = *nptr;
79 xl = *xlptr;
81 while (get_a_line(fp, line, STRLEN))
83 na = sscanf(line, "%s%s%s%s", rbuf, abuf, repbuf, dumbuf);
84 /* Check if we are reading an old format file with the number of items
85 * on the first line.
87 if (na == 1 && n == *nptr && sscanf(rbuf, "%d", &idum) == 1)
89 continue;
91 if (na != 3)
93 gmx_fatal(FARGS, "Expected a residue name and two atom names in file '%s', not '%s'", filename.c_str(), line);
96 srenew(xl, n+1);
97 xl[n].filebase = gmx_strdup(filebase);
99 /* Use wildcards... */
100 if (strcmp(rbuf, "*") != 0)
102 xl[n].res = gmx_strdup(rbuf);
104 else
106 xl[n].res = nullptr;
109 /* Replace underscores in the string by spaces */
110 while ((_ptr = strchr(abuf, '_')) != nullptr)
112 *_ptr = ' ';
115 xl[n].atom = gmx_strdup(abuf);
116 xl[n].replace = gmx_strdup(repbuf);
117 n++;
120 *nptr = n;
121 *xlptr = xl;
124 static void done_xlatom(int nxlate, t_xlate_atom *xlatom)
126 int i;
128 for (i = 0; (i < nxlate); i++)
130 sfree(xlatom[i].filebase);
131 if (xlatom[i].res != nullptr)
133 sfree(xlatom[i].res);
135 sfree(xlatom[i].atom);
136 sfree(xlatom[i].replace);
138 sfree(xlatom);
141 void rename_atoms(const char* xlfile, const char *ffdir,
142 t_atoms *atoms, t_symtab *symtab, gmx::ArrayRef<const PreprocessResidue> localPpResidue,
143 bool bResname, ResidueType *rt, bool bReorderNum,
144 bool bVerbose)
146 int nxlate, a, i, resind;
147 t_xlate_atom *xlatom;
148 char c, *rnm, atombuf[32];
149 bool bReorderedNum, bRenamed, bMatch;
150 bool bStartTerm, bEndTerm;
152 nxlate = 0;
153 xlatom = nullptr;
154 if (xlfile != nullptr)
156 gmx::FilePtr fp = gmx::openLibraryFile(xlfile);
157 get_xlatoms(xlfile, fp.get(), &nxlate, &xlatom);
159 else
161 std::vector<std::string> fns = fflib_search_file_end(ffdir, ".arn", FALSE);
162 for (const auto &filename : fns)
164 FILE * fp = fflib_open(filename);
165 get_xlatoms(filename, fp, &nxlate, &xlatom);
166 gmx_ffclose(fp);
170 for (a = 0; (a < atoms->nr); a++)
172 resind = atoms->atom[a].resind;
174 bStartTerm = (resind == 0) || atoms->resinfo[resind].chainnum != atoms->resinfo[resind-1].chainnum;
175 bEndTerm = (resind >= atoms->nres-1) || atoms->resinfo[resind].chainnum != atoms->resinfo[resind+1].chainnum;
177 if (bResname)
179 rnm = *(atoms->resinfo[resind].name);
181 else
183 rnm = *(atoms->resinfo[resind].rtp);
186 strcpy(atombuf, *(atoms->atomname[a]));
187 bReorderedNum = FALSE;
188 if (bReorderNum)
190 if (isdigit(atombuf[0]))
192 c = atombuf[0];
193 for (i = 0; (static_cast<size_t>(i) < strlen(atombuf)-1); i++)
195 atombuf[i] = atombuf[i+1];
197 atombuf[i] = c;
198 bReorderedNum = TRUE;
201 bRenamed = FALSE;
202 for (i = 0; (i < nxlate) && !bRenamed; i++)
204 /* Check if the base file name of the rtp and arn entry match */
205 if (localPpResidue.empty() ||
206 gmx::equalCaseInsensitive(localPpResidue[resind].filebase, xlatom[i].filebase))
208 /* Match the residue name */
209 bMatch = (xlatom[i].res == nullptr ||
210 (gmx_strcasecmp("protein-nterm", xlatom[i].res) == 0 &&
211 rt->namedResidueHasType(rnm, "Protein") && bStartTerm) ||
212 (gmx_strcasecmp("protein-cterm", xlatom[i].res) == 0 &&
213 rt->namedResidueHasType(rnm, "Protein") && bEndTerm) ||
214 (gmx_strcasecmp("protein", xlatom[i].res) == 0 &&
215 rt->namedResidueHasType(rnm, "Protein")) ||
216 (gmx_strcasecmp("DNA", xlatom[i].res) == 0 &&
217 rt->namedResidueHasType(rnm, "DNA")) ||
218 (gmx_strcasecmp("RNA", xlatom[i].res) == 0 &&
219 rt->namedResidueHasType(rnm, "RNA")));
220 if (!bMatch)
222 const char *ptr0 = rnm;
223 const char *ptr1 = xlatom[i].res;
224 while (ptr0[0] != '\0' && ptr1[0] != '\0' &&
225 (ptr0[0] == ptr1[0] || ptr1[0] == '?'))
227 ptr0++;
228 ptr1++;
230 bMatch = (ptr0[0] == '\0' && ptr1[0] == '\0');
232 if (bMatch && strcmp(atombuf, xlatom[i].atom) == 0)
234 /* We have a match. */
235 /* Don't free the old atomname,
236 * since it might be in the symtab.
238 const char *ptr0 = xlatom[i].replace;
239 if (bVerbose)
241 printf("Renaming atom '%s' in residue %d %s to '%s'\n",
242 *atoms->atomname[a],
243 atoms->resinfo[resind].nr,
244 *atoms->resinfo[resind].name,
245 ptr0);
247 atoms->atomname[a] = put_symtab(symtab, ptr0);
248 bRenamed = TRUE;
252 if (bReorderedNum && !bRenamed)
254 atoms->atomname[a] = put_symtab(symtab, atombuf);
258 done_xlatom(nxlate, xlatom);