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