Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / tomorse.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) 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 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <stdlib.h>
43 #include <math.h>
44 #include <ctype.h>
45 #include "typedefs.h"
46 #include "string2.h"
47 #include "grompp-impl.h"
48 #include "gromacs/fileio/futil.h"
49 #include "smalloc.h"
50 #include "toputil.h"
51 #include "gmx_fatal.h"
52 #include "gpp_atomtype.h"
53 #include "macros.h"
54
55 #include "tomorse.h"
56
57 typedef struct {
58     char *ai, *aj;
59     real  e_diss;
60 } t_2morse;
61
62 static t_2morse *read_dissociation_energies(int *n2morse)
63 {
64     FILE       *fp;
65     char        ai[32], aj[32];
66     double      e_diss;
67     const char *fn     = "edissoc.dat";
68     t_2morse   *t2m    = NULL;
69     int         maxn2m = 0, n2m = 0;
70     int         nread;
71
72     /* Open the file with dissociation energies */
73     fp = libopen(fn);
74     do
75     {
76         /* Try and read two atom names and an energy term from it */
77         nread = fscanf(fp, "%s%s%lf", ai, aj, &e_diss);
78         if (nread == 3)
79         {
80             /* If we got three terms, it probably was OK, no further checking */
81             if (n2m >= maxn2m)
82             {
83                 /* Increase memory for 16 at once, some mallocs are stupid */
84                 maxn2m += 16;
85                 srenew(t2m, maxn2m);
86             }
87             /* Copy the values */
88             t2m[n2m].ai     = strdup(ai);
89             t2m[n2m].aj     = strdup(aj);
90             t2m[n2m].e_diss = e_diss;
91             /* Increment counter */
92             n2m++;
93         }
94         /* If we did not read three items, quit reading */
95     }
96     while (nread == 3);
97     gmx_ffclose(fp);
98
99     /* Set the return values */
100     *n2morse = n2m;
101
102     return t2m;
103 }
104
105 static int nequal(char *a1, char *a2)
106 {
107     int i;
108
109     /* Count the number of (case insensitive) characters that are equal in
110      * two strings. If they are equally long their respective null characters are
111      * counted also.
112      */
113     for (i = 0; (a1[i] != '\0') && (a2[i] != '\0'); i++)
114     {
115         if (toupper(a1[i]) != toupper(a2[i]))
116         {
117             break;
118         }
119     }
120     if ((a1[i] == '\0') && (a2[i] == '\0'))
121     {
122         i++;
123     }
124
125     return i;
126 }
127
128 static real search_e_diss(int n2m, t_2morse t2m[], char *ai, char *aj)
129 {
130     int  i;
131     int  ibest = -1;
132     int  nii, njj, nbstii = 0, nbstjj = 0;
133     real ediss = 400;
134
135     /* Do a best match search for dissociation energies */
136     for (i = 0; (i < n2m); i++)
137     {
138         /* Check for a perfect match */
139         if (((gmx_strcasecmp(t2m[i].ai, ai) == 0) && (gmx_strcasecmp(t2m[i].aj, aj) == 0)) ||
140             ((gmx_strcasecmp(t2m[i].aj, ai) == 0) && (gmx_strcasecmp(t2m[i].ai, aj) == 0)))
141         {
142             ibest = i;
143             break;
144         }
145         else
146         {
147             /* Otherwise count the number of equal characters in the strings ai and aj
148              * and the ones from the file
149              */
150             nii = nequal(t2m[i].ai, ai);
151             njj = nequal(t2m[i].aj, aj);
152             if (((nii >  nbstii) && (njj >= nbstjj)) ||
153                 ((nii >= nbstii) && (njj >  nbstjj)))
154             {
155                 if ((nii > 0) && (njj > 0))
156                 {
157                     ibest  = i;
158                     nbstii = nii;
159                     nbstjj = njj;
160                 }
161             }
162             else
163             {
164                 /* Swap ai and aj (at least in counting the number of equal chars) */
165                 nii = nequal(t2m[i].ai, aj);
166                 njj = nequal(t2m[i].aj, ai);
167                 if (((nii >  nbstii) && (njj >= nbstjj)) ||
168                     ((nii >= nbstii) && (njj >  nbstjj)))
169                 {
170                     if ((nii > 0) && (njj > 0))
171                     {
172                         ibest  = i;
173                         nbstii = nii;
174                         nbstjj = njj;
175                     }
176                 }
177             }
178         }
179     }
180     /* Return the dissocation energy corresponding to the best match, if we have
181      * found one. Do some debug output anyway.
182      */
183     if (ibest == -1)
184     {
185         if (debug)
186         {
187             fprintf(debug, "MORSE: Couldn't find E_diss for bond %s - %s, using default %g\n", ai, aj, ediss);
188         }
189         return ediss;
190     }
191     else
192     {
193         if (debug)
194         {
195             fprintf(debug, "MORSE: Dissoc. E (%10.3f) for bond %4s-%4s taken from bond %4s-%4s\n",
196                     t2m[ibest].e_diss, ai, aj, t2m[ibest].ai, t2m[ibest].aj);
197         }
198         return t2m[ibest].e_diss;
199     }
200 }
201
202 void convert_harmonics(int nrmols, t_molinfo mols[], gpp_atomtype_t atype)
203 {
204     int       n2m;
205     t_2morse *t2m;
206
207     int       i, j, k, last, ni, nj;
208     int       nrharm, nrmorse, bb;
209     real      edis, kb, b0, beta;
210     gmx_bool *bRemoveHarm;
211
212     /* First get the data */
213     t2m = read_dissociation_energies(&n2m);
214     if (debug)
215     {
216         fprintf(debug, "MORSE: read %d dissoc energies\n", n2m);
217     }
218     if (n2m <= 0)
219     {
220         fprintf(stderr, "No dissocation energies read\n");
221         return;
222     }
223
224     /* For all the molecule types */
225     for (i = 0; (i < nrmols); i++)
226     {
227         /* Check how many morse and harmonic BONDSs there are, increase size of
228          * morse with the number of harmonics
229          */
230         nrmorse = mols[i].plist[F_MORSE].nr;
231
232         for (bb = 0; (bb < F_NRE); bb++)
233         {
234             if ((interaction_function[bb].flags & IF_BTYPE) && (bb != F_MORSE))
235             {
236                 nrharm  = mols[i].plist[bb].nr;
237                 pr_alloc(nrharm, &(mols[i].plist[F_MORSE]));
238                 snew(bRemoveHarm, nrharm);
239
240                 /* Now loop over the harmonics, trying to convert them */
241                 for (j = 0; (j < nrharm); j++)
242                 {
243                     ni   = mols[i].plist[bb].param[j].AI;
244                     nj   = mols[i].plist[bb].param[j].AJ;
245                     edis =
246                         search_e_diss(n2m, t2m,
247                                       get_atomtype_name(mols[i].atoms.atom[ni].type, atype),
248                                       get_atomtype_name(mols[i].atoms.atom[nj].type, atype));
249                     if (edis != 0)
250                     {
251                         bRemoveHarm[j] = TRUE;
252                         b0             = mols[i].plist[bb].param[j].c[0];
253                         kb             = mols[i].plist[bb].param[j].c[1];
254                         beta           = sqrt(kb/(2*edis));
255                         mols[i].plist[F_MORSE].param[nrmorse].a[0] = ni;
256                         mols[i].plist[F_MORSE].param[nrmorse].a[1] = nj;
257                         mols[i].plist[F_MORSE].param[nrmorse].c[0] = b0;
258                         mols[i].plist[F_MORSE].param[nrmorse].c[1] = edis;
259                         mols[i].plist[F_MORSE].param[nrmorse].c[2] = beta;
260                         nrmorse++;
261                     }
262                 }
263                 mols[i].plist[F_MORSE].nr = nrmorse;
264
265                 /* Now remove the harmonics */
266                 for (j = last = 0; (j < nrharm); j++)
267                 {
268                     if (!bRemoveHarm[j])
269                     {
270                         /* Copy it to the last position */
271                         for (k = 0; (k < MAXATOMLIST); k++)
272                         {
273                             mols[i].plist[bb].param[last].a[k] =
274                                 mols[i].plist[bb].param[j].a[k];
275                         }
276                         for (k = 0; (k < MAXFORCEPARAM); k++)
277                         {
278                             mols[i].plist[bb].param[last].c[k] =
279                                 mols[i].plist[bb].param[j].c[k];
280                         }
281                         last++;
282                     }
283                 }
284                 sfree(bRemoveHarm);
285                 fprintf(stderr, "Converted %d out of %d %s to morse bonds for mol %d\n",
286                         nrharm-last, nrharm, interaction_function[bb].name, i);
287                 mols[i].plist[bb].nr = last;
288             }
289         }
290     }
291     sfree(t2m);
292 }