changed reading hint
[gromacs/adressmacs.git] / src / gmxlib / h_db.c
blob04e535e7ccb65b4f301005407d23bf45599c92a6
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 2.0
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
16 * Please refer to:
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
23 * or e-mail to:
24 * gromacs@chem.rug.nl
26 * And Hey:
27 * Green Red Orange Magenta Azure Cyan Skyblue
29 static char *SRCID_h_db_c = "$Id$";
31 #include <string.h>
32 #include "string2.h"
33 #include "sysstuff.h"
34 #include "smalloc.h"
35 #include "futil.h"
36 #include "symtab.h"
37 #include "h_db.h"
39 /* There are 9 types of adding hydrogens, numbered from
40 * 1 thru 9. Each of these has a specific number of
41 * control atoms, that determine how the hydrogens are added.
42 * Here these number are given. Because arrays start at 0 an
43 * extra dummy for index 0 is added
45 /* 1 2 3 4 5 6 7 8 9 */
46 int ncontrol[12] = { -1, 3, 3, 3, 3, 4, 3, 1, 3, 3 };
48 int compaddh(const void *a,const void *b)
50 t_hackblock *ah,*bh;
52 ah=(t_hackblock *)a;
53 bh=(t_hackblock *)b;
54 return strcasecmp(ah->name,bh->name);
57 void read_ab(char *line,char *fn,t_hack *hack)
59 int i,n,nh,tp,ncntl;
60 char buf[80];
62 if (sscanf(line,"%d%d%n",&nh,&tp,&n) != 2)
63 fatal_error(0,"wrong format in input file %s on line\n%s\n",fn,line);
64 line+=n;
65 hack->nr=nh;
66 hack->tp=tp;
67 ncntl=ncontrol[tp];
68 for(i=0; i<ncntl; i++) {
69 if (sscanf(line,"%s%n",buf,&n) != 1)
70 fatal_error(0,"Expected %d control atoms instead of %d while reading Hydrogen Database %s on line\n%s\n",ncntl,i-1,fn,line);
71 hack->a[i]=strdup(buf);
72 line+=n;
74 for( ; i<4; i++)
75 hack->a[i]=NULL;
76 hack->oname=NULL;
77 hack->nname=NULL;
78 hack->atom=NULL;
79 hack->cgnr=NOTSET;
80 for(i=0; i<DIM; i++)
81 hack->newx[i]=NOTSET;
84 int read_h_db(char *fn,t_hackblock **ah)
86 FILE *in;
87 char hfn[STRLEN], line[STRLEN], buf[STRLEN];
88 int i, n, nab, nah;
89 t_hackblock *aah;
91 sprintf(hfn,"%s.hdb",fn);
92 in=libopen(hfn);
93 if (debug) fprintf(debug,"Hydrogen Database (%s):\n",hfn);
94 nah=0;
95 aah=NULL;
96 while (fgets(line,STRLEN,in)) {
97 sscanf(line,"%s%n",buf,&n);
98 if (debug) fprintf(debug,"%s",buf);
99 srenew(aah,++nah);
100 aah[nah-1].name=strdup(buf);
101 sscanf(line+n,"%d",&nab);
102 if (debug) fprintf(debug,"\t%d\n",nab);
103 snew(aah[nah-1].hack,nab);
104 aah[nah-1].nhack=nab;
105 for(i=0; (i<nab); i++) {
106 if (feof(in))
107 fatal_error(0, "Expected %d lines of hydrogens, found only %d "
108 "while reading Hydrogen Database %s residue %s",
109 nab, i-1, aah[nah-1].name, hfn);
110 fgets(buf, STRLEN, in);
111 read_ab(buf,hfn,&(aah[nah-1].hack[i]));
112 if (debug) print_ab(debug, &(aah[nah-1].hack[i]));
115 fclose(in);
117 /* Sort the list (necessary to be able to use bsearch */
118 qsort(aah,nah,(size_t)sizeof(**ah),compaddh);
120 *ah=aah;
121 return nah;
124 void print_ab(FILE *out,t_hack *hack)
126 int i;
128 fprintf(out,"%d\t%d",hack->nr,hack->tp);
129 for(i=0; i < ncontrol[hack->tp]; i++)
130 fprintf(out,"\t%s",hack->a[i]);
131 fprintf(out,"\n");
134 void print_h_db(FILE *out,int nh,t_hackblock ah[])
136 int i,j;
138 for(i=0; (i<nh); i++) {
139 fprintf(out,"%s\t%d\n",ah[i].name,ah[i].nhack);
140 for(j=0; (j<ah[i].nhack); j++) {
141 fprintf(out,"\t");
142 print_ab(out,&(ah[i].hack[j]));
147 t_hackblock *search_h_db(int nh,t_hackblock ah[],char *key)
149 t_hackblock ahkey,*result;
151 if (nh <= 0)
152 return NULL;
154 ahkey.name=key;
156 result=(t_hackblock *)bsearch(&ahkey,ah,nh,(size_t)sizeof(ah[0]),compaddh);
158 return result;