Merge remote-tracking branch 'origin/release-4-6' into HEAD
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / resall.c
1 /*  -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
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 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  * 
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  * 
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include "sysstuff.h"
41 #include <ctype.h>
42 #include "string2.h"
43 #include "strdb.h"
44 #include "futil.h"
45 #include "smalloc.h"
46 #include "gmx_fatal.h"
47 #include "symtab.h"
48 #include "macros.h"
49 #include "resall.h"
50 #include "pgutil.h"
51 #include "fflibutil.h"
52
53 gpp_atomtype_t read_atype(const char *ffdir,t_symtab *tab)
54 {
55     int        nfile,f;
56     char       **file;
57     FILE       *in;
58     char       buf[STRLEN],name[STRLEN];
59     double     m;
60     int        nratt=0;
61     gpp_atomtype_t at;
62     t_atom     *a;
63     t_param    *nb;
64     
65     nfile = fflib_search_file_end(ffdir,".atp",TRUE,&file);
66     at = init_atomtype();
67     snew(a,1);
68     snew(nb,1);
69     
70     for(f=0; f<nfile; f++)
71     {
72         in = fflib_open(file[f]);
73         while(!feof(in))
74         {
75             /* Skip blank or comment-only lines */
76             do
77             {
78                 fgets2(buf,STRLEN,in);
79                 if (NULL != buf)
80                 {
81                     strip_comment(buf);
82                     trim(buf);
83                 }
84             }
85             while (!feof(in) && NULL!=buf && strlen(buf)==0);
86             
87             if ((buf != NULL) && (sscanf(buf,"%s%lf",name,&m) == 2))
88             {
89                 a->m = m;
90                 add_atomtype(at,tab,a,name,nb, 0 , 0.0 , 0.0 , 0.0 , 0.0 , 0.0 , 0.0 );
91                 fprintf(stderr,"\rAtomtype %d",nratt+1);
92             }
93         }
94         ffclose(in);
95         sfree(file[f]);
96     }
97     fprintf(stderr,"\n");
98     sfree(file);
99   
100     return at;
101 }
102
103 static void print_resatoms(FILE *out,gpp_atomtype_t atype,t_restp *rtp)
104 {
105   int j,tp;
106   char *tpnm;
107   
108   /* fprintf(out,"%5s\n",rtp->resname);
109      fprintf(out,"%5d\n",rtp->natom); */
110   fprintf(out,"[ %s ]\n",rtp->resname);
111   fprintf(out," [ atoms ]\n");
112   
113   for(j=0; (j<rtp->natom); j++) {
114     tp = rtp->atom[j].type;
115     tpnm = get_atomtype_name(tp,atype);
116     if (tpnm == NULL)
117       gmx_fatal(FARGS,"Incorrect atomtype (%d)",tp);
118     fprintf(out,"%6s  %6s  %8.3f  %6d\n",
119             *(rtp->atomname[j]),tpnm,rtp->atom[j].q,rtp->cgnr[j]);
120   }
121 }
122
123 static gmx_bool read_atoms(FILE *in,char *line,
124                        t_restp *r0,t_symtab *tab,gpp_atomtype_t atype)
125 {
126   int    i,j,cg,maxentries;
127   char   buf[256],buf1[256];
128   double q;
129
130   /* Read Atoms */
131   maxentries=0;
132   r0->atom=     NULL;
133   r0->atomname= NULL;
134   r0->cgnr=     NULL;
135   i=0;
136   while (get_a_line(in,line,STRLEN) && (strchr(line,'[')==NULL)) { 
137     if (sscanf(line,"%s%s%lf%d",buf,buf1,&q,&cg) != 4)
138       return FALSE;
139     if (i>=maxentries) {
140       maxentries+=100;
141       srenew(r0->atom,     maxentries);
142       srenew(r0->atomname, maxentries);
143       srenew(r0->cgnr,     maxentries);
144     }
145     r0->atomname[i] = put_symtab(tab,buf);
146     r0->atom[i].q=q;
147     r0->cgnr[i]=cg;
148     j = get_atomtype_type(buf1,atype);
149     if (j == NOTSET)
150       gmx_fatal(FARGS,"Atom type %s (residue %s) not found in atomtype "
151                   "database",buf1,r0->resname);
152     r0->atom[i].type=j;
153     r0->atom[i].m=get_atomtype_massA(j,atype);
154     i++;
155   }
156   r0->natom=i;
157   srenew(r0->atom,i);
158   srenew(r0->atomname,i);
159   srenew(r0->cgnr,i);
160
161   return TRUE;
162 }
163
164 gmx_bool read_bondeds(int bt, FILE *in, char *line, t_restp *rtp)
165 {
166   char str[STRLEN];
167   int  j,n,ni,maxrb;
168   
169   maxrb = rtp->rb[bt].nb;
170   while (get_a_line(in,line,STRLEN) && (strchr(line,'[')==NULL)) {
171     if ( rtp->rb[bt].nb >= maxrb ) {
172       maxrb+=100;
173       srenew(rtp->rb[bt].b,maxrb);
174     }
175     n=0;
176     for(j=0; j < btsNiatoms[bt]; j++) {
177       if ( sscanf(line+n,"%s%n",str,&ni)==1 )
178         rtp->rb[bt].b[rtp->rb[bt].nb].a[j]=strdup(str);
179       else
180         return FALSE;
181       n+=ni;
182     }
183     for(  ; j < MAXATOMLIST; j++)
184       rtp->rb[bt].b[rtp->rb[bt].nb].a[j]=NULL;
185     while (isspace(line[n]))
186       n++;
187     rtrim(line+n);
188     rtp->rb[bt].b[rtp->rb[bt].nb].s=strdup(line+n);
189     rtp->rb[bt].nb++;
190   }
191   /* give back unused memory */
192   srenew(rtp->rb[bt].b,rtp->rb[bt].nb);
193   
194   return TRUE;
195 }
196
197 static void print_resbondeds(FILE *out, int bt, t_restp *rtp)
198 {
199   int i,j;
200   
201   if (rtp->rb[bt].nb) {
202     fprintf(out," [ %s ]\n",btsNames[bt]);
203     
204     for(i=0; i < rtp->rb[bt].nb; i++) {
205       for(j=0; j<btsNiatoms[bt]; j++)
206         fprintf(out,"%6s ",rtp->rb[bt].b[i].a[j]);
207       if (rtp->rb[bt].b[i].s[0])
208         fprintf(out,"    %s",rtp->rb[bt].b[i].s);
209       fprintf(out,"\n");
210     }
211   }
212 }
213
214 static void check_rtp(int nrtp,t_restp rtp[],char *libfn)
215 {
216   int i;
217
218   /* check for double entries, assuming list is already sorted */
219   for(i=1; (i<nrtp); i++) {
220     if (gmx_strcasecmp(rtp[i-1].resname,rtp[i].resname) == 0)
221       fprintf(stderr,"WARNING double entry %s in file %s\n",
222               rtp[i].resname,libfn);
223   }
224 }
225
226 static int comprtp(const void *a,const void *b)
227 {
228   t_restp *ra,*rb;
229
230   ra=(t_restp *)a;
231   rb=(t_restp *)b;
232
233   return gmx_strcasecmp(ra->resname,rb->resname);
234 }
235
236 int get_bt(char* header)
237 {
238   int i;
239
240   for(i=0; i<ebtsNR; i++)
241     if ( gmx_strcasecmp(btsNames[i],header)==0 )
242       return i;
243   return NOTSET;
244 }
245
246 void clear_t_restp(t_restp *rrtp)
247 {
248   memset((void *)rrtp, 0, sizeof(t_restp));
249 }
250
251 void read_resall(char *rrdb, int *nrtpptr, t_restp **rtp, 
252                  gpp_atomtype_t atype, t_symtab *tab,
253                  gmx_bool bAllowOverrideRTP)
254 {
255   FILE      *in;
256   char      filebase[STRLEN],*ptr,line[STRLEN],header[STRLEN];
257   int       i,nrtp,maxrtp,bt,nparam;
258   int       dum1,dum2,dum3;
259   t_restp   *rrtp;
260   gmx_bool      bNextResidue,bError;
261   int       bts[ebtsNR];
262   gmx_bool      bAlldih;
263   int       nrexcl;
264   gmx_bool      HH14;
265   gmx_bool      bRemoveDih;
266   int       firstrtp;
267
268   fflib_filename_base(rrdb,filebase,STRLEN);
269
270   in = fflib_open(rrdb);
271   
272   if (debug) {
273     fprintf(debug,"%9s %5s", "Residue", "atoms");
274     for(i=0; i<ebtsNR; i++)
275       fprintf(debug," %10s",btsNames[i]);
276     fprintf(debug,"\n");
277   }
278
279   /* these bonded parameters will overwritten be when  *
280    * there is a [ bondedtypes ] entry in the .rtp file */
281   bts[ebtsBONDS]  = 1; /* normal bonds     */
282   bts[ebtsANGLES] = 1; /* normal angles    */
283   bts[ebtsPDIHS]  = 1; /* normal dihedrals */
284   bts[ebtsIDIHS]  = 2; /* normal impropers */
285   bts[ebtsEXCLS]  = 1; /* normal exclusions */
286   bts[ebtsCMAP]   = 1; /* normal cmap torsions */
287
288   bAlldih    = FALSE;
289   nrexcl     = 3;
290   HH14       = TRUE;
291   bRemoveDih = TRUE;
292   
293   /* Column 5 & 6 aren't really bonded types, but we include
294    * them here to avoid introducing a new section:
295    * Column 5: 1 means generate all dihedrals, 0 not.
296    * Column 6: Number of bonded neighbors to exclude.
297    * Coulmn 7: Generate 1,4 interactions between pairs of hydrogens
298    * Column 8: Remove impropers over the same bond as a proper dihedral
299    */
300   
301   get_a_line(in,line,STRLEN);
302   if (!get_header(line,header))
303     gmx_fatal(FARGS,"in .rtp file at line:\n%s\n",line);
304   if (gmx_strncasecmp("bondedtypes",header,5)==0) {
305     get_a_line(in,line,STRLEN);
306     if ((nparam=sscanf(line,"%d %d %d %d %d %d %d %d",
307                        &bts[ebtsBONDS],&bts[ebtsANGLES],
308                        &bts[ebtsPDIHS],&bts[ebtsIDIHS],
309                        &dum1,&nrexcl,&dum2,&dum3)) < 4 )
310     {
311       gmx_fatal(FARGS,"need at least 4 (up to 8) parameters in .rtp file at line:\n%s\n",line);
312     }
313     bAlldih    = (dum1 != 0);
314     HH14       = (dum2 != 0);
315     bRemoveDih = (dum3 != 0);
316     get_a_line(in,line,STRLEN);
317     if(nparam<5) {
318       fprintf(stderr,"Using default: not generating all possible dihedrals\n");
319       bAlldih = FALSE;
320     }
321     if(nparam<6) {
322       fprintf(stderr,"Using default: excluding 3 bonded neighbors\n");
323       nrexcl = 3;
324     }
325     if(nparam<7) {
326       fprintf(stderr,"Using default: generating 1,4 H--H interactions\n");
327       HH14 = TRUE;
328     }
329     if(nparam<8) {
330       fprintf(stderr,"Using default: removing impropers on same bond as a proper\n");
331       bRemoveDih = TRUE;
332     }
333   } else {
334     fprintf(stderr,
335             "Reading .rtp file without '[ bondedtypes ]' directive,\n"
336             "Will proceed as if the entry\n"
337             "\n"
338             "\n[ bondedtypes ]"
339             "\n; bonds  angles  dihedrals  impropers all_dihedrals nr_exclusions HH14 remove_dih"
340             "\n   %3d     %3d        %3d        %3d           %3d           %3d   %3d    %3d"
341             "\n"
342             "was present at the beginning of %s",
343             bts[0],bts[1],bts[2],bts[3], bAlldih ? 1 : 0,nrexcl,HH14,bRemoveDih,rrdb);
344   }
345   /* We don't know the current size of rrtp, but simply realloc immediately */
346   nrtp = *nrtpptr;
347   rrtp = *rtp;
348   maxrtp = nrtp;
349   while (!feof(in)) {
350     if (nrtp >= maxrtp) {
351       maxrtp+=100;
352       srenew(rrtp,maxrtp);
353     }
354     clear_t_restp(&rrtp[nrtp]);
355     if (!get_header(line,header))
356       gmx_fatal(FARGS,"in .rtp file at line:\n%s\n",line);
357     rrtp[nrtp].resname  = strdup(header);
358     rrtp[nrtp].filebase = strdup(filebase);
359
360     /* Set the bonded types */
361     rrtp[nrtp].bAlldih    = bAlldih;
362     rrtp[nrtp].nrexcl     = nrexcl;
363     rrtp[nrtp].HH14       = HH14;
364     rrtp[nrtp].bRemoveDih = bRemoveDih;
365     for(i=0; i<ebtsNR; i++) {
366       rrtp[nrtp].rb[i].type = bts[i];
367     }
368
369     get_a_line(in,line,STRLEN);
370     bError=FALSE;
371     bNextResidue=FALSE;
372     do {
373       if (!get_header(line,header)) {
374         bError = TRUE;
375       } else {
376         bt = get_bt(header);
377         if (bt != NOTSET) {
378           /* header is an bonded directive */
379           bError = !read_bondeds(bt,in,line,&rrtp[nrtp]);
380         } else if (gmx_strncasecmp("atoms",header,5) == 0) {
381           /* header is the atoms directive */
382           bError = !read_atoms(in,line,&(rrtp[nrtp]),tab,atype);
383         } else {
384           /* else header must be a residue name */
385           bNextResidue = TRUE;
386         }
387       }
388       if (bError)
389         gmx_fatal(FARGS,"in .rtp file in residue %s at line:\n%s\n",
390                     rrtp[nrtp].resname,line);
391     } while (!feof(in) && !bNextResidue);
392
393     if (rrtp[nrtp].natom == 0)
394       gmx_fatal(FARGS,"No atoms found in .rtp file in residue %s\n",
395                   rrtp[nrtp].resname);
396     if (debug) {
397       fprintf(debug,"%3d %5s %5d",
398               nrtp+1,rrtp[nrtp].resname,rrtp[nrtp].natom);
399       for(i=0; i<ebtsNR; i++)
400         fprintf(debug," %10d",rrtp[nrtp].rb[i].nb);
401       fprintf(debug,"\n");
402     }
403     
404     firstrtp = -1;
405     for(i=0; i<nrtp; i++) {
406         if (gmx_strcasecmp(rrtp[i].resname,rrtp[nrtp].resname) == 0) {
407             firstrtp = i;
408         }
409     }
410     
411     if (firstrtp == -1) {
412         nrtp++;
413         fprintf(stderr,"\rResidue %d",nrtp);
414     } else {
415         if (firstrtp >= *nrtpptr)
416         {
417             gmx_fatal(FARGS,"Found a second entry for '%s' in '%s'",
418                       rrtp[nrtp].resname,rrdb);
419         }
420         if (bAllowOverrideRTP)
421         {
422             fprintf(stderr,"\n");
423             fprintf(stderr,"Found another rtp entry for '%s' in '%s', ignoring this entry and keeping the one from '%s.rtp'\n",
424                     rrtp[nrtp].resname,rrdb,rrtp[firstrtp].filebase);
425             /* We should free all the data for this entry.
426              * The current code gives a lot of dangling pointers.
427              */
428             clear_t_restp(&rrtp[nrtp]);
429         }
430         else
431         {
432             gmx_fatal(FARGS,"Found rtp entries for '%s' in both '%s' and '%s'. If you want the first definition to override the second one, set the -rtpo option of pdb2gmx.",rrtp[nrtp].resname,rrtp[firstrtp].filebase,rrdb);
433         }
434     }
435   }
436   ffclose(in);
437   /* give back unused memory */
438   srenew(rrtp,nrtp);
439   
440   fprintf(stderr,"\nSorting it all out...\n");
441   qsort(rrtp,nrtp,(size_t)sizeof(rrtp[0]),comprtp);
442   
443   check_rtp(nrtp,rrtp,rrdb);
444
445   *nrtpptr = nrtp;
446   *rtp     = rrtp;
447 }
448
449 void print_resall(FILE *out, int nrtp, t_restp rtp[],
450                   gpp_atomtype_t atype)
451 {
452   int i,bt;
453
454   if (nrtp == 0) {
455     return;
456   }
457
458   /* print all the ebtsNR type numbers */
459   fprintf(out,"[ bondedtypes ]\n");
460   fprintf(out,"; bonds  angles  dihedrals  impropers all_dihedrals nr_exclusions  HH14  remove_dih\n");
461   fprintf(out," %5d  %6d  %9d  %9d  %14d  %14d %14d %14d\n\n",
462           rtp[0].rb[0].type,
463           rtp[0].rb[1].type,
464           rtp[0].rb[2].type,
465           rtp[0].rb[3].type,
466           rtp[0].bAlldih,rtp[0].nrexcl,rtp[0].HH14,rtp[0].bRemoveDih);
467
468   for(i=0; i<nrtp; i++) {
469     if (rtp[i].natom > 0) {
470       print_resatoms(out,atype,&rtp[i]);
471       for(bt=0; bt<ebtsNR; bt++)
472         print_resbondeds(out,bt,&rtp[i]);
473     }
474   }
475 }
476
477 /************************************************************
478  *
479  *                  SEARCH   ROUTINES
480  * 
481  ***********************************************************/
482 static gmx_bool is_sign(char c)
483 {
484     return (c == '+' || c == '-');
485 }
486
487 /* Compares if the strins match except for a sign at the end
488  * of (only) one of the two.
489  */
490 static int neq_str_sign(const char *a1,const char *a2)
491 {
492     int l1,l2,lm;
493   
494     l1 = (int)strlen(a1);
495     l2 = (int)strlen(a2);
496     lm = min(l1,l2);
497
498     if (lm >= 1 &&
499         ((l1 == l2+1 && is_sign(a1[l1-1])) ||
500          (l2 == l1+1 && is_sign(a2[l2-1]))) &&
501         gmx_strncasecmp(a1,a2,lm) == 0)
502     {
503         return lm;
504     }
505     else
506     {
507         return 0;
508     }
509 }
510
511 char *search_rtp(const char *key,int nrtp,t_restp rtp[])
512 {
513     int  i,n,nbest,best,besti;
514     char bestbuf[STRLEN];
515
516     nbest =  0;
517     besti = -1;
518     /* We want to match at least one character */
519     best  =  1;
520     for(i=0; (i<nrtp); i++)
521     {
522         if (gmx_strcasecmp(key,rtp[i].resname) == 0)
523         {
524              besti = i;
525              nbest = 1;
526              break;
527         }
528         else
529         {
530             /* Allow a mismatch of at most a sign character (with warning) */
531             n = neq_str_sign(key,rtp[i].resname);
532             if (n >= best &&
533                 n+1 >= (int)strlen(key) &&
534                 n+1 >= (int)strlen(rtp[i].resname))
535             {
536                 if (n == best)
537                 {
538                     if (nbest == 1)
539                     {
540                         strcpy(bestbuf,rtp[besti].resname);
541                     }
542                     if (nbest >= 1)
543                     {
544                         strcat(bestbuf," or ");
545                         strcat(bestbuf,rtp[i].resname);
546                     }
547                 }
548                 else
549                 {
550                     nbest = 0;
551                 }
552                 besti = i;
553                 best  = n;
554                 nbest++;
555             }
556         }
557     }
558     if (nbest > 1)
559     {
560         gmx_fatal(FARGS,"Residue '%s' not found in residue topology database, looks a bit like %s",key,bestbuf);
561     }
562     else if (besti == -1)
563     {
564         gmx_fatal(FARGS,"Residue '%s' not found in residue topology database",key);
565     }
566     if (gmx_strcasecmp(rtp[besti].resname,key) != 0)
567     {
568         fprintf(stderr,
569                 "\nWARNING: '%s' not found in residue topology database, "
570                 "trying to use '%s'\n\n", key, rtp[besti].resname);
571     }
572   
573     return rtp[besti].resname;
574 }
575
576 t_restp *get_restp(const char *rtpname,int nrtp,t_restp rtp[])
577 {
578     int  i;
579
580     i = 0;
581     while (i < nrtp && gmx_strcasecmp(rtpname,rtp[i].resname) != 0)
582     {
583         i++;
584     }
585     if (i >= nrtp)
586     {
587         /* This should never happen, since search_rtp should have been called
588          * before calling get_restp.
589          */
590         gmx_fatal(FARGS,"Residue type '%s' not found in residue topology database",rtpname);
591     }
592   
593     return &rtp[i];
594 }