3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 /* This file is completely threadsafe - keep it that way! */
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
)
66 return strcasecmp(ah
->name
,bh
->name
);
69 void print_ab(FILE *out
,t_hack
*hack
,char *nname
)
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
]);
80 void read_ab(char *line
,const char *fn
,t_hack
*hack
)
86 ns
= sscanf(line
,"%d%d%s%s%s%s%s",&nh
,&tp
,hn
,a
[0],a
[1],a
[2],a
[3]);
88 gmx_fatal(FARGS
,"wrong format in input file %s on line\n%s\n",fn
,line
);
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
);
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
]);
105 hack
->nname
=strdup(hn
);
110 hack
->newx
[i
]=NOTSET
;
113 static void dump_h_db(const char *fn
,int nah
,t_hackblock
*ah
)
116 char buf
[STRLEN
],nname
[STRLEN
];
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]);
126 print_ab(fp
,&ah
[i
].hack
[k
],nname
);
132 static void read_h_db_file(const char *hfn
,int *nahptr
,t_hackblock
**ah
)
135 char filebase
[STRLEN
],line
[STRLEN
], buf
[STRLEN
];
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
147 in
= fflib_open(hfn
);
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",
157 if (debug
) fprintf(debug
,"%s",buf
);
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
++) {
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
]));
183 /* Sort the list (necessary to be able to use bsearch */
184 qsort(aah
,nah
,(size_t)sizeof(**ah
),compaddh
);
187 dump_h_db(hfn
,nah
,aah
);
193 int read_h_db(const char *ffdir
,t_hackblock
**ah
)
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
);
206 for(f
=0; f
<nhdbf
; f
++) {
207 read_h_db_file(hdbf
[f
],&nah
,ah
);
215 t_hackblock
*search_h_db(int nh
,t_hackblock ah
[],char *key
)
217 t_hackblock ahkey
,*result
;
224 result
=(t_hackblock
*)bsearch(&ahkey
,ah
,nh
,(size_t)sizeof(ah
[0]),compaddh
);