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, 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.
44 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/utility/futil.h"
50 #include "gromacs/utility/cstringutil.h"
51 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/fileio/gmxfio.h"
55 #include "fflibutil.h"
57 #include "gromacs/fileio/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 const char *kw_names[ekwNR] = {
65 "replace", "add", "delete"
68 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, gmx_bool bAdd,
93 char **nname, t_atom *a, gpp_atomtype_t atype, int *cgnr)
96 char buf[5][30], type[12];
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 = 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;
157 char buf[STRLEN], nname[STRLEN];
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 != NULL && tb[i].hack[j].nname != NULL)
176 else if (tb[i].hack[j].oname == NULL && tb[i].hack[j].nname != NULL)
180 else if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname == NULL)
184 else if (tb[i].hack[j].oname == NULL && tb[i].hack[j].nname == NULL)
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 != NULL && tb[i].hack[j].nname != NULL)
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 == NULL && tb[i].hack[j].nname != NULL)
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 != NULL && tb[i].hack[j].nname == NULL)
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(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, '.');
269 fprintf(debug, "Opened %s\n", fn);
276 get_a_line(in, line, STRLEN);
279 if (get_header(line, header))
281 /* this is a new block, or a new keyword */
282 kwnr = find_kw(header);
287 /* here starts a new block */
293 clear_t_hackblock(&tb[nb]);
294 tb[nb].name = strdup(header);
295 tb[nb].filebase = strdup(filebase);
302 gmx_fatal(FARGS, "reading termini database: "
303 "directive expected before line:\n%s\n"
304 "This might be a file in an old format.", line);
306 /* this is not a header, so it must be data */
309 /* this is a hack: add/rename/delete atoms */
310 /* make space for hacks */
311 if (tb[nb].nhack >= tb[nb].maxhack)
313 tb[nb].maxhack += 10;
314 srenew(tb[nb].hack, tb[nb].maxhack);
317 clear_t_hack(&(tb[nb].hack[nh]));
318 for (i = 0; i < 4; i++)
320 tb[nb].hack[nh].a[i] = NULL;
326 if (kwnr == ekwRepl || kwnr == ekwDel)
328 if (sscanf(line, "%s%n", buf, &n) != 1)
330 gmx_fatal(FARGS, "Reading Termini Database '%s': "
331 "expected atom name on line\n%s", fn, line);
333 tb[nb].hack[nh].oname = strdup(buf);
334 /* we only replace or delete one atom at a time */
335 tb[nb].hack[nh].nr = 1;
337 else if (kwnr == ekwAdd)
339 read_ab(line, fn, &(tb[nb].hack[nh]));
340 get_a_line(in, line, STRLEN);
344 gmx_fatal(FARGS, "unimplemented keyword number %d (%s:%d)",
345 kwnr, __FILE__, __LINE__);
347 if (kwnr == ekwRepl || kwnr == ekwAdd)
349 snew(tb[nb].hack[nh].atom, 1);
350 read_atom(line+n, kwnr == ekwAdd,
351 &tb[nb].hack[nh].nname, tb[nb].hack[nh].atom, atype,
352 &tb[nb].hack[nh].cgnr);
353 if (tb[nb].hack[nh].nname == NULL)
355 if (tb[nb].hack[nh].oname != NULL)
357 tb[nb].hack[nh].nname = strdup(tb[nb].hack[nh].oname);
361 gmx_fatal(FARGS, "Reading Termini Database '%s': don't know which name the new atom should have on line\n%s", fn, line);
366 else if (kwnr >= 0 && kwnr < ebtsNR)
368 /* this is bonded data: bonds, angles, dihedrals or impropers */
369 srenew(tb[nb].rb[kwnr].b, tb[nb].rb[kwnr].nb+1);
371 for (j = 0; j < btsNiatoms[kwnr]; j++)
373 if (sscanf(line+n, "%s%n", buf, &ni) == 1)
375 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].a[j] = strdup(buf);
379 gmx_fatal(FARGS, "Reading Termini Database '%s': expected %d atom names (found %d) on line\n%s", fn, btsNiatoms[kwnr], j-1, line);
383 for (; j < MAXATOMLIST; j++)
385 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].a[j] = NULL;
388 sscanf(line+n, "%s", buf);
389 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].s = strdup(buf);
390 tb[nb].rb[kwnr].nb++;
394 gmx_fatal(FARGS, "Reading Termini Database: Expecting a header at line\n"
398 get_a_line(in, line, STRLEN);
409 int read_ter_db(const char *ffdir, char ter,
410 t_hackblock **tbptr, gpp_atomtype_t atype)
417 sprintf(ext, ".%c.tdb", ter);
419 /* Search for termini database files.
420 * Do not generate an error when none are found.
422 ntdbf = fflib_search_file_end(ffdir, ext, FALSE, &tdbf);
425 for (f = 0; f < ntdbf; f++)
427 read_ter_db_file(tdbf[f], &ntb, tbptr, atype);
434 print_ter_db("new", ter, ntb, *tbptr, atype);
440 t_hackblock **filter_ter(int nrtp, t_restp rtp[],
441 int nb, t_hackblock tb[],
446 /* Since some force fields (e.g. OPLS) needs different
447 * atomtypes for different residues there could be a lot
448 * of entries in the databases for specific residues
449 * (e.g. GLY-NH3+, SER-NH3+, ALA-NH3+).
451 * To reduce the database size, we assume that a terminus specifier liker
453 * [ GLY|SER|ALA-NH3+ ]
455 * would cover all of the three residue types above.
456 * Wildcards (*,?) are not OK. Don't worry about e.g. GLU vs. GLUH since
457 * pdb2gmx only uses the first 3 letters when calling this routine.
459 * To automate this, this routines scans a list of termini
460 * for the residue name "resname" and returns an allocated list of
461 * pointers to the termini that could be applied to the
462 * residue in question. The variable pointed to by nret will
463 * contain the number of valid pointers in the list.
464 * Remember to free the list when you are done with it...
468 int i, j, n, len, none_idx;
470 char *rtpname_match, *s, *s2, *c;
473 rtpname_match = search_rtp(rtpname, nrtp, rtp);
474 restp = get_restp(rtpname_match, nrtp, rtp);
479 for (i = 0; i < nb; i++)
485 /* The residue name should appear in a tdb file with the same base name
486 * as the file containing the rtp entry.
487 * This makes termini selection for different molecule types
490 if (gmx_strcasecmp(restp->filebase, tb[i].filebase) == 0 &&
491 gmx_strncasecmp(resname, s, 3) == 0)
500 /* advance to next |-separated field */
508 while (!found && s != NULL);
511 /* All residue-specific termini have been added. See if there
512 * are some generic ones by searching for the occurence of
513 * '-' in the name prior to the last position (which indicates charge).
514 * The [ None ] alternative is special since we don't want that
515 * to be the default, so we put it last in the list we return.
518 for (i = 0; i < nb; i++)
521 /* The residue name should appear in a tdb file with the same base name
522 * as the file containing the rtp entry.
523 * This makes termini selection for different molecule types
526 if (gmx_strcasecmp(restp->filebase, tb[i].filebase) == 0)
528 if (!gmx_strcasecmp("None", s))
535 if (c == NULL || ((c-s+1) == strlen(s)))
537 /* Check that we haven't already added a residue-specific version
540 for (j = 0; j < n && strstr((*list[j]).name, s) == NULL; j++)
557 list[n] = &(tb[none_idx]);
566 t_hackblock *choose_ter(int nb, t_hackblock **tb, const char *title)
570 printf("%s\n", title);
571 for (i = 0; (i < nb); i++)
573 char *advice_string = "";
574 if (0 == gmx_wcmatch("*ZWITTERION*", (*tb[i]).name))
576 advice_string = " (only use with zwitterions containing exactly one residue)";
578 printf("%2d: %s%s\n", i, (*tb[i]).name, advice_string);
582 ret = fscanf(stdin, "%d", &sel);
584 while ((ret != 1) || (sel < 0) || (sel >= nb));