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