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