Improve use of gmxpreprocess headers
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / nm2type.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-2008, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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 "nm2type.h"
41
42 #include <cstring>
43
44 #include <algorithm>
45 #include <string>
46 #include <vector>
47
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/gmxpreprocess/fflibutil.h"
50 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
51 #include "gromacs/gmxpreprocess/gpp_nextnb.h"
52 #include "gromacs/gmxpreprocess/grompp-impl.h"
53 #include "gromacs/gmxpreprocess/notset.h"
54 #include "gromacs/gmxpreprocess/pdb2top.h"
55 #include "gromacs/gmxpreprocess/toppush.h"
56 #include "gromacs/math/utilities.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdtypes/md_enums.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
63
64 static void rd_nm2type_file(const std::string &filename, int *nnm, t_nm2type **nmp)
65 {
66     FILE         *fp;
67     bool          bCont;
68     char          libfilename[128];
69     char          format[128], f1[128];
70     char          buf[1024], elem[16], type[16], nbbuf[16], **newbuf;
71     int           i, nb, nnnm, line = 1;
72     double        qq, mm;
73     t_nm2type    *nm2t = nullptr;
74
75     fp = fflib_open(filename);
76     if (nullptr == fp)
77     {
78         gmx_fatal(FARGS, "Can not find %s in library directory", filename.c_str());
79     }
80
81     nnnm = *nnm;
82     nm2t = *nmp;
83     do
84     {
85         /* Read a line from the file */
86         bCont = (fgets2(buf, 1023, fp) != nullptr);
87
88         if (bCont)
89         {
90             /* Remove comment */
91             strip_comment(buf);
92             if (sscanf(buf, "%s%s%lf%lf%d", elem, type, &qq, &mm, &nb) == 5)
93             {
94                 /* If we can read the first four, there probably is more */
95                 srenew(nm2t, nnnm+1);
96                 snew(nm2t[nnnm].blen, nb);
97                 if (nb > 0)
98                 {
99                     snew(newbuf, nb);
100                     strcpy(format, "%*s%*s%*s%*s%*s");
101                     for (i = 0; (i < nb); i++)
102                     {
103                         /* Complicated format statement */
104                         strcpy(f1, format);
105                         strcat(f1, "%s%lf");
106                         if (sscanf(buf, f1, nbbuf, &(nm2t[nnnm].blen[i])) != 2)
107                         {
108                             gmx_fatal(FARGS, "Error on line %d of %s", line, libfilename);
109                         }
110                         newbuf[i] = gmx_strdup(nbbuf);
111                         strcat(format, "%*s%*s");
112                     }
113                 }
114                 else
115                 {
116                     newbuf = nullptr;
117                 }
118                 nm2t[nnnm].elem   = gmx_strdup(elem);
119                 nm2t[nnnm].type   = gmx_strdup(type);
120                 nm2t[nnnm].q      = qq;
121                 nm2t[nnnm].m      = mm;
122                 nm2t[nnnm].nbonds = nb;
123                 nm2t[nnnm].bond   = newbuf;
124                 nnnm++;
125             }
126             line++;
127         }
128     }
129     while (bCont);
130     gmx_ffclose(fp);
131
132     *nnm = nnnm;
133     *nmp = nm2t;
134 }
135
136 t_nm2type *rd_nm2type(const char *ffdir, int *nnm)
137 {
138     std::vector<std::string> ff  = fflib_search_file_end(ffdir, ".n2t", FALSE);
139     *nnm = 0;
140     t_nm2type               *nm = nullptr;
141     for (const auto &filename : ff)
142     {
143         rd_nm2type_file(filename, nnm, &nm);
144     }
145     return nm;
146 }
147
148 void dump_nm2type(FILE *fp, int nnm, t_nm2type nm2t[])
149 {
150     int i, j;
151
152     fprintf(fp, "; nm2type database\n");
153     for (i = 0; (i < nnm); i++)
154     {
155         fprintf(fp, "%-8s %-8s %8.4f %8.4f %-4d",
156                 nm2t[i].elem, nm2t[i].type,
157                 nm2t[i].q, nm2t[i].m, nm2t[i].nbonds);
158         for (j = 0; (j < nm2t[i].nbonds); j++)
159         {
160             fprintf(fp, " %-5s %6.4f", nm2t[i].bond[j], nm2t[i].blen[j]);
161         }
162         fprintf(fp, "\n");
163     }
164 }
165
166 enum {
167     ematchNone, ematchWild, ematchElem, ematchExact, ematchNR
168 };
169
170 static int match_str(const char *atom, const char *template_string)
171 {
172     if (!atom || !template_string)
173     {
174         return ematchNone;
175     }
176     else if (gmx_strcasecmp(atom, template_string) == 0)
177     {
178         return ematchExact;
179     }
180     else if (atom[0] == template_string[0])
181     {
182         return ematchElem;
183     }
184     else if (strcmp(template_string, "*") == 0)
185     {
186         return ematchWild;
187     }
188     else
189     {
190         return ematchNone;
191     }
192 }
193
194 int nm2type(int nnm, t_nm2type nm2t[], t_symtab *tab, t_atoms *atoms,
195             gpp_atomtype *atype, int *nbonds, t_params *bonds)
196 {
197     int      cur = 0;
198 #define prev (1-cur)
199     int      i, j, k, m, n, nresolved, nb, maxbond, ai, aj, best, im, nqual[2][ematchNR];
200     int     *bbb, *n_mask, *m_mask, **match;
201     char    *aname_i, *aname_m, *aname_n, *type;
202     double   qq, mm;
203     t_param *param;
204
205     snew(param, 1);
206     maxbond = 0;
207     for (i = 0; (i < atoms->nr); i++)
208     {
209         maxbond = std::max(maxbond, nbonds[i]);
210     }
211     if (debug)
212     {
213         fprintf(debug, "Max number of bonds per atom is %d\n", maxbond);
214     }
215     snew(bbb, maxbond);
216     snew(n_mask, maxbond);
217     snew(m_mask, maxbond);
218     snew(match, maxbond);
219     for (i = 0; (i < maxbond); i++)
220     {
221         snew(match[i], maxbond);
222     }
223
224     nresolved = 0;
225     for (i = 0; (i < atoms->nr); i++)
226     {
227         aname_i = *atoms->atomname[i];
228         nb      = 0;
229         for (j = 0; (j < bonds->nr); j++)
230         {
231             ai = bonds->param[j].ai();
232             aj = bonds->param[j].aj();
233             if (ai == i)
234             {
235                 bbb[nb++] = aj;
236             }
237             else if (aj == i)
238             {
239                 bbb[nb++] = ai;
240             }
241         }
242         if (nb != nbonds[i])
243         {
244             gmx_fatal(FARGS, "Counting number of bonds nb = %d, nbonds[%d] = %d",
245                       nb, i, nbonds[i]);
246         }
247         if (debug)
248         {
249             fprintf(debug, "%4s has bonds to", aname_i);
250             for (j = 0; (j < nb); j++)
251             {
252                 fprintf(debug, " %4s", *atoms->atomname[bbb[j]]);
253             }
254             fprintf(debug, "\n");
255         }
256         best = -1;
257         for (k = 0; (k < ematchNR); k++)
258         {
259             nqual[prev][k] = 0;
260         }
261
262         /* First check for names */
263         for (k = 0; (k < nnm); k++)
264         {
265             if (nm2t[k].nbonds == nb)
266             {
267                 im = match_str(*atoms->atomname[i], nm2t[k].elem);
268                 if (im > ematchWild)
269                 {
270                     for (j = 0; (j < ematchNR); j++)
271                     {
272                         nqual[cur][j] = 0;
273                     }
274
275                     /* Fill a matrix with matching quality */
276                     for (m = 0; (m < nb); m++)
277                     {
278                         aname_m = *atoms->atomname[bbb[m]];
279                         for (n = 0; (n < nb); n++)
280                         {
281                             aname_n     = nm2t[k].bond[n];
282                             match[m][n] = match_str(aname_m, aname_n);
283                         }
284                     }
285                     /* Now pick the best matches */
286                     for (m = 0; (m < nb); m++)
287                     {
288                         n_mask[m] = 0;
289                         m_mask[m] = 0;
290                     }
291                     for (j = ematchNR-1; (j > 0); j--)
292                     {
293                         for (m = 0; (m < nb); m++)
294                         {
295                             for (n = 0; (n < nb); n++)
296                             {
297                                 if ((n_mask[n] == 0) &&
298                                     (m_mask[m] == 0) &&
299                                     (match[m][n] == j))
300                                 {
301                                     n_mask[n] = 1;
302                                     m_mask[m] = 1;
303                                     nqual[cur][j]++;
304                                 }
305                             }
306                         }
307                     }
308                     if ((nqual[cur][ematchExact]+
309                          nqual[cur][ematchElem]+
310                          nqual[cur][ematchWild]) == nb)
311                     {
312                         if ((nqual[cur][ematchExact] > nqual[prev][ematchExact]) ||
313
314                             ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
315                              (nqual[cur][ematchElem] > nqual[prev][ematchElem])) ||
316
317                             ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
318                              (nqual[cur][ematchElem] == nqual[prev][ematchElem]) &&
319                              (nqual[cur][ematchWild] > nqual[prev][ematchWild])))
320                         {
321                             best = k;
322                             cur  = prev;
323                         }
324                     }
325                 }
326             }
327         }
328         if (best != -1)
329         {
330             int  atomnr = 0;
331             real alpha  = 0;
332
333             qq   = nm2t[best].q;
334             mm   = nm2t[best].m;
335             type = nm2t[best].type;
336
337             if ((k = get_atomtype_type(type, atype)) == NOTSET)
338             {
339                 atoms->atom[i].qB = alpha;
340                 atoms->atom[i].m  = atoms->atom[i].mB = mm;
341                 k                 = add_atomtype(atype, tab, &(atoms->atom[i]), type, param,
342                                                  atoms->atom[i].type, atomnr);
343             }
344             atoms->atom[i].type  = k;
345             atoms->atom[i].typeB = k;
346             atoms->atom[i].q     = qq;
347             atoms->atom[i].qB    = qq;
348             atoms->atom[i].m     = mm;
349             atoms->atom[i].mB    = mm;
350             nresolved++;
351         }
352         else
353         {
354             fprintf(stderr, "Can not find forcefield for atom %s-%d with %d bonds\n",
355                     *atoms->atomname[i], i+1, nb);
356         }
357     }
358     sfree(bbb);
359     sfree(n_mask);
360     sfree(m_mask);
361     sfree(param);
362
363     return nresolved;
364 }