use separate PME nodes with #nodes>18
[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     fclose(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[])
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   gmx_residuetype_t rt;
453     
454   nmissat = 0;
455
456   resind=-1;
457   bProt=FALSE;
458   bNterm=FALSE;
459   i0=0;
460   snew(*cgnr,at->nr);
461   qt=0;
462   prevcg=NOTSET;
463   curcg=0;
464   cg=-1;
465   j=NOTSET;
466   gmx_residuetype_init(&rt);
467   
468   for(i=0; (i<at->nr); i++) {
469     prevresind=resind;
470     if (at->atom[i].resind != resind) {
471       resind = at->atom[i].resind;
472       bProt = gmx_residuetype_is_protein(rt,*(at->resinfo[resind].name));
473       bNterm=bProt && (resind == 0);
474       if (resind > 0) {
475           nmissat += missing_atoms(&restp[prevresind],prevresind,at,i0,i);
476       }
477       i0=i;
478     }
479     if (at->atom[i].m == 0) {
480       if (debug)
481         fprintf(debug,"atom %d%s: curcg=%d, prevcg=%d, cg=%d\n",
482                 i+1,*(at->atomname[i]),curcg,prevcg,
483                 j==NOTSET ? NOTSET : restp[resind].cgnr[j]);
484       qt=0;
485       prevcg=cg;
486       name=*(at->atomname[i]);
487       j=search_jtype(&restp[resind],name,bNterm);
488       at->atom[i].type = restp[resind].atom[j].type;
489       at->atom[i].q    = restp[resind].atom[j].q;
490       at->atom[i].m    = get_atomtype_massA(restp[resind].atom[j].type,
491                                             atype);
492       cg = restp[resind].cgnr[j];
493       /* A charge group number -1 signals a separate charge group
494        * for this atom.
495        */
496       if ( (cg == -1) || (cg != prevcg) || (resind != prevresind) ) {
497           curcg++;
498       }
499     } else {
500       if (debug)
501         fprintf(debug,"atom %d%s: curcg=%d, qt=%g, is_int=%d\n",
502                 i+1,*(at->atomname[i]),curcg,qt,is_int(qt));
503       cg=-1;
504       if (is_int(qt)) {
505         qt=0;
506         curcg++;
507       }
508       qt+=at->atom[i].q;
509     }
510     (*cgnr)[i]=curcg;
511     at->atom[i].typeB = at->atom[i].type;
512     at->atom[i].qB    = at->atom[i].q;
513     at->atom[i].mB    = at->atom[i].m;
514   }
515   nmissat += missing_atoms(&restp[resind],resind,at,i0,i);
516
517   gmx_residuetype_destroy(rt);
518                            
519   return nmissat;
520 }
521
522 static void print_top_heavy_H(FILE *out, real mHmult)
523 {
524   if (mHmult == 2.0) 
525     fprintf(out,"; Using deuterium instead of hydrogen\n\n");
526   else if (mHmult == 4.0)
527     fprintf(out,"#define HEAVY_H\n\n");
528   else if (mHmult != 1.0)
529     fprintf(stderr,"WARNING: unsupported proton mass multiplier (%g) "
530             "in pdb2top\n",mHmult);
531 }
532
533 void print_top_comment(FILE *out,
534                        const char *filename,
535                        const char *generator,
536                        const char *ffdir,
537                        gmx_bool bITP)
538 {
539   char tmp[256]; 
540   char ffdir_parent[STRLEN];
541   char *p;
542         
543   nice_header(out,filename);
544   fprintf(out,";\tThis is a %s topology file\n;\n",bITP ? "include" : "standalone");
545   fprintf(out,";\tIt was generated using program:\n;\t%s\n;\n",
546           (NULL == generator) ? "unknown" : generator);
547   fprintf(out,";\tCommand line was:\n;\t%s\n;\n",command_line());
548
549   if(strchr(ffdir,'/')==NULL)
550   {
551       fprintf(out,";\tForce field was read from the standard Gromacs share directory.\n;\n\n");
552   }
553   else if(ffdir[0]=='.')
554   {
555       fprintf(out,";\tForce field was read from current directory or a relative path - path added.\n;\n\n");
556   }
557   else
558   {
559       strncpy(ffdir_parent,ffdir,STRLEN-1);
560       p=strrchr(ffdir_parent,'/');
561
562       *p='\0';
563       
564       fprintf(out,
565               ";\tForce field data was read from:\n"
566               ";\t%s\n"
567               ";\n"
568               ";\tNote:\n"
569               ";\tThis might be a non-standard force field location. When you use this topology, the\n"
570               ";\tforce field must either be present in the current directory, or the location\n"
571               ";\tspecified in the GMXLIB path variable or with the 'include' mdp file option.\n;\n\n",
572               ffdir_parent);
573   }
574 }
575
576 void print_top_header(FILE *out,const char *filename, 
577                       const char *title,gmx_bool bITP,const char *ffdir,real mHmult)
578 {
579     const char *p;
580     
581     print_top_comment(out,filename,title,ffdir,bITP);
582     
583     print_top_heavy_H(out, mHmult);
584     fprintf(out,"; Include forcefield parameters\n");
585
586     p=strrchr(ffdir,'/');        
587     p = (ffdir[0]=='.' || p==NULL) ? ffdir : p+1;
588
589     fprintf(out,"#include \"%s/%s\"\n\n",p,fflib_forcefield_itp());
590 }
591
592 static void print_top_posre(FILE *out,const char *pr)
593 {
594   fprintf(out,"; Include Position restraint file\n");
595   fprintf(out,"#ifdef POSRES\n");
596   fprintf(out,"#include \"%s\"\n",pr);
597   fprintf(out,"#endif\n\n");
598 }
599   
600 static void print_top_water(FILE *out,const char *ffdir,const char *water)
601 {
602   const char *p;
603   char  buf[STRLEN];
604     
605   fprintf(out,"; Include water topology\n");
606
607   p=strrchr(ffdir,'/');        
608   p = (ffdir[0]=='.' || p==NULL) ? ffdir : p+1;
609   fprintf(out,"#include \"%s/%s.itp\"\n",p,water);
610   
611   fprintf(out,"\n");
612   fprintf(out,"#ifdef POSRES_WATER\n");
613   fprintf(out,"; Position restraint for each water oxygen\n");
614   fprintf(out,"[ position_restraints ]\n");
615   fprintf(out,";%3s %5s %9s %10s %10s\n","i","funct","fcx","fcy","fcz");
616   fprintf(out,"%4d %4d %10g %10g %10g\n",1,1,1000.0,1000.0,1000.0);
617   fprintf(out,"#endif\n");
618   fprintf(out,"\n");
619
620   sprintf(buf,"%s/ions.itp",p);
621
622   if (fflib_fexist(buf))
623   {
624     fprintf(out,"; Include topology for ions\n");
625     fprintf(out,"#include \"%s\"\n",buf);
626     fprintf(out,"\n");
627   }
628 }
629
630 static void print_top_system(FILE *out, const char *title)
631 {
632   fprintf(out,"[ %s ]\n",dir2str(d_system));
633   fprintf(out,"; Name\n");
634   fprintf(out,"%s\n\n",title[0]?title:"Protein");
635 }
636
637 void print_top_mols(FILE *out,
638                     const char *title, const char *ffdir, const char *water,
639                     int nincl, char **incls, int nmol, t_mols *mols)
640 {
641   int  i;
642   char *incl;
643
644   if (nincl>0) {
645     fprintf(out,"; Include chain topologies\n");
646     for (i=0; (i<nincl); i++) {
647         incl = strrchr(incls[i],DIR_SEPARATOR);
648         if (incl == NULL) {
649             incl = incls[i];
650         } else {
651             /* Remove the path from the include name */
652             incl = incl + 1;
653         }
654       fprintf(out,"#include \"%s\"\n",incl);
655     }
656     fprintf(out,"\n");
657   }
658
659     if (water)
660     {
661       print_top_water(out,ffdir,water);
662     }
663     print_top_system(out, title);
664   
665   if (nmol) {
666     fprintf(out,"[ %s ]\n",dir2str(d_molecules));
667     fprintf(out,"; %-15s %5s\n","Compound","#mols");
668     for (i=0; (i<nmol); i++)
669       fprintf(out,"%-15s %5d\n",mols[i].name,mols[i].nr);
670   }
671 }
672
673 void write_top(FILE *out, char *pr,char *molname,
674                t_atoms *at,gmx_bool bRTPresname,
675                int bts[],t_params plist[],t_excls excls[],
676                gpp_atomtype_t atype,int *cgnr, int nrexcl)
677      /* NOTE: nrexcl is not the size of *excl! */
678 {
679   if (at && atype && cgnr) {
680     fprintf(out,"[ %s ]\n",dir2str(d_moleculetype));
681     fprintf(out,"; %-15s %5s\n","Name","nrexcl");
682     fprintf(out,"%-15s %5d\n\n",molname?molname:"Protein",nrexcl);
683     
684     print_atoms(out, atype, at, cgnr, bRTPresname);
685     print_bondeds(out,at->nr,d_bonds,      F_BONDS,    bts[ebtsBONDS], plist);
686     print_bondeds(out,at->nr,d_constraints,F_CONSTR,   0,              plist);
687     print_bondeds(out,at->nr,d_constraints,F_CONSTRNC, 0,              plist);
688     print_bondeds(out,at->nr,d_pairs,      F_LJ14,     0,              plist);
689     print_excl(out,at->nr,excls);
690     print_bondeds(out,at->nr,d_angles,     F_ANGLES,   bts[ebtsANGLES],plist);
691     print_bondeds(out,at->nr,d_dihedrals,  F_PDIHS,    bts[ebtsPDIHS], plist);
692     print_bondeds(out,at->nr,d_dihedrals,  F_IDIHS,    bts[ebtsIDIHS], plist);
693     print_bondeds(out,at->nr,d_cmap,       F_CMAP,     bts[ebtsCMAP],  plist);
694     print_bondeds(out,at->nr,d_polarization,F_POLARIZATION,   0,       plist);
695     print_bondeds(out,at->nr,d_thole_polarization,F_THOLE_POL,0,       plist);
696     print_bondeds(out,at->nr,d_vsites2,    F_VSITE2,   0,              plist);
697     print_bondeds(out,at->nr,d_vsites3,    F_VSITE3,   0,              plist);
698     print_bondeds(out,at->nr,d_vsites3,    F_VSITE3FD, 0,              plist);
699     print_bondeds(out,at->nr,d_vsites3,    F_VSITE3FAD,0,              plist);
700     print_bondeds(out,at->nr,d_vsites3,    F_VSITE3OUT,0,              plist);
701     print_bondeds(out,at->nr,d_vsites4,    F_VSITE4FD, 0,              plist);
702     print_bondeds(out,at->nr,d_vsites4,    F_VSITE4FDN, 0,             plist);
703     
704     if (pr)
705       print_top_posre(out,pr);
706   }
707 }
708
709 static atom_id search_res_atom(const char *type,int resind,
710                                int natom,t_atom at[],
711                                char ** const *aname,
712                                const char *bondtype,gmx_bool bAllowMissing)
713 {
714   int i;
715
716   for(i=0; (i<natom); i++)
717     if (at[i].resind == resind)
718       return search_atom(type,i,natom,at,aname,bondtype,bAllowMissing);
719   
720   return NO_ATID;
721 }
722
723 static void do_ssbonds(t_params *ps,int natoms,t_atom atom[],char **aname[],
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,natoms,atom,aname,
733                          "special bond",bAllowMissing);
734     aj = search_res_atom(ssbonds[i].a2,rj,natoms,atom,aname,
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 gmx_bool inter_res_bond(const t_rbonded *b)
744 {
745     return (b->AI[0] == '-' || b->AI[0] == '+' ||
746             b->AJ[0] == '-' || b->AJ[0] == '+');
747 }
748
749 static void at2bonds(t_params *psb, t_hackblock *hb,
750                      int natoms, t_atom atom[], char **aname[], 
751                      int nres, rvec x[], 
752                      real long_bond_dist, real short_bond_dist,
753                      gmx_bool bAllowMissing)
754 {
755   int     resind,i,j,k;
756   atom_id ai,aj;
757   real    dist2, long_bond_dist2, short_bond_dist2;
758   const char *ptr;
759
760   long_bond_dist2  = sqr(long_bond_dist);
761   short_bond_dist2 = sqr(short_bond_dist);
762
763   if (debug)
764     ptr = "bond";
765   else
766     ptr = "check";
767
768   fprintf(stderr,"Making bonds...\n");
769   i=0;
770   for(resind=0; (resind < nres) && (i<natoms); resind++) {
771     /* add bonds from list of bonded interactions */
772     for(j=0; j < hb[resind].rb[ebtsBONDS].nb; j++) {
773       /* Unfortunately we can not issue errors or warnings
774        * for missing atoms in bonds, as the hydrogens and terminal atoms
775        * have not been added yet.
776        */
777       ai=search_atom(hb[resind].rb[ebtsBONDS].b[j].AI,i,natoms,atom,aname,
778                      ptr,TRUE);
779       aj=search_atom(hb[resind].rb[ebtsBONDS].b[j].AJ,i,natoms,atom,aname,
780                      ptr,TRUE);
781       if (ai != NO_ATID && aj != NO_ATID) {
782           dist2 = distance2(x[ai],x[aj]);
783           if (dist2 > long_bond_dist2 )
784           {
785               fprintf(stderr,"Warning: Long Bond (%d-%d = %g nm)\n",
786                       ai+1,aj+1,sqrt(dist2));
787           }
788           else if (dist2 < short_bond_dist2 )
789           {
790               fprintf(stderr,"Warning: Short Bond (%d-%d = %g nm)\n",
791                       ai+1,aj+1,sqrt(dist2));
792           }
793           add_param(psb,ai,aj,NULL,hb[resind].rb[ebtsBONDS].b[j].s);
794       }
795     }
796     /* add bonds from list of hacks (each added atom gets a bond) */
797     while( (i<natoms) && (atom[i].resind == resind) ) {
798       for(j=0; j < hb[resind].nhack; j++)
799         if ( ( hb[resind].hack[j].tp > 0 ||
800                hb[resind].hack[j].oname==NULL ) &&
801              strcmp(hb[resind].hack[j].AI,*(aname[i])) == 0 ) {
802           switch(hb[resind].hack[j].tp) {
803           case 9:          /* COOH terminus */
804             add_param(psb,i,i+1,NULL,NULL);     /* C-O  */
805             add_param(psb,i,i+2,NULL,NULL);     /* C-OA */
806             add_param(psb,i+2,i+3,NULL,NULL);   /* OA-H */
807             break;
808           default:
809             for(k=0; (k<hb[resind].hack[j].nr); k++)
810               add_param(psb,i,i+k+1,NULL,NULL);
811           }
812         }
813       i++;
814     }
815     /* we're now at the start of the next residue */
816   }
817 }
818
819 static int pcompar(const void *a, const void *b)
820 {
821   t_param *pa,*pb;
822   int     d;
823   pa=(t_param *)a;
824   pb=(t_param *)b;
825   
826   d = pa->AI - pb->AI;
827   if (d == 0) 
828     d = pa->AJ - pb->AJ;
829   if (d == 0)
830     return strlen(pb->s) - strlen(pa->s);
831   else
832     return d;
833 }
834
835 static void clean_bonds(t_params *ps)
836 {
837   int     i,j;
838   atom_id a;
839   
840   if (ps->nr > 0) {
841     /* swap atomnumbers in bond if first larger than second: */
842     for(i=0; (i<ps->nr); i++)
843       if ( ps->param[i].AJ < ps->param[i].AI ) {
844         a = ps->param[i].AI;
845         ps->param[i].AI = ps->param[i].AJ;
846         ps->param[i].AJ = a;
847       }
848     
849     /* Sort bonds */
850     qsort(ps->param,ps->nr,(size_t)sizeof(ps->param[0]),pcompar);
851     
852     /* remove doubles, keep the first one always. */
853     j = 1;
854     for(i=1; (i<ps->nr); i++) {
855       if ((ps->param[i].AI != ps->param[j-1].AI) ||
856           (ps->param[i].AJ != ps->param[j-1].AJ) ) {
857         if (j != i) {
858           cp_param(&(ps->param[j]),&(ps->param[i]));
859         }
860         j++;
861       } 
862     }
863     fprintf(stderr,"Number of bonds was %d, now %d\n",ps->nr,j);
864     ps->nr=j;
865   }
866   else
867     fprintf(stderr,"No bonds\n");
868 }
869
870 void print_sums(t_atoms *atoms, gmx_bool bSystem)
871 {
872   double m,qtot;
873   int    i;
874   const char *where;
875   
876   if (bSystem)
877     where=" in system";
878   else
879     where="";
880   
881   m=0;
882   qtot=0;
883   for(i=0; (i<atoms->nr); i++) {
884     m+=atoms->atom[i].m;
885     qtot+=atoms->atom[i].q;
886   }
887   
888   fprintf(stderr,"Total mass%s %.3f a.m.u.\n",where,m);
889   fprintf(stderr,"Total charge%s %.3f e\n",where,qtot);
890 }
891
892 static void check_restp_type(const char *name,int t1,int t2)
893 {
894     if (t1 != t2)
895     {
896         gmx_fatal(FARGS,"Residues in one molecule have a different '%s' type: %d and %d",name,t1,t2);
897     }
898 }
899
900 static void check_restp_types(t_restp *r0,t_restp *r1)
901 {
902     int i;
903
904     check_restp_type("all dihedrals",r0->bAlldih,r1->bAlldih);
905     check_restp_type("nrexcl",r0->nrexcl,r1->nrexcl);
906     check_restp_type("HH14",r0->HH14,r1->HH14);
907     check_restp_type("remove dihedrals",r0->bRemoveDih,r1->bRemoveDih);
908
909     for(i=0; i<ebtsNR; i++)
910     {
911         check_restp_type(btsNames[i],r0->rb[i].type,r1->rb[i].type);
912     }
913 }
914
915 void add_atom_to_restp(t_restp *restp,int resnr,int at_start,const t_hack *hack)
916 {
917     char buf[STRLEN];
918     int  k;
919     const char *Hnum="123456";
920
921     /*if (debug)
922     {
923         fprintf(debug,"adding atom(s) %s to atom %s in res %d%s in rtp\n",
924                 hack->nname,
925                 *restp->atomname[at_start], resnr, restp->resname);
926                 }*/
927     strcpy(buf, hack->nname);
928     buf[strlen(buf)+1]='\0';
929     if ( hack->nr > 1 )
930     {
931         buf[strlen(buf)]='-';
932     }
933     /* make space */
934     restp->natom += hack->nr;
935     srenew(restp->atom,     restp->natom);
936     srenew(restp->atomname, restp->natom);
937     srenew(restp->cgnr,     restp->natom);
938     /* shift the rest */
939     for(k=restp->natom-1; k > at_start+hack->nr; k--)
940     {
941         restp->atom[k] =
942             restp->atom    [k - hack->nr];
943         restp->atomname[k] =
944             restp->atomname[k - hack->nr];
945         restp->cgnr[k] =
946             restp->cgnr    [k - hack->nr];
947     }
948     /* now add them */
949     for(k=0; k < hack->nr; k++)
950     {
951         /* set counter in atomname */
952         if ( hack->nr > 1 )
953         {
954             buf[strlen(buf)-1] = Hnum[k];
955         }
956         snew( restp->atomname[at_start+1+k], 1);
957         restp->atom     [at_start+1+k] = *hack->atom;
958         *restp->atomname[at_start+1+k] = strdup(buf);
959         if ( hack->cgnr != NOTSET )
960         {
961             restp->cgnr   [at_start+1+k] = hack->cgnr;
962         }
963         else
964         {
965             restp->cgnr   [at_start+1+k] = restp->cgnr[at_start];
966         }
967     }
968 }
969
970 void get_hackblocks_rtp(t_hackblock **hb, t_restp **restp, 
971                                int nrtp, t_restp rtp[],
972                                int nres, t_resinfo *resinfo, 
973                                int nterpairs,
974                                t_hackblock **ntdb, t_hackblock **ctdb,
975                                int *rn, int *rc)
976 {
977   int i, j, k, l;
978   char *key;
979   t_restp *res;
980   char buf[STRLEN];
981   const char *Hnum="123456";
982   int tern,terc;
983   gmx_bool bN,bC,bRM;
984
985   snew(*hb,nres);
986   snew(*restp,nres);
987   /* first the termini */
988   for(i=0; i<nterpairs; i++) {
989       if (rn[i] >= 0 && ntdb[i] != NULL) {
990           copy_t_hackblock(ntdb[i], &(*hb)[rn[i]]);
991       }
992       if (rc[i] >= 0 && ctdb[i] != NULL) {
993           merge_t_hackblock(ctdb[i], &(*hb)[rc[i]]);
994       }
995   }  
996
997   /* then the whole rtp */
998   for(i=0; i < nres; i++) {
999     /* Here we allow a mismatch of one character when looking for the rtp entry.
1000      * For such a mismatch there should be only one mismatching name.
1001      * This is mainly useful for small molecules such as ions.
1002      * Note that this will usually not work for protein, DNA and RNA,
1003      * since there the residue names should be listed in residuetypes.dat
1004      * and an error will have been generated earlier in the process.
1005      */
1006     key = *resinfo[i].rtp;
1007     snew(resinfo[i].rtp,1);
1008     *resinfo[i].rtp = search_rtp(key,nrtp,rtp);
1009     res = get_restp(*resinfo[i].rtp,nrtp,rtp);
1010     copy_t_restp(res, &(*restp)[i]);
1011
1012     /* Check that we do not have different bonded types in one molecule */
1013     check_restp_types(&(*restp)[0],&(*restp)[i]);
1014
1015     tern = -1;
1016     for(j=0; j<nterpairs && tern==-1; j++) {
1017         if (i == rn[j]) {
1018             tern = j;
1019         }
1020     }
1021     terc = -1;
1022     for(j=0; j<nterpairs && terc == -1; j++) {
1023         if (i == rc[j]) {
1024             terc = j;
1025         }
1026     }
1027     bRM = merge_t_bondeds(res->rb, (*hb)[i].rb,tern>=0,terc>=0);
1028
1029     if (bRM && ((tern >= 0 && ntdb[tern] == NULL) ||
1030                 (terc >= 0 && ctdb[terc] == NULL))) {
1031         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.");
1032     }
1033     if (bRM && ((tern >= 0 && ntdb[tern]->nhack == 0) ||
1034                 (terc >= 0 && ctdb[terc]->nhack == 0))) {
1035         gmx_fatal(FARGS,"There is a dangling bond at at least one of the terminal ends. Select a proper terminal entry.");
1036     }
1037   }
1038   
1039   /* now perform t_hack's on t_restp's,
1040      i.e. add's and deletes from termini database will be 
1041      added to/removed from residue topology 
1042      we'll do this on one big dirty loop, so it won't make easy reading! */
1043     for(i=0; i < nres; i++)
1044     {
1045         for(j=0; j < (*hb)[i].nhack; j++)
1046         {
1047             if ( (*hb)[i].hack[j].nr )
1048             {
1049                 /* find atom in restp */
1050                 for(l=0; l < (*restp)[i].natom; l++)
1051                     if ( ( (*hb)[i].hack[j].oname==NULL && 
1052                            strcmp((*hb)[i].hack[j].AI, *(*restp)[i].atomname[l])==0 ) ||
1053                          ( (*hb)[i].hack[j].oname!=NULL &&
1054                            strcmp((*hb)[i].hack[j].oname,*(*restp)[i].atomname[l])==0 ) )
1055                         break;
1056                 if (l == (*restp)[i].natom)
1057                 {
1058                     /* If we are doing an atom rename only, we don't need
1059                      * to generate a fatal error if the old name is not found
1060                      * in the rtp.
1061                      */
1062                     /* Deleting can happen also only on the input atoms,
1063                      * not necessarily always on the rtp entry.
1064                      */
1065                     if (!((*hb)[i].hack[j].oname != NULL &&
1066                           (*hb)[i].hack[j].nname != NULL) &&
1067                         !((*hb)[i].hack[j].oname != NULL &&
1068                           (*hb)[i].hack[j].nname == NULL))
1069                     {
1070                         gmx_fatal(FARGS,
1071                                   "atom %s not found in buiding block %d%s "
1072                                   "while combining tdb and rtp",
1073                                   (*hb)[i].hack[j].oname!=NULL ? 
1074                                   (*hb)[i].hack[j].oname : (*hb)[i].hack[j].AI, 
1075                                   i+1,*resinfo[i].rtp);
1076                     }
1077                 }
1078                 else
1079                 {
1080                     if ( (*hb)[i].hack[j].oname==NULL ) {
1081                         /* we're adding: */
1082                         add_atom_to_restp(&(*restp)[i],resinfo[i].nr,l,
1083                                           &(*hb)[i].hack[j]);
1084                     }
1085                     else
1086                     {
1087                         /* oname != NULL */
1088                         if ( (*hb)[i].hack[j].nname==NULL ) {
1089                             /* we're deleting */
1090                             if (debug) 
1091                                 fprintf(debug, "deleting atom %s from res %d%s in rtp\n",
1092                                         *(*restp)[i].atomname[l], 
1093                                         i+1,(*restp)[i].resname);
1094                             /* shift the rest */
1095                             (*restp)[i].natom--;
1096                             for(k=l; k < (*restp)[i].natom; k++) {
1097                                 (*restp)[i].atom    [k] = (*restp)[i].atom    [k+1];
1098                                 (*restp)[i].atomname[k] = (*restp)[i].atomname[k+1];
1099                                 (*restp)[i].cgnr    [k] = (*restp)[i].cgnr    [k+1];
1100                             }
1101                             /* give back space */
1102                             srenew((*restp)[i].atom,     (*restp)[i].natom);
1103                             srenew((*restp)[i].atomname, (*restp)[i].natom);
1104                             srenew((*restp)[i].cgnr,     (*restp)[i].natom);
1105                         } else { /* nname != NULL */
1106                             /* we're replacing */
1107                             if (debug) 
1108                                 fprintf(debug, "replacing atom %s by %s in res %d%s in rtp\n",
1109                                         *(*restp)[i].atomname[l], (*hb)[i].hack[j].nname, 
1110                                         i+1,(*restp)[i].resname);
1111                             snew( (*restp)[i].atomname[l], 1);
1112                             (*restp)[i].atom[l]      =       *(*hb)[i].hack[j].atom;
1113                             *(*restp)[i].atomname[l] = strdup((*hb)[i].hack[j].nname);
1114                             if ( (*hb)[i].hack[j].cgnr != NOTSET )
1115                                 (*restp)[i].cgnr   [l] = (*hb)[i].hack[j].cgnr;
1116                         }
1117                     }
1118                 }
1119             }
1120         }
1121     }
1122 }
1123
1124 static gmx_bool atomname_cmp_nr(const char *anm,t_hack *hack,int *nr)
1125 {
1126
1127     if (hack->nr == 1)
1128     {
1129         *nr = 0;
1130         
1131         return (gmx_strcasecmp(anm,hack->nname) == 0);
1132     }
1133     else
1134     {
1135         if (isdigit(anm[strlen(anm)-1]))
1136         {
1137             *nr = anm[strlen(anm)-1] - '0';
1138         }
1139         else
1140         {
1141             *nr = 0;
1142         }
1143         if (*nr <= 0 || *nr > hack->nr)
1144         {
1145             return FALSE;
1146         }
1147         else
1148         {
1149             return (strlen(anm) == strlen(hack->nname) + 1 &&
1150                     gmx_strncasecmp(anm,hack->nname,strlen(hack->nname)) == 0);
1151         }
1152     }
1153 }
1154
1155 static gmx_bool match_atomnames_with_rtp_atom(t_atoms *pdba,rvec *x,int atind,
1156                                           t_restp *rptr,t_hackblock *hbr,
1157                                           gmx_bool bVerbose)
1158 {
1159     int  resnr;
1160     int  i,j,k;
1161     char *oldnm,*newnm;
1162     int  anmnr;
1163     char *start_at,buf[STRLEN];
1164     int  start_nr;
1165     gmx_bool bReplaceReplace,bFoundInAdd;
1166     gmx_bool bDeleted;
1167
1168     oldnm = *pdba->atomname[atind];
1169     resnr = pdba->resinfo[pdba->atom[atind].resind].nr;
1170
1171     bDeleted = FALSE;
1172     for(j=0; j<hbr->nhack; j++)
1173     {
1174         if (hbr->hack[j].oname != NULL && hbr->hack[j].nname != NULL &&
1175             gmx_strcasecmp(oldnm,hbr->hack[j].oname) == 0)
1176         {
1177             /* This is a replace entry. */
1178             /* Check if we are not replacing a replaced atom. */
1179             bReplaceReplace = FALSE;
1180             for(k=0; k<hbr->nhack; k++) {
1181                 if (k != j &&
1182                     hbr->hack[k].oname != NULL && hbr->hack[k].nname != NULL &&
1183                     gmx_strcasecmp(hbr->hack[k].nname,hbr->hack[j].oname) == 0)
1184                 {
1185                     /* The replace in hack[j] replaces an atom that
1186                      * was already replaced in hack[k], we do not want
1187                      * second or higher level replaces at this stage.
1188                      */
1189                     bReplaceReplace = TRUE;
1190                 }
1191             }
1192             if (bReplaceReplace)
1193             {
1194                 /* Skip this replace. */
1195                 continue;
1196             }
1197
1198             /* This atom still has the old name, rename it */
1199             newnm = hbr->hack[j].nname;
1200             for(k=0; k<rptr->natom; k++)
1201             {
1202                 if (gmx_strcasecmp(newnm,*rptr->atomname[k]) == 0)
1203                 {
1204                     break;
1205                 }
1206             }
1207             if (k == rptr->natom)
1208             {
1209                 /* The new name is not present in the rtp.
1210                  * We need to apply the replace also to the rtp entry.
1211                  */
1212                 
1213                 /* We need to find the add hack that can add this atom
1214                  * to find out after which atom it should be added.
1215                  */
1216                 bFoundInAdd = FALSE;
1217                 for(k=0; k<hbr->nhack; k++)
1218                 {
1219                     if (hbr->hack[k].oname == NULL &&
1220                         hbr->hack[k].nname != NULL &&
1221                         atomname_cmp_nr(newnm,&hbr->hack[k],&anmnr))
1222                     {
1223                         if (anmnr <= 1)
1224                         {
1225                             start_at = hbr->hack[k].a[0];
1226                         }
1227                         else
1228                         {
1229                             sprintf(buf,"%s%d",hbr->hack[k].nname,anmnr-1);
1230                             start_at = buf;
1231                         }
1232                         for(start_nr=0; start_nr<rptr->natom; start_nr++)
1233                         {
1234                             if (gmx_strcasecmp(start_at,(*rptr->atomname[start_nr])) == 0)
1235                             {
1236                                 break;
1237                             }
1238                         }
1239                         if (start_nr == rptr->natom)
1240                         {
1241                             gmx_fatal(FARGS,"Could not find atom '%s' in residue building block '%s' to add atom '%s' to",
1242                                       start_at,rptr->resname,newnm);
1243                         }
1244                         /* We can add the atom after atom start_nr */
1245                         add_atom_to_restp(rptr,resnr,start_nr,
1246                                           &hbr->hack[j]);
1247                         
1248                         bFoundInAdd = TRUE;
1249                     }
1250                 }
1251
1252                 if (!bFoundInAdd)
1253                 {
1254                     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'",
1255                               newnm,
1256                               hbr->hack[j].oname,hbr->hack[j].nname,
1257                               rptr->resname);
1258                 }
1259             }
1260                 
1261             if (bVerbose)
1262             {
1263                 printf("Renaming atom '%s' in residue '%s' %d to '%s'\n",
1264                        oldnm,rptr->resname,resnr,newnm);
1265             }
1266             /* Rename the atom in pdba */
1267             snew(pdba->atomname[atind],1);
1268             *pdba->atomname[atind] = strdup(newnm);
1269         }
1270         else if (hbr->hack[j].oname != NULL && hbr->hack[j].nname == NULL &&
1271                  gmx_strcasecmp(oldnm,hbr->hack[j].oname) == 0)
1272         {
1273             /* This is a delete entry, check if this atom is present
1274              * in the rtp entry of this residue.
1275              */
1276             for(k=0; k<rptr->natom; k++)
1277             {
1278                 if (gmx_strcasecmp(oldnm,*rptr->atomname[k]) == 0)
1279                 {
1280                     break;
1281                 }
1282             }
1283             if (k == rptr->natom)
1284             {
1285                 /* This atom is not present in the rtp entry,
1286                  * delete is now from the input pdba.
1287                  */
1288                 if (bVerbose)
1289                 {
1290                     printf("Deleting atom '%s' in residue '%s' %d\n",
1291                            oldnm,rptr->resname,resnr);
1292                 }
1293                 /* We should free the atom name,
1294                  * but it might be used multiple times in the symtab.
1295                  * sfree(pdba->atomname[atind]);
1296                  */
1297                 for(k=atind+1; k<pdba->nr; k++)
1298                 {
1299                     pdba->atom[k-1]     = pdba->atom[k];
1300                     pdba->atomname[k-1] = pdba->atomname[k];
1301                     copy_rvec(x[k],x[k-1]);
1302                 }
1303                 pdba->nr--;
1304                 bDeleted = TRUE;
1305             }
1306         }
1307     }
1308
1309     return bDeleted;
1310 }
1311     
1312 void match_atomnames_with_rtp(t_restp restp[],t_hackblock hb[],
1313                               t_atoms *pdba,rvec *x,
1314                               gmx_bool bVerbose)
1315 {
1316     int  i,j,k;
1317     char *oldnm,*newnm;
1318     int  resnr;
1319     t_restp *rptr;
1320     t_hackblock *hbr;
1321     int  anmnr;
1322     char *start_at,buf[STRLEN];
1323     int  start_nr;
1324     gmx_bool bFoundInAdd;
1325     
1326     for(i=0; i<pdba->nr; i++)
1327     {
1328         oldnm = *pdba->atomname[i];
1329         resnr = pdba->resinfo[pdba->atom[i].resind].nr;
1330         rptr  = &restp[pdba->atom[i].resind];
1331         for(j=0; (j<rptr->natom); j++)
1332         {
1333             if (gmx_strcasecmp(oldnm,*(rptr->atomname[j])) == 0)
1334             {
1335                 break;
1336             }
1337         }
1338         if (j == rptr->natom)
1339         {
1340             /* Not found yet, check if we have to rename this atom */
1341             if (match_atomnames_with_rtp_atom(pdba,x,i,
1342                                               rptr,&(hb[pdba->atom[i].resind]),
1343                                               bVerbose))
1344             {
1345                 /* We deleted this atom, decrease the atom counter by 1. */
1346                 i--;
1347             }
1348         }
1349     }
1350 }
1351
1352 void gen_cmap(t_params *psb, t_restp *restp, int natoms, t_atom atom[], char **aname[], int nres)
1353 {
1354         int     residx,i,ii,j,k;
1355         atom_id ai,aj,ak,al,am;
1356         const char *ptr;
1357         
1358         if (debug)
1359                 ptr = "cmap";
1360         else
1361                 ptr = "check";
1362         
1363         fprintf(stderr,"Making cmap torsions...");
1364         i=0;
1365         /* End loop at nres-1, since the very last residue does not have a +N atom, and
1366          * therefore we get a valgrind invalid 4 byte read error with atom am */
1367         for(residx=0; residx<nres-1; residx++)
1368         {
1369                 /* Add CMAP terms from the list of CMAP interactions */
1370                 for(j=0;j<restp[residx].rb[ebtsCMAP].nb; j++)
1371                 {
1372                         ai=search_atom(restp[residx].rb[ebtsCMAP].b[j].a[0],i,natoms,atom,aname,
1373                                                    ptr,TRUE);
1374                         aj=search_atom(restp[residx].rb[ebtsCMAP].b[j].a[1],i,natoms,atom,aname,
1375                                                    ptr,TRUE);
1376                         ak=search_atom(restp[residx].rb[ebtsCMAP].b[j].a[2],i,natoms,atom,aname,
1377                                                    ptr,TRUE);
1378                         al=search_atom(restp[residx].rb[ebtsCMAP].b[j].a[3],i,natoms,atom,aname,
1379                                                    ptr,TRUE);
1380                         am=search_atom(restp[residx].rb[ebtsCMAP].b[j].a[4],i,natoms,atom,aname,
1381                                                    ptr,TRUE);
1382                         
1383                         /* The first and last residues no not have cmap torsions */
1384                         if(ai!=NO_ATID && aj!=NO_ATID && ak!=NO_ATID && al!=NO_ATID && am!=NO_ATID)
1385                         {
1386                                 add_cmap_param(psb,ai,aj,ak,al,am,restp[residx].rb[ebtsCMAP].b[j].s);
1387                         }
1388                 }
1389                 
1390                 if(residx<nres-1)
1391                 {
1392                         while(atom[i].resind<residx+1)
1393                         {
1394                                 i++;
1395                         }
1396                 }
1397         }
1398         
1399         /* Start the next residue */
1400 }
1401
1402 static void 
1403 scrub_charge_groups(int *cgnr, int natoms)
1404 {
1405         int i;
1406         
1407         for(i=0;i<natoms;i++)
1408         {
1409                 cgnr[i]=i+1;
1410         }
1411 }
1412
1413
1414 void pdb2top(FILE *top_file, char *posre_fn, char *molname,
1415              t_atoms *atoms, rvec **x, gpp_atomtype_t atype, t_symtab *tab,
1416              int nrtp, t_restp rtp[],
1417              t_restp *restp, t_hackblock *hb,
1418              int nterpairs,t_hackblock **ntdb, t_hackblock **ctdb,
1419              gmx_bool bAllowMissing,
1420              gmx_bool bVsites, gmx_bool bVsiteAromatics,
1421              const char *ff, const char *ffdir,
1422              real mHmult,
1423              int nssbonds, t_ssbond *ssbonds,
1424              real long_bond_dist, real short_bond_dist,
1425              gmx_bool bDeuterate, gmx_bool bChargeGroups, gmx_bool bCmap,
1426              gmx_bool bRenumRes,gmx_bool bRTPresname)
1427 {
1428     /*
1429   t_hackblock *hb;
1430   t_restp  *restp;
1431     */
1432   t_params plist[F_NRE];
1433   t_excls  *excls;
1434   t_nextnb nnb;
1435   int      *cgnr;
1436   int      *vsite_type;
1437   int      i,nmissat;
1438   int      bts[ebtsNR];
1439   
1440   init_plist(plist);
1441
1442   if (debug) {
1443     print_resall(debug, atoms->nres, restp, atype);
1444     dump_hb(debug, atoms->nres, hb);
1445   }
1446   
1447   /* Make bonds */
1448   at2bonds(&(plist[F_BONDS]), hb, 
1449            atoms->nr, atoms->atom, atoms->atomname, atoms->nres, *x, 
1450            long_bond_dist, short_bond_dist, bAllowMissing);
1451   
1452   /* specbonds: disulphide bonds & heme-his */
1453   do_ssbonds(&(plist[F_BONDS]),
1454              atoms->nr, atoms->atom, atoms->atomname, nssbonds, ssbonds,
1455              bAllowMissing);
1456   
1457   nmissat = name2type(atoms, &cgnr, atype, restp);
1458   if (nmissat) {
1459     if (bAllowMissing)
1460       fprintf(stderr,"There were %d missing atoms in molecule %s\n",
1461               nmissat,molname);
1462     else
1463       gmx_fatal(FARGS,"There were %d missing atoms in molecule %s, if you want to use this incomplete topology anyhow, use the option -missing",
1464                   nmissat,molname);
1465   }
1466   
1467   /* Cleanup bonds (sort and rm doubles) */ 
1468   clean_bonds(&(plist[F_BONDS]));
1469   
1470   snew(vsite_type,atoms->nr);
1471   for(i=0; i<atoms->nr; i++)
1472     vsite_type[i]=NOTSET;
1473   if (bVsites) {
1474     /* determine which atoms will be vsites and add dummy masses 
1475        also renumber atom numbers in plist[0..F_NRE]! */
1476     do_vsites(nrtp, rtp, atype, atoms, tab, x, plist, 
1477               &vsite_type, &cgnr, mHmult, bVsiteAromatics, ffdir);
1478   }
1479   
1480   /* Make Angles and Dihedrals */
1481   fprintf(stderr,"Generating angles, dihedrals and pairs...\n");
1482   snew(excls,atoms->nr);
1483   init_nnb(&nnb,atoms->nr,4);
1484   gen_nnb(&nnb,plist);
1485   print_nnb(&nnb,"NNB");
1486   gen_pad(&nnb,atoms,restp[0].nrexcl,restp[0].HH14,
1487           plist,excls,hb,restp[0].bAlldih,restp[0].bRemoveDih,
1488           bAllowMissing);
1489   done_nnb(&nnb);
1490   
1491     /* Make CMAP */
1492     if (TRUE == bCmap)
1493     {
1494                 gen_cmap(&(plist[F_CMAP]), restp, atoms->nr, atoms->atom, atoms->atomname, atoms->nres);
1495         if (plist[F_CMAP].nr > 0)
1496         {
1497             fprintf(stderr, "There are %4d cmap torsion pairs\n",
1498                     plist[F_CMAP].nr);
1499         }
1500     }
1501
1502   /* set mass of all remaining hydrogen atoms */
1503   if (mHmult != 1.0)
1504     do_h_mass(&(plist[F_BONDS]),vsite_type,atoms,mHmult,bDeuterate);
1505   sfree(vsite_type);
1506   
1507   /* Cleanup bonds (sort and rm doubles) */ 
1508   /* clean_bonds(&(plist[F_BONDS]));*/
1509    
1510   fprintf(stderr,
1511           "There are %4d dihedrals, %4d impropers, %4d angles\n"
1512           "          %4d pairs,     %4d bonds and  %4d virtual sites\n",
1513           plist[F_PDIHS].nr, plist[F_IDIHS].nr, plist[F_ANGLES].nr,
1514           plist[F_LJ14].nr, plist[F_BONDS].nr,
1515           plist[F_VSITE2].nr +
1516           plist[F_VSITE3].nr +
1517           plist[F_VSITE3FD].nr +
1518           plist[F_VSITE3FAD].nr +
1519           plist[F_VSITE3OUT].nr +
1520       plist[F_VSITE4FD].nr +
1521       plist[F_VSITE4FDN].nr );
1522   
1523   print_sums(atoms, FALSE);
1524   
1525   if (FALSE == bChargeGroups)
1526   {
1527           scrub_charge_groups(cgnr, atoms->nr);
1528   }
1529
1530     if (bRenumRes)
1531     {
1532         for(i=0; i<atoms->nres; i++) 
1533         {
1534             atoms->resinfo[i].nr = i + 1;
1535             atoms->resinfo[i].ic = ' ';
1536         }
1537     }
1538         
1539   if (top_file) {
1540     fprintf(stderr,"Writing topology\n");
1541     /* We can copy the bonded types from the first restp,
1542      * since the types have to be identical for all residues in one molecule.
1543      */
1544     for(i=0; i<ebtsNR; i++) {
1545         bts[i] = restp[0].rb[i].type;
1546     }
1547     write_top(top_file, posre_fn, molname,
1548               atoms, bRTPresname, 
1549               bts, plist, excls, atype, cgnr, restp[0].nrexcl);
1550   }
1551   
1552   /* cleaning up */
1553   free_t_hackblock(atoms->nres, &hb);
1554   free_t_restp(atoms->nres, &restp);
1555     
1556   /* we should clean up hb and restp here, but that is a *L*O*T* of work! */
1557   sfree(cgnr);
1558   for (i=0; i<F_NRE; i++)
1559     sfree(plist[i].param);
1560   for (i=0; i<atoms->nr; i++)
1561     sfree(excls[i].e);
1562   sfree(excls);
1563 }