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