Remove unused function generate_excls and make clean_excls static
[gromacs.git] / src / gromacs / gmxpreprocess / nm2type.cpp
blob437030c57b655622bb325f0eaed88210fa1a57b6
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
40 #include "nm2type.h"
42 #include <cstring>
44 #include <algorithm>
45 #include <string>
46 #include <vector>
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"
64 static void rd_nm2type_file(const std::string &filename, int *nnm, t_nm2type **nmp)
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;
75 fp = fflib_open(filename);
76 if (nullptr == fp)
78 gmx_fatal(FARGS, "Can not find %s in library directory", filename.c_str());
81 nnnm = *nnm;
82 nm2t = *nmp;
85 /* Read a line from the file */
86 bCont = (fgets2(buf, 1023, fp) != nullptr);
88 if (bCont)
90 /* Remove comment */
91 strip_comment(buf);
92 if (sscanf(buf, "%s%s%lf%lf%d", elem, type, &qq, &mm, &nb) == 5)
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)
99 snew(newbuf, nb);
100 strcpy(format, "%*s%*s%*s%*s%*s");
101 for (i = 0; (i < nb); i++)
103 /* Complicated format statement */
104 strcpy(f1, format);
105 strcat(f1, "%s%lf");
106 if (sscanf(buf, f1, nbbuf, &(nm2t[nnnm].blen[i])) != 2)
108 gmx_fatal(FARGS, "Error on line %d of %s", line, libfilename);
110 newbuf[i] = gmx_strdup(nbbuf);
111 strcat(format, "%*s%*s");
114 else
116 newbuf = nullptr;
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++;
126 line++;
129 while (bCont);
130 gmx_ffclose(fp);
132 *nnm = nnnm;
133 *nmp = nm2t;
136 t_nm2type *rd_nm2type(const char *ffdir, int *nnm)
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)
143 rd_nm2type_file(filename, nnm, &nm);
145 return nm;
148 void dump_nm2type(FILE *fp, int nnm, t_nm2type nm2t[])
150 int i, j;
152 fprintf(fp, "; nm2type database\n");
153 for (i = 0; (i < nnm); i++)
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++)
160 fprintf(fp, " %-5s %6.4f", nm2t[i].bond[j], nm2t[i].blen[j]);
162 fprintf(fp, "\n");
166 enum {
167 ematchNone, ematchWild, ematchElem, ematchExact, ematchNR
170 static int match_str(const char *atom, const char *template_string)
172 if (!atom || !template_string)
174 return ematchNone;
176 else if (gmx_strcasecmp(atom, template_string) == 0)
178 return ematchExact;
180 else if (atom[0] == template_string[0])
182 return ematchElem;
184 else if (strcmp(template_string, "*") == 0)
186 return ematchWild;
188 else
190 return ematchNone;
194 int nm2type(int nnm, t_nm2type nm2t[], t_symtab *tab, t_atoms *atoms,
195 PreprocessingAtomTypes *atype, int *nbonds, InteractionsOfType *bonds)
197 int cur = 0;
198 #define prev (1-cur)
199 int nresolved, nb, maxbond, nqual[2][ematchNR];
200 int *bbb, *n_mask, *m_mask, **match;
201 char *aname_i, *aname_n;
203 maxbond = 0;
204 for (int i = 0; (i < atoms->nr); i++)
206 maxbond = std::max(maxbond, nbonds[i]);
208 if (debug)
210 fprintf(debug, "Max number of bonds per atom is %d\n", maxbond);
212 snew(bbb, maxbond);
213 snew(n_mask, maxbond);
214 snew(m_mask, maxbond);
215 snew(match, maxbond);
216 for (int i = 0; (i < maxbond); i++)
218 snew(match[i], maxbond);
221 nresolved = 0;
222 for (int i = 0; (i < atoms->nr); i++)
224 aname_i = *atoms->atomname[i];
225 nb = 0;
226 for (const auto &bond : bonds->interactionTypes)
228 int ai = bond.ai();
229 int aj = bond.aj();
230 if (ai == i)
232 bbb[nb++] = aj;
234 else if (aj == i)
236 bbb[nb++] = ai;
239 if (nb != nbonds[i])
241 gmx_fatal(FARGS, "Counting number of bonds nb = %d, nbonds[%d] = %d",
242 nb, i, nbonds[i]);
244 if (debug)
246 fprintf(debug, "%4s has bonds to", aname_i);
247 for (int j = 0; (j < nb); j++)
249 fprintf(debug, " %4s", *atoms->atomname[bbb[j]]);
251 fprintf(debug, "\n");
253 int best = -1;
254 for (int k = 0; (k < ematchNR); k++)
256 nqual[prev][k] = 0;
259 /* First check for names */
260 for (int k = 0; (k < nnm); k++)
262 if (nm2t[k].nbonds == nb)
264 int im = match_str(*atoms->atomname[i], nm2t[k].elem);
265 if (im > ematchWild)
267 for (int j = 0; (j < ematchNR); j++)
269 nqual[cur][j] = 0;
272 /* Fill a matrix with matching quality */
273 for (int m = 0; (m < nb); m++)
275 const char *aname_m = *atoms->atomname[bbb[m]];
276 for (int n = 0; (n < nb); n++)
278 aname_n = nm2t[k].bond[n];
279 match[m][n] = match_str(aname_m, aname_n);
282 /* Now pick the best matches */
283 for (int m = 0; (m < nb); m++)
285 n_mask[m] = 0;
286 m_mask[m] = 0;
288 for (int j = ematchNR-1; (j > 0); j--)
290 for (int m = 0; (m < nb); m++)
292 for (int n = 0; (n < nb); n++)
294 if ((n_mask[n] == 0) &&
295 (m_mask[m] == 0) &&
296 (match[m][n] == j))
298 n_mask[n] = 1;
299 m_mask[m] = 1;
300 nqual[cur][j]++;
305 if ((nqual[cur][ematchExact]+
306 nqual[cur][ematchElem]+
307 nqual[cur][ematchWild]) == nb)
309 if ((nqual[cur][ematchExact] > nqual[prev][ematchExact]) ||
311 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
312 (nqual[cur][ematchElem] > nqual[prev][ematchElem])) ||
314 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
315 (nqual[cur][ematchElem] == nqual[prev][ematchElem]) &&
316 (nqual[cur][ematchWild] > nqual[prev][ematchWild])))
318 best = k;
319 cur = prev;
325 if (best != -1)
327 int atomnr = 0;
328 real alpha = 0;
330 double qq = nm2t[best].q;
331 double mm = nm2t[best].m;
332 const char *type = nm2t[best].type;
334 int k;
335 if ((k = atype->atomTypeFromName(type)) == NOTSET)
337 atoms->atom[i].qB = alpha;
338 atoms->atom[i].m = atoms->atom[i].mB = mm;
339 k = atype->addType(tab, atoms->atom[i], type, InteractionOfType({}, {}),
340 atoms->atom[i].type, atomnr);
342 atoms->atom[i].type = k;
343 atoms->atom[i].typeB = k;
344 atoms->atom[i].q = qq;
345 atoms->atom[i].qB = qq;
346 atoms->atom[i].m = mm;
347 atoms->atom[i].mB = mm;
348 nresolved++;
350 else
352 fprintf(stderr, "Can not find forcefield for atom %s-%d with %d bonds\n",
353 *atoms->atomname[i], i+1, nb);
356 sfree(bbb);
357 sfree(n_mask);
358 sfree(m_mask);
360 return nresolved;