Move physics.* to math/units.*
[gromacs.git] / src / gromacs / gmxpreprocess / nm2type.c
blob7f1f21a797cb20e593674b9d254a7ae321e0eb22
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, 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 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
42 #include <string.h>
44 #include "gromacs/math/utilities.h"
45 #include "macros.h"
46 #include "bondf.h"
47 #include "gromacs/utility/cstringutil.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/fileio/confio.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/math/3dview.h"
52 #include "txtdump.h"
53 #include "readinp.h"
54 #include "names.h"
55 #include "toppush.h"
56 #include "pdb2top.h"
57 #include "gpp_nextnb.h"
58 #include "gpp_atomtype.h"
59 #include "fflibutil.h"
61 #include "nm2type.h"
63 #include "gromacs/utility/fatalerror.h"
65 static void rd_nm2type_file(const char *fn, int *nnm, t_nm2type **nmp)
67 FILE *fp;
68 gmx_bool bCont;
69 char libfilename[128];
70 char format[128], f1[128];
71 char buf[1024], elem[16], type[16], nbbuf[16], **newbuf;
72 int i, nb, nnnm, line = 1;
73 double qq, mm, *blen;
74 t_nm2type *nm2t = NULL;
76 fp = fflib_open(fn);
77 if (NULL == fp)
79 gmx_fatal(FARGS, "Can not find %s in library directory", fn);
82 nnnm = *nnm;
83 nm2t = *nmp;
86 /* Read a line from the file */
87 bCont = (fgets2(buf, 1023, fp) != NULL);
89 if (bCont)
91 /* Remove comment */
92 strip_comment(buf);
93 if (sscanf(buf, "%s%s%lf%lf%d", elem, type, &qq, &mm, &nb) == 5)
95 /* If we can read the first four, there probably is more */
96 srenew(nm2t, nnnm+1);
97 snew(nm2t[nnnm].blen, nb);
98 if (nb > 0)
100 snew(newbuf, nb);
101 strcpy(format, "%*s%*s%*s%*s%*s");
102 for (i = 0; (i < nb); i++)
104 /* Complicated format statement */
105 strcpy(f1, format);
106 strcat(f1, "%s%lf");
107 if (sscanf(buf, f1, nbbuf, &(nm2t[nnnm].blen[i])) != 2)
109 gmx_fatal(FARGS, "Error on line %d of %s", line, libfilename);
111 newbuf[i] = strdup(nbbuf);
112 strcat(format, "%*s%*s");
115 else
117 newbuf = NULL;
119 nm2t[nnnm].elem = strdup(elem);
120 nm2t[nnnm].type = strdup(type);
121 nm2t[nnnm].q = qq;
122 nm2t[nnnm].m = mm;
123 nm2t[nnnm].nbonds = nb;
124 nm2t[nnnm].bond = newbuf;
125 nnnm++;
127 line++;
130 while (bCont);
131 gmx_ffclose(fp);
133 *nnm = nnnm;
134 *nmp = nm2t;
137 t_nm2type *rd_nm2type(const char *ffdir, int *nnm)
139 int nff, f;
140 char **ff;
141 t_nm2type *nm;
143 nff = fflib_search_file_end(ffdir, ".n2t", FALSE, &ff);
144 *nnm = 0;
145 nm = NULL;
146 for (f = 0; f < nff; f++)
148 rd_nm2type_file(ff[f], nnm, &nm);
149 sfree(ff[f]);
151 sfree(ff);
153 return nm;
156 void dump_nm2type(FILE *fp, int nnm, t_nm2type nm2t[])
158 int i, j;
160 fprintf(fp, "; nm2type database\n");
161 for (i = 0; (i < nnm); i++)
163 fprintf(fp, "%-8s %-8s %8.4f %8.4f %-4d",
164 nm2t[i].elem, nm2t[i].type,
165 nm2t[i].q, nm2t[i].m, nm2t[i].nbonds);
166 for (j = 0; (j < nm2t[i].nbonds); j++)
168 fprintf(fp, " %-5s %6.4f", nm2t[i].bond[j], nm2t[i].blen[j]);
170 fprintf(fp, "\n");
174 enum {
175 ematchNone, ematchWild, ematchElem, ematchExact, ematchNR
178 static int match_str(const char *atom, const char *template_string)
180 if (!atom || !template_string)
182 return ematchNone;
184 else if (gmx_strcasecmp(atom, template_string) == 0)
186 return ematchExact;
188 else if (atom[0] == template_string[0])
190 return ematchElem;
192 else if (strcmp(template_string, "*") == 0)
194 return ematchWild;
196 else
198 return ematchNone;
202 int nm2type(int nnm, t_nm2type nm2t[], t_symtab *tab, t_atoms *atoms,
203 gpp_atomtype_t atype, int *nbonds, t_params *bonds)
205 int cur = 0;
206 #define prev (1-cur)
207 int i, j, k, m, n, nresolved, nb, maxbond, ai, aj, best, im, nqual[2][ematchNR];
208 int *bbb, *n_mask, *m_mask, **match, **quality;
209 char *aname_i, *aname_m, *aname_n, *type;
210 double qq, mm;
211 t_param *param;
213 snew(param, 1);
214 maxbond = 0;
215 for (i = 0; (i < atoms->nr); i++)
217 maxbond = max(maxbond, nbonds[i]);
219 if (debug)
221 fprintf(debug, "Max number of bonds per atom is %d\n", maxbond);
223 snew(bbb, maxbond);
224 snew(n_mask, maxbond);
225 snew(m_mask, maxbond);
226 snew(match, maxbond);
227 for (i = 0; (i < maxbond); i++)
229 snew(match[i], maxbond);
232 nresolved = 0;
233 for (i = 0; (i < atoms->nr); i++)
235 aname_i = *atoms->atomname[i];
236 nb = 0;
237 for (j = 0; (j < bonds->nr); j++)
239 ai = bonds->param[j].AI;
240 aj = bonds->param[j].AJ;
241 if (ai == i)
243 bbb[nb++] = aj;
245 else if (aj == i)
247 bbb[nb++] = ai;
250 if (nb != nbonds[i])
252 gmx_fatal(FARGS, "Counting number of bonds nb = %d, nbonds[%d] = %d",
253 nb, i, nbonds[i]);
255 if (debug)
257 fprintf(debug, "%4s has bonds to", aname_i);
258 for (j = 0; (j < nb); j++)
260 fprintf(debug, " %4s", *atoms->atomname[bbb[j]]);
262 fprintf(debug, "\n");
264 best = -1;
265 for (k = 0; (k < ematchNR); k++)
267 nqual[prev][k] = 0;
270 /* First check for names */
271 for (k = 0; (k < nnm); k++)
273 if (nm2t[k].nbonds == nb)
275 im = match_str(*atoms->atomname[i], nm2t[k].elem);
276 if (im > ematchWild)
278 for (j = 0; (j < ematchNR); j++)
280 nqual[cur][j] = 0;
283 /* Fill a matrix with matching quality */
284 for (m = 0; (m < nb); m++)
286 aname_m = *atoms->atomname[bbb[m]];
287 for (n = 0; (n < nb); n++)
289 aname_n = nm2t[k].bond[n];
290 match[m][n] = match_str(aname_m, aname_n);
293 /* Now pick the best matches */
294 for (m = 0; (m < nb); m++)
296 n_mask[m] = 0;
297 m_mask[m] = 0;
299 for (j = ematchNR-1; (j > 0); j--)
301 for (m = 0; (m < nb); m++)
303 for (n = 0; (n < nb); n++)
305 if ((n_mask[n] == 0) &&
306 (m_mask[m] == 0) &&
307 (match[m][n] == j))
309 n_mask[n] = 1;
310 m_mask[m] = 1;
311 nqual[cur][j]++;
316 if ((nqual[cur][ematchExact]+
317 nqual[cur][ematchElem]+
318 nqual[cur][ematchWild]) == nb)
320 if ((nqual[cur][ematchExact] > nqual[prev][ematchExact]) ||
322 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
323 (nqual[cur][ematchElem] > nqual[prev][ematchElem])) ||
325 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
326 (nqual[cur][ematchElem] == nqual[prev][ematchElem]) &&
327 (nqual[cur][ematchWild] > nqual[prev][ematchWild])))
329 best = k;
330 cur = prev;
336 if (best != -1)
338 int atomnr = 0;
339 real alpha = 0;
341 qq = nm2t[best].q;
342 mm = nm2t[best].m;
343 type = nm2t[best].type;
345 if ((k = get_atomtype_type(type, atype)) == NOTSET)
347 atoms->atom[i].qB = alpha;
348 atoms->atom[i].m = atoms->atom[i].mB = mm;
349 k = add_atomtype(atype, tab, &(atoms->atom[i]), type, param,
350 atoms->atom[i].type, 0, 0, 0, atomnr, 0, 0);
352 atoms->atom[i].type = k;
353 atoms->atom[i].typeB = k;
354 atoms->atom[i].q = qq;
355 atoms->atom[i].qB = qq;
356 atoms->atom[i].m = mm;
357 atoms->atom[i].mB = mm;
358 nresolved++;
360 else
362 fprintf(stderr, "Can not find forcefield for atom %s-%d with %d bonds\n",
363 *atoms->atomname[i], i+1, nb);
366 sfree(bbb);
367 sfree(n_mask);
368 sfree(m_mask);
369 sfree(param);
371 return nresolved;