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