7c8a9423711a236bee35b0de6ff81da4c1e159be
[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/grompp_impl.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     } while (bCont);
128     gmx_ffclose(fp);
129
130     *nnm = nnnm;
131     *nmp = nm2t;
132 }
133
134 t_nm2type* rd_nm2type(const char* ffdir, int* nnm)
135 {
136     std::vector<std::string> ff = fflib_search_file_end(ffdir, ".n2t", FALSE);
137     *nnm                        = 0;
138     t_nm2type* nm               = nullptr;
139     for (const auto& filename : ff)
140     {
141         rd_nm2type_file(filename, nnm, &nm);
142     }
143     return nm;
144 }
145
146 void dump_nm2type(FILE* fp, int nnm, t_nm2type nm2t[])
147 {
148     int i, j;
149
150     fprintf(fp, "; nm2type database\n");
151     for (i = 0; (i < nnm); i++)
152     {
153         fprintf(fp, "%-8s %-8s %8.4f %8.4f %-4d", nm2t[i].elem, nm2t[i].type, nm2t[i].q, nm2t[i].m,
154                 nm2t[i].nbonds);
155         for (j = 0; (j < nm2t[i].nbonds); j++)
156         {
157             fprintf(fp, " %-5s %6.4f", nm2t[i].bond[j], nm2t[i].blen[j]);
158         }
159         fprintf(fp, "\n");
160     }
161 }
162
163 enum
164 {
165     ematchNone,
166     ematchWild,
167     ematchElem,
168     ematchExact,
169     ematchNR
170 };
171
172 static int match_str(const char* atom, const char* template_string)
173 {
174     if (!atom || !template_string)
175     {
176         return ematchNone;
177     }
178     else if (gmx_strcasecmp(atom, template_string) == 0)
179     {
180         return ematchExact;
181     }
182     else if (atom[0] == template_string[0])
183     {
184         return ematchElem;
185     }
186     else if (strcmp(template_string, "*") == 0)
187     {
188         return ematchWild;
189     }
190     else
191     {
192         return ematchNone;
193     }
194 }
195
196 int nm2type(int                     nnm,
197             t_nm2type               nm2t[],
198             t_symtab*               tab,
199             t_atoms*                atoms,
200             PreprocessingAtomTypes* atype,
201             int*                    nbonds,
202             InteractionsOfType*     bonds)
203 {
204     int cur = 0;
205 #define prev (1 - cur)
206     int   nresolved, nb, maxbond, nqual[2][ematchNR];
207     int * bbb, *n_mask, *m_mask, **match;
208     char *aname_i, *aname_n;
209
210     maxbond = 0;
211     for (int i = 0; (i < atoms->nr); i++)
212     {
213         maxbond = std::max(maxbond, nbonds[i]);
214     }
215     if (debug)
216     {
217         fprintf(debug, "Max number of bonds per atom is %d\n", maxbond);
218     }
219     snew(bbb, maxbond);
220     snew(n_mask, maxbond);
221     snew(m_mask, maxbond);
222     snew(match, maxbond);
223     for (int i = 0; (i < maxbond); i++)
224     {
225         snew(match[i], maxbond);
226     }
227
228     nresolved = 0;
229     for (int i = 0; (i < atoms->nr); i++)
230     {
231         aname_i = *atoms->atomname[i];
232         nb      = 0;
233         for (const auto& bond : bonds->interactionTypes)
234         {
235             int ai = bond.ai();
236             int aj = bond.aj();
237             if (ai == i)
238             {
239                 bbb[nb++] = aj;
240             }
241             else if (aj == i)
242             {
243                 bbb[nb++] = ai;
244             }
245         }
246         if (nb != nbonds[i])
247         {
248             gmx_fatal(FARGS, "Counting number of bonds nb = %d, nbonds[%d] = %d", nb, i, nbonds[i]);
249         }
250         if (debug)
251         {
252             fprintf(debug, "%4s has bonds to", aname_i);
253             for (int j = 0; (j < nb); j++)
254             {
255                 fprintf(debug, " %4s", *atoms->atomname[bbb[j]]);
256             }
257             fprintf(debug, "\n");
258         }
259         int best = -1;
260         for (int k = 0; (k < ematchNR); k++)
261         {
262             nqual[prev][k] = 0;
263         }
264
265         /* First check for names */
266         for (int k = 0; (k < nnm); k++)
267         {
268             if (nm2t[k].nbonds == nb)
269             {
270                 int im = match_str(*atoms->atomname[i], nm2t[k].elem);
271                 if (im > ematchWild)
272                 {
273                     for (int j = 0; (j < ematchNR); j++)
274                     {
275                         nqual[cur][j] = 0;
276                     }
277
278                     /* Fill a matrix with matching quality */
279                     for (int m = 0; (m < nb); m++)
280                     {
281                         const char* aname_m = *atoms->atomname[bbb[m]];
282                         for (int n = 0; (n < nb); n++)
283                         {
284                             aname_n     = nm2t[k].bond[n];
285                             match[m][n] = match_str(aname_m, aname_n);
286                         }
287                     }
288                     /* Now pick the best matches */
289                     for (int m = 0; (m < nb); m++)
290                     {
291                         n_mask[m] = 0;
292                         m_mask[m] = 0;
293                     }
294                     for (int j = ematchNR - 1; (j > 0); j--)
295                     {
296                         for (int m = 0; (m < nb); m++)
297                         {
298                             for (int n = 0; (n < nb); n++)
299                             {
300                                 if ((n_mask[n] == 0) && (m_mask[m] == 0) && (match[m][n] == j))
301                                 {
302                                     n_mask[n] = 1;
303                                     m_mask[m] = 1;
304                                     nqual[cur][j]++;
305                                 }
306                             }
307                         }
308                     }
309                     if ((nqual[cur][ematchExact] + nqual[cur][ematchElem] + 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
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             double      qq   = nm2t[best].q;
334             double      mm   = nm2t[best].m;
335             const char* type = nm2t[best].type;
336
337             int k;
338             if ((k = atype->atomTypeFromName(type)) == NOTSET)
339             {
340                 atoms->atom[i].qB = alpha;
341                 atoms->atom[i].m = atoms->atom[i].mB = mm;
342                 k = atype->addType(tab, atoms->atom[i], type, InteractionOfType({}, {}),
343                                    atoms->atom[i].type, atomnr);
344             }
345             atoms->atom[i].type  = k;
346             atoms->atom[i].typeB = k;
347             atoms->atom[i].q     = qq;
348             atoms->atom[i].qB    = qq;
349             atoms->atom[i].m     = mm;
350             atoms->atom[i].mB    = mm;
351             nresolved++;
352         }
353         else
354         {
355             fprintf(stderr, "Can not find forcefield for atom %s-%d with %d bonds\n",
356                     *atoms->atomname[i], i + 1, nb);
357         }
358     }
359     sfree(bbb);
360     sfree(n_mask);
361     sfree(m_mask);
362
363     return nresolved;
364 }