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