21bb2a94b76aa87da48b09002dfcbddeea0ddb34
[alexxy/gromacs.git] / src / gromacs / gmxlib / atomprop.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <ctype.h>
40 #include "sysstuff.h"
41 #include "smalloc.h"
42 #include "string2.h"
43 #include "futil.h"
44 #include "maths.h"
45 #include "gmx_fatal.h"
46 #include "atomprop.h"
47 #include "maths.h"
48 #include "macros.h"
49 #include "index.h"
50 #include "strdb.h"
51 #include "copyrite.h"
52 #include "statutil.h"
53
54 typedef struct {
55     gmx_bool    bSet;
56     int         nprop, maxprop;
57     char       *db;
58     double      def;
59     char      **atomnm;
60     char      **resnm;
61     gmx_bool   *bAvail;
62     real       *value;
63 } aprop_t;
64
65 typedef struct gmx_atomprop {
66     gmx_bool          bWarned, bWarnVDW;
67     aprop_t           prop[epropNR];
68     gmx_residuetype_t restype;
69 } gmx_atomprop;
70
71
72
73 /* NOTFOUND should be smallest, others larger in increasing priority */
74 enum {
75     NOTFOUND = -4, WILDCARD, WILDPROT, PROTEIN
76 };
77
78 /* return number of matching characters,
79    or NOTFOUND if not at least all characters in char *database match */
80 static int dbcmp_len(char *search, char *database)
81 {
82     int i;
83
84     i = 0;
85     while (search[i] && database[i] && (search[i] == database[i]) )
86     {
87         i++;
88     }
89
90     if (database[i])
91     {
92         i = NOTFOUND;
93     }
94     return i;
95 }
96
97 static int get_prop_index(aprop_t *ap, gmx_residuetype_t restype,
98                           char *resnm, char *atomnm,
99                           gmx_bool *bExact)
100 {
101     int      i, j = NOTFOUND;
102     long int alen, rlen;
103     long int malen, mrlen;
104     gmx_bool bProtein, bProtWild;
105
106     bProtein  = gmx_residuetype_is_protein(restype, resnm);
107     bProtWild = (strcmp(resnm, "AAA") == 0);
108     malen     = NOTFOUND;
109     mrlen     = NOTFOUND;
110     for (i = 0; (i < ap->nprop); i++)
111     {
112         rlen = dbcmp_len(resnm, ap->resnm[i]);
113         if (rlen == NOTFOUND)
114         {
115             if ( (strcmp(ap->resnm[i], "*") == 0) ||
116                  (strcmp(ap->resnm[i], "???") == 0) )
117             {
118                 rlen = WILDCARD;
119             }
120             else if (strcmp(ap->resnm[i], "AAA") == 0)
121             {
122                 rlen = WILDPROT;
123             }
124         }
125         alen = dbcmp_len(atomnm, ap->atomnm[i]);
126         if ( (alen > NOTFOUND) && (rlen > NOTFOUND))
127         {
128             if ( ( (alen > malen) && (rlen >= mrlen)) ||
129                  ( (rlen > mrlen) && (alen >= malen) ) )
130             {
131                 malen = alen;
132                 mrlen = rlen;
133                 j     = i;
134             }
135         }
136     }
137
138     *bExact = ((malen == (long int)strlen(atomnm)) &&
139                ((mrlen == (long int)strlen(resnm)) ||
140                 ((mrlen == WILDPROT) && bProtWild) ||
141                 ((mrlen == WILDCARD) && !bProtein && !bProtWild)));
142
143     if (debug)
144     {
145         fprintf(debug, "searching residue: %4s atom: %4s\n", resnm, atomnm);
146         if (j == NOTFOUND)
147         {
148             fprintf(debug, " not successful\n");
149         }
150         else
151         {
152             fprintf(debug, " match: %4s %4s\n", ap->resnm[j], ap->atomnm[j]);
153         }
154     }
155     return j;
156 }
157
158 static void add_prop(aprop_t *ap, gmx_residuetype_t restype,
159                      char *resnm, char *atomnm,
160                      real p, int line)
161 {
162     int      i, j;
163     gmx_bool bExact;
164
165     j = get_prop_index(ap, restype, resnm, atomnm, &bExact);
166
167     if (!bExact)
168     {
169         if (ap->nprop >= ap->maxprop)
170         {
171             ap->maxprop += 10;
172             srenew(ap->resnm, ap->maxprop);
173             srenew(ap->atomnm, ap->maxprop);
174             srenew(ap->value, ap->maxprop);
175             srenew(ap->bAvail, ap->maxprop);
176             for (i = ap->nprop; (i < ap->maxprop); i++)
177             {
178                 ap->atomnm[i] = NULL;
179                 ap->resnm[i]  = NULL;
180                 ap->value[i]  = 0;
181                 ap->bAvail[i] = FALSE;
182             }
183         }
184         ap->atomnm[ap->nprop] = strdup(atomnm);
185         ap->resnm[ap->nprop]  = strdup(resnm);
186         j                     = ap->nprop;
187         ap->nprop++;
188     }
189     if (ap->bAvail[j])
190     {
191         if (ap->value[j] == p)
192         {
193             fprintf(stderr, "Warning double identical entries for %s %s %g on line %d in file %s\n",
194                     resnm, atomnm, p, line, ap->db);
195         }
196         else
197         {
198             fprintf(stderr, "Warning double different entries %s %s %g and %g on line %d in file %s\n"
199                     "Using last entry (%g)\n",
200                     resnm, atomnm, p, ap->value[j], line, ap->db, p);
201             ap->value[j] = p;
202         }
203     }
204     else
205     {
206         ap->bAvail[j] = TRUE;
207         ap->value[j]  = p;
208     }
209 }
210
211 static void read_prop(gmx_atomprop_t aps, int eprop, double factor)
212 {
213     gmx_atomprop *ap2 = (gmx_atomprop*) aps;
214     FILE         *fp;
215     char          line[STRLEN], resnm[32], atomnm[32];
216     double        pp;
217     int           line_no;
218     aprop_t      *ap;
219
220     ap = &ap2->prop[eprop];
221
222     fp      = libopen(ap->db);
223     line_no = 0;
224     while (get_a_line(fp, line, STRLEN))
225     {
226         line_no++;
227         if (sscanf(line, "%s %s %lf", resnm, atomnm, &pp) == 3)
228         {
229             pp *= factor;
230             add_prop(ap, aps->restype, resnm, atomnm, pp, line_no);
231         }
232         else
233         {
234             fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n",
235                     ap->db, line_no);
236         }
237     }
238
239     /* for libraries we can use the low-level close routines */
240     ffclose(fp);
241
242     ap->bSet = TRUE;
243 }
244
245 static void set_prop(gmx_atomprop_t aps, int eprop)
246 {
247     gmx_atomprop *ap2           = (gmx_atomprop*) aps;
248     const char   *fns[epropNR]  = { "atommass.dat", "vdwradii.dat", "dgsolv.dat", "electroneg.dat", "elements.dat" };
249     double        fac[epropNR]  = { 1.0,    1.0,  418.4, 1.0, 1.0 };
250     double        def[epropNR]  = { 12.011, 0.14, 0.0, 2.2, -1 };
251     aprop_t      *ap;
252
253     ap = &ap2->prop[eprop];
254     if (!ap->bSet)
255     {
256         ap->db  = strdup(fns[eprop]);
257         ap->def = def[eprop];
258         read_prop(aps, eprop, fac[eprop]);
259
260         if (debug)
261         {
262             fprintf(debug, "Entries in %s: %d\n", ap->db, ap->nprop);
263         }
264
265         if ( ( (!aps->bWarned) && (eprop == epropMass) ) || (eprop == epropVDW))
266         {
267             printf("\n"
268                    "WARNING: Masses and atomic (Van der Waals) radii will be guessed\n"
269                    "         based on residue and atom names, since they could not be\n"
270                    "         definitively assigned from the information in your input\n"
271                    "         files. These guessed numbers might deviate from the mass\n"
272                    "         and radius of the atom type. Please check the output\n"
273                    "         files if necessary.\n\n");
274             aps->bWarned = TRUE;
275         }
276     }
277 }
278
279 gmx_atomprop_t gmx_atomprop_init(void)
280 {
281     gmx_atomprop *aps;
282     int           p;
283
284     snew(aps, 1);
285
286     gmx_residuetype_init(&aps->restype);
287     aps->bWarned  = FALSE;
288     aps->bWarnVDW = FALSE;
289
290     return (gmx_atomprop_t)aps;
291 }
292
293 static void destroy_prop(aprop_t *ap)
294 {
295     int i;
296
297     if (ap->bSet)
298     {
299         sfree(ap->db);
300
301         for (i = 0; i < ap->nprop; i++)
302         {
303             sfree(ap->atomnm[i]);
304             sfree(ap->resnm[i]);
305         }
306         sfree(ap->atomnm);
307         sfree(ap->resnm);
308         sfree(ap->bAvail);
309         sfree(ap->value);
310     }
311 }
312
313 void gmx_atomprop_destroy(gmx_atomprop_t aps)
314 {
315     gmx_atomprop *ap = (gmx_atomprop*) aps;
316     int           p;
317
318     if (aps == NULL)
319     {
320         printf("\nWARNING: gmx_atomprop_destroy called with a NULL pointer\n\n");
321         return;
322     }
323
324     for (p = 0; p < epropNR; p++)
325     {
326         destroy_prop(&ap->prop[p]);
327     }
328
329     gmx_residuetype_destroy(ap->restype);
330
331     sfree(ap);
332 }
333
334 static void vdw_warning(FILE *fp)
335 {
336     if (NULL != fp)
337     {
338         fprintf(fp, "NOTE: From version 5.0 %s uses the Van der Waals radii\n",
339                 ShortProgram());
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");
343     }
344 }
345
346 gmx_bool gmx_atomprop_query(gmx_atomprop_t aps,
347                             int eprop, const char *resnm, const char *atomnm,
348                             real *value)
349 {
350     gmx_atomprop *ap = (gmx_atomprop*) aps;
351     size_t        i;
352     int           j;
353 #define MAXQ 32
354     char          atomname[MAXQ], resname[MAXQ];
355     gmx_bool      bExact;
356
357     set_prop(aps, eprop);
358     if ((strlen(atomnm) > MAXQ-1) || (strlen(resnm) > MAXQ-1))
359     {
360         if (debug)
361         {
362             fprintf(debug, "WARNING: will only compare first %d characters\n",
363                     MAXQ-1);
364         }
365     }
366     if (isdigit(atomnm[0]))
367     {
368         /* put digit after atomname */
369         for (i = 1; (i < min(MAXQ-1, strlen(atomnm))); i++)
370         {
371             atomname[i-1] = atomnm[i];
372         }
373         atomname[i-1] = atomnm[0];
374         atomname[i]   = '\0';
375     }
376     else
377     {
378         strncpy(atomname, atomnm, MAXQ-1);
379     }
380     strncpy(resname, resnm, MAXQ-1);
381
382     j = get_prop_index(&(ap->prop[eprop]), ap->restype, resname,
383                        atomname, &bExact);
384
385     if (eprop == epropVDW && !ap->bWarnVDW)
386     {
387         vdw_warning(stdout);
388         ap->bWarnVDW = TRUE;
389     }
390     if (j >= 0)
391     {
392         *value = ap->prop[eprop].value[j];
393         return TRUE;
394     }
395     else
396     {
397         *value = ap->prop[eprop].def;
398         return FALSE;
399     }
400 }
401
402 char *gmx_atomprop_element(gmx_atomprop_t aps, int atomnumber)
403 {
404     gmx_atomprop *ap = (gmx_atomprop*) aps;
405     int           i;
406
407     set_prop(aps, epropElement);
408     for (i = 0; (i < ap->prop[epropElement].nprop); i++)
409     {
410         if (gmx_nint(ap->prop[epropElement].value[i]) == atomnumber)
411         {
412             return ap->prop[epropElement].atomnm[i];
413         }
414     }
415     return NULL;
416 }
417
418 int gmx_atomprop_atomnumber(gmx_atomprop_t aps, const char *elem)
419 {
420     gmx_atomprop *ap = (gmx_atomprop*) aps;
421     int           i;
422
423     set_prop(aps, epropElement);
424     for (i = 0; (i < ap->prop[epropElement].nprop); i++)
425     {
426         if (gmx_strcasecmp(ap->prop[epropElement].atomnm[i], elem) == 0)
427         {
428             return gmx_nint(ap->prop[epropElement].value[i]);
429         }
430     }
431     return NOTSET;
432 }