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