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,2015,2017,2018,2019, 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! */
47 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
48 #include "gromacs/gmxpreprocess/grompp_impl.h"
49 #include "gromacs/gmxpreprocess/toputil.h"
50 #include "gromacs/topology/ifunc.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/smalloc.h"
61 static t_2morse *read_dissociation_energies(int *n2morse)
65 const char *fn = "edissoc.dat";
66 t_2morse *t2m = nullptr;
67 int maxn2m = 0, n2m = 0;
70 /* Open the file with dissociation energies */
71 gmx::FilePtr fp = gmx::openLibraryFile(fn);
74 /* Try and read two atom names and an energy term from it */
75 nread = fscanf(fp.get(), "%s%s%lf", ai, aj, &e_diss);
78 /* If we got three terms, it probably was OK, no further checking */
81 /* Increase memory for 16 at once, some mallocs are stupid */
86 t2m[n2m].ai = gmx_strdup(ai);
87 t2m[n2m].aj = gmx_strdup(aj);
88 t2m[n2m].e_diss = e_diss;
89 /* Increment counter */
92 /* 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
186 return t2m[ibest].e_diss;
190 void convert_harmonics(int nrmols, t_molinfo mols[], gpp_atomtype *atype)
195 int i, j, k, last, ni, nj;
196 int nrharm, nrmorse, bb;
197 real edis, kb, b0, beta;
200 /* First get the data */
201 t2m = read_dissociation_energies(&n2m);
204 fprintf(stderr, "No dissocation energies read\n");
208 /* For all the molecule types */
209 for (i = 0; (i < nrmols); i++)
211 /* Check how many morse and harmonic BONDSs there are, increase size of
212 * morse with the number of harmonics
214 nrmorse = mols[i].plist[F_MORSE].nr;
216 for (bb = 0; (bb < F_NRE); bb++)
218 if ((interaction_function[bb].flags & IF_BTYPE) && (bb != F_MORSE))
220 nrharm = mols[i].plist[bb].nr;
221 pr_alloc(nrharm, &(mols[i].plist[F_MORSE]));
222 snew(bRemoveHarm, nrharm);
224 /* Now loop over the harmonics, trying to convert them */
225 for (j = 0; (j < nrharm); j++)
227 ni = mols[i].plist[bb].param[j].ai();
228 nj = mols[i].plist[bb].param[j].aj();
230 search_e_diss(n2m, t2m,
231 get_atomtype_name(mols[i].atoms.atom[ni].type, atype),
232 get_atomtype_name(mols[i].atoms.atom[nj].type, atype));
235 bRemoveHarm[j] = TRUE;
236 b0 = mols[i].plist[bb].param[j].c[0];
237 kb = mols[i].plist[bb].param[j].c[1];
238 beta = std::sqrt(kb/(2*edis));
239 mols[i].plist[F_MORSE].param[nrmorse].a[0] = ni;
240 mols[i].plist[F_MORSE].param[nrmorse].a[1] = nj;
241 mols[i].plist[F_MORSE].param[nrmorse].c[0] = b0;
242 mols[i].plist[F_MORSE].param[nrmorse].c[1] = edis;
243 mols[i].plist[F_MORSE].param[nrmorse].c[2] = beta;
247 mols[i].plist[F_MORSE].nr = nrmorse;
249 /* Now remove the harmonics */
250 for (j = last = 0; (j < nrharm); j++)
254 /* Copy it to the last position */
255 for (k = 0; (k < MAXATOMLIST); k++)
257 mols[i].plist[bb].param[last].a[k] =
258 mols[i].plist[bb].param[j].a[k];
260 for (k = 0; (k < MAXFORCEPARAM); k++)
262 mols[i].plist[bb].param[last].c[k] =
263 mols[i].plist[bb].param[j].c[k];
269 fprintf(stderr, "Converted %d out of %d %s to morse bonds for mol %d\n",
270 nrharm-last, nrharm, interaction_function[bb].name, i);
271 mols[i].plist[bb].nr = last;