708b18a824f361e0eac00af162bc664810175273
[alexxy/gromacs.git] / src / programs / gmx / tomorse.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <stdlib.h>
41 #include <math.h>
42 #include <ctype.h>
43 #include "typedefs.h"
44 #include "string2.h"
45 #include "grompp.h"
46 #include "futil.h"
47 #include "smalloc.h"
48 #include "toputil.h"
49 #include "gmx_fatal.h"
50 #include "gpp_atomtype.h"
51 #include "macros.h"
52
53 #include "tomorse.h"
54
55 typedef struct {
56     char *ai, *aj;
57     real  e_diss;
58 } t_2morse;
59
60 static t_2morse *read_dissociation_energies(int *n2morse)
61 {
62     FILE       *fp;
63     char        ai[32], aj[32];
64     double      e_diss;
65     const char *fn     = "edissoc.dat";
66     t_2morse   *t2m    = NULL;
67     int         maxn2m = 0, n2m = 0;
68     int         nread;
69
70     /* Open the file with dissociation energies */
71     fp = libopen(fn);
72     do
73     {
74         /* Try and read two atom names and an energy term from it */
75         nread = fscanf(fp, "%s%s%lf", ai, aj, &e_diss);
76         if (nread == 3)
77         {
78             /* If we got three terms, it probably was OK, no further checking */
79             if (n2m >= maxn2m)
80             {
81                 /* Increase memory for 16 at once, some mallocs are stupid */
82                 maxn2m += 16;
83                 srenew(t2m, maxn2m);
84             }
85             /* Copy the values */
86             t2m[n2m].ai     = strdup(ai);
87             t2m[n2m].aj     = strdup(aj);
88             t2m[n2m].e_diss = e_diss;
89             /* Increment counter */
90             n2m++;
91         }
92         /* If we did not read three items, quit reading */
93     }
94     while (nread == 3);
95     ffclose(fp);
96
97     /* Set the return values */
98     *n2morse = n2m;
99
100     return t2m;
101 }
102
103 static int nequal(char *a1, char *a2)
104 {
105     int i;
106
107     /* Count the number of (case insensitive) characters that are equal in
108      * two strings. If they are equally long their respective null characters are
109      * counted also.
110      */
111     for (i = 0; (a1[i] != '\0') && (a2[i] != '\0'); i++)
112     {
113         if (toupper(a1[i]) != toupper(a2[i]))
114         {
115             break;
116         }
117     }
118     if ((a1[i] == '\0') && (a2[i] == '\0'))
119     {
120         i++;
121     }
122
123     return i;
124 }
125
126 static real search_e_diss(int n2m, t_2morse t2m[], char *ai, char *aj)
127 {
128     int  i;
129     int  ibest = -1;
130     int  nii, njj, nbstii = 0, nbstjj = 0;
131     real ediss = 400;
132
133     /* Do a best match search for dissociation energies */
134     for (i = 0; (i < n2m); i++)
135     {
136         /* Check for a perfect match */
137         if (((gmx_strcasecmp(t2m[i].ai, ai) == 0) && (gmx_strcasecmp(t2m[i].aj, aj) == 0)) ||
138             ((gmx_strcasecmp(t2m[i].aj, ai) == 0) && (gmx_strcasecmp(t2m[i].ai, aj) == 0)))
139         {
140             ibest = i;
141             break;
142         }
143         else
144         {
145             /* Otherwise count the number of equal characters in the strings ai and aj
146              * and the ones from the file
147              */
148             nii = nequal(t2m[i].ai, ai);
149             njj = nequal(t2m[i].aj, aj);
150             if (((nii >  nbstii) && (njj >= nbstjj)) ||
151                 ((nii >= nbstii) && (njj >  nbstjj)))
152             {
153                 if ((nii > 0) && (njj > 0))
154                 {
155                     ibest  = i;
156                     nbstii = nii;
157                     nbstjj = njj;
158                 }
159             }
160             else
161             {
162                 /* Swap ai and aj (at least in counting the number of equal chars) */
163                 nii = nequal(t2m[i].ai, aj);
164                 njj = nequal(t2m[i].aj, ai);
165                 if (((nii >  nbstii) && (njj >= nbstjj)) ||
166                     ((nii >= nbstii) && (njj >  nbstjj)))
167                 {
168                     if ((nii > 0) && (njj > 0))
169                     {
170                         ibest  = i;
171                         nbstii = nii;
172                         nbstjj = njj;
173                     }
174                 }
175             }
176         }
177     }
178     /* Return the dissocation energy corresponding to the best match, if we have
179      * found one. Do some debug output anyway.
180      */
181     if (ibest == -1)
182     {
183         if (debug)
184         {
185             fprintf(debug, "MORSE: Couldn't find E_diss for bond %s - %s, using default %g\n", ai, aj, ediss);
186         }
187         return ediss;
188     }
189     else
190     {
191         if (debug)
192         {
193             fprintf(debug, "MORSE: Dissoc. E (%10.3f) for bond %4s-%4s taken from bond %4s-%4s\n",
194                     t2m[ibest].e_diss, ai, aj, t2m[ibest].ai, t2m[ibest].aj);
195         }
196         return t2m[ibest].e_diss;
197     }
198 }
199
200 void convert_harmonics(int nrmols, t_molinfo mols[], gpp_atomtype_t atype)
201 {
202     int       n2m;
203     t_2morse *t2m;
204
205     int       i, j, k, last, ni, nj;
206     int       nrharm, nrmorse, bb;
207     real      edis, kb, b0, beta;
208     gmx_bool *bRemoveHarm;
209
210     /* First get the data */
211     t2m = read_dissociation_energies(&n2m);
212     if (debug)
213     {
214         fprintf(debug, "MORSE: read %d dissoc energies\n", n2m);
215     }
216     if (n2m <= 0)
217     {
218         fprintf(stderr, "No dissocation energies read\n");
219         return;
220     }
221
222     /* For all the molecule types */
223     for (i = 0; (i < nrmols); i++)
224     {
225         /* Check how many morse and harmonic BONDSs there are, increase size of
226          * morse with the number of harmonics
227          */
228         nrmorse = mols[i].plist[F_MORSE].nr;
229
230         for (bb = 0; (bb < F_NRE); bb++)
231         {
232             if ((interaction_function[bb].flags & IF_BTYPE) && (bb != F_MORSE))
233             {
234                 nrharm  = mols[i].plist[bb].nr;
235                 pr_alloc(nrharm, &(mols[i].plist[F_MORSE]));
236                 snew(bRemoveHarm, nrharm);
237
238                 /* Now loop over the harmonics, trying to convert them */
239                 for (j = 0; (j < nrharm); j++)
240                 {
241                     ni   = mols[i].plist[bb].param[j].AI;
242                     nj   = mols[i].plist[bb].param[j].AJ;
243                     edis =
244                         search_e_diss(n2m, t2m,
245                                       get_atomtype_name(mols[i].atoms.atom[ni].type, atype),
246                                       get_atomtype_name(mols[i].atoms.atom[nj].type, atype));
247                     if (edis != 0)
248                     {
249                         bRemoveHarm[j] = TRUE;
250                         b0             = mols[i].plist[bb].param[j].c[0];
251                         kb             = mols[i].plist[bb].param[j].c[1];
252                         beta           = sqrt(kb/(2*edis));
253                         mols[i].plist[F_MORSE].param[nrmorse].a[0] = ni;
254                         mols[i].plist[F_MORSE].param[nrmorse].a[1] = nj;
255                         mols[i].plist[F_MORSE].param[nrmorse].c[0] = b0;
256                         mols[i].plist[F_MORSE].param[nrmorse].c[1] = edis;
257                         mols[i].plist[F_MORSE].param[nrmorse].c[2] = beta;
258                         nrmorse++;
259                     }
260                 }
261                 mols[i].plist[F_MORSE].nr = nrmorse;
262
263                 /* Now remove the harmonics */
264                 for (j = last = 0; (j < nrharm); j++)
265                 {
266                     if (!bRemoveHarm[j])
267                     {
268                         /* Copy it to the last position */
269                         for (k = 0; (k < MAXATOMLIST); k++)
270                         {
271                             mols[i].plist[bb].param[last].a[k] =
272                                 mols[i].plist[bb].param[j].a[k];
273                         }
274                         for (k = 0; (k < MAXFORCEPARAM); k++)
275                         {
276                             mols[i].plist[bb].param[last].c[k] =
277                                 mols[i].plist[bb].param[j].c[k];
278                         }
279                         last++;
280                     }
281                 }
282                 sfree(bRemoveHarm);
283                 fprintf(stderr, "Converted %d out of %d %s to morse bonds for mol %d\n",
284                         nrharm-last, nrharm, interaction_function[bb].name, i);
285                 mols[i].plist[bb].nr = last;
286             }
287         }
288     }
289     sfree(t2m);
290 }