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