8b79fb4fef91d3fb8dbf968afed9ffb9639e63dc
[alexxy/gromacs.git] / src / kernel / nm2type.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2008, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /* This file is completely threadsafe - keep it that way! */
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include "assert.h"
44 #include "maths.h"
45 #include "macros.h"
46 #include "copyrite.h"
47 #include "bondf.h"
48 #include "string2.h"
49 #include "smalloc.h"
50 #include "sysstuff.h"
51 #include "confio.h"
52 #include "physics.h"
53 #include "statutil.h"
54 #include "vec.h"
55 #include "random.h"
56 #include "3dview.h"
57 #include "txtdump.h"
58 #include "readinp.h"
59 #include "names.h"
60 #include "toppush.h"
61 #include "pdb2top.h"
62 #include "gpp_nextnb.h"
63 #include "gpp_atomtype.h"
64 #include "g_x2top.h"
65 #include "fflibutil.h"
66
67 static void rd_nm2type_file(const char *fn,int *nnm,t_nm2type **nmp)
68 {
69   FILE      *fp;
70   gmx_bool      bCont;
71   char      libfilename[128];
72   char      format[128],f1[128];
73   char      buf[1024],elem[16],type[16],nbbuf[16],**newbuf;
74   int       i,nb,nnnm,line=1;
75   double    qq,mm,*blen;
76   t_nm2type *nm2t=NULL;
77   
78   fp = fflib_open(fn);
79   if (NULL == fp)
80     gmx_fatal(FARGS,"Can not find %s in library directory",fn);
81     
82   nnnm = *nnm;
83   nm2t = *nmp;
84   do {
85     /* Read a line from the file */
86     bCont = (fgets2(buf,1023,fp) != NULL); 
87     
88     if (bCont) {
89       /* Remove comment */
90       strip_comment(buf);
91       if (sscanf(buf,"%s%s%lf%lf%d",elem,type,&qq,&mm,&nb) == 5) {
92         /* If we can read the first four, there probably is more */
93         srenew(nm2t,nnnm+1);
94         snew(nm2t[nnnm].blen,nb);
95         if (nb > 0) {
96           snew(newbuf,nb);
97           strcpy(format,"%*s%*s%*s%*s%*s");
98           for(i=0; (i<nb); i++) {
99             /* Complicated format statement */
100             strcpy(f1,format);
101             strcat(f1,"%s%lf");
102             if (sscanf(buf,f1,nbbuf,&(nm2t[nnnm].blen[i])) != 2)
103               gmx_fatal(FARGS,"Error on line %d of %s",line,libfilename);
104             newbuf[i] = strdup(nbbuf);
105             strcat(format,"%*s%*s");
106           }
107         }
108         else
109           newbuf = NULL;
110         nm2t[nnnm].elem   = strdup(elem);
111         nm2t[nnnm].type   = strdup(type);
112         nm2t[nnnm].q      = qq;
113         nm2t[nnnm].m      = mm;
114         nm2t[nnnm].nbonds = nb;
115         nm2t[nnnm].bond   = newbuf;
116         nnnm++;
117       }
118       line++;
119     }
120   } while(bCont);
121   ffclose(fp);
122   
123   *nnm = nnnm;
124   *nmp = nm2t;
125 }
126
127 t_nm2type *rd_nm2type(const char *ffdir,int *nnm)
128 {
129   int  nff,f;
130   char **ff;
131   t_nm2type *nm;
132
133   nff = fflib_search_file_end(ffdir,".n2t",FALSE,&ff);
134   *nnm = 0;
135   nm   = NULL;
136   for(f=0; f<nff; f++) {
137     rd_nm2type_file(ff[f],nnm,&nm);
138     sfree(ff[f]);
139   }
140   sfree(ff);
141
142   return nm;
143 }
144
145 void dump_nm2type(FILE *fp,int nnm,t_nm2type nm2t[])
146 {
147   int i,j;
148   
149   fprintf(fp,"; nm2type database\n");
150   for(i=0; (i<nnm); i++) {
151     fprintf(fp,"%-8s %-8s %8.4f %8.4f %-4d",
152             nm2t[i].elem,nm2t[i].type,
153             nm2t[i].q,nm2t[i].m,nm2t[i].nbonds);
154     for(j=0; (j<nm2t[i].nbonds); j++)
155       fprintf(fp," %-5s %6.4f",nm2t[i].bond[j],nm2t[i].blen[j]);
156     fprintf(fp,"\n");
157   }
158 }
159
160 enum { ematchNone, ematchWild, ematchElem, ematchExact, ematchNR };
161
162 static int match_str(const char *atom,const char *template_string)
163 {
164   if (!atom || !template_string)
165     return ematchNone;
166   else if (gmx_strcasecmp(atom,template_string) == 0) 
167     return ematchExact;
168   else if (atom[0] == template_string[0])
169     return ematchElem;
170   else if (strcmp(template_string,"*") == 0) 
171     return ematchWild;
172   else
173     return ematchNone;
174 }
175
176 int nm2type(int nnm,t_nm2type nm2t[],t_symtab *tab,t_atoms *atoms,
177             gpp_atomtype_t atype,int *nbonds,t_params *bonds)
178 {
179   int cur = 0;
180 #define prev (1-cur)
181   int i,j,k,m,n,nresolved,nb,maxbond,ai,aj,best,im,nqual[2][ematchNR];
182   int *bbb,*n_mask,*m_mask,**match,**quality;
183   char *aname_i,*aname_m,*aname_n,*type;
184   double qq,mm;
185   t_param *param;
186       
187   snew(param,1);
188   maxbond = 0;
189   for(i=0; (i<atoms->nr); i++) 
190     maxbond = max(maxbond,nbonds[i]);
191   if (debug)
192     fprintf(debug,"Max number of bonds per atom is %d\n",maxbond);
193   snew(bbb,maxbond);
194   snew(n_mask,maxbond);
195   snew(m_mask,maxbond);
196   snew(match,maxbond);
197   for(i=0; (i<maxbond); i++)
198     snew(match[i],maxbond);
199     
200   nresolved = 0;
201   for(i=0; (i<atoms->nr); i++) {
202     aname_i = *atoms->atomname[i];
203     nb = 0;
204     for(j=0; (j<bonds->nr); j++) {
205       ai = bonds->param[j].AI;
206       aj = bonds->param[j].AJ;
207       if (ai == i)
208         bbb[nb++] = aj;
209       else if (aj == i)
210         bbb[nb++] = ai;
211     }
212     if (nb != nbonds[i])
213       gmx_fatal(FARGS,"Counting number of bonds nb = %d, nbonds[%d] = %d",
214                 nb,i,nbonds[i]);
215     if(debug) {
216       fprintf(debug,"%4s has bonds to",aname_i);
217       for(j=0; (j<nb); j++)
218         fprintf(debug," %4s",*atoms->atomname[bbb[j]]);
219       fprintf(debug,"\n");
220     }    
221     best = -1;
222     for(k=0; (k<ematchNR); k++) 
223       nqual[prev][k] = 0;
224     
225     /* First check for names */  
226     for(k=0; (k<nnm); k++) {
227       if (nm2t[k].nbonds == nb) {
228         im = match_str(*atoms->atomname[i],nm2t[k].elem);
229         if (im > ematchWild) {
230           for(j=0; (j<ematchNR); j++) 
231             nqual[cur][j] = 0;
232
233           /* Fill a matrix with matching quality */
234           for(m=0; (m<nb); m++) {
235             aname_m = *atoms->atomname[bbb[m]];
236             for(n=0; (n<nb); n++) {
237               aname_n = nm2t[k].bond[n];
238               match[m][n] = match_str(aname_m,aname_n);
239             }
240           }
241           /* Now pick the best matches */
242           for(m=0; (m<nb); m++) {
243             n_mask[m] = 0;
244             m_mask[m] = 0;
245           }
246           for(j=ematchNR-1; (j>0); j--) {
247             for(m=0; (m<nb); m++) {
248               for(n=0; (n<nb); n++) {
249                 if ((n_mask[n] == 0) && 
250                     (m_mask[m] == 0) &&
251                     (match[m][n] == j)) {
252                   n_mask[n] = 1;
253                   m_mask[m] = 1;
254                   nqual[cur][j]++;
255                 }
256               }
257             }
258           }
259           if ((nqual[cur][ematchExact]+
260                nqual[cur][ematchElem]+
261                nqual[cur][ematchWild]) == nb) {
262             if ((nqual[cur][ematchExact] > nqual[prev][ematchExact]) ||
263                 
264                 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
265                  (nqual[cur][ematchElem] > nqual[prev][ematchElem])) ||
266                 
267                 ((nqual[cur][ematchExact] == nqual[prev][ematchExact]) &&
268                  (nqual[cur][ematchElem] == nqual[prev][ematchElem]) &&
269                  (nqual[cur][ematchWild] > nqual[prev][ematchWild]))) {
270               best = k;
271               cur  = prev;
272             }
273           }
274         }
275       }
276     }
277     if (best != -1) {
278       int  atomnr=0;
279       real alpha=0;
280         
281       qq   = nm2t[best].q;
282       mm   = nm2t[best].m;
283       type = nm2t[best].type;
284       
285       if ((k = get_atomtype_type(type,atype)) == NOTSET) {
286         atoms->atom[i].qB = alpha;
287         atoms->atom[i].m = atoms->atom[i].mB = mm;
288         k = add_atomtype(atype,tab,&(atoms->atom[i]),type,param,
289                          atoms->atom[i].type,0,0,0,atomnr,0,0);
290       }
291       atoms->atom[i].type  = k;
292       atoms->atom[i].typeB = k;
293       atoms->atom[i].q  = qq;
294       atoms->atom[i].qB = qq;
295       atoms->atom[i].m  = mm;
296       atoms->atom[i].mB = mm;
297       nresolved++;
298     }
299     else {
300       fprintf(stderr,"Can not find forcefield for atom %s-%d with %d bonds\n",
301               *atoms->atomname[i],i+1,nb);
302     }
303   }
304   sfree(bbb);
305   sfree(n_mask);
306   sfree(m_mask);
307   sfree(param);
308   
309   return nresolved;
310 }
311