Remove unnecessary config.h includes
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / gpp_atomtype.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-2004, The GROMACS development team.
6  * Copyright (c) 2011,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 #include "gmxpre.h"
38
39 #include <math.h>
40 #include <string.h>
41
42 #include "topdirs.h"
43 #include "toputil.h"
44 #include "gromacs/legacyheaders/txtdump.h"
45 #include "gpp_atomtype.h"
46
47 #include "gromacs/topology/symtab.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/smalloc.h"
50
51 typedef struct gpp_atomtype {
52     int              nr;           /* The number of atomtypes           */
53     t_atom          *atom;         /* Array of atoms                    */
54     char          ***atomname;     /* Names of the atomtypes            */
55     t_param         *nb;           /* Nonbonded force default params    */
56     int             *bondatomtype; /* The bond_atomtype for each atomtype  */
57     real            *radius;       /* Radius for GBSA stuff                */
58     real            *vol;          /* Effective volume for GBSA            */
59     real            *surftens;     /* Surface tension with water, for GBSA */
60     real            *gb_radius;    /* Radius for Still model               */
61     real            *S_hct;        /* Overlap factor for HCT model         */
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 NULL;
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, gpp_atomtype_t 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_radius(int nt, gpp_atomtype_t ga)
167 {
168     if ((nt < 0) || (nt >= ga->nr))
169     {
170         return NOTSET;
171     }
172
173     return ga->radius[nt];
174 }
175
176 real get_atomtype_vol(int nt, gpp_atomtype_t ga)
177 {
178     if ((nt < 0) || (nt >= ga->nr))
179     {
180         return NOTSET;
181     }
182
183     return ga->vol[nt];
184 }
185
186 real get_atomtype_surftens(int nt, gpp_atomtype_t ga)
187 {
188     if ((nt < 0) || (nt >= ga->nr))
189     {
190         return NOTSET;
191     }
192
193     return ga->surftens[nt];
194 }
195
196 real get_atomtype_gb_radius(int nt, gpp_atomtype_t ga)
197 {
198     if ((nt < 0) || (nt >= ga->nr))
199     {
200         return NOTSET;
201     }
202
203     return ga->gb_radius[nt];
204 }
205
206 real get_atomtype_S_hct(int nt, gpp_atomtype_t ga)
207 {
208     if ((nt < 0) || (nt >= ga->nr))
209     {
210         return NOTSET;
211     }
212
213     return ga->S_hct[nt];
214 }
215
216 real get_atomtype_nbparam(int nt, int param, gpp_atomtype_t ga)
217 {
218     if ((nt < 0) || (nt >= ga->nr))
219     {
220         return NOTSET;
221     }
222     if ((param < 0) || (param >= MAXFORCEPARAM))
223     {
224         return NOTSET;
225     }
226     return ga->nb[nt].c[param];
227 }
228
229 gpp_atomtype_t init_atomtype(void)
230 {
231     gpp_atomtype_t ga;
232
233     snew(ga, 1);
234
235     ga->nr           = 0;
236     ga->atom         = NULL;
237     ga->atomname     = NULL;
238     ga->nb           = NULL;
239     ga->bondatomtype = NULL;
240     ga->radius       = NULL;
241     ga->vol          = NULL;
242     ga->surftens     = NULL;
243     ga->atomnumber   = NULL;
244     ga->gb_radius    = NULL;
245     ga->S_hct        = NULL;
246
247     return ga;
248 }
249
250 int
251 set_atomtype_gbparam(gpp_atomtype_t ga, int i,
252                      real radius, real vol, real surftens,
253                      real gb_radius, real S_hct)
254 {
255     if ( (i < 0) || (i >= ga->nr))
256     {
257         return NOTSET;
258     }
259
260     ga->radius[i]    = radius;
261     ga->vol[i]       = vol;
262     ga->surftens[i]  = surftens;
263     ga->gb_radius[i] = gb_radius;
264     ga->S_hct[i]     = S_hct;
265
266     return i;
267 }
268
269
270 int set_atomtype(int nt, gpp_atomtype_t ga, t_symtab *tab,
271                  t_atom *a, const char *name, t_param *nb,
272                  int bondatomtype,
273                  real radius, real vol, real surftens, int atomnumber,
274                  real gb_radius, real S_hct)
275 {
276     if ((nt < 0) || (nt >= ga->nr))
277     {
278         return NOTSET;
279     }
280
281     ga->atom[nt]         = *a;
282     ga->atomname[nt]     = put_symtab(tab, name);
283     ga->nb[nt]           = *nb;
284     ga->bondatomtype[nt] = bondatomtype;
285     ga->radius[nt]       = radius;
286     ga->vol[nt]          = vol;
287     ga->surftens[nt]     = surftens;
288     ga->atomnumber[nt]   = atomnumber;
289     ga->gb_radius[nt]    = gb_radius;
290     ga->S_hct[nt]        = S_hct;
291
292     return nt;
293 }
294
295 int add_atomtype(gpp_atomtype_t ga, t_symtab *tab,
296                  t_atom *a, const char *name, t_param *nb,
297                  int bondatomtype,
298                  real radius, real vol, real surftens, real atomnumber,
299                  real gb_radius, real S_hct)
300 {
301     int i;
302
303     for (i = 0; (i < ga->nr); i++)
304     {
305         if (strcmp(*ga->atomname[i], name) == 0)
306         {
307             if (NULL != debug)
308             {
309                 fprintf(debug, "Trying to add atomtype %s again. Skipping it.\n", name);
310             }
311             break;
312         }
313     }
314     if (i == ga->nr)
315     {
316         ga->nr++;
317         srenew(ga->atom, ga->nr);
318         srenew(ga->atomname, ga->nr);
319         srenew(ga->nb, ga->nr);
320         srenew(ga->bondatomtype, ga->nr);
321         srenew(ga->radius, ga->nr);
322         srenew(ga->vol, ga->nr);
323         srenew(ga->surftens, ga->nr);
324         srenew(ga->atomnumber, ga->nr);
325         srenew(ga->gb_radius, ga->nr);
326         srenew(ga->S_hct, ga->nr);
327
328         return set_atomtype(ga->nr-1, ga, tab, a, name, nb, bondatomtype, radius,
329                             vol, surftens, atomnumber, gb_radius, S_hct);
330     }
331     else
332     {
333         return i;
334     }
335 }
336
337 void print_at (FILE * out, gpp_atomtype_t ga)
338 {
339     int         i;
340     t_atom     *atom = ga->atom;
341     t_param    *nb   = ga->nb;
342
343     fprintf (out, "[ %s ]\n", dir2str(d_atomtypes));
344     fprintf (out, "; %6s  %8s  %8s  %8s  %12s  %12s\n",
345              "type", "mass", "charge", "particle", "c6", "c12");
346     for (i = 0; (i < ga->nr); i++)
347     {
348         fprintf(out, "%8s  %8.3f  %8.3f  %8s  %12e  %12e\n",
349                 *(ga->atomname[i]), atom[i].m, atom[i].q, "A",
350                 nb[i].C0, nb[i].C1);
351     }
352
353     fprintf (out, "\n");
354 }
355
356 void done_atomtype(gpp_atomtype_t ga)
357 {
358     sfree(ga->atom);
359     sfree(ga->atomname);
360     sfree(ga->nb);
361     sfree(ga->bondatomtype);
362     sfree(ga->radius);
363     sfree(ga->vol);
364     sfree(ga->gb_radius);
365     sfree(ga->S_hct);
366     sfree(ga->surftens);
367     sfree(ga->atomnumber);
368     ga->nr = 0;
369     sfree(ga);
370 }
371
372 static int search_atomtypes(gpp_atomtype_t ga, int *n, int typelist[],
373                             int thistype,
374                             t_param param[], int ftype)
375 {
376     int      i, nn, nrfp, j, k, ntype, tli;
377     gmx_bool bFound = FALSE;
378
379     nn    = *n;
380     nrfp  = NRFP(ftype);
381     ntype = get_atomtype_ntypes(ga);
382
383     for (i = 0; (i < nn); i++)
384     {
385         if (typelist[i] == thistype)
386         {
387             /* This type number has already been added */
388             break;
389         }
390
391         /* Otherwise, check if the parameters are identical to any previously added type */
392
393         bFound = TRUE;
394         for (j = 0; j < ntype && bFound; j++)
395         {
396             /* Check nonbonded parameters */
397             for (k = 0; k < nrfp && bFound; k++)
398             {
399                 bFound = (param[ntype*typelist[i]+j].c[k] == param[ntype*thistype+j].c[k]);
400             }
401
402             /* Check radius, volume, surftens */
403             tli    = typelist[i];
404             bFound = bFound &&
405                 (get_atomtype_radius(tli, ga) == get_atomtype_radius(thistype, ga)) &&
406                 (get_atomtype_vol(tli, ga) == get_atomtype_vol(thistype, ga)) &&
407                 (get_atomtype_surftens(tli, ga) == get_atomtype_surftens(thistype, ga)) &&
408                 (get_atomtype_atomnumber(tli, ga) == get_atomtype_atomnumber(thistype, ga)) &&
409                 (get_atomtype_gb_radius(tli, ga) == get_atomtype_gb_radius(thistype, ga)) &&
410                 (get_atomtype_S_hct(tli, ga) == get_atomtype_S_hct(thistype, ga));
411         }
412         if (bFound)
413         {
414             break;
415         }
416     }
417
418     if (i == nn)
419     {
420         if (debug)
421         {
422             fprintf(debug, "Renumbering atomtype %d to %d\n", thistype, nn);
423         }
424         if (nn == ntype)
425         {
426             gmx_fatal(FARGS, "Atomtype horror n = %d, %s, %d", nn, __FILE__, __LINE__);
427         }
428         typelist[nn] = thistype;
429         nn++;
430     }
431     *n = nn;
432
433     return i;
434 }
435
436 void renum_atype(t_params plist[], gmx_mtop_t *mtop,
437                  int *wall_atomtype,
438                  gpp_atomtype_t ga, gmx_bool bVerbose)
439 {
440     int         i, j, k, l, molt, mi, mj, nat, nrfp, ftype, ntype;
441     t_atoms    *atoms;
442     t_param    *nbsnew;
443     int        *typelist;
444     real       *new_radius;
445     real       *new_vol;
446     real       *new_surftens;
447     real       *new_gb_radius;
448     real       *new_S_hct;
449     int        *new_atomnumber;
450     char     ***new_atomname;
451
452     ntype = get_atomtype_ntypes(ga);
453     snew(typelist, ntype);
454
455     if (bVerbose)
456     {
457         fprintf(stderr, "renumbering atomtypes...\n");
458     }
459
460     /* Since the bonded interactions have been assigned now,
461      * we want to reduce the number of atom types by merging
462      * ones with identical nonbonded interactions, in addition
463      * to removing unused ones.
464      *
465      * With Generalized-Born electrostatics, or implicit solvent
466      * we also check that the atomtype radius, effective_volume
467      * and surface tension match.
468      *
469      * With QM/MM we also check that the atom numbers match
470      */
471
472     /* Get nonbonded interaction type */
473     if (plist[F_LJ].nr > 0)
474     {
475         ftype = F_LJ;
476     }
477     else
478     {
479         ftype = F_BHAM;
480     }
481
482     /* Renumber atomtypes by first making a list of which ones are actually used.
483      * We provide the list of nonbonded parameters so search_atomtypes
484      * can determine if two types should be merged.
485      */
486     nat = 0;
487     for (molt = 0; molt < mtop->nmoltype; molt++)
488     {
489         atoms = &mtop->moltype[molt].atoms;
490         for (i = 0; (i < atoms->nr); i++)
491         {
492             atoms->atom[i].type =
493                 search_atomtypes(ga, &nat, typelist, atoms->atom[i].type,
494                                  plist[ftype].param, ftype);
495             atoms->atom[i].typeB =
496                 search_atomtypes(ga, &nat, typelist, atoms->atom[i].typeB,
497                                  plist[ftype].param, ftype);
498         }
499     }
500
501     for (i = 0; i < 2; i++)
502     {
503         if (wall_atomtype[i] >= 0)
504         {
505             wall_atomtype[i] = search_atomtypes(ga, &nat, typelist, wall_atomtype[i],
506                                                 plist[ftype].param, ftype);
507         }
508     }
509
510     snew(new_radius, nat);
511     snew(new_vol, nat);
512     snew(new_surftens, nat);
513     snew(new_atomnumber, nat);
514     snew(new_gb_radius, nat);
515     snew(new_S_hct, nat);
516     snew(new_atomname, nat);
517
518     /* We now have a list of unique atomtypes in typelist */
519
520     if (debug)
521     {
522         pr_ivec(debug, 0, "typelist", typelist, nat, TRUE);
523     }
524
525     /* Renumber nlist */
526     nbsnew = NULL;
527     snew(nbsnew, plist[ftype].nr);
528
529     nrfp  = NRFP(ftype);
530
531     for (i = k = 0; (i < nat); i++)
532     {
533         mi = typelist[i];
534         for (j = 0; (j < nat); j++, k++)
535         {
536             mj = typelist[j];
537             for (l = 0; (l < nrfp); l++)
538             {
539                 nbsnew[k].c[l] = plist[ftype].param[ntype*mi+mj].c[l];
540             }
541         }
542         new_radius[i]     = get_atomtype_radius(mi, ga);
543         new_vol[i]        = get_atomtype_vol(mi, ga);
544         new_surftens[i]   = get_atomtype_surftens(mi, ga);
545         new_atomnumber[i] = get_atomtype_atomnumber(mi, ga);
546         new_gb_radius[i]  = get_atomtype_gb_radius(mi, ga);
547         new_S_hct[i]      = get_atomtype_S_hct(mi, ga);
548         new_atomname[i]   = ga->atomname[mi];
549     }
550
551     for (i = 0; (i < nat*nat); i++)
552     {
553         for (l = 0; (l < nrfp); l++)
554         {
555             plist[ftype].param[i].c[l] = nbsnew[i].c[l];
556         }
557     }
558     plist[ftype].nr     = i;
559     mtop->ffparams.atnr = nat;
560
561     sfree(ga->radius);
562     sfree(ga->vol);
563     sfree(ga->surftens);
564     sfree(ga->atomnumber);
565     sfree(ga->gb_radius);
566     sfree(ga->S_hct);
567     /* Dangling atomname pointers ? */
568     sfree(ga->atomname);
569
570     ga->radius     = new_radius;
571     ga->vol        = new_vol;
572     ga->surftens   = new_surftens;
573     ga->atomnumber = new_atomnumber;
574     ga->gb_radius  = new_gb_radius;
575     ga->S_hct      = new_S_hct;
576     ga->atomname   = new_atomname;
577
578     ga->nr = nat;
579
580     sfree(nbsnew);
581     sfree(typelist);
582 }
583
584 void copy_atomtype_atomtypes(gpp_atomtype_t ga, t_atomtypes *atomtypes)
585 {
586     int i, ntype;
587
588     /* Copy the atomtype data to the topology atomtype list */
589     ntype         = get_atomtype_ntypes(ga);
590     atomtypes->nr = ntype;
591     snew(atomtypes->radius, ntype);
592     snew(atomtypes->vol, ntype);
593     snew(atomtypes->surftens, ntype);
594     snew(atomtypes->atomnumber, ntype);
595     snew(atomtypes->gb_radius, ntype);
596     snew(atomtypes->S_hct, ntype);
597
598     for (i = 0; i < ntype; i++)
599     {
600         atomtypes->radius[i]     = ga->radius[i];
601         atomtypes->vol[i]        = ga->vol[i];
602         atomtypes->surftens[i]   = ga->surftens[i];
603         atomtypes->atomnumber[i] = ga->atomnumber[i];
604         atomtypes->gb_radius[i]  = ga->gb_radius[i];
605         atomtypes->S_hct[i]      = ga->S_hct[i];
606     }
607 }