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, 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/fileio/gmxfio.h"
48 #include "gromacs/gmxpreprocess/fflibutil.h"
49 #include "gromacs/gmxpreprocess/h_db.h"
50 #include "gromacs/gmxpreprocess/notset.h"
51 #include "gromacs/gmxpreprocess/resall.h"
52 #include "gromacs/gmxpreprocess/toputil.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/utility/smalloc.h"
57 #include "gromacs/utility/strdb.h"
59 /* use bonded types definitions in hackblock.h */
60 #define ekwRepl (ebtsNR+1)
61 #define ekwAdd (ebtsNR+2)
62 #define ekwDel (ebtsNR+3)
64 static const char *kw_names[ekwNR] = {
65 "replace", "add", "delete"
68 static int find_kw(char *keyw)
72 for (i = 0; i < ebtsNR; i++)
74 if (gmx_strcasecmp(btsNames[i], keyw) == 0)
79 for (i = 0; i < ekwNR; i++)
81 if (gmx_strcasecmp(kw_names[i], keyw) == 0)
83 return ebtsNR + 1 + i;
90 #define FATAL() gmx_fatal(FARGS, "Reading Termini Database: not enough items on line\n%s", line)
92 static void read_atom(char *line, bool bAdd,
93 char **nname, t_atom *a, gpp_atomtype_t atype, int *cgnr)
99 /* This code is messy, because of support for different formats:
100 * for replace: [new name] <atom type> <m> <q> [cgnr (old format)]
101 * for add: <atom type> <m> <q> [cgnr]
103 nr = sscanf(line, "%s %s %s %s %s", buf[0], buf[1], buf[2], buf[3], buf[4]);
105 /* Here there an ambiguity due to the old replace format with cgnr,
106 * which was read for years, but ignored in the rest of the code.
107 * We have to assume that the atom type does not start with a digit
108 * to make a line with 4 entries uniquely interpretable.
110 if (!bAdd && nr == 4 && isdigit(buf[1][0]))
115 if (nr < 3 || nr > 4)
117 gmx_fatal(FARGS, "Reading Termini Database: expected %d or %d items of atom data in stead of %d on line\n%s", 3, 4, nr, line);
124 *nname = gmx_strdup(buf[i++]);
131 a->type = get_atomtype_type(buf[i++], atype);
132 sscanf(buf[i++], "%lf", &m);
134 sscanf(buf[i++], "%lf", &q);
138 sscanf(buf[i++], "%d", cgnr);
146 static void print_atom(FILE *out, t_atom *a, gpp_atomtype_t atype)
148 fprintf(out, "\t%s\t%g\t%g\n",
149 get_atomtype_name(a->type, atype), a->m, a->q);
152 static void print_ter_db(const char *ff, char C, int nb, t_hackblock tb[],
153 gpp_atomtype_t atype)
156 int i, j, k, bt, nrepl, nadd, ndel;
159 sprintf(buf, "%s-%c.tdb", ff, C);
160 out = gmx_fio_fopen(buf, "w");
162 for (i = 0; (i < nb); i++)
164 fprintf(out, "[ %s ]\n", tb[i].name);
170 for (j = 0; j < tb[i].nhack; j++)
172 if (tb[i].hack[j].oname != nullptr && tb[i].hack[j].nname != nullptr)
176 else if (tb[i].hack[j].oname == nullptr && tb[i].hack[j].nname != nullptr)
180 else if (tb[i].hack[j].oname != nullptr && tb[i].hack[j].nname == nullptr)
184 else if (tb[i].hack[j].oname == nullptr && tb[i].hack[j].nname == nullptr)
186 gmx_fatal(FARGS, "invalid hack (%s) in termini database", tb[i].name);
191 fprintf(out, "[ %s ]\n", kw_names[ekwRepl-ebtsNR-1]);
192 for (j = 0; j < tb[i].nhack; j++)
194 if (tb[i].hack[j].oname != nullptr && tb[i].hack[j].nname != nullptr)
196 fprintf(out, "%s\t", tb[i].hack[j].oname);
197 print_atom(out, tb[i].hack[j].atom, atype);
203 fprintf(out, "[ %s ]\n", kw_names[ekwAdd-ebtsNR-1]);
204 for (j = 0; j < tb[i].nhack; j++)
206 if (tb[i].hack[j].oname == nullptr && tb[i].hack[j].nname != nullptr)
208 print_ab(out, &(tb[i].hack[j]), tb[i].hack[j].nname);
209 print_atom(out, tb[i].hack[j].atom, atype);
215 fprintf(out, "[ %s ]\n", kw_names[ekwDel-ebtsNR-1]);
216 for (j = 0; j < tb[i].nhack; j++)
218 if (tb[i].hack[j].oname != nullptr && tb[i].hack[j].nname == nullptr)
220 fprintf(out, "%s\n", tb[i].hack[j].oname);
224 for (bt = 0; bt < ebtsNR; bt++)
228 fprintf(out, "[ %s ]\n", btsNames[bt]);
229 for (j = 0; j < tb[i].rb[bt].nb; j++)
231 for (k = 0; k < btsNiatoms[bt]; k++)
233 fprintf(out, "%s%s", k ? "\t" : "", tb[i].rb[bt].b[j].a[k]);
235 if (tb[i].rb[bt].b[j].s)
237 fprintf(out, "\t%s", tb[i].rb[bt].b[j].s);
248 static void read_ter_db_file(const char *fn,
249 int *ntbptr, t_hackblock **tbptr,
250 gpp_atomtype_t atype)
252 char filebase[STRLEN], *ptr;
254 char header[STRLEN], buf[STRLEN], line[STRLEN];
256 int i, j, n, ni, kwnr, nb, maxnb, nh;
258 fflib_filename_base(fn, filebase, STRLEN);
259 /* Remove the C/N termini extension */
260 ptr = strrchr(filebase, '.');
272 get_a_line(in, line, STRLEN);
275 if (get_header(line, header))
277 /* this is a new block, or a new keyword */
278 kwnr = find_kw(header);
283 /* here starts a new block */
289 clear_t_hackblock(&tb[nb]);
290 tb[nb].name = gmx_strdup(header);
291 tb[nb].filebase = gmx_strdup(filebase);
298 gmx_fatal(FARGS, "reading termini database: "
299 "directive expected before line:\n%s\n"
300 "This might be a file in an old format.", line);
302 /* this is not a header, so it must be data */
305 /* this is a hack: add/rename/delete atoms */
306 /* make space for hacks */
307 if (tb[nb].nhack >= tb[nb].maxhack)
309 tb[nb].maxhack += 10;
310 srenew(tb[nb].hack, tb[nb].maxhack);
313 clear_t_hack(&(tb[nb].hack[nh]));
314 for (i = 0; i < 4; i++)
316 tb[nb].hack[nh].a[i] = nullptr;
322 if (kwnr == ekwRepl || kwnr == ekwDel)
324 if (sscanf(line, "%s%n", buf, &n) != 1)
326 gmx_fatal(FARGS, "Reading Termini Database '%s': "
327 "expected atom name on line\n%s", fn, line);
329 tb[nb].hack[nh].oname = gmx_strdup(buf);
330 /* we only replace or delete one atom at a time */
331 tb[nb].hack[nh].nr = 1;
333 else if (kwnr == ekwAdd)
335 read_ab(line, fn, &(tb[nb].hack[nh]));
336 get_a_line(in, line, STRLEN);
340 gmx_fatal(FARGS, "unimplemented keyword number %d (%s:%d)",
341 kwnr, __FILE__, __LINE__);
343 if (kwnr == ekwRepl || kwnr == ekwAdd)
345 snew(tb[nb].hack[nh].atom, 1);
346 read_atom(line+n, kwnr == ekwAdd,
347 &tb[nb].hack[nh].nname, tb[nb].hack[nh].atom, atype,
348 &tb[nb].hack[nh].cgnr);
349 if (tb[nb].hack[nh].nname == nullptr)
351 if (tb[nb].hack[nh].oname != nullptr)
353 tb[nb].hack[nh].nname = gmx_strdup(tb[nb].hack[nh].oname);
357 gmx_fatal(FARGS, "Reading Termini Database '%s': don't know which name the new atom should have on line\n%s", fn, line);
362 else if (kwnr >= 0 && kwnr < ebtsNR)
364 /* this is bonded data: bonds, angles, dihedrals or impropers */
365 srenew(tb[nb].rb[kwnr].b, tb[nb].rb[kwnr].nb+1);
367 for (j = 0; j < btsNiatoms[kwnr]; j++)
369 if (sscanf(line+n, "%s%n", buf, &ni) == 1)
371 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].a[j] = gmx_strdup(buf);
375 gmx_fatal(FARGS, "Reading Termini Database '%s': expected %d atom names (found %d) on line\n%s", fn, btsNiatoms[kwnr], j-1, line);
379 for (; j < MAXATOMLIST; j++)
381 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].a[j] = nullptr;
384 sscanf(line+n, "%s", buf);
385 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].s = gmx_strdup(buf);
386 tb[nb].rb[kwnr].nb++;
390 gmx_fatal(FARGS, "Reading Termini Database: Expecting a header at line\n"
394 get_a_line(in, line, STRLEN);
405 int read_ter_db(const char *ffdir, char ter,
406 t_hackblock **tbptr, gpp_atomtype_t atype)
411 sprintf(ext, ".%c.tdb", ter);
413 /* Search for termini database files.
414 * Do not generate an error when none are found.
416 std::vector<std::string> tdbf = fflib_search_file_end(ffdir, ext, FALSE);
419 for (const auto &filename : tdbf)
421 read_ter_db_file(filename.c_str(), &ntb, tbptr, atype);
426 print_ter_db("new", ter, ntb, *tbptr, atype);
432 t_hackblock **filter_ter(int nb, t_hackblock tb[],
436 // TODO Four years later, no force fields have ever used this, so decide status of this feature
437 /* Since some force fields (e.g. OPLS) needs different
438 * atomtypes for different residues there could be a lot
439 * of entries in the databases for specific residues
440 * (e.g. GLY-NH3+, SER-NH3+, ALA-NH3+).
442 * To reduce the database size, we assume that a terminus specifier liker
444 * [ GLY|SER|ALA-NH3+ ]
446 * would cover all of the three residue types above.
447 * Wildcards (*,?) are not OK. Don't worry about e.g. GLU vs. GLUH since
448 * pdb2gmx only uses the first 3 letters when calling this routine.
450 * To automate this, this routines scans a list of termini
451 * for the residue name "resname" and returns an allocated list of
452 * pointers to the termini that could be applied to the
453 * residue in question. The variable pointed to by nret will
454 * contain the number of valid pointers in the list.
455 * Remember to free the list when you are done with it...
458 int i, j, n, none_idx;
466 for (i = 0; i < nb; i++)
472 if (gmx_strncasecmp(resname, s, 3) == 0)
481 /* advance to next |-separated field */
489 while (!found && s != nullptr);
492 /* All residue-specific termini have been added. We might have to fall
493 * back on generic termini, which are characterized by not having
494 * '-' in the name prior to the last position (which indicates charge).
495 * The [ None ] alternative is special since we don't want that
496 * to be the default, so we put it last in the list we return.
499 for (i = 0; i < nb; i++)
502 if (!gmx_strcasecmp("None", s))
508 /* Time to see if there's a generic terminus that matches.
509 Is there a hyphen? */
510 char *c = strchr(s, '-');
512 /* A conjunction hyphen normally indicates a residue-specific
513 terminus, which is named like "GLY-COOH". A generic terminus
514 won't have a hyphen. */
515 bool bFoundAnyHyphen = (c != nullptr);
516 /* '-' as the last character indicates charge, so if that's
517 the only one found e.g. "COO-", then it was not a conjunction
518 hyphen, so this is a generic terminus */
519 bool bOnlyFoundChargeHyphen = (bFoundAnyHyphen &&
521 /* Thus, "GLY-COO-" is not recognized as a generic terminus. */
522 bool bFoundGenericTerminus = !bFoundAnyHyphen || bOnlyFoundChargeHyphen;
523 if (bFoundGenericTerminus)
525 /* Check that we haven't already added a residue-specific version
528 for (j = 0; j < n && strstr((*list[j]).name, s) == nullptr; j++)
544 list[n] = &(tb[none_idx]);
553 t_hackblock *choose_ter(int nb, t_hackblock **tb, const char *title)
557 printf("%s\n", title);
558 for (i = 0; (i < nb); i++)
560 bool bIsZwitterion = (0 == gmx_wcmatch("*ZWITTERION*", (*tb[i]).name));
561 printf("%2d: %s%s\n", i, (*tb[i]).name,
562 bIsZwitterion ? " (only use with zwitterions containing exactly one residue)" : "");
566 ret = fscanf(stdin, "%d", &sel);
568 while ((ret != 1) || (sel < 0) || (sel >= nb));