Split txtdump.*, part 1
[gromacs.git] / src / gromacs / gmxpreprocess / nm2type.cpp
blob3f14d74c55d3a6c434404d2724c47ce98c373581
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, 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 <string.h>
44 #include <algorithm>
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/gmxlib/readinp.h"
48 #include "gromacs/gmxpreprocess/fflibutil.h"
49 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
50 #include "gromacs/gmxpreprocess/gpp_nextnb.h"
51 #include "gromacs/gmxpreprocess/notset.h"
52 #include "gromacs/gmxpreprocess/pdb2top.h"
53 #include "gromacs/gmxpreprocess/toppush.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
62 static void rd_nm2type_file(const char *fn, int *nnm, t_nm2type **nmp)
64 FILE *fp;
65 gmx_bool bCont;
66 char libfilename[128];
67 char format[128], f1[128];
68 char buf[1024], elem[16], type[16], nbbuf[16], **newbuf;
69 int i, nb, nnnm, line = 1;
70 double qq, mm;
71 t_nm2type *nm2t = NULL;
73 fp = fflib_open(fn);
74 if (NULL == fp)
76 gmx_fatal(FARGS, "Can not find %s in library directory", fn);
79 nnnm = *nnm;
80 nm2t = *nmp;
83 /* Read a line from the file */
84 bCont = (fgets2(buf, 1023, fp) != NULL);
86 if (bCont)
88 /* Remove comment */
89 strip_comment(buf);
90 if (sscanf(buf, "%s%s%lf%lf%d", elem, type, &qq, &mm, &nb) == 5)
92 /* If we can read the first four, there probably is more */
93 srenew(nm2t, nnnm+1);
94 snew(nm2t[nnnm].blen, nb);
95 if (nb > 0)
97 snew(newbuf, nb);
98 strcpy(format, "%*s%*s%*s%*s%*s");
99 for (i = 0; (i < nb); i++)
101 /* Complicated format statement */
102 strcpy(f1, format);
103 strcat(f1, "%s%lf");
104 if (sscanf(buf, f1, nbbuf, &(nm2t[nnnm].blen[i])) != 2)
106 gmx_fatal(FARGS, "Error on line %d of %s", line, libfilename);
108 newbuf[i] = gmx_strdup(nbbuf);
109 strcat(format, "%*s%*s");
112 else
114 newbuf = NULL;
116 nm2t[nnnm].elem = gmx_strdup(elem);
117 nm2t[nnnm].type = gmx_strdup(type);
118 nm2t[nnnm].q = qq;
119 nm2t[nnnm].m = mm;
120 nm2t[nnnm].nbonds = nb;
121 nm2t[nnnm].bond = newbuf;
122 nnnm++;
124 line++;
127 while (bCont);
128 gmx_ffclose(fp);
130 *nnm = nnnm;
131 *nmp = nm2t;
134 t_nm2type *rd_nm2type(const char *ffdir, int *nnm)
136 int nff, f;
137 char **ff;
138 t_nm2type *nm;
140 nff = fflib_search_file_end(ffdir, ".n2t", FALSE, &ff);
141 *nnm = 0;
142 nm = NULL;
143 for (f = 0; f < nff; f++)
145 rd_nm2type_file(ff[f], nnm, &nm);
146 sfree(ff[f]);
148 sfree(ff);
150 return nm;
153 void dump_nm2type(FILE *fp, int nnm, t_nm2type nm2t[])
155 int i, j;
157 fprintf(fp, "; nm2type database\n");
158 for (i = 0; (i < nnm); i++)
160 fprintf(fp, "%-8s %-8s %8.4f %8.4f %-4d",
161 nm2t[i].elem, nm2t[i].type,
162 nm2t[i].q, nm2t[i].m, nm2t[i].nbonds);
163 for (j = 0; (j < nm2t[i].nbonds); j++)
165 fprintf(fp, " %-5s %6.4f", nm2t[i].bond[j], nm2t[i].blen[j]);
167 fprintf(fp, "\n");
171 enum {
172 ematchNone, ematchWild, ematchElem, ematchExact, ematchNR
175 static int match_str(const char *atom, const char *template_string)
177 if (!atom || !template_string)
179 return ematchNone;
181 else if (gmx_strcasecmp(atom, template_string) == 0)
183 return ematchExact;
185 else if (atom[0] == template_string[0])
187 return ematchElem;
189 else if (strcmp(template_string, "*") == 0)
191 return ematchWild;
193 else
195 return ematchNone;
199 int nm2type(int nnm, t_nm2type nm2t[], struct t_symtab *tab, t_atoms *atoms,
200 gpp_atomtype_t atype, int *nbonds, t_params *bonds)
202 int cur = 0;
203 #define prev (1-cur)
204 int i, j, k, m, n, nresolved, nb, maxbond, ai, aj, best, im, nqual[2][ematchNR];
205 int *bbb, *n_mask, *m_mask, **match;
206 char *aname_i, *aname_m, *aname_n, *type;
207 double qq, mm;
208 t_param *param;
210 snew(param, 1);
211 maxbond = 0;
212 for (i = 0; (i < atoms->nr); i++)
214 maxbond = std::max(maxbond, nbonds[i]);
216 if (debug)
218 fprintf(debug, "Max number of bonds per atom is %d\n", maxbond);
220 snew(bbb, maxbond);
221 snew(n_mask, maxbond);
222 snew(m_mask, maxbond);
223 snew(match, maxbond);
224 for (i = 0; (i < maxbond); i++)
226 snew(match[i], maxbond);
229 nresolved = 0;
230 for (i = 0; (i < atoms->nr); i++)
232 aname_i = *atoms->atomname[i];
233 nb = 0;
234 for (j = 0; (j < bonds->nr); j++)
236 ai = bonds->param[j].ai();
237 aj = bonds->param[j].aj();
238 if (ai == i)
240 bbb[nb++] = aj;
242 else if (aj == i)
244 bbb[nb++] = ai;
247 if (nb != nbonds[i])
249 gmx_fatal(FARGS, "Counting number of bonds nb = %d, nbonds[%d] = %d",
250 nb, i, nbonds[i]);
252 if (debug)
254 fprintf(debug, "%4s has bonds to", aname_i);
255 for (j = 0; (j < nb); j++)
257 fprintf(debug, " %4s", *atoms->atomname[bbb[j]]);
259 fprintf(debug, "\n");
261 best = -1;
262 for (k = 0; (k < ematchNR); k++)
264 nqual[prev][k] = 0;
267 /* First check for names */
268 for (k = 0; (k < nnm); k++)
270 if (nm2t[k].nbonds == nb)
272 im = match_str(*atoms->atomname[i], nm2t[k].elem);
273 if (im > ematchWild)
275 for (j = 0; (j < ematchNR); j++)
277 nqual[cur][j] = 0;
280 /* Fill a matrix with matching quality */
281 for (m = 0; (m < nb); m++)
283 aname_m = *atoms->atomname[bbb[m]];
284 for (n = 0; (n < nb); n++)
286 aname_n = nm2t[k].bond[n];
287 match[m][n] = match_str(aname_m, aname_n);
290 /* Now pick the best matches */
291 for (m = 0; (m < nb); m++)
293 n_mask[m] = 0;
294 m_mask[m] = 0;
296 for (j = ematchNR-1; (j > 0); j--)
298 for (m = 0; (m < nb); m++)
300 for (n = 0; (n < nb); n++)
302 if ((n_mask[n] == 0) &&
303 (m_mask[m] == 0) &&
304 (match[m][n] == j))
306 n_mask[n] = 1;
307 m_mask[m] = 1;
308 nqual[cur][j]++;
313 if ((nqual[cur][ematchExact]+
314 nqual[cur][ematchElem]+
315 nqual[cur][ematchWild]) == nb)
317 if ((nqual[cur][ematchExact] > nqual[prev][ematchExact]) ||
319 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
320 (nqual[cur][ematchElem] > nqual[prev][ematchElem])) ||
322 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
323 (nqual[cur][ematchElem] == nqual[prev][ematchElem]) &&
324 (nqual[cur][ematchWild] > nqual[prev][ematchWild])))
326 best = k;
327 cur = prev;
333 if (best != -1)
335 int atomnr = 0;
336 real alpha = 0;
338 qq = nm2t[best].q;
339 mm = nm2t[best].m;
340 type = nm2t[best].type;
342 if ((k = get_atomtype_type(type, atype)) == NOTSET)
344 atoms->atom[i].qB = alpha;
345 atoms->atom[i].m = atoms->atom[i].mB = mm;
346 k = add_atomtype(atype, tab, &(atoms->atom[i]), type, param,
347 atoms->atom[i].type, 0, 0, 0, atomnr, 0, 0);
349 atoms->atom[i].type = k;
350 atoms->atom[i].typeB = k;
351 atoms->atom[i].q = qq;
352 atoms->atom[i].qB = qq;
353 atoms->atom[i].m = mm;
354 atoms->atom[i].mB = mm;
355 nresolved++;
357 else
359 fprintf(stderr, "Can not find forcefield for atom %s-%d with %d bonds\n",
360 *atoms->atomname[i], i+1, nb);
363 sfree(bbb);
364 sfree(n_mask);
365 sfree(m_mask);
366 sfree(param);
368 return nresolved;