71da9aece88da71844c2c2e3266763ae4fde770c
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / h_db.cpp
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,2015,2016,2017,2018, 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 <cstdlib>
43 #include <cstring>
44
45 #include <string>
46 #include <vector>
47
48 #include "gromacs/gmxpreprocess/fflibutil.h"
49 #include "gromacs/gmxpreprocess/notset.h"
50 #include "gromacs/utility/arraysize.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/smalloc.h"
55 #include "gromacs/utility/stringutil.h"
56
57 /* Number of control atoms for each 'add' type.
58  *
59  * There are 11 types of adding hydrogens, numbered from 1 thru
60  * 11. Each of these has a specific number of control atoms, that
61  * determine how the hydrogens are added.  Here these number are
62  * given. Because arrays start at 0, an extra dummy for index 0 is
63  * added.
64  */
65 const int ncontrol[] = { -1, 3, 3, 3, 3, 4, 3, 1, 3, 3, 1, 1 };
66 #define maxcontrol asize(ncontrol)
67
68 int compaddh(const void *a, const void *b)
69 {
70     const t_hackblock *ah, *bh;
71
72     ah = static_cast<const t_hackblock *>(a);
73     bh = static_cast<const t_hackblock *>(b);
74     return gmx_strcasecmp(ah->name, bh->name);
75 }
76
77 void print_ab(FILE *out, t_hack *hack, char *nname)
78 {
79     int i;
80
81     fprintf(out, "%d\t%d\t%s", hack->nr, hack->tp, nname);
82     for (i = 0; (i < hack->nctl); i++)
83     {
84         fprintf(out, "\t%s", hack->a[i]);
85     }
86     fprintf(out, "\n");
87 }
88
89
90 void read_ab(char *line, const char *fn, t_hack *hack)
91 {
92     int  i, nh, tp, ns;
93     char a[4][12];
94     char hn[32];
95
96     ns = sscanf(line, "%d%d%s%s%s%s%s", &nh, &tp, hn, a[0], a[1], a[2], a[3]);
97     if (ns < 4)
98     {
99         gmx_fatal(FARGS, "wrong format in input file %s on line\n%s\n", fn, line);
100     }
101
102     hack->nr = nh;
103     hack->tp = tp;
104     if ((tp < 1) || (tp >= maxcontrol))
105     {
106         gmx_fatal(FARGS, "Error in hdb file %s:\nH-type should be in 1-%d. Offending line:\n%s", fn, maxcontrol-1, line);
107     }
108
109     hack->nctl = ns - 3;
110     if ((hack->nctl != ncontrol[hack->tp]) && (ncontrol[hack->tp] != -1))
111     {
112         gmx_fatal(FARGS, "Error in hdb file %s:\nWrong number of control atoms (%d instead of %d) on line:\n%s\n", fn, hack->nctl, ncontrol[hack->tp], line);
113     }
114     for (i = 0; (i < hack->nctl); i++)
115     {
116         hack->a[i] = gmx_strdup(a[i]);
117     }
118     for (; i < 4; i++)
119     {
120         hack->a[i] = nullptr;
121     }
122     hack->oname = nullptr;
123     hack->nname = gmx_strdup(hn);
124     hack->atom  = nullptr;
125     hack->cgnr  = NOTSET;
126     hack->bXSet = FALSE;
127     for (i = 0; i < DIM; i++)
128     {
129         hack->newx[i] = NOTSET;
130     }
131 }
132
133 static void read_h_db_file(const char *hfn, int *nahptr, t_hackblock **ah)
134 {
135     FILE        *in;
136     char         filebase[STRLEN], line[STRLEN], buf[STRLEN];
137     int          i, n, nab, nah;
138     t_hackblock *aah;
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         // Skip lines that are only whitespace
153         if (gmx::countWords(line) == 0)
154         {
155             continue;
156         }
157         if (sscanf(line, "%s%n", buf, &n) != 1)
158         {
159             fprintf(stderr, "Error in hdb file: nah = %d\nline = '%s'\n",
160                     nah, line);
161             break;
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             snew(aah[nah].hack, nab);
171             aah[nah].nhack = nab;
172             for (i = 0; (i < nab); i++)
173             {
174                 if (feof(in))
175                 {
176                     gmx_fatal(FARGS, "Expected %d lines of hydrogens, found only %d "
177                               "while reading Hydrogen Database %s residue %s",
178                               nab, i-1, aah[nah].name, hfn);
179                 }
180                 if (nullptr == fgets(buf, STRLEN, in))
181                 {
182                     gmx_fatal(FARGS, "Error reading from file %s", hfn);
183                 }
184                 read_ab(buf, hfn, &(aah[nah].hack[i]));
185             }
186         }
187         nah++;
188     }
189     gmx_ffclose(in);
190
191     if (nah > 0)
192     {
193         /* Sort the list (necessary to be able to use bsearch */
194         qsort(aah, nah, static_cast<size_t>(sizeof(**ah)), compaddh);
195     }
196
197     *nahptr = nah;
198     *ah     = aah;
199 }
200
201 int read_h_db(const char *ffdir, t_hackblock **ah)
202 {
203     /* Read the hydrogen database file(s).
204      * Do not generate an error when no files are found.
205      */
206
207     std::vector<std::string> hdbf = fflib_search_file_end(ffdir, ".hdb", FALSE);
208     int nah                       = 0;
209     *ah   = nullptr;
210     for (const auto &filename : hdbf)
211     {
212         read_h_db_file(filename.c_str(), &nah, ah);
213     }
214     return nah;
215 }
216
217 t_hackblock *search_h_db(int nh, t_hackblock ah[], char *key)
218 {
219     t_hackblock ahkey, *result;
220
221     if (nh <= 0)
222     {
223         return nullptr;
224     }
225
226     ahkey.name = key;
227
228     result = static_cast<t_hackblock *>(bsearch(&ahkey, ah, nh, static_cast<size_t>(sizeof(ah[0])), compaddh));
229
230     return result;
231 }