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.
42 #include "gromacs/utility/smalloc.h"
43 #include "gromacs/legacyheaders/typedefs.h"
44 #include "gromacs/utility/futil.h"
47 #include "gromacs/utility/cstringutil.h"
48 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/fileio/gmxfio.h"
52 #include "fflibutil.h"
54 #include "gromacs/fileio/strdb.h"
56 /* use bonded types definitions in hackblock.h */
57 #define ekwRepl ebtsNR+1
58 #define ekwAdd ebtsNR+2
59 #define ekwDel ebtsNR+3
61 const char *kw_names[ekwNR] = {
62 "replace", "add", "delete"
65 int find_kw(char *keyw)
69 for (i = 0; i < ebtsNR; i++)
71 if (gmx_strcasecmp(btsNames[i], keyw) == 0)
76 for (i = 0; i < ekwNR; i++)
78 if (gmx_strcasecmp(kw_names[i], keyw) == 0)
80 return ebtsNR + 1 + i;
87 #define FATAL() gmx_fatal(FARGS, "Reading Termini Database: not enough items on line\n%s", line)
89 static void read_atom(char *line, gmx_bool bAdd,
90 char **nname, t_atom *a, gpp_atomtype_t atype, int *cgnr)
93 char buf[5][30], type[12];
96 /* This code is messy, because of support for different formats:
97 * for replace: [new name] <atom type> <m> <q> [cgnr (old format)]
98 * for add: <atom type> <m> <q> [cgnr]
100 nr = sscanf(line, "%s %s %s %s %s", buf[0], buf[1], buf[2], buf[3], buf[4]);
102 /* Here there an ambiguity due to the old replace format with cgnr,
103 * which was read for years, but ignored in the rest of the code.
104 * We have to assume that the atom type does not start with a digit
105 * to make a line with 4 entries uniquely interpretable.
107 if (!bAdd && nr == 4 && isdigit(buf[1][0]))
112 if (nr < 3 || nr > 4)
114 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);
121 *nname = gmx_strdup(buf[i++]);
128 a->type = get_atomtype_type(buf[i++], atype);
129 sscanf(buf[i++], "%lf", &m);
131 sscanf(buf[i++], "%lf", &q);
135 sscanf(buf[i++], "%d", cgnr);
143 static void print_atom(FILE *out, t_atom *a, gpp_atomtype_t atype)
145 fprintf(out, "\t%s\t%g\t%g\n",
146 get_atomtype_name(a->type, atype), a->m, a->q);
149 static void print_ter_db(const char *ff, char C, int nb, t_hackblock tb[],
150 gpp_atomtype_t atype)
153 int i, j, k, bt, nrepl, nadd, ndel;
154 char buf[STRLEN], nname[STRLEN];
156 sprintf(buf, "%s-%c.tdb", ff, C);
157 out = gmx_fio_fopen(buf, "w");
159 for (i = 0; (i < nb); i++)
161 fprintf(out, "[ %s ]\n", tb[i].name);
167 for (j = 0; j < tb[i].nhack; j++)
169 if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname != NULL)
173 else if (tb[i].hack[j].oname == NULL && tb[i].hack[j].nname != NULL)
177 else if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname == NULL)
181 else if (tb[i].hack[j].oname == NULL && tb[i].hack[j].nname == NULL)
183 gmx_fatal(FARGS, "invalid hack (%s) in termini database", tb[i].name);
188 fprintf(out, "[ %s ]\n", kw_names[ekwRepl-ebtsNR-1]);
189 for (j = 0; j < tb[i].nhack; j++)
191 if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname != NULL)
193 fprintf(out, "%s\t", tb[i].hack[j].oname);
194 print_atom(out, tb[i].hack[j].atom, atype);
200 fprintf(out, "[ %s ]\n", kw_names[ekwAdd-ebtsNR-1]);
201 for (j = 0; j < tb[i].nhack; j++)
203 if (tb[i].hack[j].oname == NULL && tb[i].hack[j].nname != NULL)
205 print_ab(out, &(tb[i].hack[j]), tb[i].hack[j].nname);
206 print_atom(out, tb[i].hack[j].atom, atype);
212 fprintf(out, "[ %s ]\n", kw_names[ekwDel-ebtsNR-1]);
213 for (j = 0; j < tb[i].nhack; j++)
215 if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname == NULL)
217 fprintf(out, "%s\n", tb[i].hack[j].oname);
221 for (bt = 0; bt < ebtsNR; bt++)
225 fprintf(out, "[ %s ]\n", btsNames[bt]);
226 for (j = 0; j < tb[i].rb[bt].nb; j++)
228 for (k = 0; k < btsNiatoms[bt]; k++)
230 fprintf(out, "%s%s", k ? "\t" : "", tb[i].rb[bt].b[j].a[k]);
232 if (tb[i].rb[bt].b[j].s)
234 fprintf(out, "\t%s", tb[i].rb[bt].b[j].s);
245 static void read_ter_db_file(char *fn,
246 int *ntbptr, t_hackblock **tbptr,
247 gpp_atomtype_t atype)
249 char filebase[STRLEN], *ptr;
251 char header[STRLEN], buf[STRLEN], line[STRLEN];
253 int i, j, n, ni, kwnr, nb, maxnb, nh;
255 fflib_filename_base(fn, filebase, STRLEN);
256 /* Remove the C/N termini extension */
257 ptr = strrchr(filebase, '.');
266 fprintf(debug, "Opened %s\n", fn);
273 get_a_line(in, line, STRLEN);
276 if (get_header(line, header))
278 /* this is a new block, or a new keyword */
279 kwnr = find_kw(header);
284 /* here starts a new block */
290 clear_t_hackblock(&tb[nb]);
291 tb[nb].name = gmx_strdup(header);
292 tb[nb].filebase = gmx_strdup(filebase);
299 gmx_fatal(FARGS, "reading termini database: "
300 "directive expected before line:\n%s\n"
301 "This might be a file in an old format.", line);
303 /* this is not a header, so it must be data */
306 /* this is a hack: add/rename/delete atoms */
307 /* make space for hacks */
308 if (tb[nb].nhack >= tb[nb].maxhack)
310 tb[nb].maxhack += 10;
311 srenew(tb[nb].hack, tb[nb].maxhack);
314 clear_t_hack(&(tb[nb].hack[nh]));
315 for (i = 0; i < 4; i++)
317 tb[nb].hack[nh].a[i] = NULL;
323 if (kwnr == ekwRepl || kwnr == ekwDel)
325 if (sscanf(line, "%s%n", buf, &n) != 1)
327 gmx_fatal(FARGS, "Reading Termini Database '%s': "
328 "expected atom name on line\n%s", fn, line);
330 tb[nb].hack[nh].oname = gmx_strdup(buf);
331 /* we only replace or delete one atom at a time */
332 tb[nb].hack[nh].nr = 1;
334 else if (kwnr == ekwAdd)
336 read_ab(line, fn, &(tb[nb].hack[nh]));
337 get_a_line(in, line, STRLEN);
341 gmx_fatal(FARGS, "unimplemented keyword number %d (%s:%d)",
342 kwnr, __FILE__, __LINE__);
344 if (kwnr == ekwRepl || kwnr == ekwAdd)
346 snew(tb[nb].hack[nh].atom, 1);
347 read_atom(line+n, kwnr == ekwAdd,
348 &tb[nb].hack[nh].nname, tb[nb].hack[nh].atom, atype,
349 &tb[nb].hack[nh].cgnr);
350 if (tb[nb].hack[nh].nname == NULL)
352 if (tb[nb].hack[nh].oname != NULL)
354 tb[nb].hack[nh].nname = gmx_strdup(tb[nb].hack[nh].oname);
358 gmx_fatal(FARGS, "Reading Termini Database '%s': don't know which name the new atom should have on line\n%s", fn, line);
363 else if (kwnr >= 0 && kwnr < ebtsNR)
365 /* this is bonded data: bonds, angles, dihedrals or impropers */
366 srenew(tb[nb].rb[kwnr].b, tb[nb].rb[kwnr].nb+1);
368 for (j = 0; j < btsNiatoms[kwnr]; j++)
370 if (sscanf(line+n, "%s%n", buf, &ni) == 1)
372 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].a[j] = gmx_strdup(buf);
376 gmx_fatal(FARGS, "Reading Termini Database '%s': expected %d atom names (found %d) on line\n%s", fn, btsNiatoms[kwnr], j-1, line);
380 for (; j < MAXATOMLIST; j++)
382 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].a[j] = NULL;
385 sscanf(line+n, "%s", buf);
386 tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].s = gmx_strdup(buf);
387 tb[nb].rb[kwnr].nb++;
391 gmx_fatal(FARGS, "Reading Termini Database: Expecting a header at line\n"
395 get_a_line(in, line, STRLEN);
406 int read_ter_db(const char *ffdir, char ter,
407 t_hackblock **tbptr, gpp_atomtype_t atype)
414 sprintf(ext, ".%c.tdb", ter);
416 /* Search for termini database files.
417 * Do not generate an error when none are found.
419 ntdbf = fflib_search_file_end(ffdir, ext, FALSE, &tdbf);
422 for (f = 0; f < ntdbf; f++)
424 read_ter_db_file(tdbf[f], &ntb, tbptr, atype);
431 print_ter_db("new", ter, ntb, *tbptr, atype);
437 t_hackblock **filter_ter(int nrtp, t_restp rtp[],
438 int nb, t_hackblock tb[],
443 /* Since some force fields (e.g. OPLS) needs different
444 * atomtypes for different residues there could be a lot
445 * of entries in the databases for specific residues
446 * (e.g. GLY-NH3+, SER-NH3+, ALA-NH3+).
448 * To reduce the database size, we assume that a terminus specifier liker
450 * [ GLY|SER|ALA-NH3+ ]
452 * would cover all of the three residue types above.
453 * Wildcards (*,?) are not OK. Don't worry about e.g. GLU vs. GLUH since
454 * pdb2gmx only uses the first 3 letters when calling this routine.
456 * To automate this, this routines scans a list of termini
457 * for the residue name "resname" and returns an allocated list of
458 * pointers to the termini that could be applied to the
459 * residue in question. The variable pointed to by nret will
460 * contain the number of valid pointers in the list.
461 * Remember to free the list when you are done with it...
465 int i, j, n, len, none_idx;
467 char *rtpname_match, *s, *s2, *c;
470 rtpname_match = search_rtp(rtpname, nrtp, rtp);
471 restp = get_restp(rtpname_match, nrtp, rtp);
476 for (i = 0; i < nb; i++)
482 /* The residue name should appear in a tdb file with the same base name
483 * as the file containing the rtp entry.
484 * This makes termini selection for different molecule types
487 if (gmx_strcasecmp(restp->filebase, tb[i].filebase) == 0 &&
488 gmx_strncasecmp(resname, s, 3) == 0)
497 /* advance to next |-separated field */
505 while (!found && s != NULL);
508 /* All residue-specific termini have been added. See if there
509 * are some generic ones by searching for the occurence of
510 * '-' in the name prior to the last position (which indicates charge).
511 * The [ None ] alternative is special since we don't want that
512 * to be the default, so we put it last in the list we return.
515 for (i = 0; i < nb; i++)
518 /* The residue name should appear in a tdb file with the same base name
519 * as the file containing the rtp entry.
520 * This makes termini selection for different molecule types
523 if (gmx_strcasecmp(restp->filebase, tb[i].filebase) == 0)
525 if (!gmx_strcasecmp("None", s))
532 if (c == NULL || ((c-s+1) == strlen(s)))
534 /* Check that we haven't already added a residue-specific version
537 for (j = 0; j < n && strstr((*list[j]).name, s) == NULL; j++)
554 list[n] = &(tb[none_idx]);
563 t_hackblock *choose_ter(int nb, t_hackblock **tb, const char *title)
567 printf("%s\n", title);
568 for (i = 0; (i < nb); i++)
570 char *advice_string = "";
571 if (0 == gmx_wcmatch("*ZWITTERION*", (*tb[i]).name))
573 advice_string = " (only use with zwitterions containing exactly one residue)";
575 printf("%2d: %s%s\n", i, (*tb[i]).name, advice_string);
579 ret = fscanf(stdin, "%d", &sel);
581 while ((ret != 1) || (sel < 0) || (sel >= nb));