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, 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! */
42 #include "gromacs/math/utilities.h"
48 #include "gromacs/fileio/confio.h"
51 #include "gromacs/math/3dview.h"
57 #include "gpp_nextnb.h"
58 #include "gpp_atomtype.h"
59 #include "fflibutil.h"
63 static void rd_nm2type_file(const char *fn, int *nnm, t_nm2type **nmp)
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;
72 t_nm2type *nm2t = NULL;
77 gmx_fatal(FARGS, "Can not find %s in library directory", fn);
84 /* Read a line from the file */
85 bCont = (fgets2(buf, 1023, fp) != NULL);
91 if (sscanf(buf, "%s%s%lf%lf%d", elem, type, &qq, &mm, &nb) == 5)
93 /* If we can read the first four, there probably is more */
95 snew(nm2t[nnnm].blen, nb);
99 strcpy(format, "%*s%*s%*s%*s%*s");
100 for (i = 0; (i < nb); i++)
102 /* Complicated format statement */
105 if (sscanf(buf, f1, nbbuf, &(nm2t[nnnm].blen[i])) != 2)
107 gmx_fatal(FARGS, "Error on line %d of %s", line, libfilename);
109 newbuf[i] = strdup(nbbuf);
110 strcat(format, "%*s%*s");
117 nm2t[nnnm].elem = strdup(elem);
118 nm2t[nnnm].type = strdup(type);
121 nm2t[nnnm].nbonds = nb;
122 nm2t[nnnm].bond = newbuf;
135 t_nm2type *rd_nm2type(const char *ffdir, int *nnm)
141 nff = fflib_search_file_end(ffdir, ".n2t", FALSE, &ff);
144 for (f = 0; f < nff; f++)
146 rd_nm2type_file(ff[f], nnm, &nm);
154 void dump_nm2type(FILE *fp, int nnm, t_nm2type nm2t[])
158 fprintf(fp, "; nm2type database\n");
159 for (i = 0; (i < nnm); i++)
161 fprintf(fp, "%-8s %-8s %8.4f %8.4f %-4d",
162 nm2t[i].elem, nm2t[i].type,
163 nm2t[i].q, nm2t[i].m, nm2t[i].nbonds);
164 for (j = 0; (j < nm2t[i].nbonds); j++)
166 fprintf(fp, " %-5s %6.4f", nm2t[i].bond[j], nm2t[i].blen[j]);
173 ematchNone, ematchWild, ematchElem, ematchExact, ematchNR
176 static int match_str(const char *atom, const char *template_string)
178 if (!atom || !template_string)
182 else if (gmx_strcasecmp(atom, template_string) == 0)
186 else if (atom[0] == template_string[0])
190 else if (strcmp(template_string, "*") == 0)
200 int nm2type(int nnm, t_nm2type nm2t[], t_symtab *tab, t_atoms *atoms,
201 gpp_atomtype_t atype, int *nbonds, t_params *bonds)
205 int i, j, k, m, n, nresolved, nb, maxbond, ai, aj, best, im, nqual[2][ematchNR];
206 int *bbb, *n_mask, *m_mask, **match, **quality;
207 char *aname_i, *aname_m, *aname_n, *type;
213 for (i = 0; (i < atoms->nr); i++)
215 maxbond = max(maxbond, nbonds[i]);
219 fprintf(debug, "Max number of bonds per atom is %d\n", maxbond);
222 snew(n_mask, maxbond);
223 snew(m_mask, maxbond);
224 snew(match, maxbond);
225 for (i = 0; (i < maxbond); i++)
227 snew(match[i], maxbond);
231 for (i = 0; (i < atoms->nr); i++)
233 aname_i = *atoms->atomname[i];
235 for (j = 0; (j < bonds->nr); j++)
237 ai = bonds->param[j].AI;
238 aj = bonds->param[j].AJ;
250 gmx_fatal(FARGS, "Counting number of bonds nb = %d, nbonds[%d] = %d",
255 fprintf(debug, "%4s has bonds to", aname_i);
256 for (j = 0; (j < nb); j++)
258 fprintf(debug, " %4s", *atoms->atomname[bbb[j]]);
260 fprintf(debug, "\n");
263 for (k = 0; (k < ematchNR); k++)
268 /* First check for names */
269 for (k = 0; (k < nnm); k++)
271 if (nm2t[k].nbonds == nb)
273 im = match_str(*atoms->atomname[i], nm2t[k].elem);
276 for (j = 0; (j < ematchNR); j++)
281 /* Fill a matrix with matching quality */
282 for (m = 0; (m < nb); m++)
284 aname_m = *atoms->atomname[bbb[m]];
285 for (n = 0; (n < nb); n++)
287 aname_n = nm2t[k].bond[n];
288 match[m][n] = match_str(aname_m, aname_n);
291 /* Now pick the best matches */
292 for (m = 0; (m < nb); m++)
297 for (j = ematchNR-1; (j > 0); j--)
299 for (m = 0; (m < nb); m++)
301 for (n = 0; (n < nb); n++)
303 if ((n_mask[n] == 0) &&
314 if ((nqual[cur][ematchExact]+
315 nqual[cur][ematchElem]+
316 nqual[cur][ematchWild]) == nb)
318 if ((nqual[cur][ematchExact] > nqual[prev][ematchExact]) ||
320 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
321 (nqual[cur][ematchElem] > nqual[prev][ematchElem])) ||
323 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
324 (nqual[cur][ematchElem] == nqual[prev][ematchElem]) &&
325 (nqual[cur][ematchWild] > nqual[prev][ematchWild])))
341 type = nm2t[best].type;
343 if ((k = get_atomtype_type(type, atype)) == NOTSET)
345 atoms->atom[i].qB = alpha;
346 atoms->atom[i].m = atoms->atom[i].mB = mm;
347 k = add_atomtype(atype, tab, &(atoms->atom[i]), type, param,
348 atoms->atom[i].type, 0, 0, 0, atomnr, 0, 0);
350 atoms->atom[i].type = k;
351 atoms->atom[i].typeB = k;
352 atoms->atom[i].q = qq;
353 atoms->atom[i].qB = qq;
354 atoms->atom[i].m = mm;
355 atoms->atom[i].mB = mm;
360 fprintf(stderr, "Can not find forcefield for atom %s-%d with %d bonds\n",
361 *atoms->atomname[i], i+1, nb);