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