6b1b536ca98942c8f52e0e70f6ca46ba65b685ed
[alexxy/gromacs.git] / src / kernel / nm2type.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.3.3
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-2008, 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 Machine for Chemical Simulation
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 "assert.h"
41 #include "maths.h"
42 #include "macros.h"
43 #include "copyrite.h"
44 #include "bondf.h"
45 #include "string2.h"
46 #include "smalloc.h"
47 #include "sysstuff.h"
48 #include "confio.h"
49 #include "physics.h"
50 #include "statutil.h"
51 #include "vec.h"
52 #include "random.h"
53 #include "3dview.h"
54 #include "txtdump.h"
55 #include "readinp.h"
56 #include "names.h"
57 #include "toppush.h"
58 #include "pdb2top.h"
59 #include "gpp_nextnb.h"
60 #include "gpp_atomtype.h"
61 #include "g_x2top.h"
62 #include "fflibutil.h"
63
64 static void rd_nm2type_file(const char *fn,int *nnm,t_nm2type **nmp)
65 {
66   FILE      *fp;
67   gmx_bool      bCont;
68   char      libfilename[128];
69   char      format[128],f1[128];
70   char      buf[1024],elem[16],type[16],nbbuf[16],**newbuf;
71   int       i,nb,nnnm,line=1;
72   double    qq,mm,*blen;
73   t_nm2type *nm2t=NULL;
74   
75   fp = fflib_open(fn);
76   if (NULL == fp)
77     gmx_fatal(FARGS,"Can not find %s in library directory",fn);
78     
79   nnnm = *nnm;
80   nm2t = *nmp;
81   do {
82     /* Read a line from the file */
83     bCont = (fgets2(buf,1023,fp) != NULL); 
84     
85     if (bCont) {
86       /* Remove comment */
87       strip_comment(buf);
88       if (sscanf(buf,"%s%s%lf%lf%d",elem,type,&qq,&mm,&nb) == 5) {
89         /* If we can read the first four, there probably is more */
90         srenew(nm2t,nnnm+1);
91         snew(nm2t[nnnm].blen,nb);
92         if (nb > 0) {
93           snew(newbuf,nb);
94           strcpy(format,"%*s%*s%*s%*s%*s");
95           for(i=0; (i<nb); i++) {
96             /* Complicated format statement */
97             strcpy(f1,format);
98             strcat(f1,"%s%lf");
99             if (sscanf(buf,f1,nbbuf,&(nm2t[nnnm].blen[i])) != 2)
100               gmx_fatal(FARGS,"Error on line %d of %s",line,libfilename);
101             newbuf[i] = strdup(nbbuf);
102             strcat(format,"%*s%*s");
103           }
104         }
105         else
106           newbuf = NULL;
107         nm2t[nnnm].elem   = strdup(elem);
108         nm2t[nnnm].type   = strdup(type);
109         nm2t[nnnm].q      = qq;
110         nm2t[nnnm].m      = mm;
111         nm2t[nnnm].nbonds = nb;
112         nm2t[nnnm].bond   = newbuf;
113         nnnm++;
114       }
115       line++;
116     }
117   } while(bCont);
118   ffclose(fp);
119   
120   *nnm = nnnm;
121   *nmp = nm2t;
122 }
123
124 t_nm2type *rd_nm2type(const char *ffdir,int *nnm)
125 {
126   int  nff,f;
127   char **ff;
128   t_nm2type *nm;
129
130   nff = fflib_search_file_end(ffdir,".n2t",FALSE,&ff);
131   *nnm = 0;
132   nm   = NULL;
133   for(f=0; f<nff; f++) {
134     rd_nm2type_file(ff[f],nnm,&nm);
135     sfree(ff[f]);
136   }
137   sfree(ff);
138
139   return nm;
140 }
141
142 void dump_nm2type(FILE *fp,int nnm,t_nm2type nm2t[])
143 {
144   int i,j;
145   
146   fprintf(fp,"; nm2type database\n");
147   for(i=0; (i<nnm); i++) {
148     fprintf(fp,"%-8s %-8s %8.4f %8.4f %-4d",
149             nm2t[i].elem,nm2t[i].type,
150             nm2t[i].q,nm2t[i].m,nm2t[i].nbonds);
151     for(j=0; (j<nm2t[i].nbonds); j++)
152       fprintf(fp," %-5s %6.4f",nm2t[i].bond[j],nm2t[i].blen[j]);
153     fprintf(fp,"\n");
154   }
155 }
156
157 enum { ematchNone, ematchWild, ematchElem, ematchExact, ematchNR };
158
159 static int match_str(const char *atom,const char *template_string)
160 {
161   if (!atom || !template_string)
162     return ematchNone;
163   else if (gmx_strcasecmp(atom,template_string) == 0) 
164     return ematchExact;
165   else if (atom[0] == template_string[0])
166     return ematchElem;
167   else if (strcmp(template_string,"*") == 0) 
168     return ematchWild;
169   else
170     return ematchNone;
171 }
172
173 int nm2type(int nnm,t_nm2type nm2t[],t_symtab *tab,t_atoms *atoms,
174             gpp_atomtype_t atype,int *nbonds,t_params *bonds)
175 {
176   int cur = 0;
177 #define prev (1-cur)
178   int i,j,k,m,n,nresolved,nb,maxbond,ai,aj,best,im,nqual[2][ematchNR];
179   int *bbb,*n_mask,*m_mask,**match,**quality;
180   char *aname_i,*aname_m,*aname_n,*type;
181   double qq,mm;
182   t_param *param;
183       
184   snew(param,1);
185   maxbond = 0;
186   for(i=0; (i<atoms->nr); i++) 
187     maxbond = max(maxbond,nbonds[i]);
188   if (debug)
189     fprintf(debug,"Max number of bonds per atom is %d\n",maxbond);
190   snew(bbb,maxbond);
191   snew(n_mask,maxbond);
192   snew(m_mask,maxbond);
193   snew(match,maxbond);
194   for(i=0; (i<maxbond); i++)
195     snew(match[i],maxbond);
196     
197   nresolved = 0;
198   for(i=0; (i<atoms->nr); i++) {
199     aname_i = *atoms->atomname[i];
200     nb = 0;
201     for(j=0; (j<bonds->nr); j++) {
202       ai = bonds->param[j].AI;
203       aj = bonds->param[j].AJ;
204       if (ai == i)
205         bbb[nb++] = aj;
206       else if (aj == i)
207         bbb[nb++] = ai;
208     }
209     if (nb != nbonds[i])
210       gmx_fatal(FARGS,"Counting number of bonds nb = %d, nbonds[%d] = %d",
211                 nb,i,nbonds[i]);
212     if(debug) {
213       fprintf(debug,"%4s has bonds to",aname_i);
214       for(j=0; (j<nb); j++)
215         fprintf(debug," %4s",*atoms->atomname[bbb[j]]);
216       fprintf(debug,"\n");
217     }    
218     best = -1;
219     for(k=0; (k<ematchNR); k++) 
220       nqual[prev][k] = 0;
221     
222     /* First check for names */  
223     for(k=0; (k<nnm); k++) {
224       if (nm2t[k].nbonds == nb) {
225         im = match_str(*atoms->atomname[i],nm2t[k].elem);
226         if (im > ematchWild) {
227           for(j=0; (j<ematchNR); j++) 
228             nqual[cur][j] = 0;
229
230           /* Fill a matrix with matching quality */
231           for(m=0; (m<nb); m++) {
232             aname_m = *atoms->atomname[bbb[m]];
233             for(n=0; (n<nb); n++) {
234               aname_n = nm2t[k].bond[n];
235               match[m][n] = match_str(aname_m,aname_n);
236             }
237           }
238           /* Now pick the best matches */
239           for(m=0; (m<nb); m++) {
240             n_mask[m] = 0;
241             m_mask[m] = 0;
242           }
243           for(j=ematchNR-1; (j>0); j--) {
244             for(m=0; (m<nb); m++) {
245               for(n=0; (n<nb); n++) {
246                 if ((n_mask[n] == 0) && 
247                     (m_mask[m] == 0) &&
248                     (match[m][n] == j)) {
249                   n_mask[n] = 1;
250                   m_mask[m] = 1;
251                   nqual[cur][j]++;
252                 }
253               }
254             }
255           }
256           if ((nqual[cur][ematchExact]+
257                nqual[cur][ematchElem]+
258                nqual[cur][ematchWild]) == nb) {
259             if ((nqual[cur][ematchExact] > nqual[prev][ematchExact]) ||
260                 
261                 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
262                  (nqual[cur][ematchElem] > nqual[prev][ematchElem])) ||
263                 
264                 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
265                  (nqual[cur][ematchElem] == nqual[prev][ematchElem]) &&
266                  (nqual[cur][ematchWild] > nqual[prev][ematchWild]))) {
267               best = k;
268               cur  = prev;
269             }
270           }
271         }
272       }
273     }
274     if (best != -1) {
275       int  atomnr=0;
276       real alpha=0;
277         
278       qq   = nm2t[best].q;
279       mm   = nm2t[best].m;
280       type = nm2t[best].type;
281       
282       if ((k = get_atomtype_type(type,atype)) == NOTSET) {
283         atoms->atom[i].qB = alpha;
284         atoms->atom[i].m = atoms->atom[i].mB = mm;
285         k = add_atomtype(atype,tab,&(atoms->atom[i]),type,param,
286                          atoms->atom[i].type,0,0,0,atomnr,0,0);
287       }
288       atoms->atom[i].type  = k;
289       atoms->atom[i].typeB = k;
290       atoms->atom[i].q  = qq;
291       atoms->atom[i].qB = qq;
292       atoms->atom[i].m  = mm;
293       atoms->atom[i].mB = mm;
294       nresolved++;
295     }
296     else {
297       fprintf(stderr,"Can not find forcefield for atom %s-%d with %d bonds\n",
298               *atoms->atomname[i],i+1,nb);
299     }
300   }
301   sfree(bbb);
302   sfree(n_mask);
303   sfree(m_mask);
304   sfree(param);
305   
306   return nresolved;
307 }
308