Merge remote-tracking branch 'gerrit/release-4-5-patches' into release-4-6
[alexxy/gromacs.git] / src / kernel / pdb2top.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 <stdio.h>
41 #include <math.h>
42 #include <ctype.h>
43
44 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
45 #include <direct.h>
46 #include <io.h>
47 #endif
48
49
50 #include "vec.h"
51 #include "copyrite.h"
52 #include "smalloc.h"
53 #include "macros.h"
54 #include "symtab.h"
55 #include "futil.h"
56 #include "statutil.h"
57 #include "gmx_fatal.h"
58 #include "pdb2top.h"
59 #include "gpp_nextnb.h"
60 #include "topdirs.h"
61 #include "toputil.h"
62 #include "h_db.h"
63 #include "pgutil.h"
64 #include "resall.h"
65 #include "topio.h"
66 #include "string2.h"
67 #include "physics.h"
68 #include "pdbio.h"
69 #include "gen_ad.h"
70 #include "filenm.h"
71 #include "index.h"
72 #include "gen_vsite.h"
73 #include "add_par.h"
74 #include "toputil.h"
75 #include "fflibutil.h"
76 #include "strdb.h"
77
78 /* this must correspond to enum in pdb2top.h */
79 const char *hh[ehisNR]   = { "HISD", "HISE", "HISH", "HIS1" };
80
81 static int missing_atoms(t_restp *rp, int resind,t_atoms *at, int i0, int i)
82 {
83     int  j,k,nmiss;
84     char *name;
85     gmx_bool bFound, bRet;
86     
87     nmiss = 0;
88     for (j=0; j<rp->natom; j++)
89     {
90         name=*(rp->atomname[j]);
91         bFound=FALSE;
92         for (k=i0; k<i; k++) 
93         {
94             bFound = (bFound || !gmx_strcasecmp(*(at->atomname[k]),name));
95         }
96         if (!bFound)
97         {
98             nmiss++;
99             fprintf(stderr,"\nWARNING: "
100                     "atom %s is missing in residue %s %d in the pdb file\n",
101                     name,*(at->resinfo[resind].name),at->resinfo[resind].nr);
102             if (name[0]=='H' || name[0]=='h')
103             {
104                 fprintf(stderr,"         You might need to add atom %s to the hydrogen database of building block %s\n"
105                         "         in the file %s.hdb (see the manual)\n",
106                         name,*(at->resinfo[resind].rtp),rp->filebase);
107             }
108             fprintf(stderr,"\n");
109         }
110     }
111   
112     return nmiss;
113 }
114
115 gmx_bool is_int(double x)
116 {
117   const double tol = 1e-4;
118   int   ix;
119   
120   if (x < 0)
121     x=-x;
122   ix=gmx_nint(x);
123   
124   return (fabs(x-ix) < tol);
125 }
126
127 static void swap_strings(char **s,int i,int j)
128 {
129     char *tmp;
130
131     tmp  = s[i];
132     s[i] = s[j];
133     s[j] = tmp;
134 }
135
136 void
137 choose_ff(const char *ffsel,
138           char *forcefield, int ff_maxlen,
139           char *ffdir, int ffdir_maxlen)
140 {
141     int  nff;
142     char **ffdirs,**ffs,**ffs_dir,*ptr;
143     int  i,j,sel,cwdsel,nfound;
144     char buf[STRLEN],**desc;
145     FILE *fp;
146     char *pret;
147     
148     nff = fflib_search_file_in_dirend(fflib_forcefield_itp(),
149                                       fflib_forcefield_dir_ext(),
150                                       &ffdirs);
151
152     if (nff == 0)
153     {
154         gmx_fatal(FARGS,"No force fields found (files with name '%s' in subdirectories ending on '%s')",
155                   fflib_forcefield_itp(),fflib_forcefield_dir_ext());
156     }
157
158     /* Replace with unix path separators */
159     if(DIR_SEPARATOR!='/')
160     {
161         for(i=0;i<nff;i++)
162         {
163             while( (ptr=strchr(ffdirs[i],DIR_SEPARATOR))!=NULL )
164             {
165                 *ptr='/';
166             }
167         }
168     }
169     
170     /* Store the force field names in ffs */
171     snew(ffs,nff);
172     snew(ffs_dir,nff);
173     for(i=0; i<nff; i++)
174     {
175         /* Remove the path from the ffdir name - use our unix standard here! */
176         ptr = strrchr(ffdirs[i],'/');
177         if (ptr == NULL)
178         {
179             ffs[i] = strdup(ffdirs[i]);
180             ffs_dir[i] = low_gmxlibfn(ffdirs[i],FALSE,FALSE);
181             if (ffs_dir[i] == NULL)
182             {
183                 gmx_fatal(FARGS,"Can no longer find file '%s'",ffdirs[i]);
184             }
185         }
186         else
187         {
188             ffs[i] = strdup(ptr+1);
189             ffs_dir[i] = strdup(ffdirs[i]);
190         }
191         ffs_dir[i][strlen(ffs_dir[i])-strlen(ffs[i])-1] = '\0';
192         /* Remove the extension from the ffdir name */
193         ffs[i][strlen(ffs[i])-strlen(fflib_forcefield_dir_ext())] = '\0';
194     }
195
196     if (ffsel != NULL)
197     {
198         sel     = -1;
199         cwdsel  = -1;
200         nfound  = 0;
201         for(i=0; i<nff; i++)
202         {
203             if ( strcmp(ffs[i],ffsel)==0 )
204             {
205                 /* Matching ff name */
206                 sel = i;
207                 nfound++;
208                 
209                 if( strncmp(ffs_dir[i],".",1)==0 )
210                 {
211                     cwdsel = i;
212                 }
213             }
214         }
215         
216         if(cwdsel != -1)
217         {
218             sel = cwdsel;
219         }
220         
221         if(nfound>1)
222         {
223             if(cwdsel!=-1)
224             {
225                 fprintf(stderr,
226                         "Force field '%s' occurs in %d places. pdb2gmx is using the one in the\n"
227                         "current directory. Use interactive selection (not the -ff option) if\n"
228                         "you would prefer a different one.\n",ffsel,nfound);
229             }
230             else
231             {
232                 gmx_fatal(FARGS,
233                           "Force field '%s' occurs in %d places, but not in the current directory.\n"
234                           "Run without the -ff switch and select the force field interactively.",ffsel,nfound);
235             }
236         }
237         else if (nfound==0)
238         {
239             gmx_fatal(FARGS,"Could not find force field '%s' in current directory, install tree or GMXDATA path.",ffsel);
240         }
241     }
242     else if (nff > 1)
243     {
244         snew(desc,nff);
245         for(i=0; (i<nff); i++)
246         {
247             sprintf(buf,"%s%c%s%s%c%s",
248                     ffs_dir[i],DIR_SEPARATOR,
249                     ffs[i],fflib_forcefield_dir_ext(),DIR_SEPARATOR,
250                     fflib_forcefield_doc());
251             if (gmx_fexist(buf))
252             {
253                 /* We don't use fflib_open, because we don't want printf's */
254                 fp = ffopen(buf,"r");
255                 snew(desc[i],STRLEN);
256                 get_a_line(fp,desc[i],STRLEN);
257                 ffclose(fp);
258             }
259             else
260             {
261                 desc[i] = strdup(ffs[i]);
262             }
263         }
264         /* Order force fields from the same dir alphabetically
265          * and put deprecated force fields at the end.
266          */
267         for(i=0; (i<nff); i++)
268         {
269             for(j=i+1; (j<nff); j++)
270             {
271                 if (strcmp(ffs_dir[i],ffs_dir[j]) == 0 &&
272                     ((desc[i][0] == '[' && desc[j][0] != '[') ||
273                      ((desc[i][0] == '[' || desc[j][0] != '[') &&
274                       gmx_strcasecmp(desc[i],desc[j]) > 0)))
275                 {
276                     swap_strings(ffdirs,i,j);
277                     swap_strings(ffs   ,i,j);
278                     swap_strings(desc  ,i,j);
279                 }
280             }
281         }
282
283         printf("\nSelect the Force Field:\n");
284         for(i=0; (i<nff); i++)
285         {
286             if (i == 0 || strcmp(ffs_dir[i-1],ffs_dir[i]) != 0)
287             {
288                 if( strcmp(ffs_dir[i],".")==0 )
289                 {
290                     printf("From current directory:\n");
291                 }
292                 else
293                 {
294                     printf("From '%s':\n",ffs_dir[i]);
295                 }
296             }
297             printf("%2d: %s\n",i+1,desc[i]);
298             sfree(desc[i]);
299         }
300         sfree(desc);
301
302         do
303         {
304             pret = fgets(buf,STRLEN,stdin);
305             
306             if (pret != NULL)
307             {
308                 sscanf(buf,"%d",&sel);
309                 sel--;
310             }
311         }
312         while ( pret==NULL || (sel < 0) || (sel >= nff));
313
314         /* Check for a current limitation of the fflib code.
315          * It will always read from the first ff directory in the list.
316          * This check assumes that the order of ffs matches the order
317          * in which fflib_open searches ff library files.
318          */
319         for(i=0; i<sel; i++)
320         {
321             if (strcmp(ffs[i],ffs[sel]) == 0)
322             {
323                 gmx_fatal(FARGS,"Can only select the first of multiple force field entries with directory name '%s%s' in the list. If you want to use the next entry, run pdb2gmx in a different directory or rename or move the force field directory present in the current working directory.",
324                           ffs[sel],fflib_forcefield_dir_ext());
325             }
326         }
327     }
328     else
329     {
330         sel = 0;
331     }
332
333     if (strlen(ffs[sel]) >= (size_t)ff_maxlen)
334     {
335         gmx_fatal(FARGS,"Length of force field name (%d) >= maxlen (%d)",
336                   strlen(ffs[sel]),ff_maxlen);
337     }
338     strcpy(forcefield,ffs[sel]);
339
340     if (strlen(ffdirs[sel]) >= (size_t)ffdir_maxlen)
341     {
342         gmx_fatal(FARGS,"Length of force field dir (%d) >= maxlen (%d)",
343                   strlen(ffdirs[sel]),ffdir_maxlen);
344     }
345     strcpy(ffdir,ffdirs[sel]);
346
347     for(i=0; (i<nff); i++)
348     {
349         sfree(ffdirs[i]);
350         sfree(ffs[i]);
351         sfree(ffs_dir[i]);
352     }
353     sfree(ffdirs);
354     sfree(ffs);
355     sfree(ffs_dir);
356 }
357
358 void choose_watermodel(const char *wmsel,const char *ffdir,
359                        char **watermodel)
360 {
361     const char *fn_watermodels="watermodels.dat";
362     char fn_list[STRLEN];
363     FILE *fp;
364     char buf[STRLEN];
365     int  nwm,sel,i;
366     char **model;
367     char *pret;
368
369     if (strcmp(wmsel,"none") == 0)
370     {
371         *watermodel = NULL;
372         
373         return;
374     }
375     else if (strcmp(wmsel,"select") != 0)
376     {
377         *watermodel = strdup(wmsel);
378
379         return;
380     }
381
382     sprintf(fn_list,"%s%c%s",ffdir,DIR_SEPARATOR,fn_watermodels);
383     
384     if (!fflib_fexist(fn_list))
385     {
386         fprintf(stderr,"No file '%s' found, will not include a water model\n",
387                 fn_watermodels);
388         *watermodel = NULL;
389         
390         return;
391     }
392
393     fp = fflib_open(fn_list);
394     printf("\nSelect the Water Model:\n");
395     nwm = 0;
396     model = NULL;
397     while (get_a_line(fp,buf,STRLEN))
398     {
399         srenew(model,nwm+1);
400         snew(model[nwm],STRLEN);
401         sscanf(buf,"%s%n",model[nwm],&i);
402         if (i > 0)
403         {
404             ltrim(buf+i);
405             fprintf(stderr,"%2d: %s\n",nwm+1,buf+i);
406             nwm++;
407         }
408         else
409         {
410             sfree(model[nwm]);
411         }
412     }
413     ffclose(fp);
414     fprintf(stderr,"%2d: %s\n",nwm+1,"None");
415
416     do
417     {
418         pret = fgets(buf,STRLEN,stdin);
419         
420         if (pret != NULL)
421         {
422             sscanf(buf,"%d",&sel);
423             sel--;
424         }
425     }
426     while (pret == NULL || sel < 0 || sel > nwm);
427
428     if (sel == nwm)
429     {
430         *watermodel = NULL;
431     }
432     else
433     {
434         *watermodel = strdup(model[sel]);
435     }
436
437     for(i=0; i<nwm; i++)
438     {
439         sfree(model[i]);
440     }
441     sfree(model);
442 }
443
444 static int name2type(t_atoms *at, int **cgnr, gpp_atomtype_t atype, 
445                      t_restp restp[], gmx_residuetype_t rt)
446 {
447   int     i,j,prevresind,resind,i0,prevcg,cg,curcg;
448   char    *name;
449   gmx_bool    bProt, bNterm;
450   double  qt;
451   int     nmissat;
452     
453   nmissat = 0;
454
455   resind=-1;
456   bProt=FALSE;
457   bNterm=FALSE;
458   i0=0;
459   snew(*cgnr,at->nr);
460   qt=0;
461   prevcg=NOTSET;
462   curcg=0;
463   cg=-1;
464   j=NOTSET;
465   
466   for(i=0; (i<at->nr); i++) {
467     prevresind=resind;
468     if (at->atom[i].resind != resind) {
469       resind = at->atom[i].resind;
470       bProt = gmx_residuetype_is_protein(rt,*(at->resinfo[resind].name));
471       bNterm=bProt && (resind == 0);
472       if (resind > 0) {
473           nmissat += missing_atoms(&restp[prevresind],prevresind,at,i0,i);
474       }
475       i0=i;
476     }
477     if (at->atom[i].m == 0) {
478       if (debug)
479         fprintf(debug,"atom %d%s: curcg=%d, prevcg=%d, cg=%d\n",
480                 i+1,*(at->atomname[i]),curcg,prevcg,
481                 j==NOTSET ? NOTSET : restp[resind].cgnr[j]);
482       qt=0;
483       prevcg=cg;
484       name=*(at->atomname[i]);
485       j=search_jtype(&restp[resind],name,bNterm);
486       at->atom[i].type = restp[resind].atom[j].type;
487       at->atom[i].q    = restp[resind].atom[j].q;
488       at->atom[i].m    = get_atomtype_massA(restp[resind].atom[j].type,
489                                             atype);
490       cg = restp[resind].cgnr[j];
491       /* A charge group number -1 signals a separate charge group
492        * for this atom.
493        */
494       if ( (cg == -1) || (cg != prevcg) || (resind != prevresind) ) {
495           curcg++;
496       }
497     } else {
498       if (debug)
499         fprintf(debug,"atom %d%s: curcg=%d, qt=%g, is_int=%d\n",
500                 i+1,*(at->atomname[i]),curcg,qt,is_int(qt));
501       cg=-1;
502       if (is_int(qt)) {
503         qt=0;
504         curcg++;
505       }
506       qt+=at->atom[i].q;
507     }
508     (*cgnr)[i]=curcg;
509     at->atom[i].typeB = at->atom[i].type;
510     at->atom[i].qB    = at->atom[i].q;
511     at->atom[i].mB    = at->atom[i].m;
512   }
513   nmissat += missing_atoms(&restp[resind],resind,at,i0,i);
514
515   return nmissat;
516 }
517
518 static void print_top_heavy_H(FILE *out, real mHmult)
519 {
520   if (mHmult == 2.0) 
521     fprintf(out,"; Using deuterium instead of hydrogen\n\n");
522   else if (mHmult == 4.0)
523     fprintf(out,"#define HEAVY_H\n\n");
524   else if (mHmult != 1.0)
525     fprintf(stderr,"WARNING: unsupported proton mass multiplier (%g) "
526             "in pdb2top\n",mHmult);
527 }
528
529 void print_top_comment(FILE *out,
530                        const char *filename,
531                        const char *generator,
532                        const char *ffdir,
533                        gmx_bool bITP)
534 {
535   char tmp[256]; 
536   char ffdir_parent[STRLEN];
537   char *p;
538         
539   nice_header(out,filename);
540   fprintf(out,";\tThis is a %s topology file\n;\n",bITP ? "include" : "standalone");
541   fprintf(out,";\tIt was generated using program:\n;\t%s\n;\n",
542           (NULL == generator) ? "unknown" : generator);
543   fprintf(out,";\tCommand line was:\n;\t%s\n;\n",command_line());
544
545   if(strchr(ffdir,'/')==NULL)
546   {
547       fprintf(out,";\tForce field was read from the standard Gromacs share directory.\n;\n\n");
548   }
549   else if(ffdir[0]=='.')
550   {
551       fprintf(out,";\tForce field was read from current directory or a relative path - path added.\n;\n\n");
552   }
553   else
554   {
555       strncpy(ffdir_parent,ffdir,STRLEN-1);
556       ffdir_parent[STRLEN-1]='\0'; /*make sure it is 0-terminated even for long string*/
557       p=strrchr(ffdir_parent,'/');
558
559       *p='\0';
560       
561       fprintf(out,
562               ";\tForce field data was read from:\n"
563               ";\t%s\n"
564               ";\n"
565               ";\tNote:\n"
566               ";\tThis might be a non-standard force field location. When you use this topology, the\n"
567               ";\tforce field must either be present in the current directory, or the location\n"
568               ";\tspecified in the GMXLIB path variable or with the 'include' mdp file option.\n;\n\n",
569               ffdir_parent);
570   }
571 }
572
573 void print_top_header(FILE *out,const char *filename, 
574                       const char *title,gmx_bool bITP,const char *ffdir,real mHmult)
575 {
576     const char *p;
577     
578     print_top_comment(out,filename,title,ffdir,bITP);
579     
580     print_top_heavy_H(out, mHmult);
581     fprintf(out,"; Include forcefield parameters\n");
582
583     p=strrchr(ffdir,'/');        
584     p = (ffdir[0]=='.' || p==NULL) ? ffdir : p+1;
585
586     fprintf(out,"#include \"%s/%s\"\n\n",p,fflib_forcefield_itp());
587 }
588
589 static void print_top_posre(FILE *out,const char *pr)
590 {
591   fprintf(out,"; Include Position restraint file\n");
592   fprintf(out,"#ifdef POSRES\n");
593   fprintf(out,"#include \"%s\"\n",pr);
594   fprintf(out,"#endif\n\n");
595 }
596   
597 static void print_top_water(FILE *out,const char *ffdir,const char *water)
598 {
599   const char *p;
600   char  buf[STRLEN];
601     
602   fprintf(out,"; Include water topology\n");
603
604   p=strrchr(ffdir,'/');        
605   p = (ffdir[0]=='.' || p==NULL) ? ffdir : p+1;
606   fprintf(out,"#include \"%s/%s.itp\"\n",p,water);
607   
608   fprintf(out,"\n");
609   fprintf(out,"#ifdef POSRES_WATER\n");
610   fprintf(out,"; Position restraint for each water oxygen\n");
611   fprintf(out,"[ position_restraints ]\n");
612   fprintf(out,";%3s %5s %9s %10s %10s\n","i","funct","fcx","fcy","fcz");
613   fprintf(out,"%4d %4d %10g %10g %10g\n",1,1,1000.0,1000.0,1000.0);
614   fprintf(out,"#endif\n");
615   fprintf(out,"\n");
616
617   sprintf(buf,"%s/ions.itp",p);
618
619   if (fflib_fexist(buf))
620   {
621     fprintf(out,"; Include topology for ions\n");
622     fprintf(out,"#include \"%s\"\n",buf);
623     fprintf(out,"\n");
624   }
625 }
626
627 static void print_top_system(FILE *out, const char *title)
628 {
629   fprintf(out,"[ %s ]\n",dir2str(d_system));
630   fprintf(out,"; Name\n");
631   fprintf(out,"%s\n\n",title[0]?title:"Protein");
632 }
633
634 void print_top_mols(FILE *out,
635                     const char *title, const char *ffdir, const char *water,
636                     int nincl, char **incls, int nmol, t_mols *mols)
637 {
638   int  i;
639   char *incl;
640
641   if (nincl>0) {
642     fprintf(out,"; Include chain topologies\n");
643     for (i=0; (i<nincl); i++) {
644         incl = strrchr(incls[i],DIR_SEPARATOR);
645         if (incl == NULL) {
646             incl = incls[i];
647         } else {
648             /* Remove the path from the include name */
649             incl = incl + 1;
650         }
651       fprintf(out,"#include \"%s\"\n",incl);
652     }
653     fprintf(out,"\n");
654   }
655
656     if (water)
657     {
658       print_top_water(out,ffdir,water);
659     }
660     print_top_system(out, title);
661   
662   if (nmol) {
663     fprintf(out,"[ %s ]\n",dir2str(d_molecules));
664     fprintf(out,"; %-15s %5s\n","Compound","#mols");
665     for (i=0; (i<nmol); i++)
666       fprintf(out,"%-15s %5d\n",mols[i].name,mols[i].nr);
667   }
668 }
669
670 void write_top(FILE *out, char *pr,char *molname,
671                t_atoms *at,gmx_bool bRTPresname,
672                int bts[],t_params plist[],t_excls excls[],
673                gpp_atomtype_t atype,int *cgnr, int nrexcl)
674      /* NOTE: nrexcl is not the size of *excl! */
675 {
676   if (at && atype && cgnr) {
677     fprintf(out,"[ %s ]\n",dir2str(d_moleculetype));
678     fprintf(out,"; %-15s %5s\n","Name","nrexcl");
679     fprintf(out,"%-15s %5d\n\n",molname?molname:"Protein",nrexcl);
680     
681     print_atoms(out, atype, at, cgnr, bRTPresname);
682     print_bondeds(out,at->nr,d_bonds,      F_BONDS,    bts[ebtsBONDS], plist);
683     print_bondeds(out,at->nr,d_constraints,F_CONSTR,   0,              plist);
684     print_bondeds(out,at->nr,d_constraints,F_CONSTRNC, 0,              plist);
685     print_bondeds(out,at->nr,d_pairs,      F_LJ14,     0,              plist);
686     print_excl(out,at->nr,excls);
687     print_bondeds(out,at->nr,d_angles,     F_ANGLES,   bts[ebtsANGLES],plist);
688     print_bondeds(out,at->nr,d_dihedrals,  F_PDIHS,    bts[ebtsPDIHS], plist);
689     print_bondeds(out,at->nr,d_dihedrals,  F_IDIHS,    bts[ebtsIDIHS], plist);
690     print_bondeds(out,at->nr,d_cmap,       F_CMAP,     bts[ebtsCMAP],  plist);
691     print_bondeds(out,at->nr,d_polarization,F_POLARIZATION,   0,       plist);
692     print_bondeds(out,at->nr,d_thole_polarization,F_THOLE_POL,0,       plist);
693     print_bondeds(out,at->nr,d_vsites2,    F_VSITE2,   0,              plist);
694     print_bondeds(out,at->nr,d_vsites3,    F_VSITE3,   0,              plist);
695     print_bondeds(out,at->nr,d_vsites3,    F_VSITE3FD, 0,              plist);
696     print_bondeds(out,at->nr,d_vsites3,    F_VSITE3FAD,0,              plist);
697     print_bondeds(out,at->nr,d_vsites3,    F_VSITE3OUT,0,              plist);
698     print_bondeds(out,at->nr,d_vsites4,    F_VSITE4FD, 0,              plist);
699     print_bondeds(out,at->nr,d_vsites4,    F_VSITE4FDN, 0,             plist);
700     
701     if (pr)
702       print_top_posre(out,pr);
703   }
704 }
705
706 static atom_id search_res_atom(const char *type,int resind,
707                    t_atoms *atoms,
708                                const char *bondtype,gmx_bool bAllowMissing)
709 {
710   int i;
711
712   for(i=0; (i<atoms->nr); i++)
713   {
714     if (atoms->atom[i].resind == resind)
715     {
716       return search_atom(type,i,atoms,bondtype,bAllowMissing);
717     }
718   }
719   
720   return NO_ATID;
721 }
722
723 static void do_ssbonds(t_params *ps,t_atoms *atoms,
724                        int nssbonds,t_ssbond *ssbonds,gmx_bool bAllowMissing)
725 {
726   int     i,ri,rj;
727   atom_id ai,aj;
728   
729   for(i=0; (i<nssbonds); i++) {
730     ri = ssbonds[i].res1;
731     rj = ssbonds[i].res2;
732     ai = search_res_atom(ssbonds[i].a1,ri,atoms,
733                          "special bond",bAllowMissing);
734     aj = search_res_atom(ssbonds[i].a2,rj,atoms,
735                          "special bond",bAllowMissing);
736     if ((ai == NO_ATID) || (aj == NO_ATID))
737       gmx_fatal(FARGS,"Trying to make impossible special bond (%s-%s)!",
738                   ssbonds[i].a1,ssbonds[i].a2);
739     add_param(ps,ai,aj,NULL,NULL);
740   }
741 }
742
743 static void at2bonds(t_params *psb, t_hackblock *hb,
744                      t_atoms *atoms,
745                      rvec x[],
746                      real long_bond_dist, real short_bond_dist,
747                      gmx_bool bAllowMissing)
748 {
749   int     resind,i,j,k;
750   atom_id ai,aj;
751   real    dist2, long_bond_dist2, short_bond_dist2;
752   const char *ptr;
753
754   long_bond_dist2  = sqr(long_bond_dist);
755   short_bond_dist2 = sqr(short_bond_dist);
756
757   if (debug)
758     ptr = "bond";
759   else
760     ptr = "check";
761
762   fprintf(stderr,"Making bonds...\n");
763   i=0;
764   for(resind=0; (resind < atoms->nres) && (i<atoms->nr); resind++) {
765     /* add bonds from list of bonded interactions */
766     for(j=0; j < hb[resind].rb[ebtsBONDS].nb; j++) {
767       /* Unfortunately we can not issue errors or warnings
768        * for missing atoms in bonds, as the hydrogens and terminal atoms
769        * have not been added yet.
770        */
771       ai=search_atom(hb[resind].rb[ebtsBONDS].b[j].AI,i,atoms,
772                      ptr,TRUE);
773       aj=search_atom(hb[resind].rb[ebtsBONDS].b[j].AJ,i,atoms,
774                      ptr,TRUE);
775       if (ai != NO_ATID && aj != NO_ATID) {
776           dist2 = distance2(x[ai],x[aj]);
777           if (dist2 > long_bond_dist2 )
778           {
779               fprintf(stderr,"Warning: Long Bond (%d-%d = %g nm)\n",
780                       ai+1,aj+1,sqrt(dist2));
781           }
782           else if (dist2 < short_bond_dist2 )
783           {
784               fprintf(stderr,"Warning: Short Bond (%d-%d = %g nm)\n",
785                       ai+1,aj+1,sqrt(dist2));
786           }
787           add_param(psb,ai,aj,NULL,hb[resind].rb[ebtsBONDS].b[j].s);
788       }
789     }
790     /* add bonds from list of hacks (each added atom gets a bond) */
791     while( (i<atoms->nr) && (atoms->atom[i].resind == resind) ) {
792       for(j=0; j < hb[resind].nhack; j++)
793         if ( ( hb[resind].hack[j].tp > 0 ||
794                hb[resind].hack[j].oname==NULL ) &&
795              strcmp(hb[resind].hack[j].AI,*(atoms->atomname[i])) == 0 ) {
796           switch(hb[resind].hack[j].tp) {
797           case 9:          /* COOH terminus */
798             add_param(psb,i,i+1,NULL,NULL);     /* C-O  */
799             add_param(psb,i,i+2,NULL,NULL);     /* C-OA */
800             add_param(psb,i+2,i+3,NULL,NULL);   /* OA-H */
801             break;
802           default:
803             for(k=0; (k<hb[resind].hack[j].nr); k++)
804               add_param(psb,i,i+k+1,NULL,NULL);
805           }
806         }
807       i++;
808     }
809     /* we're now at the start of the next residue */
810   }
811 }
812
813 static int pcompar(const void *a, const void *b)
814 {
815   t_param *pa,*pb;
816   int     d;
817   pa=(t_param *)a;
818   pb=(t_param *)b;
819   
820   d = pa->AI - pb->AI;
821   if (d == 0) 
822     d = pa->AJ - pb->AJ;
823   if (d == 0)
824     return strlen(pb->s) - strlen(pa->s);
825   else
826     return d;
827 }
828
829 static void clean_bonds(t_params *ps)
830 {
831   int     i,j;
832   atom_id a;
833   
834   if (ps->nr > 0) {
835     /* swap atomnumbers in bond if first larger than second: */
836     for(i=0; (i<ps->nr); i++)
837       if ( ps->param[i].AJ < ps->param[i].AI ) {
838         a = ps->param[i].AI;
839         ps->param[i].AI = ps->param[i].AJ;
840         ps->param[i].AJ = a;
841       }
842     
843     /* Sort bonds */
844     qsort(ps->param,ps->nr,(size_t)sizeof(ps->param[0]),pcompar);
845     
846     /* remove doubles, keep the first one always. */
847     j = 1;
848     for(i=1; (i<ps->nr); i++) {
849       if ((ps->param[i].AI != ps->param[j-1].AI) ||
850           (ps->param[i].AJ != ps->param[j-1].AJ) ) {
851         if (j != i) {
852           cp_param(&(ps->param[j]),&(ps->param[i]));
853         }
854         j++;
855       } 
856     }
857     fprintf(stderr,"Number of bonds was %d, now %d\n",ps->nr,j);
858     ps->nr=j;
859   }
860   else
861     fprintf(stderr,"No bonds\n");
862 }
863
864 void print_sums(t_atoms *atoms, gmx_bool bSystem)
865 {
866   double m,qtot;
867   int    i;
868   const char *where;
869   
870   if (bSystem)
871     where=" in system";
872   else
873     where="";
874   
875   m=0;
876   qtot=0;
877   for(i=0; (i<atoms->nr); i++) {
878     m+=atoms->atom[i].m;
879     qtot+=atoms->atom[i].q;
880   }
881   
882   fprintf(stderr,"Total mass%s %.3f a.m.u.\n",where,m);
883   fprintf(stderr,"Total charge%s %.3f e\n",where,qtot);
884 }
885
886 static void check_restp_type(const char *name,int t1,int t2)
887 {
888     if (t1 != t2)
889     {
890         gmx_fatal(FARGS,"Residues in one molecule have a different '%s' type: %d and %d",name,t1,t2);
891     }
892 }
893
894 static void check_restp_types(t_restp *r0,t_restp *r1)
895 {
896     int i;
897
898     check_restp_type("all dihedrals",r0->bAlldih,r1->bAlldih);
899     check_restp_type("nrexcl",r0->nrexcl,r1->nrexcl);
900     check_restp_type("HH14",r0->HH14,r1->HH14);
901     check_restp_type("remove dihedrals",r0->bRemoveDih,r1->bRemoveDih);
902
903     for(i=0; i<ebtsNR; i++)
904     {
905         check_restp_type(btsNames[i],r0->rb[i].type,r1->rb[i].type);
906     }
907 }
908
909 void add_atom_to_restp(t_restp *restp,int resnr,int at_start,const t_hack *hack)
910 {
911     char buf[STRLEN];
912     int  k;
913     const char *Hnum="123456";
914
915     /*if (debug)
916     {
917         fprintf(debug,"adding atom(s) %s to atom %s in res %d%s in rtp\n",
918                 hack->nname,
919                 *restp->atomname[at_start], resnr, restp->resname);
920                 }*/
921     strcpy(buf, hack->nname);
922     buf[strlen(buf)+1]='\0';
923     if ( hack->nr > 1 )
924     {
925         buf[strlen(buf)]='-';
926     }
927     /* make space */
928     restp->natom += hack->nr;
929     srenew(restp->atom,     restp->natom);
930     srenew(restp->atomname, restp->natom);
931     srenew(restp->cgnr,     restp->natom);
932     /* shift the rest */
933     for(k=restp->natom-1; k > at_start+hack->nr; k--)
934     {
935         restp->atom[k] =
936             restp->atom    [k - hack->nr];
937         restp->atomname[k] =
938             restp->atomname[k - hack->nr];
939         restp->cgnr[k] =
940             restp->cgnr    [k - hack->nr];
941     }
942     /* now add them */
943     for(k=0; k < hack->nr; k++)
944     {
945         /* set counter in atomname */
946         if ( hack->nr > 1 )
947         {
948             buf[strlen(buf)-1] = Hnum[k];
949         }
950         snew( restp->atomname[at_start+1+k], 1);
951         restp->atom     [at_start+1+k] = *hack->atom;
952         *restp->atomname[at_start+1+k] = strdup(buf);
953         if ( hack->cgnr != NOTSET )
954         {
955             restp->cgnr   [at_start+1+k] = hack->cgnr;
956         }
957         else
958         {
959             restp->cgnr   [at_start+1+k] = restp->cgnr[at_start];
960         }
961     }
962 }
963
964 void get_hackblocks_rtp(t_hackblock **hb, t_restp **restp, 
965                                int nrtp, t_restp rtp[],
966                                int nres, t_resinfo *resinfo, 
967                                int nterpairs,
968                                t_hackblock **ntdb, t_hackblock **ctdb,
969                                int *rn, int *rc)
970 {
971   int i, j, k, l;
972   char *key;
973   t_restp *res;
974   char buf[STRLEN];
975   const char *Hnum="123456";
976   int tern,terc;
977   gmx_bool bN,bC,bRM;
978
979   snew(*hb,nres);
980   snew(*restp,nres);
981   /* first the termini */
982   for(i=0; i<nterpairs; i++) {
983       if (rn[i] >= 0 && ntdb[i] != NULL) {
984           copy_t_hackblock(ntdb[i], &(*hb)[rn[i]]);
985       }
986       if (rc[i] >= 0 && ctdb[i] != NULL) {
987           merge_t_hackblock(ctdb[i], &(*hb)[rc[i]]);
988       }
989   }  
990
991   /* then the whole rtp */
992   for(i=0; i < nres; i++) {
993     /* Here we allow a mismatch of one character when looking for the rtp entry.
994      * For such a mismatch there should be only one mismatching name.
995      * This is mainly useful for small molecules such as ions.
996      * Note that this will usually not work for protein, DNA and RNA,
997      * since there the residue names should be listed in residuetypes.dat
998      * and an error will have been generated earlier in the process.
999      */
1000     key = *resinfo[i].rtp;
1001     snew(resinfo[i].rtp,1);
1002     *resinfo[i].rtp = search_rtp(key,nrtp,rtp);
1003     res = get_restp(*resinfo[i].rtp,nrtp,rtp);
1004     copy_t_restp(res, &(*restp)[i]);
1005
1006     /* Check that we do not have different bonded types in one molecule */
1007     check_restp_types(&(*restp)[0],&(*restp)[i]);
1008
1009     tern = -1;
1010     for(j=0; j<nterpairs && tern==-1; j++) {
1011         if (i == rn[j]) {
1012             tern = j;
1013         }
1014     }
1015     terc = -1;
1016     for(j=0; j<nterpairs && terc == -1; j++) {
1017         if (i == rc[j]) {
1018             terc = j;
1019         }
1020     }
1021     bRM = merge_t_bondeds(res->rb, (*hb)[i].rb,tern>=0,terc>=0);
1022
1023     if (bRM && ((tern >= 0 && ntdb[tern] == NULL) ||
1024                 (terc >= 0 && ctdb[terc] == NULL))) {
1025         gmx_fatal(FARGS,"There is a dangling bond at at least one of the terminal ends and the force field does not provide terminal entries or files. Edit a .n.tdb and/or .c.tdb file.");
1026     }
1027     if (bRM && ((tern >= 0 && ntdb[tern]->nhack == 0) ||
1028                 (terc >= 0 && ctdb[terc]->nhack == 0))) {
1029         gmx_fatal(FARGS,"There is a dangling bond at at least one of the terminal ends. Select a proper terminal entry.");
1030     }
1031   }
1032   
1033   /* now perform t_hack's on t_restp's,
1034      i.e. add's and deletes from termini database will be 
1035      added to/removed from residue topology 
1036      we'll do this on one big dirty loop, so it won't make easy reading! */
1037     for(i=0; i < nres; i++)
1038     {
1039         for(j=0; j < (*hb)[i].nhack; j++)
1040         {
1041             if ( (*hb)[i].hack[j].nr )
1042             {
1043                 /* find atom in restp */
1044                 for(l=0; l < (*restp)[i].natom; l++)
1045                     if ( ( (*hb)[i].hack[j].oname==NULL && 
1046                            strcmp((*hb)[i].hack[j].AI, *(*restp)[i].atomname[l])==0 ) ||
1047                          ( (*hb)[i].hack[j].oname!=NULL &&
1048                            strcmp((*hb)[i].hack[j].oname,*(*restp)[i].atomname[l])==0 ) )
1049                         break;
1050                 if (l == (*restp)[i].natom)
1051                 {
1052                     /* If we are doing an atom rename only, we don't need
1053                      * to generate a fatal error if the old name is not found
1054                      * in the rtp.
1055                      */
1056                     /* Deleting can happen also only on the input atoms,
1057                      * not necessarily always on the rtp entry.
1058                      */
1059                     if (!((*hb)[i].hack[j].oname != NULL &&
1060                           (*hb)[i].hack[j].nname != NULL) &&
1061                         !((*hb)[i].hack[j].oname != NULL &&
1062                           (*hb)[i].hack[j].nname == NULL))
1063                     {
1064                         gmx_fatal(FARGS,
1065                                   "atom %s not found in buiding block %d%s "
1066                                   "while combining tdb and rtp",
1067                                   (*hb)[i].hack[j].oname!=NULL ? 
1068                                   (*hb)[i].hack[j].oname : (*hb)[i].hack[j].AI, 
1069                                   i+1,*resinfo[i].rtp);
1070                     }
1071                 }
1072                 else
1073                 {
1074                     if ( (*hb)[i].hack[j].oname==NULL ) {
1075                         /* we're adding: */
1076                         add_atom_to_restp(&(*restp)[i],resinfo[i].nr,l,
1077                                           &(*hb)[i].hack[j]);
1078                     }
1079                     else
1080                     {
1081                         /* oname != NULL */
1082                         if ( (*hb)[i].hack[j].nname==NULL ) {
1083                             /* we're deleting */
1084                             if (debug) 
1085                                 fprintf(debug, "deleting atom %s from res %d%s in rtp\n",
1086                                         *(*restp)[i].atomname[l], 
1087                                         i+1,(*restp)[i].resname);
1088                             /* shift the rest */
1089                             (*restp)[i].natom--;
1090                             for(k=l; k < (*restp)[i].natom; k++) {
1091                                 (*restp)[i].atom    [k] = (*restp)[i].atom    [k+1];
1092                                 (*restp)[i].atomname[k] = (*restp)[i].atomname[k+1];
1093                                 (*restp)[i].cgnr    [k] = (*restp)[i].cgnr    [k+1];
1094                             }
1095                             /* give back space */
1096                             srenew((*restp)[i].atom,     (*restp)[i].natom);
1097                             srenew((*restp)[i].atomname, (*restp)[i].natom);
1098                             srenew((*restp)[i].cgnr,     (*restp)[i].natom);
1099                         } else { /* nname != NULL */
1100                             /* we're replacing */
1101                             if (debug) 
1102                                 fprintf(debug, "replacing atom %s by %s in res %d%s in rtp\n",
1103                                         *(*restp)[i].atomname[l], (*hb)[i].hack[j].nname, 
1104                                         i+1,(*restp)[i].resname);
1105                             snew( (*restp)[i].atomname[l], 1);
1106                             (*restp)[i].atom[l]      =       *(*hb)[i].hack[j].atom;
1107                             *(*restp)[i].atomname[l] = strdup((*hb)[i].hack[j].nname);
1108                             if ( (*hb)[i].hack[j].cgnr != NOTSET )
1109                                 (*restp)[i].cgnr   [l] = (*hb)[i].hack[j].cgnr;
1110                         }
1111                     }
1112                 }
1113             }
1114         }
1115     }
1116 }
1117
1118 static gmx_bool atomname_cmp_nr(const char *anm,t_hack *hack,int *nr)
1119 {
1120
1121     if (hack->nr == 1)
1122     {
1123         *nr = 0;
1124         
1125         return (gmx_strcasecmp(anm,hack->nname) == 0);
1126     }
1127     else
1128     {
1129         if (isdigit(anm[strlen(anm)-1]))
1130         {
1131             *nr = anm[strlen(anm)-1] - '0';
1132         }
1133         else
1134         {
1135             *nr = 0;
1136         }
1137         if (*nr <= 0 || *nr > hack->nr)
1138         {
1139             return FALSE;
1140         }
1141         else
1142         {
1143             return (strlen(anm) == strlen(hack->nname) + 1 &&
1144                     gmx_strncasecmp(anm,hack->nname,strlen(hack->nname)) == 0);
1145         }
1146     }
1147 }
1148
1149 static gmx_bool match_atomnames_with_rtp_atom(t_atoms *pdba,rvec *x,int atind,
1150                                           t_restp *rptr,t_hackblock *hbr,
1151                                           gmx_bool bVerbose)
1152 {
1153     int  resnr;
1154     int  i,j,k;
1155     char *oldnm,*newnm;
1156     int  anmnr;
1157     char *start_at,buf[STRLEN];
1158     int  start_nr;
1159     gmx_bool bReplaceReplace,bFoundInAdd;
1160     gmx_bool bDeleted;
1161
1162     oldnm = *pdba->atomname[atind];
1163     resnr = pdba->resinfo[pdba->atom[atind].resind].nr;
1164
1165     bDeleted = FALSE;
1166     for(j=0; j<hbr->nhack; j++)
1167     {
1168         if (hbr->hack[j].oname != NULL && hbr->hack[j].nname != NULL &&
1169             gmx_strcasecmp(oldnm,hbr->hack[j].oname) == 0)
1170         {
1171             /* This is a replace entry. */
1172             /* Check if we are not replacing a replaced atom. */
1173             bReplaceReplace = FALSE;
1174             for(k=0; k<hbr->nhack; k++) {
1175                 if (k != j &&
1176                     hbr->hack[k].oname != NULL && hbr->hack[k].nname != NULL &&
1177                     gmx_strcasecmp(hbr->hack[k].nname,hbr->hack[j].oname) == 0)
1178                 {
1179                     /* The replace in hack[j] replaces an atom that
1180                      * was already replaced in hack[k], we do not want
1181                      * second or higher level replaces at this stage.
1182                      */
1183                     bReplaceReplace = TRUE;
1184                 }
1185             }
1186             if (bReplaceReplace)
1187             {
1188                 /* Skip this replace. */
1189                 continue;
1190             }
1191
1192             /* This atom still has the old name, rename it */
1193             newnm = hbr->hack[j].nname;
1194             for(k=0; k<rptr->natom; k++)
1195             {
1196                 if (gmx_strcasecmp(newnm,*rptr->atomname[k]) == 0)
1197                 {
1198                     break;
1199                 }
1200             }
1201             if (k == rptr->natom)
1202             {
1203                 /* The new name is not present in the rtp.
1204                  * We need to apply the replace also to the rtp entry.
1205                  */
1206                 
1207                 /* We need to find the add hack that can add this atom
1208                  * to find out after which atom it should be added.
1209                  */
1210                 bFoundInAdd = FALSE;
1211                 for(k=0; k<hbr->nhack; k++)
1212                 {
1213                     if (hbr->hack[k].oname == NULL &&
1214                         hbr->hack[k].nname != NULL &&
1215                         atomname_cmp_nr(newnm,&hbr->hack[k],&anmnr))
1216                     {
1217                         if (anmnr <= 1)
1218                         {
1219                             start_at = hbr->hack[k].a[0];
1220                         }
1221                         else
1222                         {
1223                             sprintf(buf,"%s%d",hbr->hack[k].nname,anmnr-1);
1224                             start_at = buf;
1225                         }
1226                         for(start_nr=0; start_nr<rptr->natom; start_nr++)
1227                         {
1228                             if (gmx_strcasecmp(start_at,(*rptr->atomname[start_nr])) == 0)
1229                             {
1230                                 break;
1231                             }
1232                         }
1233                         if (start_nr == rptr->natom)
1234                         {
1235                             gmx_fatal(FARGS,"Could not find atom '%s' in residue building block '%s' to add atom '%s' to",
1236                                       start_at,rptr->resname,newnm);
1237                         }
1238                         /* We can add the atom after atom start_nr */
1239                         add_atom_to_restp(rptr,resnr,start_nr,
1240                                           &hbr->hack[j]);
1241                         
1242                         bFoundInAdd = TRUE;
1243                     }
1244                 }
1245
1246                 if (!bFoundInAdd)
1247                 {
1248                     gmx_fatal(FARGS,"Could not find an 'add' entry for atom named '%s' corresponding to the 'replace' entry from atom name '%s' to '%s' for tdb or hdb database of residue type '%s'",
1249                               newnm,
1250                               hbr->hack[j].oname,hbr->hack[j].nname,
1251                               rptr->resname);
1252                 }
1253             }
1254                 
1255             if (bVerbose)
1256             {
1257                 printf("Renaming atom '%s' in residue '%s' %d to '%s'\n",
1258                        oldnm,rptr->resname,resnr,newnm);
1259             }
1260             /* Rename the atom in pdba */
1261             snew(pdba->atomname[atind],1);
1262             *pdba->atomname[atind] = strdup(newnm);
1263         }
1264         else if (hbr->hack[j].oname != NULL && hbr->hack[j].nname == NULL &&
1265                  gmx_strcasecmp(oldnm,hbr->hack[j].oname) == 0)
1266         {
1267             /* This is a delete entry, check if this atom is present
1268              * in the rtp entry of this residue.
1269              */
1270             for(k=0; k<rptr->natom; k++)
1271             {
1272                 if (gmx_strcasecmp(oldnm,*rptr->atomname[k]) == 0)
1273                 {
1274                     break;
1275                 }
1276             }
1277             if (k == rptr->natom)
1278             {
1279                 /* This atom is not present in the rtp entry,
1280                  * delete is now from the input pdba.
1281                  */
1282                 if (bVerbose)
1283                 {
1284                     printf("Deleting atom '%s' in residue '%s' %d\n",
1285                            oldnm,rptr->resname,resnr);
1286                 }
1287                 /* We should free the atom name,
1288                  * but it might be used multiple times in the symtab.
1289                  * sfree(pdba->atomname[atind]);
1290                  */
1291                 for(k=atind+1; k<pdba->nr; k++)
1292                 {
1293                     pdba->atom[k-1]     = pdba->atom[k];
1294                     pdba->atomname[k-1] = pdba->atomname[k];
1295                     copy_rvec(x[k],x[k-1]);
1296                 }
1297                 pdba->nr--;
1298                 bDeleted = TRUE;
1299             }
1300         }
1301     }
1302
1303     return bDeleted;
1304 }
1305     
1306 void match_atomnames_with_rtp(t_restp restp[],t_hackblock hb[],
1307                               t_atoms *pdba,rvec *x,
1308                               gmx_bool bVerbose)
1309 {
1310     int  i,j,k;
1311     char *oldnm,*newnm;
1312     int  resnr;
1313     t_restp *rptr;
1314     t_hackblock *hbr;
1315     int  anmnr;
1316     char *start_at,buf[STRLEN];
1317     int  start_nr;
1318     gmx_bool bFoundInAdd;
1319     
1320     for(i=0; i<pdba->nr; i++)
1321     {
1322         oldnm = *pdba->atomname[i];
1323         resnr = pdba->resinfo[pdba->atom[i].resind].nr;
1324         rptr  = &restp[pdba->atom[i].resind];
1325         for(j=0; (j<rptr->natom); j++)
1326         {
1327             if (gmx_strcasecmp(oldnm,*(rptr->atomname[j])) == 0)
1328             {
1329                 break;
1330             }
1331         }
1332         if (j == rptr->natom)
1333         {
1334             /* Not found yet, check if we have to rename this atom */
1335             if (match_atomnames_with_rtp_atom(pdba,x,i,
1336                                               rptr,&(hb[pdba->atom[i].resind]),
1337                                               bVerbose))
1338             {
1339                 /* We deleted this atom, decrease the atom counter by 1. */
1340                 i--;
1341             }
1342         }
1343     }
1344 }
1345
1346 #define NUM_CMAP_ATOMS 5
1347 static void gen_cmap(t_params *psb, t_restp *restp, t_atoms *atoms, gmx_residuetype_t rt)
1348 {
1349     int residx,i,j,k;
1350     const char *ptr;
1351     t_resinfo *resinfo = atoms->resinfo;
1352     int nres = atoms->nres;
1353     gmx_bool bAddCMAP;
1354     atom_id cmap_atomid[NUM_CMAP_ATOMS];
1355     int cmap_chainnum=-1, this_residue_index;
1356
1357         if (debug)
1358                 ptr = "cmap";
1359         else
1360                 ptr = "check";
1361         
1362         fprintf(stderr,"Making cmap torsions...");
1363         i=0;
1364         /* End loop at nres-1, since the very last residue does not have a +N atom, and
1365          * therefore we get a valgrind invalid 4 byte read error with atom am */
1366         for(residx=0; residx<nres-1; residx++)
1367         {
1368                 /* Add CMAP terms from the list of CMAP interactions */
1369                 for(j=0;j<restp[residx].rb[ebtsCMAP].nb; j++)
1370                 {
1371             bAddCMAP = TRUE;
1372             /* Loop over atoms in a candidate CMAP interaction and
1373              * check that they exist, are from the same chain and are
1374              * from residues labelled as protein. */
1375             for(k = 0; k < NUM_CMAP_ATOMS && bAddCMAP; k++)
1376             {
1377                 cmap_atomid[k] = search_atom(restp[residx].rb[ebtsCMAP].b[j].a[k],
1378                                              i,atoms,ptr,TRUE);
1379                 bAddCMAP = bAddCMAP && (cmap_atomid[k] != NO_ATID);
1380                 if (!bAddCMAP)
1381                 {
1382                     /* This break is necessary, because cmap_atomid[k]
1383                      * == NO_ATID cannot be safely used as an index
1384                      * into the atom array. */
1385                     break;
1386                 }
1387                 this_residue_index = atoms->atom[cmap_atomid[k]].resind;
1388                 if (0 == k)
1389                 {
1390                     cmap_chainnum = resinfo[this_residue_index].chainnum;
1391                 }
1392                 else
1393                 {
1394                     /* Does the residue for this atom have the same
1395                      * chain number as the residues for previous
1396                      * atoms? */
1397                     bAddCMAP = bAddCMAP &&
1398                         cmap_chainnum == resinfo[this_residue_index].chainnum;
1399                 }
1400                 bAddCMAP = bAddCMAP && gmx_residuetype_is_protein(rt,*(resinfo[this_residue_index].name));
1401             }
1402
1403             if(bAddCMAP)
1404             {
1405                 add_cmap_param(psb,cmap_atomid[0],cmap_atomid[1],cmap_atomid[2],cmap_atomid[3],cmap_atomid[4],restp[residx].rb[ebtsCMAP].b[j].s);
1406                         }
1407                 }
1408                 
1409                 if(residx<nres-1)
1410                 {
1411                         while(atoms->atom[i].resind<residx+1)
1412                         {
1413                                 i++;
1414                         }
1415                 }
1416         }
1417         
1418         /* Start the next residue */
1419 }
1420
1421 static void 
1422 scrub_charge_groups(int *cgnr, int natoms)
1423 {
1424         int i;
1425         
1426         for(i=0;i<natoms;i++)
1427         {
1428                 cgnr[i]=i+1;
1429         }
1430 }
1431
1432
1433 void pdb2top(FILE *top_file, char *posre_fn, char *molname,
1434              t_atoms *atoms, rvec **x, gpp_atomtype_t atype, t_symtab *tab,
1435              int nrtp, t_restp rtp[],
1436              t_restp *restp, t_hackblock *hb,
1437              int nterpairs,t_hackblock **ntdb, t_hackblock **ctdb,
1438              gmx_bool bAllowMissing,
1439              gmx_bool bVsites, gmx_bool bVsiteAromatics,
1440              const char *ff, const char *ffdir,
1441              real mHmult,
1442              int nssbonds, t_ssbond *ssbonds,
1443              real long_bond_dist, real short_bond_dist,
1444              gmx_bool bDeuterate, gmx_bool bChargeGroups, gmx_bool bCmap,
1445              gmx_bool bRenumRes,gmx_bool bRTPresname)
1446 {
1447     /*
1448   t_hackblock *hb;
1449   t_restp  *restp;
1450     */
1451   t_params plist[F_NRE];
1452   t_excls  *excls;
1453   t_nextnb nnb;
1454   int      *cgnr;
1455   int      *vsite_type;
1456   int      i,nmissat;
1457   int      bts[ebtsNR];
1458   gmx_residuetype_t rt;
1459   
1460   init_plist(plist);
1461   gmx_residuetype_init(&rt);
1462
1463   if (debug) {
1464     print_resall(debug, atoms->nres, restp, atype);
1465     dump_hb(debug, atoms->nres, hb);
1466   }
1467   
1468   /* Make bonds */
1469   at2bonds(&(plist[F_BONDS]), hb, 
1470            atoms, *x,
1471            long_bond_dist, short_bond_dist, bAllowMissing);
1472   
1473   /* specbonds: disulphide bonds & heme-his */
1474   do_ssbonds(&(plist[F_BONDS]),
1475              atoms, nssbonds, ssbonds,
1476              bAllowMissing);
1477   
1478   nmissat = name2type(atoms, &cgnr, atype, restp, rt);
1479   if (nmissat) {
1480     if (bAllowMissing)
1481       fprintf(stderr,"There were %d missing atoms in molecule %s\n",
1482               nmissat,molname);
1483     else
1484       gmx_fatal(FARGS,"There were %d missing atoms in molecule %s, if you want to use this incomplete topology anyhow, use the option -missing",
1485                   nmissat,molname);
1486   }
1487   
1488   /* Cleanup bonds (sort and rm doubles) */ 
1489   clean_bonds(&(plist[F_BONDS]));
1490   
1491   snew(vsite_type,atoms->nr);
1492   for(i=0; i<atoms->nr; i++)
1493     vsite_type[i]=NOTSET;
1494   if (bVsites) {
1495     /* determine which atoms will be vsites and add dummy masses 
1496        also renumber atom numbers in plist[0..F_NRE]! */
1497     do_vsites(nrtp, rtp, atype, atoms, tab, x, plist, 
1498               &vsite_type, &cgnr, mHmult, bVsiteAromatics, ffdir);
1499   }
1500   
1501   /* Make Angles and Dihedrals */
1502   fprintf(stderr,"Generating angles, dihedrals and pairs...\n");
1503   snew(excls,atoms->nr);
1504   init_nnb(&nnb,atoms->nr,4);
1505   gen_nnb(&nnb,plist);
1506   print_nnb(&nnb,"NNB");
1507   gen_pad(&nnb,atoms,restp[0].nrexcl,restp[0].HH14,
1508           plist,excls,hb,restp[0].bAlldih,restp[0].bRemoveDih,
1509           bAllowMissing);
1510   done_nnb(&nnb);
1511   
1512     /* Make CMAP */
1513     if (TRUE == bCmap)
1514     {
1515         gen_cmap(&(plist[F_CMAP]), restp, atoms, rt);
1516         if (plist[F_CMAP].nr > 0)
1517         {
1518             fprintf(stderr, "There are %4d cmap torsion pairs\n",
1519                     plist[F_CMAP].nr);
1520         }
1521     }
1522
1523   /* set mass of all remaining hydrogen atoms */
1524   if (mHmult != 1.0)
1525     do_h_mass(&(plist[F_BONDS]),vsite_type,atoms,mHmult,bDeuterate);
1526   sfree(vsite_type);
1527   
1528   /* Cleanup bonds (sort and rm doubles) */ 
1529   /* clean_bonds(&(plist[F_BONDS]));*/
1530    
1531   fprintf(stderr,
1532           "There are %4d dihedrals, %4d impropers, %4d angles\n"
1533           "          %4d pairs,     %4d bonds and  %4d virtual sites\n",
1534           plist[F_PDIHS].nr, plist[F_IDIHS].nr, plist[F_ANGLES].nr,
1535           plist[F_LJ14].nr, plist[F_BONDS].nr,
1536           plist[F_VSITE2].nr +
1537           plist[F_VSITE3].nr +
1538           plist[F_VSITE3FD].nr +
1539           plist[F_VSITE3FAD].nr +
1540           plist[F_VSITE3OUT].nr +
1541       plist[F_VSITE4FD].nr +
1542       plist[F_VSITE4FDN].nr );
1543   
1544   print_sums(atoms, FALSE);
1545   
1546   if (FALSE == bChargeGroups)
1547   {
1548           scrub_charge_groups(cgnr, atoms->nr);
1549   }
1550
1551     if (bRenumRes)
1552     {
1553         for(i=0; i<atoms->nres; i++) 
1554         {
1555             atoms->resinfo[i].nr = i + 1;
1556             atoms->resinfo[i].ic = ' ';
1557         }
1558     }
1559         
1560   if (top_file) {
1561     fprintf(stderr,"Writing topology\n");
1562     /* We can copy the bonded types from the first restp,
1563      * since the types have to be identical for all residues in one molecule.
1564      */
1565     for(i=0; i<ebtsNR; i++) {
1566         bts[i] = restp[0].rb[i].type;
1567     }
1568     write_top(top_file, posre_fn, molname,
1569               atoms, bRTPresname, 
1570               bts, plist, excls, atype, cgnr, restp[0].nrexcl);
1571   }
1572   
1573   /* cleaning up */
1574   free_t_hackblock(atoms->nres, &hb);
1575   free_t_restp(atoms->nres, &restp);
1576   gmx_residuetype_destroy(rt);
1577     
1578   /* we should clean up hb and restp here, but that is a *L*O*T* of work! */
1579   sfree(cgnr);
1580   for (i=0; i<F_NRE; i++)
1581     sfree(plist[i].param);
1582   for (i=0; i<atoms->nr; i++)
1583     sfree(excls[i].e);
1584   sfree(excls);
1585 }