Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / nm2type.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "config.h"
39
40 #include <string.h>
41
42 #include "gromacs/math/utilities.h"
43 #include "macros.h"
44 #include "bondf.h"
45 #include "gromacs/utility/cstringutil.h"
46 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/math/vec.h"
49 #include "txtdump.h"
50 #include "readinp.h"
51 #include "names.h"
52 #include "toppush.h"
53 #include "pdb2top.h"
54 #include "gpp_nextnb.h"
55 #include "gpp_atomtype.h"
56 #include "fflibutil.h"
57
58 #include "nm2type.h"
59
60 #include "gromacs/utility/fatalerror.h"
61
62 static void rd_nm2type_file(const char *fn, int *nnm, t_nm2type **nmp)
63 {
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, *blen;
71     t_nm2type    *nm2t = NULL;
72
73     fp = fflib_open(fn);
74     if (NULL == fp)
75     {
76         gmx_fatal(FARGS, "Can not find %s in library directory", fn);
77     }
78
79     nnnm = *nnm;
80     nm2t = *nmp;
81     do
82     {
83         /* Read a line from the file */
84         bCont = (fgets2(buf, 1023, fp) != NULL);
85
86         if (bCont)
87         {
88             /* Remove comment */
89             strip_comment(buf);
90             if (sscanf(buf, "%s%s%lf%lf%d", elem, type, &qq, &mm, &nb) == 5)
91             {
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)
96                 {
97                     snew(newbuf, nb);
98                     strcpy(format, "%*s%*s%*s%*s%*s");
99                     for (i = 0; (i < nb); i++)
100                     {
101                         /* Complicated format statement */
102                         strcpy(f1, format);
103                         strcat(f1, "%s%lf");
104                         if (sscanf(buf, f1, nbbuf, &(nm2t[nnnm].blen[i])) != 2)
105                         {
106                             gmx_fatal(FARGS, "Error on line %d of %s", line, libfilename);
107                         }
108                         newbuf[i] = strdup(nbbuf);
109                         strcat(format, "%*s%*s");
110                     }
111                 }
112                 else
113                 {
114                     newbuf = NULL;
115                 }
116                 nm2t[nnnm].elem   = strdup(elem);
117                 nm2t[nnnm].type   = strdup(type);
118                 nm2t[nnnm].q      = qq;
119                 nm2t[nnnm].m      = mm;
120                 nm2t[nnnm].nbonds = nb;
121                 nm2t[nnnm].bond   = newbuf;
122                 nnnm++;
123             }
124             line++;
125         }
126     }
127     while (bCont);
128     gmx_ffclose(fp);
129
130     *nnm = nnnm;
131     *nmp = nm2t;
132 }
133
134 t_nm2type *rd_nm2type(const char *ffdir, int *nnm)
135 {
136     int        nff, f;
137     char     **ff;
138     t_nm2type *nm;
139
140     nff  = fflib_search_file_end(ffdir, ".n2t", FALSE, &ff);
141     *nnm = 0;
142     nm   = NULL;
143     for (f = 0; f < nff; f++)
144     {
145         rd_nm2type_file(ff[f], nnm, &nm);
146         sfree(ff[f]);
147     }
148     sfree(ff);
149
150     return nm;
151 }
152
153 void dump_nm2type(FILE *fp, int nnm, t_nm2type nm2t[])
154 {
155     int i, j;
156
157     fprintf(fp, "; nm2type database\n");
158     for (i = 0; (i < nnm); i++)
159     {
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++)
164         {
165             fprintf(fp, " %-5s %6.4f", nm2t[i].bond[j], nm2t[i].blen[j]);
166         }
167         fprintf(fp, "\n");
168     }
169 }
170
171 enum {
172     ematchNone, ematchWild, ematchElem, ematchExact, ematchNR
173 };
174
175 static int match_str(const char *atom, const char *template_string)
176 {
177     if (!atom || !template_string)
178     {
179         return ematchNone;
180     }
181     else if (gmx_strcasecmp(atom, template_string) == 0)
182     {
183         return ematchExact;
184     }
185     else if (atom[0] == template_string[0])
186     {
187         return ematchElem;
188     }
189     else if (strcmp(template_string, "*") == 0)
190     {
191         return ematchWild;
192     }
193     else
194     {
195         return ematchNone;
196     }
197 }
198
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)
201 {
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, **quality;
206     char    *aname_i, *aname_m, *aname_n, *type;
207     double   qq, mm;
208     t_param *param;
209
210     snew(param, 1);
211     maxbond = 0;
212     for (i = 0; (i < atoms->nr); i++)
213     {
214         maxbond = max(maxbond, nbonds[i]);
215     }
216     if (debug)
217     {
218         fprintf(debug, "Max number of bonds per atom is %d\n", maxbond);
219     }
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++)
225     {
226         snew(match[i], maxbond);
227     }
228
229     nresolved = 0;
230     for (i = 0; (i < atoms->nr); i++)
231     {
232         aname_i = *atoms->atomname[i];
233         nb      = 0;
234         for (j = 0; (j < bonds->nr); j++)
235         {
236             ai = bonds->param[j].AI;
237             aj = bonds->param[j].AJ;
238             if (ai == i)
239             {
240                 bbb[nb++] = aj;
241             }
242             else if (aj == i)
243             {
244                 bbb[nb++] = ai;
245             }
246         }
247         if (nb != nbonds[i])
248         {
249             gmx_fatal(FARGS, "Counting number of bonds nb = %d, nbonds[%d] = %d",
250                       nb, i, nbonds[i]);
251         }
252         if (debug)
253         {
254             fprintf(debug, "%4s has bonds to", aname_i);
255             for (j = 0; (j < nb); j++)
256             {
257                 fprintf(debug, " %4s", *atoms->atomname[bbb[j]]);
258             }
259             fprintf(debug, "\n");
260         }
261         best = -1;
262         for (k = 0; (k < ematchNR); k++)
263         {
264             nqual[prev][k] = 0;
265         }
266
267         /* First check for names */
268         for (k = 0; (k < nnm); k++)
269         {
270             if (nm2t[k].nbonds == nb)
271             {
272                 im = match_str(*atoms->atomname[i], nm2t[k].elem);
273                 if (im > ematchWild)
274                 {
275                     for (j = 0; (j < ematchNR); j++)
276                     {
277                         nqual[cur][j] = 0;
278                     }
279
280                     /* Fill a matrix with matching quality */
281                     for (m = 0; (m < nb); m++)
282                     {
283                         aname_m = *atoms->atomname[bbb[m]];
284                         for (n = 0; (n < nb); n++)
285                         {
286                             aname_n     = nm2t[k].bond[n];
287                             match[m][n] = match_str(aname_m, aname_n);
288                         }
289                     }
290                     /* Now pick the best matches */
291                     for (m = 0; (m < nb); m++)
292                     {
293                         n_mask[m] = 0;
294                         m_mask[m] = 0;
295                     }
296                     for (j = ematchNR-1; (j > 0); j--)
297                     {
298                         for (m = 0; (m < nb); m++)
299                         {
300                             for (n = 0; (n < nb); n++)
301                             {
302                                 if ((n_mask[n] == 0) &&
303                                     (m_mask[m] == 0) &&
304                                     (match[m][n] == j))
305                                 {
306                                     n_mask[n] = 1;
307                                     m_mask[m] = 1;
308                                     nqual[cur][j]++;
309                                 }
310                             }
311                         }
312                     }
313                     if ((nqual[cur][ematchExact]+
314                          nqual[cur][ematchElem]+
315                          nqual[cur][ematchWild]) == nb)
316                     {
317                         if ((nqual[cur][ematchExact] > nqual[prev][ematchExact]) ||
318
319                             ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
320                              (nqual[cur][ematchElem] > nqual[prev][ematchElem])) ||
321
322                             ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
323                              (nqual[cur][ematchElem] == nqual[prev][ematchElem]) &&
324                              (nqual[cur][ematchWild] > nqual[prev][ematchWild])))
325                         {
326                             best = k;
327                             cur  = prev;
328                         }
329                     }
330                 }
331             }
332         }
333         if (best != -1)
334         {
335             int  atomnr = 0;
336             real alpha  = 0;
337
338             qq   = nm2t[best].q;
339             mm   = nm2t[best].m;
340             type = nm2t[best].type;
341
342             if ((k = get_atomtype_type(type, atype)) == NOTSET)
343             {
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);
348             }
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++;
356         }
357         else
358         {
359             fprintf(stderr, "Can not find forcefield for atom %s-%d with %d bonds\n",
360                     *atoms->atomname[i], i+1, nb);
361         }
362     }
363     sfree(bbb);
364     sfree(n_mask);
365     sfree(m_mask);
366     sfree(param);
367
368     return nresolved;
369 }