Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / h_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 /* This file is completely threadsafe - keep it that way! */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <stdlib.h>
43 #include <string.h>
44
45 #include "gromacs/utility/cstringutil.h"
46 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/utility/futil.h"
48 #include "symtab.h"
49 #include "h_db.h"
50 #include "gromacs/fileio/gmxfio.h"
51 #include "fflibutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "macros.h"
54
55 /* Number of control atoms for each 'add' type.
56  *
57  * There are 11 types of adding hydrogens, numbered from 1 thru
58  * 11. Each of these has a specific number of control atoms, that
59  * determine how the hydrogens are added.  Here these number are
60  * given. Because arrays start at 0, an extra dummy for index 0 is
61  * added.
62  */
63 const int ncontrol[] = { -1, 3, 3, 3, 3, 4, 3, 1, 3, 3, 1, 1 };
64 #define maxcontrol asize(ncontrol)
65
66 int compaddh(const void *a, const void *b)
67 {
68     t_hackblock *ah, *bh;
69
70     ah = (t_hackblock *)a;
71     bh = (t_hackblock *)b;
72     return gmx_strcasecmp(ah->name, bh->name);
73 }
74
75 void print_ab(FILE *out, t_hack *hack, char *nname)
76 {
77     int i;
78
79     fprintf(out, "%d\t%d\t%s", hack->nr, hack->tp, nname);
80     for (i = 0; (i < hack->nctl); i++)
81     {
82         fprintf(out, "\t%s", hack->a[i]);
83     }
84     fprintf(out, "\n");
85 }
86
87
88 void read_ab(char *line, const char *fn, t_hack *hack)
89 {
90     int  i, nh, tp, ns;
91     char a[4][12];
92     char hn[32];
93
94     ns = sscanf(line, "%d%d%s%s%s%s%s", &nh, &tp, hn, a[0], a[1], a[2], a[3]);
95     if (ns < 4)
96     {
97         gmx_fatal(FARGS, "wrong format in input file %s on line\n%s\n", fn, line);
98     }
99
100     hack->nr = nh;
101     hack->tp = tp;
102     if ((tp < 1) || (tp >= maxcontrol))
103     {
104         gmx_fatal(FARGS, "Error in hdb file %s:\nH-type should be in 1-%d. Offending line:\n%s", fn, maxcontrol-1, line);
105     }
106
107     hack->nctl = ns - 3;
108     if ((hack->nctl != ncontrol[hack->tp]) && (ncontrol[hack->tp] != -1))
109     {
110         gmx_fatal(FARGS, "Error in hdb file %s:\nWrong number of control atoms (%d iso %d) on line:\n%s\n", fn, hack->nctl, ncontrol[hack->tp], line);
111     }
112     for (i = 0; (i < hack->nctl); i++)
113     {
114         hack->a[i] = strdup(a[i]);
115     }
116     for (; i < 4; i++)
117     {
118         hack->a[i] = NULL;
119     }
120     hack->oname = NULL;
121     hack->nname = strdup(hn);
122     hack->atom  = NULL;
123     hack->cgnr  = NOTSET;
124     hack->bXSet = FALSE;
125     for (i = 0; i < DIM; i++)
126     {
127         hack->newx[i] = NOTSET;
128     }
129 }
130
131 static void read_h_db_file(const char *hfn, int *nahptr, t_hackblock **ah)
132 {
133     FILE        *in;
134     char         filebase[STRLEN], line[STRLEN], buf[STRLEN];
135     int          i, n, nab, nah;
136     t_hackblock *aah;
137
138     if (debug)
139     {
140         fprintf(debug, "Hydrogen Database (%s):\n", hfn);
141     }
142
143     fflib_filename_base(hfn, filebase, STRLEN);
144     /* Currently filebase is read and set, but not used.
145      * hdb entries from any hdb file and be applied to rtp entries
146      * in any rtp file.
147      */
148
149     in = fflib_open(hfn);
150
151     nah = *nahptr;
152     aah = *ah;
153     while (fgets2(line, STRLEN-1, in))
154     {
155         if (sscanf(line, "%s%n", buf, &n) != 1)
156         {
157             fprintf(stderr, "Error in hdb file: nah = %d\nline = '%s'\n",
158                     nah, line);
159             break;
160         }
161         if (debug)
162         {
163             fprintf(debug, "%s", buf);
164         }
165         srenew(aah, nah+1);
166         clear_t_hackblock(&aah[nah]);
167         aah[nah].name     = strdup(buf);
168         aah[nah].filebase = strdup(filebase);
169
170         if (sscanf(line+n, "%d", &nab) == 1)
171         {
172             if (debug)
173             {
174                 fprintf(debug, "  %d\n", nab);
175             }
176             snew(aah[nah].hack, nab);
177             aah[nah].nhack = nab;
178             for (i = 0; (i < nab); i++)
179             {
180                 if (feof(in))
181                 {
182                     gmx_fatal(FARGS, "Expected %d lines of hydrogens, found only %d "
183                               "while reading Hydrogen Database %s residue %s",
184                               nab, i-1, aah[nah].name, hfn);
185                 }
186                 if (NULL == fgets(buf, STRLEN, in))
187                 {
188                     gmx_fatal(FARGS, "Error reading from file %s", hfn);
189                 }
190                 read_ab(buf, hfn, &(aah[nah].hack[i]));
191             }
192         }
193         nah++;
194     }
195     gmx_ffclose(in);
196
197     /* Sort the list (necessary to be able to use bsearch */
198     qsort(aah, nah, (size_t)sizeof(**ah), compaddh);
199
200     /*
201        if (debug)
202        dump_h_db(hfn,nah,aah);
203      */
204
205     *nahptr = nah;
206     *ah     = aah;
207 }
208
209 int read_h_db(const char *ffdir, t_hackblock **ah)
210 {
211     int    nhdbf, f;
212     char **hdbf;
213     int    nah;
214     FILE  *fp;
215
216     /* Read the hydrogen database file(s).
217      * Do not generate an error when no files are found.
218      */
219     nhdbf = fflib_search_file_end(ffdir, ".hdb", FALSE, &hdbf);
220     nah   = 0;
221     *ah   = NULL;
222     for (f = 0; f < nhdbf; f++)
223     {
224         read_h_db_file(hdbf[f], &nah, ah);
225         sfree(hdbf[f]);
226     }
227     sfree(hdbf);
228
229     return nah;
230 }
231
232 t_hackblock *search_h_db(int nh, t_hackblock ah[], char *key)
233 {
234     t_hackblock ahkey, *result;
235
236     if (nh <= 0)
237     {
238         return NULL;
239     }
240
241     ahkey.name = key;
242
243     result = (t_hackblock *)bsearch(&ahkey, ah, nh, (size_t)sizeof(ah[0]), compaddh);
244
245     return result;
246 }