Split lines with many copyright years
[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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
40
41 #include "nm2type.h"
42
43 #include <cstring>
44
45 #include <algorithm>
46 #include <string>
47 #include <vector>
48
49 #include "gromacs/fileio/confio.h"
50 #include "gromacs/gmxpreprocess/fflibutil.h"
51 #include "gromacs/gmxpreprocess/gpp_atomtype.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     } 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", nm2t[i].elem, nm2t[i].type, nm2t[i].q, nm2t[i].m,
155                 nm2t[i].nbonds);
156         for (j = 0; (j < nm2t[i].nbonds); j++)
157         {
158             fprintf(fp, " %-5s %6.4f", nm2t[i].bond[j], nm2t[i].blen[j]);
159         }
160         fprintf(fp, "\n");
161     }
162 }
163
164 enum
165 {
166     ematchNone,
167     ematchWild,
168     ematchElem,
169     ematchExact,
170     ematchNR
171 };
172
173 static int match_str(const char* atom, const char* template_string)
174 {
175     if (!atom || !template_string)
176     {
177         return ematchNone;
178     }
179     else if (gmx_strcasecmp(atom, template_string) == 0)
180     {
181         return ematchExact;
182     }
183     else if (atom[0] == template_string[0])
184     {
185         return ematchElem;
186     }
187     else if (strcmp(template_string, "*") == 0)
188     {
189         return ematchWild;
190     }
191     else
192     {
193         return ematchNone;
194     }
195 }
196
197 int nm2type(int                     nnm,
198             t_nm2type               nm2t[],
199             t_symtab*               tab,
200             t_atoms*                atoms,
201             PreprocessingAtomTypes* atype,
202             int*                    nbonds,
203             InteractionsOfType*     bonds)
204 {
205     int cur = 0;
206 #define prev (1 - cur)
207     int   nresolved, nb, maxbond, nqual[2][ematchNR];
208     int * bbb, *n_mask, *m_mask, **match;
209     char *aname_i, *aname_n;
210
211     maxbond = 0;
212     for (int i = 0; (i < atoms->nr); i++)
213     {
214         maxbond = std::max(maxbond, nbonds[i]);
215     }
216     if (debug)
217     {
218         fprintf(debug, "Max number of bonds per atom is %d\n", maxbond);
219     }
220     snew(bbb, maxbond);
221     snew(n_mask, maxbond);
222     snew(m_mask, maxbond);
223     snew(match, maxbond);
224     for (int i = 0; (i < maxbond); i++)
225     {
226         snew(match[i], maxbond);
227     }
228
229     nresolved = 0;
230     for (int i = 0; (i < atoms->nr); i++)
231     {
232         aname_i = *atoms->atomname[i];
233         nb      = 0;
234         for (const auto& bond : bonds->interactionTypes)
235         {
236             int ai = bond.ai();
237             int aj = bond.aj();
238             if (ai == i)
239             {
240                 bbb[nb++] = aj;
241             }
242             else if (aj == i)
243             {
244                 bbb[nb++] = ai;
245             }
246         }
247         if (nb != nbonds[i])
248         {
249             gmx_fatal(FARGS, "Counting number of bonds nb = %d, nbonds[%d] = %d", nb, i, nbonds[i]);
250         }
251         if (debug)
252         {
253             fprintf(debug, "%4s has bonds to", aname_i);
254             for (int j = 0; (j < nb); j++)
255             {
256                 fprintf(debug, " %4s", *atoms->atomname[bbb[j]]);
257             }
258             fprintf(debug, "\n");
259         }
260         int best = -1;
261         for (int k = 0; (k < ematchNR); k++)
262         {
263             nqual[prev][k] = 0;
264         }
265
266         /* First check for names */
267         for (int k = 0; (k < nnm); k++)
268         {
269             if (nm2t[k].nbonds == nb)
270             {
271                 int im = match_str(*atoms->atomname[i], nm2t[k].elem);
272                 if (im > ematchWild)
273                 {
274                     for (int j = 0; (j < ematchNR); j++)
275                     {
276                         nqual[cur][j] = 0;
277                     }
278
279                     /* Fill a matrix with matching quality */
280                     for (int m = 0; (m < nb); m++)
281                     {
282                         const char* aname_m = *atoms->atomname[bbb[m]];
283                         for (int n = 0; (n < nb); n++)
284                         {
285                             aname_n     = nm2t[k].bond[n];
286                             match[m][n] = match_str(aname_m, aname_n);
287                         }
288                     }
289                     /* Now pick the best matches */
290                     for (int m = 0; (m < nb); m++)
291                     {
292                         n_mask[m] = 0;
293                         m_mask[m] = 0;
294                     }
295                     for (int j = ematchNR - 1; (j > 0); j--)
296                     {
297                         for (int m = 0; (m < nb); m++)
298                         {
299                             for (int n = 0; (n < nb); n++)
300                             {
301                                 if ((n_mask[n] == 0) && (m_mask[m] == 0) && (match[m][n] == j))
302                                 {
303                                     n_mask[n] = 1;
304                                     m_mask[m] = 1;
305                                     nqual[cur][j]++;
306                                 }
307                             }
308                         }
309                     }
310                     if ((nqual[cur][ematchExact] + nqual[cur][ematchElem] + 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
318                             ((nqual[cur][ematchExact] == nqual[prev][ematchExact])
319                              && (nqual[cur][ematchElem] == nqual[prev][ematchElem])
320                              && (nqual[cur][ematchWild] > nqual[prev][ematchWild])))
321                         {
322                             best = k;
323                             cur  = prev;
324                         }
325                     }
326                 }
327             }
328         }
329         if (best != -1)
330         {
331             int  atomnr = 0;
332             real alpha  = 0;
333
334             double      qq   = nm2t[best].q;
335             double      mm   = nm2t[best].m;
336             const char* type = nm2t[best].type;
337
338             int k;
339             if ((k = atype->atomTypeFromName(type)) == NOTSET)
340             {
341                 atoms->atom[i].qB = alpha;
342                 atoms->atom[i].m = atoms->atom[i].mB = mm;
343                 k = atype->addType(tab, atoms->atom[i], type, InteractionOfType({}, {}),
344                                    atoms->atom[i].type, atomnr);
345             }
346             atoms->atom[i].type  = k;
347             atoms->atom[i].typeB = k;
348             atoms->atom[i].q     = qq;
349             atoms->atom[i].qB    = qq;
350             atoms->atom[i].m     = mm;
351             atoms->atom[i].mB    = mm;
352             nresolved++;
353         }
354         else
355         {
356             fprintf(stderr, "Can not find forcefield for atom %s-%d with %d bonds\n",
357                     *atoms->atomname[i], i + 1, nb);
358         }
359     }
360     sfree(bbb);
361     sfree(n_mask);
362     sfree(m_mask);
363
364     return nresolved;
365 }