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"
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "gromacs/utility/futil.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/fileio/gmxfio.h"
54 #include "fflibutil.h"
56 #include "gromacs/fileio/strdb.h"
58 /* use bonded types definitions in hackblock.h */
59 #define ekwRepl ebtsNR+1
60 #define ekwAdd ebtsNR+2
61 #define ekwDel ebtsNR+3
63 const char *kw_names[ekwNR] = {
64 "replace", "add", "delete"
67 int find_kw(char *keyw)
71 for (i = 0; i < ebtsNR; i++)
73 if (gmx_strcasecmp(btsNames[i], keyw) == 0)
78 for (i = 0; i < ekwNR; i++)
80 if (gmx_strcasecmp(kw_names[i], keyw) == 0)
82 return ebtsNR + 1 + i;
89 #define FATAL() gmx_fatal(FARGS, "Reading Termini Database: not enough items on line\n%s", line)
91 static void read_atom(char *line, gmx_bool bAdd,
92 char **nname, t_atom *a, gpp_atomtype_t atype, int *cgnr)
95 char buf[5][30], type[12];
98 /* This code is messy, because of support for different formats:
99 * for replace: [new name] <atom type> <m> <q> [cgnr (old format)]
100 * for add: <atom type> <m> <q> [cgnr]
102 nr = sscanf(line, "%s %s %s %s %s", buf[0], buf[1], buf[2], buf[3], buf[4]);
104 /* Here there an ambiguity due to the old replace format with cgnr,
105 * which was read for years, but ignored in the rest of the code.
106 * We have to assume that the atom type does not start with a digit
107 * to make a line with 4 entries uniquely interpretable.
109 if (!bAdd && nr == 4 && isdigit(buf[1][0]))
114 if (nr < 3 || nr > 4)
116 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);
123 *nname = gmx_strdup(buf[i++]);
130 a->type = get_atomtype_type(buf[i++], atype);
131 sscanf(buf[i++], "%lf", &m);
133 sscanf(buf[i++], "%lf", &q);
137 sscanf(buf[i++], "%d", cgnr);
145 static void print_atom(FILE *out, t_atom *a, gpp_atomtype_t atype)
147 fprintf(out, "\t%s\t%g\t%g\n",
148 get_atomtype_name(a->type, atype), a->m, a->q);
151 static void print_ter_db(const char *ff, char C, int nb, t_hackblock tb[],
152 gpp_atomtype_t atype)
155 int i, j, k, bt, nrepl, nadd, ndel;
156 char buf[STRLEN], nname[STRLEN];
158 sprintf(buf, "%s-%c.tdb", ff, C);
159 out = gmx_fio_fopen(buf, "w");
161 for (i = 0; (i < nb); i++)
163 fprintf(out, "[ %s ]\n", tb[i].name);
169 for (j = 0; j < tb[i].nhack; j++)
171 if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname != NULL)
175 else if (tb[i].hack[j].oname == NULL && tb[i].hack[j].nname != NULL)
179 else if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname == NULL)
183 else if (tb[i].hack[j].oname == NULL && tb[i].hack[j].nname == NULL)
185 gmx_fatal(FARGS, "invalid hack (%s) in termini database", tb[i].name);
190 fprintf(out, "[ %s ]\n", kw_names[ekwRepl-ebtsNR-1]);
191 for (j = 0; j < tb[i].nhack; j++)
193 if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname != NULL)
195 fprintf(out, "%s\t", tb[i].hack[j].oname);
196 print_atom(out, tb[i].hack[j].atom, atype);
202 fprintf(out, "[ %s ]\n", kw_names[ekwAdd-ebtsNR-1]);
203 for (j = 0; j < tb[i].nhack; j++)
205 if (tb[i].hack[j].oname == NULL && tb[i].hack[j].nname != NULL)
207 print_ab(out, &(tb[i].hack[j]), tb[i].hack[j].nname);
208 print_atom(out, tb[i].hack[j].atom, atype);
214 fprintf(out, "[ %s ]\n", kw_names[ekwDel-ebtsNR-1]);
215 for (j = 0; j < tb[i].nhack; j++)
217 if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname == NULL)
219 fprintf(out, "%s\n", tb[i].hack[j].oname);
223 for (bt = 0; bt < ebtsNR; bt++)
227 fprintf(out, "[ %s ]\n", btsNames[bt]);
228 for (j = 0; j < tb[i].rb[bt].nb; j++)
230 for (k = 0; k < btsNiatoms[bt]; k++)
232 fprintf(out, "%s%s", k ? "\t" : "", tb[i].rb[bt].b[j].a[k]);
234 if (tb[i].rb[bt].b[j].s)
236 fprintf(out, "\t%s", tb[i].rb[bt].b[j].s);
247 static void read_ter_db_file(char *fn,
248 int *ntbptr, t_hackblock **tbptr,
249 gpp_atomtype_t atype)
251 char filebase[STRLEN], *ptr;
253 char header[STRLEN], buf[STRLEN], line[STRLEN];
255 int i, j, n, ni, kwnr, nb, maxnb, nh;
257 fflib_filename_base(fn, filebase, STRLEN);
258 /* Remove the C/N termini extension */
259 ptr = strrchr(filebase, '.');
268 fprintf(debug, "Opened %s\n", fn);
275 get_a_line(in, line, STRLEN);
278 if (get_header(line, header))
280 /* this is a new block, or a new keyword */
281 kwnr = find_kw(header);
286 /* here starts a new block */
292 clear_t_hackblock(&tb[nb]);
293 tb[nb].name = gmx_strdup(header);
294 tb[nb].filebase = gmx_strdup(filebase);
301 gmx_fatal(FARGS, "reading termini database: "
302 "directive expected before line:\n%s\n"
303 "This might be a file in an old format.", line);
305 /* this is not a header, so it must be data */
308 /* this is a hack: add/rename/delete atoms */
309 /* make space for hacks */
310 if (tb[nb].nhack >= tb[nb].maxhack)
312 tb[nb].maxhack += 10;
313 srenew(tb[nb].hack, tb[nb].maxhack);
316 clear_t_hack(&(tb[nb].hack[nh]));
317 for (i = 0; i < 4; i++)
319 tb[nb].hack[nh].a[i] = NULL;
325 if (kwnr == ekwRepl || kwnr == ekwDel)
327 if (sscanf(line, "%s%n", buf, &n) != 1)
329 gmx_fatal(FARGS, "Reading Termini Database '%s': "
330 "expected atom name on line\n%s", fn, line);
332 tb[nb].hack[nh].oname = gmx_strdup(buf);
333 /* we only replace or delete one atom at a time */
334 tb[nb].hack[nh].nr = 1;
336 else if (kwnr == ekwAdd)
338 read_ab(line, fn, &(tb[nb].hack[nh]));
339 get_a_line(in, line, STRLEN);
343 gmx_fatal(FARGS, "unimplemented keyword number %d (%s:%d)",
344 kwnr, __FILE__, __LINE__);
346 if (kwnr == ekwRepl || kwnr == ekwAdd)
348 snew(tb[nb].hack[nh].atom, 1);
349 read_atom(line+n, kwnr == ekwAdd,
350 &tb[nb].hack[nh].nname, tb[nb].hack[nh].atom, atype,
351 &tb[nb].hack[nh].cgnr);
352 if (tb[nb].hack[nh].nname == NULL)
354 if (tb[nb].hack[nh].oname != NULL)
356 tb[nb].hack[nh].nname = gmx_strdup(tb[nb].hack[nh].oname);
360 gmx_fatal(FARGS, "Reading Termini Database '%s': don't know which name the new atom should have on line\n%s", fn, line);
365 else if (kwnr >= 0 && kwnr < ebtsNR)
367 /* this is bonded data: bonds, angles, dihedrals or impropers */
368 srenew(tb[nb].rb[kwnr].b, tb[nb].rb[kwnr].nb+1);
370 for (j = 0; j < btsNiatoms[kwnr]; j++)
372 if (sscanf(line+n, "%s%n", buf, &ni) == 1)
374 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].a[j] = gmx_strdup(buf);
378 gmx_fatal(FARGS, "Reading Termini Database '%s': expected %d atom names (found %d) on line\n%s", fn, btsNiatoms[kwnr], j-1, line);
382 for (; j < MAXATOMLIST; j++)
384 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].a[j] = NULL;
387 sscanf(line+n, "%s", buf);
388 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].s = gmx_strdup(buf);
389 tb[nb].rb[kwnr].nb++;
393 gmx_fatal(FARGS, "Reading Termini Database: Expecting a header at line\n"
397 get_a_line(in, line, STRLEN);
408 int read_ter_db(const char *ffdir, char ter,
409 t_hackblock **tbptr, gpp_atomtype_t atype)
416 sprintf(ext, ".%c.tdb", ter);
418 /* Search for termini database files.
419 * Do not generate an error when none are found.
421 ntdbf = fflib_search_file_end(ffdir, ext, FALSE, &tdbf);
424 for (f = 0; f < ntdbf; f++)
426 read_ter_db_file(tdbf[f], &ntb, tbptr, atype);
433 print_ter_db("new", ter, ntb, *tbptr, atype);
439 t_hackblock **filter_ter(int nrtp, t_restp rtp[],
440 int nb, t_hackblock tb[],
445 /* Since some force fields (e.g. OPLS) needs different
446 * atomtypes for different residues there could be a lot
447 * of entries in the databases for specific residues
448 * (e.g. GLY-NH3+, SER-NH3+, ALA-NH3+).
450 * To reduce the database size, we assume that a terminus specifier liker
452 * [ GLY|SER|ALA-NH3+ ]
454 * would cover all of the three residue types above.
455 * Wildcards (*,?) are not OK. Don't worry about e.g. GLU vs. GLUH since
456 * pdb2gmx only uses the first 3 letters when calling this routine.
458 * To automate this, this routines scans a list of termini
459 * for the residue name "resname" and returns an allocated list of
460 * pointers to the termini that could be applied to the
461 * residue in question. The variable pointed to by nret will
462 * contain the number of valid pointers in the list.
463 * Remember to free the list when you are done with it...
467 int i, j, n, len, none_idx;
469 char *rtpname_match, *s, *s2, *c;
472 rtpname_match = search_rtp(rtpname, nrtp, rtp);
473 restp = get_restp(rtpname_match, nrtp, rtp);
478 for (i = 0; i < nb; i++)
484 /* The residue name should appear in a tdb file with the same base name
485 * as the file containing the rtp entry.
486 * This makes termini selection for different molecule types
489 if (gmx_strcasecmp(restp->filebase, tb[i].filebase) == 0 &&
490 gmx_strncasecmp(resname, s, 3) == 0)
499 /* advance to next |-separated field */
507 while (!found && s != NULL);
510 /* All residue-specific termini have been added. See if there
511 * are some generic ones by searching for the occurence of
512 * '-' in the name prior to the last position (which indicates charge).
513 * The [ None ] alternative is special since we don't want that
514 * to be the default, so we put it last in the list we return.
517 for (i = 0; i < nb; i++)
520 /* The residue name should appear in a tdb file with the same base name
521 * as the file containing the rtp entry.
522 * This makes termini selection for different molecule types
525 if (gmx_strcasecmp(restp->filebase, tb[i].filebase) == 0)
527 if (!gmx_strcasecmp("None", s))
534 if (c == NULL || ((c-s+1) == strlen(s)))
536 /* Check that we haven't already added a residue-specific version
539 for (j = 0; j < n && strstr((*list[j]).name, s) == NULL; j++)
556 list[n] = &(tb[none_idx]);
565 t_hackblock *choose_ter(int nb, t_hackblock **tb, const char *title)
569 printf("%s\n", title);
570 for (i = 0; (i < nb); i++)
572 char *advice_string = "";
573 if (0 == gmx_wcmatch("*ZWITTERION*", (*tb[i]).name))
575 advice_string = " (only use with zwitterions containing exactly one residue)";
577 printf("%2d: %s%s\n", i, (*tb[i]).name, advice_string);
581 ret = fscanf(stdin, "%d", &sel);
583 while ((ret != 1) || (sel < 0) || (sel >= nb));