3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 /* This file is completely threadsafe - keep it that way! */
49 #include "gmx_fatal.h"
50 #include "gpp_atomtype.h"
51 #include "gpp_tomorse.h"
59 static t_2morse *read_dissociation_energies(int *n2morse)
64 const char *fn = "edissoc.dat";
66 int maxn2m = 0, n2m = 0;
69 /* Open the file with dissociation energies */
73 /* Try and read two atom names and an energy term from it */
74 nread = fscanf(fp, "%s%s%lf", ai, aj, &e_diss);
77 /* If we got three terms, it probably was OK, no further checking */
80 /* Increase memory for 16 at once, some mallocs are stupid */
85 t2m[n2m].ai = strdup(ai);
86 t2m[n2m].aj = strdup(aj);
87 t2m[n2m].e_diss = e_diss;
88 /* Increment counter */
91 /* If we did not read three items, quit reading */
96 /* Set the return values */
102 static int nequal(char *a1, char *a2)
106 /* Count the number of (case insensitive) characters that are equal in
107 * two strings. If they are equally long their respective null characters are
110 for (i = 0; (a1[i] != '\0') && (a2[i] != '\0'); i++)
112 if (toupper(a1[i]) != toupper(a2[i]))
117 if ((a1[i] == '\0') && (a2[i] == '\0'))
125 static real search_e_diss(int n2m, t_2morse t2m[], char *ai, char *aj)
129 int nii, njj, nbstii = 0, nbstjj = 0;
132 /* Do a best match search for dissociation energies */
133 for (i = 0; (i < n2m); i++)
135 /* Check for a perfect match */
136 if (((gmx_strcasecmp(t2m[i].ai, ai) == 0) && (gmx_strcasecmp(t2m[i].aj, aj) == 0)) ||
137 ((gmx_strcasecmp(t2m[i].aj, ai) == 0) && (gmx_strcasecmp(t2m[i].ai, aj) == 0)))
144 /* Otherwise count the number of equal characters in the strings ai and aj
145 * and the ones from the file
147 nii = nequal(t2m[i].ai, ai);
148 njj = nequal(t2m[i].aj, aj);
149 if (((nii > nbstii) && (njj >= nbstjj)) ||
150 ((nii >= nbstii) && (njj > nbstjj)))
152 if ((nii > 0) && (njj > 0))
161 /* Swap ai and aj (at least in counting the number of equal chars) */
162 nii = nequal(t2m[i].ai, aj);
163 njj = nequal(t2m[i].aj, ai);
164 if (((nii > nbstii) && (njj >= nbstjj)) ||
165 ((nii >= nbstii) && (njj > nbstjj)))
167 if ((nii > 0) && (njj > 0))
177 /* Return the dissocation energy corresponding to the best match, if we have
178 * found one. Do some debug output anyway.
184 fprintf(debug, "MORSE: Couldn't find E_diss for bond %s - %s, using default %g\n", ai, aj, ediss);
192 fprintf(debug, "MORSE: Dissoc. E (%10.3f) for bond %4s-%4s taken from bond %4s-%4s\n",
193 t2m[ibest].e_diss, ai, aj, t2m[ibest].ai, t2m[ibest].aj);
195 return t2m[ibest].e_diss;
199 void convert_harmonics(int nrmols, t_molinfo mols[], gpp_atomtype_t atype)
204 int i, j, k, last, ni, nj;
205 int nrharm, nrmorse, bb;
206 real edis, kb, b0, beta;
207 gmx_bool *bRemoveHarm;
209 /* First get the data */
210 t2m = read_dissociation_energies(&n2m);
213 fprintf(debug, "MORSE: read %d dissoc energies\n", n2m);
217 fprintf(stderr, "No dissocation energies read\n");
221 /* For all the molecule types */
222 for (i = 0; (i < nrmols); i++)
224 /* Check how many morse and harmonic BONDSs there are, increase size of
225 * morse with the number of harmonics
227 nrmorse = mols[i].plist[F_MORSE].nr;
229 for (bb = 0; (bb < F_NRE); bb++)
231 if ((interaction_function[bb].flags & IF_BTYPE) && (bb != F_MORSE))
233 nrharm = mols[i].plist[bb].nr;
234 pr_alloc(nrharm, &(mols[i].plist[F_MORSE]));
235 snew(bRemoveHarm, nrharm);
237 /* Now loop over the harmonics, trying to convert them */
238 for (j = 0; (j < nrharm); j++)
240 ni = mols[i].plist[bb].param[j].AI;
241 nj = mols[i].plist[bb].param[j].AJ;
243 search_e_diss(n2m, t2m,
244 get_atomtype_name(mols[i].atoms.atom[ni].type, atype),
245 get_atomtype_name(mols[i].atoms.atom[nj].type, atype));
248 bRemoveHarm[j] = TRUE;
249 b0 = mols[i].plist[bb].param[j].c[0];
250 kb = mols[i].plist[bb].param[j].c[1];
251 beta = sqrt(kb/(2*edis));
252 mols[i].plist[F_MORSE].param[nrmorse].a[0] = ni;
253 mols[i].plist[F_MORSE].param[nrmorse].a[1] = nj;
254 mols[i].plist[F_MORSE].param[nrmorse].c[0] = b0;
255 mols[i].plist[F_MORSE].param[nrmorse].c[1] = edis;
256 mols[i].plist[F_MORSE].param[nrmorse].c[2] = beta;
260 mols[i].plist[F_MORSE].nr = nrmorse;
262 /* Now remove the harmonics */
263 for (j = last = 0; (j < nrharm); j++)
267 /* Copy it to the last position */
268 for (k = 0; (k < MAXATOMLIST); k++)
270 mols[i].plist[bb].param[last].a[k] =
271 mols[i].plist[bb].param[j].a[k];
273 for (k = 0; (k < MAXFORCEPARAM); k++)
275 mols[i].plist[bb].param[last].c[k] =
276 mols[i].plist[bb].param[j].c[k];
282 fprintf(stderr, "Converted %d out of %d %s to morse bonds for mol %d\n",
283 nrharm-last, nrharm, interaction_function[bb].name, i);
284 mols[i].plist[bb].nr = last;