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