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