Partial commit of the project to remove all static variables.
[gromacs.git] / src / kernel / tomorse.c
blobf2bd2c39cbd15e0fc2a24c44fd89d7222fd58051
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.1
11 * Copyright (c) 1991-2001, University of Groningen, The Netherlands
12 * This program is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU General Public License
14 * as published by the Free Software Foundation; either version 2
15 * of the License, or (at your option) any later version.
17 * If you want to redistribute modifications, please consider that
18 * scientific software is very special. Version control is crucial -
19 * bugs must be traceable. We will be happy to consider code for
20 * inclusion in the official distribution, but derived work must not
21 * be called official GROMACS. Details are found in the README & COPYING
22 * files - if they are missing, get the official version at www.gromacs.org.
24 * To help us fund GROMACS development, we humbly ask that you cite
25 * the papers on the package - you can find them in the top README file.
27 * For more info, check our website at http://www.gromacs.org
29 * And Hey:
30 * Gromacs Runs One Microsecond At Cannonball Speeds
33 /* This file is completely threadsafe - keep it that way! */
34 #include <stdlib.h>
35 #include <math.h>
36 #include <ctype.h>
37 #include "typedefs.h"
38 #include "string2.h"
39 #include "grompp.h"
40 #include "futil.h"
41 #include "smalloc.h"
42 #include "fatal.h"
44 typedef struct {
45 char *ai,*aj;
46 real e_diss;
47 } t_2morse;
49 static t_2morse *read_dissociation_energies(int *n2morse)
51 FILE *fp;
52 char ai[32],aj[32];
53 double e_diss;
54 char *fn="edissoc.dat";
55 t_2morse *t2m=NULL;
56 int maxn2m=0,n2m=0;
57 int nread;
59 /* Open the file with dissociation energies */
60 fp = libopen(fn);
61 do {
62 /* Try and read two atom names and an energy term from it */
63 nread = fscanf(fp,"%s%s%lf",ai,aj,&e_diss);
64 if (nread == 3) {
65 /* If we got three terms, it probably was OK, no further checking */
66 if (n2m >= maxn2m) {
67 /* Increase memory for 16 at once, some mallocs are stupid */
68 maxn2m += 16;
69 srenew(t2m,maxn2m);
71 /* Copy the values */
72 t2m[n2m].ai = strdup(ai);
73 t2m[n2m].aj = strdup(aj);
74 t2m[n2m].e_diss = e_diss;
75 /* Increment counter */
76 n2m++;
78 /* If we did not read three items, quit reading */
79 } while (nread == 3);
80 ffclose(fp);
82 /* Set the return values */
83 *n2morse = n2m;
85 return t2m;
88 static int nequal(char *a1,char *a2)
90 int i;
92 /* Count the number of (case insensitive) characters that are equal in
93 * two strings. If they are equally long their respective null characters are
94 * counted also.
96 for(i=0; (a1[i] != '\0') && (a2[i] != '\0'); i++)
97 if (toupper(a1[i]) != toupper(a2[i]))
98 break;
99 if ((a1[i] == '\0') && (a2[i] == '\0'))
100 i++;
102 return i;
105 static real search_e_diss(int n2m,t_2morse t2m[],char *ai,char *aj)
107 int i;
108 int ibest=-1;
109 int nii,njj,nbstii=0,nbstjj=0;
110 real ediss = 400;
112 /* Do a best match search for dissociation energies */
113 for(i=0; (i<n2m); i++) {
114 /* Check for a perfect match */
115 if (((strcasecmp(t2m[i].ai,ai) == 0) && (strcasecmp(t2m[i].aj,aj) == 0)) ||
116 ((strcasecmp(t2m[i].aj,ai) == 0) && (strcasecmp(t2m[i].ai,aj) == 0))) {
117 ibest = i;
118 break;
120 else {
121 /* Otherwise count the number of equal characters in the strings ai and aj
122 * and the ones from the file
124 nii = nequal(t2m[i].ai,ai);
125 njj = nequal(t2m[i].aj,aj);
126 if (((nii > nbstii) && (njj >= nbstjj)) ||
127 ((nii >= nbstii) && (njj > nbstjj))) {
128 if ((nii > 0) && (njj > 0)) {
129 ibest = i;
130 nbstii = nii;
131 nbstjj = njj;
134 else {
135 /* Swap ai and aj (at least in counting the number of equal chars) */
136 nii = nequal(t2m[i].ai,aj);
137 njj = nequal(t2m[i].aj,ai);
138 if (((nii > nbstii) && (njj >= nbstjj)) ||
139 ((nii >= nbstii) && (njj > nbstjj))) {
140 if ((nii > 0) && (njj > 0)) {
141 ibest = i;
142 nbstii = nii;
143 nbstjj = njj;
149 /* Return the dissocation energy corresponding to the best match, if we have
150 * found one. Do some debug output anyway.
152 if (ibest == -1) {
153 if (debug)
154 fprintf(debug,"MORSE: Couldn't find E_diss for bond %s - %s, using default %g\n",ai,aj,ediss);
155 return ediss;
157 else {
158 if (debug)
159 fprintf(debug,"MORSE: Dissoc. E (%10.3f) for bond %4s-%4s taken from bond %4s-%4s\n",
160 t2m[ibest].e_diss,ai,aj,t2m[ibest].ai,t2m[ibest].aj);
161 return t2m[ibest].e_diss;
165 void convert_harmonics(int nrmols,t_molinfo mols[],t_atomtype *atype)
167 int n2m;
168 t_2morse *t2m;
170 int i,j,k,last,ni,nj;
171 int nrharm,nrmorse,bb;
172 real edis,kb,b0,beta;
173 bool *bRemoveHarm;
175 /* First get the data */
176 t2m = read_dissociation_energies(&n2m);
177 if (debug)
178 fprintf(debug,"MORSE: read %d dissoc energies\n",n2m);
179 if (n2m <= 0) {
180 fprintf(stderr,"No dissocation energies read\n");
181 return;
184 /* For all the molecule types */
185 for(i=0; (i<nrmols); i++) {
186 /* Check how many morse and harmonic BONDSs there are, increase size of
187 * morse with the number of harmonics
189 nrmorse = mols[i].plist[F_MORSE].nr;
191 for(bb=0; (bb < F_NRE); bb++) {
192 if ((interaction_function[bb].flags & IF_BTYPE) && (bb != F_MORSE)) {
193 nrharm = mols[i].plist[bb].nr;
194 srenew(mols[i].plist[F_MORSE].param,nrmorse+nrharm);
195 snew(bRemoveHarm,nrharm);
197 /* Now loop over the harmonics, trying to convert them */
198 for(j=0; (j<nrharm); j++) {
199 ni = mols[i].plist[bb].param[j].AI;
200 nj = mols[i].plist[bb].param[j].AJ;
201 edis = search_e_diss(n2m,t2m,
202 *atype->atomname[mols[i].atoms.atom[ni].type],
203 *atype->atomname[mols[i].atoms.atom[nj].type]);
204 if (edis != 0) {
205 bRemoveHarm[j] = TRUE;
206 b0 = mols[i].plist[bb].param[j].c[0];
207 kb = mols[i].plist[bb].param[j].c[1];
208 beta = sqrt(kb/(2*edis));
209 mols[i].plist[F_MORSE].param[nrmorse].a[0] = ni;
210 mols[i].plist[F_MORSE].param[nrmorse].a[1] = nj;
211 mols[i].plist[F_MORSE].param[nrmorse].c[0] = b0;
212 mols[i].plist[F_MORSE].param[nrmorse].c[1] = edis;
213 mols[i].plist[F_MORSE].param[nrmorse].c[2] = beta;
214 nrmorse++;
217 mols[i].plist[F_MORSE].nr = nrmorse;
219 /* Now remove the harmonics */
220 for(j=last=0; (j<nrharm); j++) {
221 if (!bRemoveHarm[j]) {
222 /* Copy it to the last position */
223 for(k=0; (k<MAXATOMLIST); k++)
224 mols[i].plist[bb].param[last].a[k] =
225 mols[i].plist[bb].param[j].a[k];
226 for(k=0; (k<MAXFORCEPARAM); k++)
227 mols[i].plist[bb].param[last].c[k] =
228 mols[i].plist[bb].param[j].c[k];
229 last++;
232 sfree(bRemoveHarm);
233 fprintf(stderr,"Converted %d out of %d %s to morse bonds for mol %d\n",
234 nrharm-last,nrharm,interaction_function[bb].name,i);
235 mols[i].plist[bb].nr = last;
239 sfree(t2m);