e55ec6d115109dcf30231b4e7b7f65f79b6aadc9
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / ter_db.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #include <ctype.h>
42 #include <string.h>
43
44 #include "gromacs/utility/smalloc.h"
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "gromacs/utility/futil.h"
47 #include "resall.h"
48 #include "h_db.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "ter_db.h"
52 #include "toputil.h"
53 #include "gromacs/fileio/gmxfio.h"
54 #include "fflibutil.h"
55
56 #include "gromacs/fileio/strdb.h"
57
58 /* use bonded types definitions in hackblock.h */
59 #define ekwRepl ebtsNR+1
60 #define ekwAdd  ebtsNR+2
61 #define ekwDel  ebtsNR+3
62 #define ekwNR   3
63 const char *kw_names[ekwNR] = {
64     "replace", "add", "delete"
65 };
66
67 int find_kw(char *keyw)
68 {
69     int i;
70
71     for (i = 0; i < ebtsNR; i++)
72     {
73         if (gmx_strcasecmp(btsNames[i], keyw) == 0)
74         {
75             return i;
76         }
77     }
78     for (i = 0; i < ekwNR; i++)
79     {
80         if (gmx_strcasecmp(kw_names[i], keyw) == 0)
81         {
82             return ebtsNR + 1 + i;
83         }
84     }
85
86     return NOTSET;
87 }
88
89 #define FATAL() gmx_fatal(FARGS, "Reading Termini Database: not enough items on line\n%s", line)
90
91 static void read_atom(char *line, gmx_bool bAdd,
92                       char **nname, t_atom *a, gpp_atomtype_t atype, int *cgnr)
93 {
94     int    nr, i;
95     char   buf[5][30], type[12];
96     double m, q;
97
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]
101      */
102     nr = sscanf(line, "%s %s %s %s %s", buf[0], buf[1], buf[2], buf[3], buf[4]);
103
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.
108      */
109     if (!bAdd && nr == 4 && isdigit(buf[1][0]))
110     {
111         nr = 3;
112     }
113
114     if (nr < 3 || nr > 4)
115     {
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);
117     }
118     i = 0;
119     if (!bAdd)
120     {
121         if (nr == 4)
122         {
123             *nname = gmx_strdup(buf[i++]);
124         }
125         else
126         {
127             *nname = NULL;
128         }
129     }
130     a->type = get_atomtype_type(buf[i++], atype);
131     sscanf(buf[i++], "%lf", &m);
132     a->m = m;
133     sscanf(buf[i++], "%lf", &q);
134     a->q = q;
135     if (bAdd && nr == 4)
136     {
137         sscanf(buf[i++], "%d", cgnr);
138     }
139     else
140     {
141         *cgnr = NOTSET;
142     }
143 }
144
145 static void print_atom(FILE *out, t_atom *a, gpp_atomtype_t atype)
146 {
147     fprintf(out, "\t%s\t%g\t%g\n",
148             get_atomtype_name(a->type, atype), a->m, a->q);
149 }
150
151 static void print_ter_db(const char *ff, char C, int nb, t_hackblock tb[],
152                          gpp_atomtype_t atype)
153 {
154     FILE *out;
155     int   i, j, k, bt, nrepl, nadd, ndel;
156     char  buf[STRLEN], nname[STRLEN];
157
158     sprintf(buf, "%s-%c.tdb", ff, C);
159     out = gmx_fio_fopen(buf, "w");
160
161     for (i = 0; (i < nb); i++)
162     {
163         fprintf(out, "[ %s ]\n", tb[i].name);
164
165         /* first count: */
166         nrepl = 0;
167         nadd  = 0;
168         ndel  = 0;
169         for (j = 0; j < tb[i].nhack; j++)
170         {
171             if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname != NULL)
172             {
173                 nrepl++;
174             }
175             else if (tb[i].hack[j].oname == NULL && tb[i].hack[j].nname != NULL)
176             {
177                 nadd++;
178             }
179             else if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname == NULL)
180             {
181                 ndel++;
182             }
183             else if (tb[i].hack[j].oname == NULL && tb[i].hack[j].nname == NULL)
184             {
185                 gmx_fatal(FARGS, "invalid hack (%s) in termini database", tb[i].name);
186             }
187         }
188         if (nrepl)
189         {
190             fprintf(out, "[ %s ]\n", kw_names[ekwRepl-ebtsNR-1]);
191             for (j = 0; j < tb[i].nhack; j++)
192             {
193                 if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname != NULL)
194                 {
195                     fprintf(out, "%s\t", tb[i].hack[j].oname);
196                     print_atom(out, tb[i].hack[j].atom, atype);
197                 }
198             }
199         }
200         if (nadd)
201         {
202             fprintf(out, "[ %s ]\n", kw_names[ekwAdd-ebtsNR-1]);
203             for (j = 0; j < tb[i].nhack; j++)
204             {
205                 if (tb[i].hack[j].oname == NULL && tb[i].hack[j].nname != NULL)
206                 {
207                     print_ab(out, &(tb[i].hack[j]), tb[i].hack[j].nname);
208                     print_atom(out, tb[i].hack[j].atom, atype);
209                 }
210             }
211         }
212         if (ndel)
213         {
214             fprintf(out, "[ %s ]\n", kw_names[ekwDel-ebtsNR-1]);
215             for (j = 0; j < tb[i].nhack; j++)
216             {
217                 if (tb[i].hack[j].oname != NULL && tb[i].hack[j].nname == NULL)
218                 {
219                     fprintf(out, "%s\n", tb[i].hack[j].oname);
220                 }
221             }
222         }
223         for (bt = 0; bt < ebtsNR; bt++)
224         {
225             if (tb[i].rb[bt].nb)
226             {
227                 fprintf(out, "[ %s ]\n", btsNames[bt]);
228                 for (j = 0; j < tb[i].rb[bt].nb; j++)
229                 {
230                     for (k = 0; k < btsNiatoms[bt]; k++)
231                     {
232                         fprintf(out, "%s%s", k ? "\t" : "", tb[i].rb[bt].b[j].a[k]);
233                     }
234                     if (tb[i].rb[bt].b[j].s)
235                     {
236                         fprintf(out, "\t%s", tb[i].rb[bt].b[j].s);
237                     }
238                     fprintf(out, "\n");
239                 }
240             }
241         }
242         fprintf(out, "\n");
243     }
244     gmx_fio_fclose(out);
245 }
246
247 static void read_ter_db_file(char *fn,
248                              int *ntbptr, t_hackblock **tbptr,
249                              gpp_atomtype_t atype)
250 {
251     char         filebase[STRLEN], *ptr;
252     FILE        *in;
253     char         header[STRLEN], buf[STRLEN], line[STRLEN];
254     t_hackblock *tb;
255     int          i, j, n, ni, kwnr, nb, maxnb, nh;
256
257     fflib_filename_base(fn, filebase, STRLEN);
258     /* Remove the C/N termini extension */
259     ptr = strrchr(filebase, '.');
260     if (ptr != NULL)
261     {
262         ptr[0] = '\0';
263     }
264
265     in = fflib_open(fn);
266     if (debug)
267     {
268         fprintf(debug, "Opened %s\n", fn);
269     }
270
271     tb    = *tbptr;
272     nb    = *ntbptr - 1;
273     maxnb = 0;
274     kwnr  = NOTSET;
275     get_a_line(in, line, STRLEN);
276     while (!feof(in))
277     {
278         if (get_header(line, header))
279         {
280             /* this is a new block, or a new keyword */
281             kwnr = find_kw(header);
282
283             if (kwnr == NOTSET)
284             {
285                 nb++;
286                 /* here starts a new block */
287                 if (nb >= maxnb)
288                 {
289                     maxnb = nb + 100;
290                     srenew(tb, maxnb);
291                 }
292                 clear_t_hackblock(&tb[nb]);
293                 tb[nb].name     = gmx_strdup(header);
294                 tb[nb].filebase = gmx_strdup(filebase);
295             }
296         }
297         else
298         {
299             if (nb < 0)
300             {
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);
304             }
305             /* this is not a header, so it must be data */
306             if (kwnr >= ebtsNR)
307             {
308                 /* this is a hack: add/rename/delete atoms */
309                 /* make space for hacks */
310                 if (tb[nb].nhack >= tb[nb].maxhack)
311                 {
312                     tb[nb].maxhack += 10;
313                     srenew(tb[nb].hack, tb[nb].maxhack);
314                 }
315                 nh = tb[nb].nhack;
316                 clear_t_hack(&(tb[nb].hack[nh]));
317                 for (i = 0; i < 4; i++)
318                 {
319                     tb[nb].hack[nh].a[i] = NULL;
320                 }
321                 tb[nb].nhack++;
322
323                 /* get data */
324                 n = 0;
325                 if (kwnr == ekwRepl || kwnr == ekwDel)
326                 {
327                     if (sscanf(line, "%s%n", buf, &n) != 1)
328                     {
329                         gmx_fatal(FARGS, "Reading Termini Database '%s': "
330                                   "expected atom name on line\n%s", fn, line);
331                     }
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;
335                 }
336                 else if (kwnr == ekwAdd)
337                 {
338                     read_ab(line, fn, &(tb[nb].hack[nh]));
339                     get_a_line(in, line, STRLEN);
340                 }
341                 else
342                 {
343                     gmx_fatal(FARGS, "unimplemented keyword number %d (%s:%d)",
344                               kwnr, __FILE__, __LINE__);
345                 }
346                 if (kwnr == ekwRepl || kwnr == ekwAdd)
347                 {
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)
353                     {
354                         if (tb[nb].hack[nh].oname != NULL)
355                         {
356                             tb[nb].hack[nh].nname = gmx_strdup(tb[nb].hack[nh].oname);
357                         }
358                         else
359                         {
360                             gmx_fatal(FARGS, "Reading Termini Database '%s': don't know which name the new atom should have on line\n%s", fn, line);
361                         }
362                     }
363                 }
364             }
365             else if (kwnr >= 0 && kwnr < ebtsNR)
366             {
367                 /* this is bonded data: bonds, angles, dihedrals or impropers */
368                 srenew(tb[nb].rb[kwnr].b, tb[nb].rb[kwnr].nb+1);
369                 n = 0;
370                 for (j = 0; j < btsNiatoms[kwnr]; j++)
371                 {
372                     if (sscanf(line+n, "%s%n", buf, &ni) == 1)
373                     {
374                         tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].a[j] = gmx_strdup(buf);
375                     }
376                     else
377                     {
378                         gmx_fatal(FARGS, "Reading Termini Database '%s': expected %d atom names (found %d) on line\n%s", fn, btsNiatoms[kwnr], j-1, line);
379                     }
380                     n += ni;
381                 }
382                 for (; j < MAXATOMLIST; j++)
383                 {
384                     tb[nb].rb[kwnr].b[tb[nb].rb[kwnr].nb].a[j] = NULL;
385                 }
386                 strcpy(buf, "");
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++;
390             }
391             else
392             {
393                 gmx_fatal(FARGS, "Reading Termini Database: Expecting a header at line\n"
394                           "%s", line);
395             }
396         }
397         get_a_line(in, line, STRLEN);
398     }
399     nb++;
400     srenew(tb, nb);
401
402     gmx_ffclose(in);
403
404     *ntbptr = nb;
405     *tbptr  = tb;
406 }
407
408 int read_ter_db(const char *ffdir, char ter,
409                 t_hackblock **tbptr, gpp_atomtype_t atype)
410 {
411     char   ext[STRLEN];
412     int    ntdbf, f;
413     char **tdbf;
414     int    ntb;
415
416     sprintf(ext, ".%c.tdb", ter);
417
418     /* Search for termini database files.
419      * Do not generate an error when none are found.
420      */
421     ntdbf  = fflib_search_file_end(ffdir, ext, FALSE, &tdbf);
422     ntb    = 0;
423     *tbptr = NULL;
424     for (f = 0; f < ntdbf; f++)
425     {
426         read_ter_db_file(tdbf[f], &ntb, tbptr, atype);
427         sfree(tdbf[f]);
428     }
429     sfree(tdbf);
430
431     if (debug)
432     {
433         print_ter_db("new", ter, ntb, *tbptr, atype);
434     }
435
436     return ntb;
437 }
438
439 t_hackblock **filter_ter(int nrtp, t_restp rtp[],
440                          int nb, t_hackblock tb[],
441                          const char *resname,
442                          const char *rtpname,
443                          int *nret)
444 {
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+).
449      *
450      * To reduce the database size, we assume that a terminus specifier liker
451      *
452      * [ GLY|SER|ALA-NH3+ ]
453      *
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.
457      *
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...
464      */
465
466     t_restp     *   restp;
467     int             i, j, n, len, none_idx;
468     gmx_bool        found;
469     char           *rtpname_match, *s, *s2, *c;
470     t_hackblock   **list;
471
472     rtpname_match = search_rtp(rtpname, nrtp, rtp);
473     restp         = get_restp(rtpname_match, nrtp, rtp);
474
475     n    = 0;
476     list = NULL;
477
478     for (i = 0; i < nb; i++)
479     {
480         s     = tb[i].name;
481         found = FALSE;
482         do
483         {
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
487              * much cleaner.
488              */
489             if (gmx_strcasecmp(restp->filebase, tb[i].filebase) == 0 &&
490                 gmx_strncasecmp(resname, s, 3) == 0)
491             {
492                 found = TRUE;
493                 srenew(list, n+1);
494                 list[n] = &(tb[i]);
495                 n++;
496             }
497             else
498             {
499                 /* advance to next |-separated field */
500                 s = strchr(s, '|');
501                 if (s != NULL)
502                 {
503                     s++;
504                 }
505             }
506         }
507         while (!found && s != NULL);
508     }
509
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.
515      */
516     none_idx = -1;
517     for (i = 0; i < nb; i++)
518     {
519         s = tb[i].name;
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
523          * much cleaner.
524          */
525         if (gmx_strcasecmp(restp->filebase, tb[i].filebase) == 0)
526         {
527             if (!gmx_strcasecmp("None", s))
528             {
529                 none_idx = i;
530             }
531             else
532             {
533                 c = strchr(s, '-');
534                 if (c == NULL || ((c-s+1) == strlen(s)))
535                 {
536                     /* Check that we haven't already added a residue-specific version
537                      * of this terminus.
538                      */
539                     for (j = 0; j < n && strstr((*list[j]).name, s) == NULL; j++)
540                     {
541                         ;
542                     }
543                     if (j == n)
544                     {
545                         srenew(list, n+1);
546                         list[n] = &(tb[i]);
547                         n++;
548                     }
549                 }
550             }
551         }
552     }
553     if (none_idx >= 0)
554     {
555         srenew(list, n+1);
556         list[n] = &(tb[none_idx]);
557         n++;
558     }
559
560     *nret = n;
561     return list;
562 }
563
564
565 t_hackblock *choose_ter(int nb, t_hackblock **tb, const char *title)
566 {
567     int i, sel, ret;
568
569     printf("%s\n", title);
570     for (i = 0; (i < nb); i++)
571     {
572         char *advice_string = "";
573         if (0 == gmx_wcmatch("*ZWITTERION*", (*tb[i]).name))
574         {
575             advice_string = " (only use with zwitterions containing exactly one residue)";
576         }
577         printf("%2d: %s%s\n", i, (*tb[i]).name, advice_string);
578     }
579     do
580     {
581         ret = fscanf(stdin, "%d", &sel);
582     }
583     while ((ret != 1) || (sel < 0) || (sel >= nb));
584
585     return tb[sel];
586 }