minor fixes in ditribution files
[gromacs/qmmm-gamess-us.git] / src / kernel / h_db.c
blobf82ae3a547d74722c77223d03eb1392ce1d2e47c
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
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"
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
57 const int ncontrol[] = { -1, 3, 3, 3, 3, 4, 3, 1, 3, 3, 1, 1 };
58 #define maxcontrol asize(ncontrol)
60 int compaddh(const void *a,const void *b)
62 t_hackblock *ah,*bh;
64 ah=(t_hackblock *)a;
65 bh=(t_hackblock *)b;
66 return strcasecmp(ah->name,bh->name);
69 void print_ab(FILE *out,t_hack *hack,char *nname)
71 int i;
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");
80 void read_ab(char *line,const char *fn,t_hack *hack)
82 int i,nh,tp,ns;
83 char a[4][12];
84 char hn[32];
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);
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);
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]);
101 for( ; i<4; i++) {
102 hack->a[i]=NULL;
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;
113 static void dump_h_db(const char *fn,int nah,t_hackblock *ah)
115 FILE *fp;
116 char buf[STRLEN],nname[STRLEN];
117 int i,j,k;
119 sprintf(buf,"%s_new.hdb",fn);
120 fp = gmx_fio_fopen(buf,"w");
121 for(i=0; (i<nah); i++) {
122 fprintf(fp,"%-8s%-8d\n",ah[i].name,ah[i].nhack);
123 for(k=0; (k<ah[i].nhack); k++) {
124 strcpy(nname,ah[i].hack[k].a[0]);
125 nname[0] = 'H';
126 print_ab(fp,&ah[i].hack[k],nname);
129 gmx_fio_fclose(fp);
132 static void read_h_db_file(const char *hfn,int *nahptr,t_hackblock **ah)
134 FILE *in;
135 char filebase[STRLEN],line[STRLEN], buf[STRLEN];
136 int i, n, nab, nah;
137 t_hackblock *aah;
139 if (debug) fprintf(debug,"Hydrogen Database (%s):\n",hfn);
141 fflib_filename_base(hfn,filebase,STRLEN);
142 /* Currently filebase is read and set, but not used.
143 * hdb entries from any hdb file and be applied to rtp entries
144 * in any rtp file.
147 in = fflib_open(hfn);
149 nah = *nahptr;
150 aah = *ah;
151 while (fgets2(line,STRLEN-1,in)) {
152 if (sscanf(line,"%s%n",buf,&n) != 1) {
153 fprintf(stderr,"Error in hdb file: nah = %d\nline = '%s'\n",
154 nah,line);
155 break;
157 if (debug) fprintf(debug,"%s",buf);
158 srenew(aah,nah+1);
159 clear_t_hackblock(&aah[nah]);
160 aah[nah].name = strdup(buf);
161 aah[nah].filebase = strdup(filebase);
163 if (sscanf(line+n,"%d",&nab) == 1) {
164 if (debug) fprintf(debug," %d\n",nab);
165 snew(aah[nah].hack,nab);
166 aah[nah].nhack = nab;
167 for(i=0; (i<nab); i++) {
168 if (feof(in))
169 gmx_fatal(FARGS, "Expected %d lines of hydrogens, found only %d "
170 "while reading Hydrogen Database %s residue %s",
171 nab, i-1, aah[nah].name, hfn);
172 if(NULL==fgets(buf, STRLEN, in))
174 gmx_fatal(FARGS,"Error reading from file %s",hfn);
176 read_ab(buf,hfn,&(aah[nah].hack[i]));
179 nah++;
181 ffclose(in);
183 /* Sort the list (necessary to be able to use bsearch */
184 qsort(aah,nah,(size_t)sizeof(**ah),compaddh);
186 if (debug)
187 dump_h_db(hfn,nah,aah);
189 *nahptr = nah;
190 *ah = aah;
193 int read_h_db(const char *ffdir,t_hackblock **ah)
195 int nhdbf,f;
196 char **hdbf;
197 int nah;
198 FILE *fp;
200 /* Read the hydrogen database file(s).
201 * Do not generate an error when no files are found.
203 nhdbf = fflib_search_file_end(ffdir,".hdb",FALSE,&hdbf);
204 nah = 0;
205 *ah = NULL;
206 for(f=0; f<nhdbf; f++) {
207 read_h_db_file(hdbf[f],&nah,ah);
208 sfree(hdbf[f]);
210 sfree(hdbf);
212 return nah;
215 t_hackblock *search_h_db(int nh,t_hackblock ah[],char *key)
217 t_hackblock ahkey,*result;
219 if (nh <= 0)
220 return NULL;
222 ahkey.name=key;
224 result=(t_hackblock *)bsearch(&ahkey,ah,nh,(size_t)sizeof(ah[0]),compaddh);
226 return result;