Add AWH biasing module + tests
[gromacs.git] / src / gromacs / gmxpreprocess / tomorse.cpp
blob97523622ec38fb414f814b8c3562ce937627a969
1 /*
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, 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! */
38 #include "gmxpre.h"
40 #include "tomorse.h"
42 #include <ctype.h>
43 #include <stdlib.h>
44 #include <string.h>
46 #include <cmath>
48 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
49 #include "gromacs/gmxpreprocess/grompp-impl.h"
50 #include "gromacs/gmxpreprocess/toputil.h"
51 #include "gromacs/topology/ifunc.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
57 typedef struct {
58 char *ai, *aj;
59 real e_diss;
60 } t_2morse;
62 static t_2morse *read_dissociation_energies(int *n2morse)
64 FILE *fp;
65 char ai[32], aj[32];
66 double e_diss;
67 const char *fn = "edissoc.dat";
68 t_2morse *t2m = nullptr;
69 int maxn2m = 0, n2m = 0;
70 int nread;
72 /* Open the file with dissociation energies */
73 fp = libopen(fn);
76 /* Try and read two atom names and an energy term from it */
77 nread = fscanf(fp, "%s%s%lf", ai, aj, &e_diss);
78 if (nread == 3)
80 /* If we got three terms, it probably was OK, no further checking */
81 if (n2m >= maxn2m)
83 /* Increase memory for 16 at once, some mallocs are stupid */
84 maxn2m += 16;
85 srenew(t2m, maxn2m);
87 /* Copy the values */
88 t2m[n2m].ai = gmx_strdup(ai);
89 t2m[n2m].aj = gmx_strdup(aj);
90 t2m[n2m].e_diss = e_diss;
91 /* Increment counter */
92 n2m++;
94 /* If we did not read three items, quit reading */
96 while (nread == 3);
97 gmx_ffclose(fp);
99 /* Set the return values */
100 *n2morse = n2m;
102 return t2m;
105 static int nequal(char *a1, char *a2)
107 int i;
109 /* Count the number of (case insensitive) characters that are equal in
110 * two strings. If they are equally long their respective null characters are
111 * counted also.
113 for (i = 0; (a1[i] != '\0') && (a2[i] != '\0'); i++)
115 if (toupper(a1[i]) != toupper(a2[i]))
117 break;
120 if ((a1[i] == '\0') && (a2[i] == '\0'))
122 i++;
125 return i;
128 static real search_e_diss(int n2m, t_2morse t2m[], char *ai, char *aj)
130 int i;
131 int ibest = -1;
132 int nii, njj, nbstii = 0, nbstjj = 0;
133 real ediss = 400;
135 /* Do a best match search for dissociation energies */
136 for (i = 0; (i < n2m); i++)
138 /* Check for a perfect match */
139 if (((gmx_strcasecmp(t2m[i].ai, ai) == 0) && (gmx_strcasecmp(t2m[i].aj, aj) == 0)) ||
140 ((gmx_strcasecmp(t2m[i].aj, ai) == 0) && (gmx_strcasecmp(t2m[i].ai, aj) == 0)))
142 ibest = i;
143 break;
145 else
147 /* Otherwise count the number of equal characters in the strings ai and aj
148 * and the ones from the file
150 nii = nequal(t2m[i].ai, ai);
151 njj = nequal(t2m[i].aj, aj);
152 if (((nii > nbstii) && (njj >= nbstjj)) ||
153 ((nii >= nbstii) && (njj > nbstjj)))
155 if ((nii > 0) && (njj > 0))
157 ibest = i;
158 nbstii = nii;
159 nbstjj = njj;
162 else
164 /* Swap ai and aj (at least in counting the number of equal chars) */
165 nii = nequal(t2m[i].ai, aj);
166 njj = nequal(t2m[i].aj, ai);
167 if (((nii > nbstii) && (njj >= nbstjj)) ||
168 ((nii >= nbstii) && (njj > nbstjj)))
170 if ((nii > 0) && (njj > 0))
172 ibest = i;
173 nbstii = nii;
174 nbstjj = njj;
180 /* Return the dissocation energy corresponding to the best match, if we have
181 * found one. Do some debug output anyway.
183 if (ibest == -1)
185 if (debug)
187 fprintf(debug, "MORSE: Couldn't find E_diss for bond %s - %s, using default %g\n", ai, aj, ediss);
189 return ediss;
191 else
193 if (debug)
195 fprintf(debug, "MORSE: Dissoc. E (%10.3f) for bond %4s-%4s taken from bond %4s-%4s\n",
196 t2m[ibest].e_diss, ai, aj, t2m[ibest].ai, t2m[ibest].aj);
198 return t2m[ibest].e_diss;
202 void convert_harmonics(int nrmols, t_molinfo mols[], gpp_atomtype_t atype)
204 int n2m;
205 t_2morse *t2m;
207 int i, j, k, last, ni, nj;
208 int nrharm, nrmorse, bb;
209 real edis, kb, b0, beta;
210 gmx_bool *bRemoveHarm;
212 /* First get the data */
213 t2m = read_dissociation_energies(&n2m);
214 if (debug)
216 fprintf(debug, "MORSE: read %d dissoc energies\n", n2m);
218 if (n2m <= 0)
220 fprintf(stderr, "No dissocation energies read\n");
221 return;
224 /* For all the molecule types */
225 for (i = 0; (i < nrmols); i++)
227 /* Check how many morse and harmonic BONDSs there are, increase size of
228 * morse with the number of harmonics
230 nrmorse = mols[i].plist[F_MORSE].nr;
232 for (bb = 0; (bb < F_NRE); bb++)
234 if ((interaction_function[bb].flags & IF_BTYPE) && (bb != F_MORSE))
236 nrharm = mols[i].plist[bb].nr;
237 pr_alloc(nrharm, &(mols[i].plist[F_MORSE]));
238 snew(bRemoveHarm, nrharm);
240 /* Now loop over the harmonics, trying to convert them */
241 for (j = 0; (j < nrharm); j++)
243 ni = mols[i].plist[bb].param[j].ai();
244 nj = mols[i].plist[bb].param[j].aj();
245 edis =
246 search_e_diss(n2m, t2m,
247 get_atomtype_name(mols[i].atoms.atom[ni].type, atype),
248 get_atomtype_name(mols[i].atoms.atom[nj].type, atype));
249 if (edis != 0)
251 bRemoveHarm[j] = TRUE;
252 b0 = mols[i].plist[bb].param[j].c[0];
253 kb = mols[i].plist[bb].param[j].c[1];
254 beta = std::sqrt(kb/(2*edis));
255 mols[i].plist[F_MORSE].param[nrmorse].a[0] = ni;
256 mols[i].plist[F_MORSE].param[nrmorse].a[1] = nj;
257 mols[i].plist[F_MORSE].param[nrmorse].c[0] = b0;
258 mols[i].plist[F_MORSE].param[nrmorse].c[1] = edis;
259 mols[i].plist[F_MORSE].param[nrmorse].c[2] = beta;
260 nrmorse++;
263 mols[i].plist[F_MORSE].nr = nrmorse;
265 /* Now remove the harmonics */
266 for (j = last = 0; (j < nrharm); j++)
268 if (!bRemoveHarm[j])
270 /* Copy it to the last position */
271 for (k = 0; (k < MAXATOMLIST); k++)
273 mols[i].plist[bb].param[last].a[k] =
274 mols[i].plist[bb].param[j].a[k];
276 for (k = 0; (k < MAXFORCEPARAM); k++)
278 mols[i].plist[bb].param[last].c[k] =
279 mols[i].plist[bb].param[j].c[k];
281 last++;
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;
291 sfree(t2m);