2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
48 #include "gromacs/utility/cstringutil.h"
49 #include "grompp-impl.h"
50 #include "gromacs/fileio/futil.h"
51 #include "gromacs/utility/smalloc.h"
53 #include "gmx_fatal.h"
54 #include "gpp_atomtype.h"
64 static t_2morse *read_dissociation_energies(int *n2morse)
69 const char *fn = "edissoc.dat";
71 int maxn2m = 0, n2m = 0;
74 /* Open the file with dissociation energies */
78 /* Try and read two atom names and an energy term from it */
79 nread = fscanf(fp, "%s%s%lf", ai, aj, &e_diss);
82 /* If we got three terms, it probably was OK, no further checking */
85 /* Increase memory for 16 at once, some mallocs are stupid */
90 t2m[n2m].ai = strdup(ai);
91 t2m[n2m].aj = strdup(aj);
92 t2m[n2m].e_diss = e_diss;
93 /* Increment counter */
96 /* If we did not read three items, quit reading */
101 /* Set the return values */
107 static int nequal(char *a1, char *a2)
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
115 for (i = 0; (a1[i] != '\0') && (a2[i] != '\0'); i++)
117 if (toupper(a1[i]) != toupper(a2[i]))
122 if ((a1[i] == '\0') && (a2[i] == '\0'))
130 static real search_e_diss(int n2m, t_2morse t2m[], char *ai, char *aj)
134 int nii, njj, nbstii = 0, nbstjj = 0;
137 /* Do a best match search for dissociation energies */
138 for (i = 0; (i < n2m); i++)
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)))
149 /* Otherwise count the number of equal characters in the strings ai and aj
150 * and the ones from the file
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)))
157 if ((nii > 0) && (njj > 0))
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)))
172 if ((nii > 0) && (njj > 0))
182 /* Return the dissocation energy corresponding to the best match, if we have
183 * found one. Do some debug output anyway.
189 fprintf(debug, "MORSE: Couldn't find E_diss for bond %s - %s, using default %g\n", ai, aj, ediss);
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);
200 return t2m[ibest].e_diss;
204 void convert_harmonics(int nrmols, t_molinfo mols[], gpp_atomtype_t atype)
209 int i, j, k, last, ni, nj;
210 int nrharm, nrmorse, bb;
211 real edis, kb, b0, beta;
212 gmx_bool *bRemoveHarm;
214 /* First get the data */
215 t2m = read_dissociation_energies(&n2m);
218 fprintf(debug, "MORSE: read %d dissoc energies\n", n2m);
222 fprintf(stderr, "No dissocation energies read\n");
226 /* For all the molecule types */
227 for (i = 0; (i < nrmols); i++)
229 /* Check how many morse and harmonic BONDSs there are, increase size of
230 * morse with the number of harmonics
232 nrmorse = mols[i].plist[F_MORSE].nr;
234 for (bb = 0; (bb < F_NRE); bb++)
236 if ((interaction_function[bb].flags & IF_BTYPE) && (bb != F_MORSE))
238 nrharm = mols[i].plist[bb].nr;
239 pr_alloc(nrharm, &(mols[i].plist[F_MORSE]));
240 snew(bRemoveHarm, nrharm);
242 /* Now loop over the harmonics, trying to convert them */
243 for (j = 0; (j < nrharm); j++)
245 ni = mols[i].plist[bb].param[j].AI;
246 nj = mols[i].plist[bb].param[j].AJ;
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));
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;
265 mols[i].plist[F_MORSE].nr = nrmorse;
267 /* Now remove the harmonics */
268 for (j = last = 0; (j < nrharm); j++)
272 /* Copy it to the last position */
273 for (k = 0; (k < MAXATOMLIST); k++)
275 mols[i].plist[bb].param[last].a[k] =
276 mols[i].plist[bb].param[j].a[k];
278 for (k = 0; (k < MAXFORCEPARAM); k++)
280 mols[i].plist[bb].param[last].c[k] =
281 mols[i].plist[bb].param[j].c[k];
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;