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 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.
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.
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.
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.
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.
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.
38 /* This file is completely threadsafe - keep it that way! */
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"
64 static void rd_nm2type_file(const std::string& filename, int* nnm, t_nm2type** nmp)
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;
73 t_nm2type* nm2t = nullptr;
75 fp = fflib_open(filename);
78 gmx_fatal(FARGS, "Can not find %s in library directory", filename.c_str());
85 /* Read a line from the file */
86 bCont = (fgets2(buf, 1023, fp) != nullptr);
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);
100 strcpy(format, "%*s%*s%*s%*s%*s");
101 for (i = 0; (i < nb); i++)
103 /* Complicated format statement */
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");
118 nm2t[nnnm].elem = gmx_strdup(elem);
119 nm2t[nnnm].type = gmx_strdup(type);
122 nm2t[nnnm].nbonds = nb;
123 nm2t[nnnm].bond = newbuf;
135 t_nm2type* rd_nm2type(const char* ffdir, int* nnm)
137 std::vector<std::string> ff = fflib_search_file_end(ffdir, ".n2t", FALSE);
139 t_nm2type* nm = nullptr;
140 for (const auto& filename : ff)
142 rd_nm2type_file(filename, nnm, &nm);
147 void dump_nm2type(FILE* fp, int nnm, t_nm2type nm2t[])
151 fprintf(fp, "; nm2type database\n");
152 for (i = 0; (i < nnm); i++)
154 fprintf(fp, "%-8s %-8s %8.4f %8.4f %-4d", nm2t[i].elem, nm2t[i].type, nm2t[i].q, nm2t[i].m,
156 for (j = 0; (j < nm2t[i].nbonds); j++)
158 fprintf(fp, " %-5s %6.4f", nm2t[i].bond[j], nm2t[i].blen[j]);
173 static int match_str(const char* atom, const char* template_string)
175 if (!atom || !template_string)
179 else if (gmx_strcasecmp(atom, template_string) == 0)
183 else if (atom[0] == template_string[0])
187 else if (strcmp(template_string, "*") == 0)
201 PreprocessingAtomTypes* atype,
203 InteractionsOfType* bonds)
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;
212 for (int i = 0; (i < atoms->nr); i++)
214 maxbond = std::max(maxbond, nbonds[i]);
218 fprintf(debug, "Max number of bonds per atom is %d\n", maxbond);
221 snew(n_mask, maxbond);
222 snew(m_mask, maxbond);
223 snew(match, maxbond);
224 for (int i = 0; (i < maxbond); i++)
226 snew(match[i], maxbond);
230 for (int i = 0; (i < atoms->nr); i++)
232 aname_i = *atoms->atomname[i];
234 for (const auto& bond : bonds->interactionTypes)
249 gmx_fatal(FARGS, "Counting number of bonds nb = %d, nbonds[%d] = %d", nb, i, nbonds[i]);
253 fprintf(debug, "%4s has bonds to", aname_i);
254 for (int j = 0; (j < nb); j++)
256 fprintf(debug, " %4s", *atoms->atomname[bbb[j]]);
258 fprintf(debug, "\n");
261 for (int k = 0; (k < ematchNR); k++)
266 /* First check for names */
267 for (int k = 0; (k < nnm); k++)
269 if (nm2t[k].nbonds == nb)
271 int im = match_str(*atoms->atomname[i], nm2t[k].elem);
274 for (int j = 0; (j < ematchNR); j++)
279 /* Fill a matrix with matching quality */
280 for (int m = 0; (m < nb); m++)
282 const char* aname_m = *atoms->atomname[bbb[m]];
283 for (int n = 0; (n < nb); n++)
285 aname_n = nm2t[k].bond[n];
286 match[m][n] = match_str(aname_m, aname_n);
289 /* Now pick the best matches */
290 for (int m = 0; (m < nb); m++)
295 for (int j = ematchNR - 1; (j > 0); j--)
297 for (int m = 0; (m < nb); m++)
299 for (int n = 0; (n < nb); n++)
301 if ((n_mask[n] == 0) && (m_mask[m] == 0) && (match[m][n] == j))
310 if ((nqual[cur][ematchExact] + nqual[cur][ematchElem] + nqual[cur][ematchWild]) == nb)
312 if ((nqual[cur][ematchExact] > nqual[prev][ematchExact]) ||
314 ((nqual[cur][ematchExact] == nqual[prev][ematchExact])
315 && (nqual[cur][ematchElem] > nqual[prev][ematchElem]))
318 ((nqual[cur][ematchExact] == nqual[prev][ematchExact])
319 && (nqual[cur][ematchElem] == nqual[prev][ematchElem])
320 && (nqual[cur][ematchWild] > nqual[prev][ematchWild])))
334 double qq = nm2t[best].q;
335 double mm = nm2t[best].m;
336 const char* type = nm2t[best].type;
339 if ((k = atype->atomTypeFromName(type)) == NOTSET)
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);
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;
356 fprintf(stderr, "Can not find forcefield for atom %s-%d with %d bonds\n",
357 *atoms->atomname[i], i + 1, nb);