38f1d7a742fd264e4093c260ca97e75885f17310
[alexxy/gromacs.git] / src / kernel / h_db.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  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <string.h>
41 #include "string2.h"
42 #include "sysstuff.h"
43 #include "smalloc.h"
44 #include "futil.h"
45 #include "symtab.h"
46 #include "h_db.h"
47 #include "gmxfio.h"
48 #include "fflibutil.h"
49 #include "gmx_fatal.h"
50
51 /* There are 11 types of adding hydrogens, numbered from
52  * 1 thru 11. Each of these has a specific number of
53  * control atoms, that determine how the hydrogens are added.
54  * Here these number are given. Because arrays start at 0 an
55  * extra dummy for index 0 is added 
56  */
57 const int ncontrol[] = { -1, 3, 3, 3, 3, 4, 3, 1, 3, 3, 1, 1 };
58 #define maxcontrol asize(ncontrol)
59
60 int compaddh(const void *a,const void *b)
61 {
62   t_hackblock *ah,*bh;
63
64   ah=(t_hackblock *)a;
65   bh=(t_hackblock *)b;
66   return gmx_strcasecmp(ah->name,bh->name);
67 }
68
69 void print_ab(FILE *out,t_hack *hack,char *nname)
70 {
71   int i;
72
73   fprintf(out,"%d\t%d\t%s",hack->nr,hack->tp,nname);
74   for(i=0; (i < hack->nctl); i++)
75     fprintf(out,"\t%s",hack->a[i]);
76   fprintf(out,"\n");
77 }
78
79
80 void read_ab(char *line,const char *fn,t_hack *hack)
81 {
82   int  i,nh,tp,ns;
83   char a[4][12];
84   char hn[32];
85   
86   ns = sscanf(line,"%d%d%s%s%s%s%s",&nh,&tp,hn,a[0],a[1],a[2],a[3]);
87   if (ns < 4)
88     gmx_fatal(FARGS,"wrong format in input file %s on line\n%s\n",fn,line);
89   
90   hack->nr=nh;
91   hack->tp=tp;
92   if ((tp < 1) || (tp >= maxcontrol)) 
93     gmx_fatal(FARGS,"Error in hdb file %s:\nH-type should be in 1-%d. Offending line:\n%s",fn,maxcontrol-1,line);
94   
95   hack->nctl = ns - 3;
96   if ((hack->nctl != ncontrol[hack->tp]) && (ncontrol[hack->tp] != -1))
97     gmx_fatal(FARGS,"Error in hdb file %s:\nWrong number of control atoms (%d iso %d) on line:\n%s\n",fn,hack->nctl,ncontrol[hack->tp],line);
98   for(i=0; (i<hack->nctl); i++) {
99     hack->a[i]=strdup(a[i]);
100   }
101   for(   ; i<4; i++) {
102     hack->a[i]=NULL;
103   }
104   hack->oname=NULL;
105   hack->nname=strdup(hn);
106   hack->atom=NULL;
107   hack->cgnr=NOTSET;
108   hack->bXSet=FALSE;
109   for(i=0; i<DIM; i++)
110     hack->newx[i]=NOTSET;
111 }
112
113 static void read_h_db_file(const char *hfn,int *nahptr,t_hackblock **ah)
114 {       
115   FILE   *in;
116   char   filebase[STRLEN],line[STRLEN], buf[STRLEN];
117   int    i, n, nab, nah;
118   t_hackblock *aah;
119
120   if (debug) fprintf(debug,"Hydrogen Database (%s):\n",hfn);
121
122   fflib_filename_base(hfn,filebase,STRLEN);
123   /* Currently filebase is read and set, but not used.
124    * hdb entries from any hdb file and be applied to rtp entries
125    * in any rtp file.
126    */
127
128   in = fflib_open(hfn);
129
130   nah = *nahptr;
131   aah = *ah;
132   while (fgets2(line,STRLEN-1,in)) {
133     if (sscanf(line,"%s%n",buf,&n) != 1) {
134       fprintf(stderr,"Error in hdb file: nah = %d\nline = '%s'\n",
135               nah,line);
136       break;
137     }
138     if (debug) fprintf(debug,"%s",buf);
139     srenew(aah,nah+1);
140     clear_t_hackblock(&aah[nah]);
141     aah[nah].name     = strdup(buf);
142     aah[nah].filebase = strdup(filebase);
143     
144     if (sscanf(line+n,"%d",&nab) == 1) {
145       if (debug) fprintf(debug,"  %d\n",nab);
146       snew(aah[nah].hack,nab);
147       aah[nah].nhack = nab;
148       for(i=0; (i<nab); i++) {
149         if (feof(in))
150           gmx_fatal(FARGS, "Expected %d lines of hydrogens, found only %d "
151                       "while reading Hydrogen Database %s residue %s",
152                       nab, i-1, aah[nah].name, hfn);
153         if(NULL==fgets(buf, STRLEN, in))
154         {
155           gmx_fatal(FARGS,"Error reading from file %s",hfn);
156         }
157         read_ab(buf,hfn,&(aah[nah].hack[i]));
158       }
159     }
160     nah++;
161   }
162   ffclose(in);
163   
164   /* Sort the list (necessary to be able to use bsearch */
165   qsort(aah,nah,(size_t)sizeof(**ah),compaddh);
166
167   /*
168   if (debug)
169     dump_h_db(hfn,nah,aah);
170   */
171   
172   *nahptr = nah;
173   *ah     = aah;
174 }
175
176 int read_h_db(const char *ffdir,t_hackblock **ah)
177 {
178   int  nhdbf,f;
179   char **hdbf;
180   int  nah;
181   FILE *fp;
182
183   /* Read the hydrogen database file(s).
184    * Do not generate an error when no files are found.
185    */
186   nhdbf = fflib_search_file_end(ffdir,".hdb",FALSE,&hdbf);
187   nah = 0;
188   *ah = NULL;
189   for(f=0; f<nhdbf; f++) {
190     read_h_db_file(hdbf[f],&nah,ah);
191     sfree(hdbf[f]);
192   }
193   sfree(hdbf);
194
195   return nah;
196 }
197
198 t_hackblock *search_h_db(int nh,t_hackblock ah[],char *key)
199 {
200   t_hackblock ahkey,*result;
201
202   if (nh <= 0)
203     return NULL;
204   
205   ahkey.name=key;
206
207   result=(t_hackblock *)bsearch(&ahkey,ah,nh,(size_t)sizeof(ah[0]),compaddh);
208   
209   return result;
210 }