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