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