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.
39 #include "gromacs/topology/atomprop.h"
46 #include "gromacs/fileio/strdb.h"
47 #include "gromacs/legacyheaders/copyrite.h"
48 #include "gromacs/math/utilities.h"
49 #include "gromacs/topology/residuetypes.h"
50 #include "gromacs/utility/cstringutil.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/futil.h"
53 #include "gromacs/utility/smalloc.h"
66 typedef struct gmx_atomprop {
67 gmx_bool bWarned, bWarnVDW;
68 aprop_t prop[epropNR];
69 gmx_residuetype_t *restype;
74 /* NOTFOUND should be smallest, others larger in increasing priority */
76 NOTFOUND = -4, WILDCARD, WILDPROT, PROTEIN
79 /* return number of matching characters,
80 or NOTFOUND if not at least all characters in char *database match */
81 static int dbcmp_len(char *search, char *database)
86 while (search[i] && database[i] && (search[i] == database[i]) )
98 static int get_prop_index(aprop_t *ap, gmx_residuetype_t *restype,
99 char *resnm, char *atomnm,
104 long int malen, mrlen;
105 gmx_bool bProtein, bProtWild;
107 bProtein = gmx_residuetype_is_protein(restype, resnm);
108 bProtWild = (strcmp(resnm, "AAA") == 0);
111 for (i = 0; (i < ap->nprop); i++)
113 rlen = dbcmp_len(resnm, ap->resnm[i]);
114 if (rlen == NOTFOUND)
116 if ( (strcmp(ap->resnm[i], "*") == 0) ||
117 (strcmp(ap->resnm[i], "???") == 0) )
121 else if (strcmp(ap->resnm[i], "AAA") == 0)
126 alen = dbcmp_len(atomnm, ap->atomnm[i]);
127 if ( (alen > NOTFOUND) && (rlen > NOTFOUND))
129 if ( ( (alen > malen) && (rlen >= mrlen)) ||
130 ( (rlen > mrlen) && (alen >= malen) ) )
139 *bExact = ((malen == (long int)strlen(atomnm)) &&
140 ((mrlen == (long int)strlen(resnm)) ||
141 ((mrlen == WILDPROT) && bProtWild) ||
142 ((mrlen == WILDCARD) && !bProtein && !bProtWild)));
146 fprintf(debug, "searching residue: %4s atom: %4s\n", resnm, atomnm);
149 fprintf(debug, " not successful\n");
153 fprintf(debug, " match: %4s %4s\n", ap->resnm[j], ap->atomnm[j]);
159 static void add_prop(aprop_t *ap, gmx_residuetype_t *restype,
160 char *resnm, char *atomnm,
166 j = get_prop_index(ap, restype, resnm, atomnm, &bExact);
170 if (ap->nprop >= ap->maxprop)
173 srenew(ap->resnm, ap->maxprop);
174 srenew(ap->atomnm, ap->maxprop);
175 srenew(ap->value, ap->maxprop);
176 srenew(ap->bAvail, ap->maxprop);
177 for (i = ap->nprop; (i < ap->maxprop); i++)
179 ap->atomnm[i] = NULL;
182 ap->bAvail[i] = FALSE;
185 ap->atomnm[ap->nprop] = gmx_strdup(atomnm);
186 ap->resnm[ap->nprop] = gmx_strdup(resnm);
192 if (ap->value[j] == p)
194 fprintf(stderr, "Warning double identical entries for %s %s %g on line %d in file %s\n",
195 resnm, atomnm, p, line, ap->db);
199 fprintf(stderr, "Warning double different entries %s %s %g and %g on line %d in file %s\n"
200 "Using last entry (%g)\n",
201 resnm, atomnm, p, ap->value[j], line, ap->db, p);
207 ap->bAvail[j] = TRUE;
212 static void read_prop(gmx_atomprop_t aps, int eprop, double factor)
214 gmx_atomprop *ap2 = (gmx_atomprop*) aps;
216 char line[STRLEN], resnm[32], atomnm[32];
221 ap = &ap2->prop[eprop];
223 fp = libopen(ap->db);
225 while (get_a_line(fp, line, STRLEN))
228 if (sscanf(line, "%31s %31s %20lf", resnm, atomnm, &pp) == 3)
231 add_prop(ap, aps->restype, resnm, atomnm, pp, line_no);
235 fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n",
240 /* for libraries we can use the low-level close routines */
246 static void set_prop(gmx_atomprop_t aps, int eprop)
248 gmx_atomprop *ap2 = (gmx_atomprop*) aps;
249 const char *fns[epropNR] = { "atommass.dat", "vdwradii.dat", "dgsolv.dat", "electroneg.dat", "elements.dat" };
250 double fac[epropNR] = { 1.0, 1.0, 418.4, 1.0, 1.0 };
251 double def[epropNR] = { 12.011, 0.14, 0.0, 2.2, -1 };
254 ap = &ap2->prop[eprop];
257 ap->db = gmx_strdup(fns[eprop]);
258 ap->def = def[eprop];
259 read_prop(aps, eprop, fac[eprop]);
263 fprintf(debug, "Entries in %s: %d\n", ap->db, ap->nprop);
266 if ( ( (!aps->bWarned) && (eprop == epropMass) ) || (eprop == epropVDW))
269 "WARNING: Masses and atomic (Van der Waals) radii will be guessed\n"
270 " based on residue and atom names, since they could not be\n"
271 " definitively assigned from the information in your input\n"
272 " files. These guessed numbers might deviate from the mass\n"
273 " and radius of the atom type. Please check the output\n"
274 " files if necessary.\n\n");
280 gmx_atomprop_t gmx_atomprop_init(void)
286 gmx_residuetype_init(&aps->restype);
287 aps->bWarned = FALSE;
288 aps->bWarnVDW = FALSE;
290 return (gmx_atomprop_t)aps;
293 static void destroy_prop(aprop_t *ap)
301 for (i = 0; i < ap->nprop; i++)
303 sfree(ap->atomnm[i]);
313 void gmx_atomprop_destroy(gmx_atomprop_t aps)
315 gmx_atomprop *ap = (gmx_atomprop*) aps;
320 printf("\nWARNING: gmx_atomprop_destroy called with a NULL pointer\n\n");
324 for (p = 0; p < epropNR; p++)
326 destroy_prop(&ap->prop[p]);
329 gmx_residuetype_destroy(ap->restype);
334 static void vdw_warning(FILE *fp)
338 fprintf(fp, "NOTE: From version 5.0 %s uses the Van der Waals radii\n",
340 fprintf(fp, "from the source below. This means the results may be different\n");
341 fprintf(fp, "compared to previous GROMACS versions.\n");
342 please_cite(fp, "Bondi1964a");
346 gmx_bool gmx_atomprop_query(gmx_atomprop_t aps,
347 int eprop, const char *resnm, const char *atomnm,
350 gmx_atomprop *ap = (gmx_atomprop*) aps;
353 char atomname[MAXQ], resname[MAXQ];
356 set_prop(aps, eprop);
357 if ((strlen(atomnm) > MAXQ-1) || (strlen(resnm) > MAXQ-1))
361 fprintf(debug, "WARNING: will only compare first %d characters\n",
365 if (isdigit(atomnm[0]))
368 /* put digit after atomname */
369 for (i = 1; i < MAXQ-1 && atomnm[i] != '\0'; i++)
371 atomname[i-1] = atomnm[i];
373 atomname[i-1] = atomnm[0];
378 strncpy(atomname, atomnm, MAXQ-1);
380 strncpy(resname, resnm, MAXQ-1);
382 j = get_prop_index(&(ap->prop[eprop]), ap->restype, resname,
385 if (eprop == epropVDW && !ap->bWarnVDW)
392 *value = ap->prop[eprop].value[j];
397 *value = ap->prop[eprop].def;
402 char *gmx_atomprop_element(gmx_atomprop_t aps, int atomnumber)
404 gmx_atomprop *ap = (gmx_atomprop*) aps;
407 set_prop(aps, epropElement);
408 for (i = 0; (i < ap->prop[epropElement].nprop); i++)
410 if (gmx_nint(ap->prop[epropElement].value[i]) == atomnumber)
412 return ap->prop[epropElement].atomnm[i];
418 int gmx_atomprop_atomnumber(gmx_atomprop_t aps, const char *elem)
420 gmx_atomprop *ap = (gmx_atomprop*) aps;
423 set_prop(aps, epropElement);
424 for (i = 0; (i < ap->prop[epropElement].nprop); i++)
426 if (gmx_strcasecmp(ap->prop[epropElement].atomnm[i], elem) == 0)
428 return gmx_nint(ap->prop[epropElement].value[i]);