remove most debug flags from gmxpreprocess
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / gpp_atomtype.cpp
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-2004, The GROMACS development team.
6  * Copyright (c) 2011,2014,2015,2017,2018, 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 #include "gmxpre.h"
38
39 #include "gpp_atomtype.h"
40
41 #include <string.h>
42
43 #include <climits>
44 #include <cmath>
45
46 #include "gromacs/gmxpreprocess/notset.h"
47 #include "gromacs/gmxpreprocess/topdirs.h"
48 #include "gromacs/gmxpreprocess/toputil.h"
49 #include "gromacs/math/vecdump.h"
50 #include "gromacs/topology/ifunc.h"
51 #include "gromacs/topology/symtab.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/smalloc.h"
55
56 typedef struct gpp_atomtype {
57     int              nr;           /* The number of atomtypes           */
58     t_atom          *atom;         /* Array of atoms                    */
59     char          ***atomname;     /* Names of the atomtypes            */
60     t_param         *nb;           /* Nonbonded force default params    */
61     int             *bondatomtype; /* The bond_atomtype for each atomtype  */
62     int             *atomnumber;   /* Atomic number, used for QM/MM        */
63 } t_gpp_atomtype;
64
65 int get_atomtype_type(const char *str, gpp_atomtype_t ga)
66 {
67     int i;
68
69     /* Atom types are always case sensitive */
70     for (i = 0; (i < ga->nr); i++)
71     {
72         if (strcmp(str, *(ga->atomname[i])) == 0)
73         {
74             return i;
75         }
76     }
77
78     return NOTSET;
79 }
80
81 int get_atomtype_ntypes(gpp_atomtype_t ga)
82 {
83     return ga->nr;
84 }
85
86 char *get_atomtype_name(int nt, gpp_atomtype_t ga)
87 {
88     if ((nt < 0) || (nt >= ga->nr))
89     {
90         return nullptr;
91     }
92
93     return *(ga->atomname[nt]);
94 }
95
96 real get_atomtype_massA(int nt, gpp_atomtype_t ga)
97 {
98     if ((nt < 0) || (nt >= ga->nr))
99     {
100         return NOTSET;
101     }
102
103     return ga->atom[nt].m;
104 }
105
106 real get_atomtype_massB(int nt, gpp_atomtype_t ga)
107 {
108     if ((nt < 0) || (nt >= ga->nr))
109     {
110         return NOTSET;
111     }
112
113     return ga->atom[nt].mB;
114 }
115
116 real get_atomtype_qA(int nt, gpp_atomtype_t ga)
117 {
118     if ((nt < 0) || (nt >= ga->nr))
119     {
120         return NOTSET;
121     }
122
123     return ga->atom[nt].q;
124 }
125
126 real get_atomtype_qB(int nt, gpp_atomtype_t ga)
127 {
128     if ((nt < 0) || (nt >= ga->nr))
129     {
130         return NOTSET;
131     }
132
133     return ga->atom[nt].qB;
134 }
135
136 int get_atomtype_ptype(int nt, gpp_atomtype_t ga)
137 {
138     if ((nt < 0) || (nt >= ga->nr))
139     {
140         return NOTSET;
141     }
142
143     return ga->atom[nt].ptype;
144 }
145
146 int get_atomtype_batype(int nt, const gpp_atomtype* ga)
147 {
148     if ((nt < 0) || (nt >= ga->nr))
149     {
150         return NOTSET;
151     }
152
153     return ga->bondatomtype[nt];
154 }
155
156 int get_atomtype_atomnumber(int nt, gpp_atomtype_t ga)
157 {
158     if ((nt < 0) || (nt >= ga->nr))
159     {
160         return NOTSET;
161     }
162
163     return ga->atomnumber[nt];
164 }
165
166 real get_atomtype_nbparam(int nt, int param, gpp_atomtype_t ga)
167 {
168     if ((nt < 0) || (nt >= ga->nr))
169     {
170         return NOTSET;
171     }
172     if ((param < 0) || (param >= MAXFORCEPARAM))
173     {
174         return NOTSET;
175     }
176     return ga->nb[nt].c[param];
177 }
178
179 gpp_atomtype_t init_atomtype(void)
180 {
181     gpp_atomtype_t ga;
182
183     snew(ga, 1);
184
185     ga->nr           = 0;
186     ga->atom         = nullptr;
187     ga->atomname     = nullptr;
188     ga->nb           = nullptr;
189     ga->bondatomtype = nullptr;
190     ga->atomnumber   = nullptr;
191
192     return ga;
193 }
194
195 int set_atomtype(int nt, gpp_atomtype_t ga, t_symtab *tab,
196                  t_atom *a, const char *name, t_param *nb,
197                  int bondatomtype, int atomnumber)
198 {
199     if ((nt < 0) || (nt >= ga->nr))
200     {
201         return NOTSET;
202     }
203
204     ga->atom[nt]         = *a;
205     ga->atomname[nt]     = put_symtab(tab, name);
206     ga->nb[nt]           = *nb;
207     ga->bondatomtype[nt] = bondatomtype;
208     ga->atomnumber[nt]   = atomnumber;
209
210     return nt;
211 }
212
213 int add_atomtype(gpp_atomtype_t ga, t_symtab *tab,
214                  t_atom *a, const char *name, t_param *nb,
215                  int bondatomtype, int atomnumber)
216 {
217     int i;
218
219     for (i = 0; (i < ga->nr); i++)
220     {
221         if (strcmp(*ga->atomname[i], name) == 0)
222         {
223             break;
224         }
225     }
226     if (i == ga->nr)
227     {
228         ga->nr++;
229         srenew(ga->atom, ga->nr);
230         srenew(ga->atomname, ga->nr);
231         srenew(ga->nb, ga->nr);
232         srenew(ga->bondatomtype, ga->nr);
233         srenew(ga->atomnumber, ga->nr);
234
235         return set_atomtype(ga->nr-1, ga, tab, a, name, nb, bondatomtype, atomnumber);
236     }
237     else
238     {
239         return i;
240     }
241 }
242
243 void print_at (FILE * out, gpp_atomtype_t ga)
244 {
245     int         i;
246     t_atom     *atom = ga->atom;
247     t_param    *nb   = ga->nb;
248
249     fprintf (out, "[ %s ]\n", dir2str(d_atomtypes));
250     fprintf (out, "; %6s  %8s  %8s  %8s  %12s  %12s\n",
251              "type", "mass", "charge", "particle", "c6", "c12");
252     for (i = 0; (i < ga->nr); i++)
253     {
254         fprintf(out, "%8s  %8.3f  %8.3f  %8s  %12e  %12e\n",
255                 *(ga->atomname[i]), atom[i].m, atom[i].q, "A",
256                 nb[i].c0(), nb[i].c1());
257     }
258
259     fprintf (out, "\n");
260 }
261
262 void done_atomtype(gpp_atomtype_t ga)
263 {
264     sfree(ga->atom);
265     sfree(ga->atomname);
266     sfree(ga->nb);
267     sfree(ga->bondatomtype);
268     sfree(ga->atomnumber);
269     ga->nr = 0;
270     sfree(ga);
271 }
272
273 static int search_atomtypes(gpp_atomtype_t ga, int *n, int typelist[],
274                             int thistype,
275                             t_param param[], int ftype)
276 {
277     int      i, nn, nrfp, j, k, ntype, tli;
278     bool     bFound = FALSE;
279
280     nn    = *n;
281     nrfp  = NRFP(ftype);
282     ntype = get_atomtype_ntypes(ga);
283
284     for (i = 0; (i < nn); i++)
285     {
286         if (typelist[i] == thistype)
287         {
288             /* This type number has already been added */
289             break;
290         }
291
292         /* Otherwise, check if the parameters are identical to any previously added type */
293
294         bFound = TRUE;
295         for (j = 0; j < ntype && bFound; j++)
296         {
297             /* Check nonbonded parameters */
298             for (k = 0; k < nrfp && bFound; k++)
299             {
300                 bFound = (param[ntype*typelist[i]+j].c[k] == param[ntype*thistype+j].c[k]);
301             }
302
303             /* Check atomnumber */
304             tli    = typelist[i];
305             bFound = bFound &&
306                 (get_atomtype_atomnumber(tli, ga) == get_atomtype_atomnumber(thistype, ga));
307         }
308         if (bFound)
309         {
310             break;
311         }
312     }
313
314     if (i == nn)
315     {
316         if (nn == ntype)
317         {
318             gmx_fatal(FARGS, "Atomtype horror n = %d, %s, %d", nn, __FILE__, __LINE__);
319         }
320         typelist[nn] = thistype;
321         nn++;
322     }
323     *n = nn;
324
325     return i;
326 }
327
328 void renum_atype(t_params plist[], gmx_mtop_t *mtop,
329                  int *wall_atomtype,
330                  gpp_atomtype_t ga, bool bVerbose)
331 {
332     int         i, j, k, l, mi, mj, nat, nrfp, ftype, ntype;
333     t_atoms    *atoms;
334     t_param    *nbsnew;
335     int        *typelist;
336     int        *new_atomnumber;
337     char     ***new_atomname;
338
339     ntype = get_atomtype_ntypes(ga);
340     snew(typelist, ntype);
341
342     if (bVerbose)
343     {
344         fprintf(stderr, "renumbering atomtypes...\n");
345     }
346
347     /* Since the bonded interactions have been assigned now,
348      * we want to reduce the number of atom types by merging
349      * ones with identical nonbonded interactions, in addition
350      * to removing unused ones.
351      *
352      * With QM/MM we also check that the atom numbers match
353      */
354
355     /* Get nonbonded interaction type */
356     if (plist[F_LJ].nr > 0)
357     {
358         ftype = F_LJ;
359     }
360     else
361     {
362         ftype = F_BHAM;
363     }
364
365     /* Renumber atomtypes by first making a list of which ones are actually used.
366      * We provide the list of nonbonded parameters so search_atomtypes
367      * can determine if two types should be merged.
368      */
369     nat = 0;
370     for (gmx_moltype_t &moltype : mtop->moltype)
371     {
372         atoms = &moltype.atoms;
373         for (i = 0; (i < atoms->nr); i++)
374         {
375             atoms->atom[i].type =
376                 search_atomtypes(ga, &nat, typelist, atoms->atom[i].type,
377                                  plist[ftype].param, ftype);
378             atoms->atom[i].typeB =
379                 search_atomtypes(ga, &nat, typelist, atoms->atom[i].typeB,
380                                  plist[ftype].param, ftype);
381         }
382     }
383
384     for (i = 0; i < 2; i++)
385     {
386         if (wall_atomtype[i] >= 0)
387         {
388             wall_atomtype[i] = search_atomtypes(ga, &nat, typelist, wall_atomtype[i],
389                                                 plist[ftype].param, ftype);
390         }
391     }
392
393     snew(new_atomnumber, nat);
394     snew(new_atomname, nat);
395     /* We now have a list of unique atomtypes in typelist */
396
397     /* Renumber nlist */
398     nbsnew = nullptr;
399     snew(nbsnew, plist[ftype].nr);
400
401     nrfp  = NRFP(ftype);
402
403     for (i = k = 0; (i < nat); i++)
404     {
405         mi = typelist[i];
406         for (j = 0; (j < nat); j++, k++)
407         {
408             mj = typelist[j];
409             for (l = 0; (l < nrfp); l++)
410             {
411                 nbsnew[k].c[l] = plist[ftype].param[ntype*mi+mj].c[l];
412             }
413         }
414         new_atomnumber[i] = get_atomtype_atomnumber(mi, ga);
415         new_atomname[i]   = ga->atomname[mi];
416     }
417
418     for (i = 0; (i < nat*nat); i++)
419     {
420         for (l = 0; (l < nrfp); l++)
421         {
422             plist[ftype].param[i].c[l] = nbsnew[i].c[l];
423         }
424     }
425     plist[ftype].nr     = i;
426     mtop->ffparams.atnr = nat;
427
428     sfree(ga->atomnumber);
429     /* Dangling atomname pointers ? */
430     sfree(ga->atomname);
431
432     ga->atomnumber = new_atomnumber;
433     ga->atomname   = new_atomname;
434
435     ga->nr = nat;
436
437     sfree(nbsnew);
438     sfree(typelist);
439 }
440
441 void copy_atomtype_atomtypes(gpp_atomtype_t ga, t_atomtypes *atomtypes)
442 {
443     int i, ntype;
444
445     /* Copy the atomtype data to the topology atomtype list */
446     ntype         = get_atomtype_ntypes(ga);
447     atomtypes->nr = ntype;
448     snew(atomtypes->atomnumber, ntype);
449
450     for (i = 0; i < ntype; i++)
451     {
452         atomtypes->atomnumber[i] = ga->atomnumber[i];
453     }
454 }