Merge branch 'master' of git://git.gromacs.org/gromacs
[gromacs/adressmacs.git] / src / kernel / pdb2top.c
blob1384f0a8555c70146d246b16b44a52497951f270
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <stdio.h>
41 #include <math.h>
42 #include <ctype.h>
44 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
45 #include <direct.h>
46 #include <io.h>
47 #endif
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"
78 /* this must correspond to enum in pdb2top.h */
79 const char *hh[ehisNR] = { "HISD", "HISE", "HISH", "HIS1" };
81 static int missing_atoms(t_restp *rp, int resind,t_atoms *at, int i0, int i)
83 int j,k,nmiss;
84 char *name;
85 gmx_bool bFound, bRet;
87 nmiss = 0;
88 for (j=0; j<rp->natom; j++)
90 name=*(rp->atomname[j]);
91 bFound=FALSE;
92 for (k=i0; k<i; k++)
94 bFound = (bFound || !gmx_strcasecmp(*(at->atomname[k]),name));
96 if (!bFound)
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')
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);
108 fprintf(stderr,"\n");
112 return nmiss;
115 gmx_bool is_int(double x)
117 const double tol = 1e-4;
118 int ix;
120 if (x < 0)
121 x=-x;
122 ix=gmx_nint(x);
124 return (fabs(x-ix) < tol);
127 static void swap_strings(char **s,int i,int j)
129 char *tmp;
131 tmp = s[i];
132 s[i] = s[j];
133 s[j] = tmp;
136 void
137 choose_ff(const char *ffsel,
138 char *forcefield, int ff_maxlen,
139 char *ffdir, int ffdir_maxlen)
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;
148 nff = fflib_search_file_in_dirend(fflib_forcefield_itp(),
149 fflib_forcefield_dir_ext(),
150 &ffdirs);
152 if (nff == 0)
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());
158 /* Replace with unix path separators */
159 if(DIR_SEPARATOR!='/')
161 for(i=0;i<nff;i++)
163 while( (ptr=strchr(ffdirs[i],DIR_SEPARATOR))!=NULL )
165 *ptr='/';
170 /* Store the force field names in ffs */
171 snew(ffs,nff);
172 snew(ffs_dir,nff);
173 for(i=0; i<nff; i++)
175 /* Remove the path from the ffdir name - use our unix standard here! */
176 ptr = strrchr(ffdirs[i],'/');
177 if (ptr == NULL)
179 ffs[i] = strdup(ffdirs[i]);
180 ffs_dir[i] = low_gmxlibfn(ffdirs[i],FALSE,FALSE);
181 if (ffs_dir[i] == NULL)
183 gmx_fatal(FARGS,"Can no longer find file '%s'",ffdirs[i]);
186 else
188 ffs[i] = strdup(ptr+1);
189 ffs_dir[i] = strdup(ffdirs[i]);
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';
196 if (ffsel != NULL)
198 sel = -1;
199 cwdsel = -1;
200 nfound = 0;
201 for(i=0; i<nff; i++)
203 if ( strcmp(ffs[i],ffsel)==0 )
205 /* Matching ff name */
206 sel = i;
207 nfound++;
209 if( strncmp(ffs_dir[i],".",1)==0 )
211 cwdsel = i;
216 if(cwdsel != -1)
218 sel = cwdsel;
221 if(nfound>1)
223 if(cwdsel!=-1)
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);
230 else
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);
237 else if (nfound==0)
239 gmx_fatal(FARGS,"Could not find force field '%s' in current directory, install tree or GMXDATA path.",ffsel);
242 else if (nff > 1)
244 snew(desc,nff);
245 for(i=0; (i<nff); i++)
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))
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);
259 else
261 desc[i] = strdup(ffs[i]);
264 /* Order force fields from the same dir alphabetically
265 * and put deprecated force fields at the end.
267 for(i=0; (i<nff); i++)
269 for(j=i+1; (j<nff); j++)
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)))
276 swap_strings(ffdirs,i,j);
277 swap_strings(ffs ,i,j);
278 swap_strings(desc ,i,j);
283 printf("\nSelect the Force Field:\n");
284 for(i=0; (i<nff); i++)
286 if (i == 0 || strcmp(ffs_dir[i-1],ffs_dir[i]) != 0)
288 if( strcmp(ffs_dir[i],".")==0 )
290 printf("From current directory:\n");
292 else
294 printf("From '%s':\n",ffs_dir[i]);
297 printf("%2d: %s\n",i+1,desc[i]);
298 sfree(desc[i]);
300 sfree(desc);
304 pret = fgets(buf,STRLEN,stdin);
306 if (pret != NULL)
308 sscanf(buf,"%d",&sel);
309 sel--;
312 while ( pret==NULL || (sel < 0) || (sel >= nff));
314 else
316 sel = 0;
319 if (strlen(ffs[sel]) >= (size_t)ff_maxlen)
321 gmx_fatal(FARGS,"Length of force field name (%d) >= maxlen (%d)",
322 strlen(ffs[sel]),ff_maxlen);
324 strcpy(forcefield,ffs[sel]);
326 if (strlen(ffdirs[sel]) >= (size_t)ffdir_maxlen)
328 gmx_fatal(FARGS,"Length of force field dir (%d) >= maxlen (%d)",
329 strlen(ffdirs[sel]),ffdir_maxlen);
331 strcpy(ffdir,ffdirs[sel]);
333 for(i=0; (i<nff); i++)
335 sfree(ffdirs[i]);
336 sfree(ffs[i]);
337 sfree(ffs_dir[i]);
339 sfree(ffdirs);
340 sfree(ffs);
341 sfree(ffs_dir);
344 void choose_watermodel(const char *wmsel,const char *ffdir,
345 char **watermodel)
347 const char *fn_watermodels="watermodels.dat";
348 char fn_list[STRLEN];
349 FILE *fp;
350 char buf[STRLEN];
351 int nwm,sel,i;
352 char **model;
353 char *pret;
355 if (strcmp(wmsel,"none") == 0)
357 *watermodel = NULL;
359 return;
361 else if (strcmp(wmsel,"select") != 0)
363 *watermodel = strdup(wmsel);
365 return;
368 sprintf(fn_list,"%s%c%s",ffdir,DIR_SEPARATOR,fn_watermodels);
370 if (!fflib_fexist(fn_list))
372 fprintf(stderr,"No file '%s' found, will not include a water model\n",
373 fn_watermodels);
374 *watermodel = NULL;
376 return;
379 fp = fflib_open(fn_list);
380 printf("\nSelect the Water Model:\n");
381 nwm = 0;
382 model = NULL;
383 while (get_a_line(fp,buf,STRLEN))
385 srenew(model,nwm+1);
386 snew(model[nwm],STRLEN);
387 sscanf(buf,"%s%n",model[nwm],&i);
388 if (i > 0)
390 ltrim(buf+i);
391 fprintf(stderr,"%2d: %s\n",nwm+1,buf+i);
392 nwm++;
394 else
396 sfree(model[nwm]);
399 fclose(fp);
400 fprintf(stderr,"%2d: %s\n",nwm+1,"None");
404 pret = fgets(buf,STRLEN,stdin);
406 if (pret != NULL)
408 sscanf(buf,"%d",&sel);
409 sel--;
412 while (pret == NULL || sel < 0 || sel > nwm);
414 if (sel == nwm)
416 *watermodel = NULL;
418 else
420 *watermodel = strdup(model[sel]);
423 for(i=0; i<nwm; i++)
425 sfree(model[i]);
427 sfree(model);
430 static int name2type(t_atoms *at, int **cgnr, gpp_atomtype_t atype,
431 t_restp restp[])
433 int i,j,prevresind,resind,i0,prevcg,cg,curcg;
434 char *name;
435 gmx_bool bProt, bNterm;
436 double qt;
437 int nmissat;
438 gmx_residuetype_t rt;
440 nmissat = 0;
442 resind=-1;
443 bProt=FALSE;
444 bNterm=FALSE;
445 i0=0;
446 snew(*cgnr,at->nr);
447 qt=0;
448 prevcg=NOTSET;
449 curcg=0;
450 cg=-1;
451 j=NOTSET;
452 gmx_residuetype_init(&rt);
454 for(i=0; (i<at->nr); i++) {
455 prevresind=resind;
456 if (at->atom[i].resind != resind) {
457 resind = at->atom[i].resind;
458 bProt = gmx_residuetype_is_protein(rt,*(at->resinfo[resind].name));
459 bNterm=bProt && (resind == 0);
460 if (resind > 0) {
461 nmissat += missing_atoms(&restp[prevresind],prevresind,at,i0,i);
463 i0=i;
465 if (at->atom[i].m == 0) {
466 if (debug)
467 fprintf(debug,"atom %d%s: curcg=%d, prevcg=%d, cg=%d\n",
468 i+1,*(at->atomname[i]),curcg,prevcg,
469 j==NOTSET ? NOTSET : restp[resind].cgnr[j]);
470 qt=0;
471 prevcg=cg;
472 name=*(at->atomname[i]);
473 j=search_jtype(&restp[resind],name,bNterm);
474 at->atom[i].type = restp[resind].atom[j].type;
475 at->atom[i].q = restp[resind].atom[j].q;
476 at->atom[i].m = get_atomtype_massA(restp[resind].atom[j].type,
477 atype);
478 cg = restp[resind].cgnr[j];
479 /* A charge group number -1 signals a separate charge group
480 * for this atom.
482 if ( (cg == -1) || (cg != prevcg) || (resind != prevresind) ) {
483 curcg++;
485 } else {
486 if (debug)
487 fprintf(debug,"atom %d%s: curcg=%d, qt=%g, is_int=%d\n",
488 i+1,*(at->atomname[i]),curcg,qt,is_int(qt));
489 cg=-1;
490 if (is_int(qt)) {
491 qt=0;
492 curcg++;
494 qt+=at->atom[i].q;
496 (*cgnr)[i]=curcg;
497 at->atom[i].typeB = at->atom[i].type;
498 at->atom[i].qB = at->atom[i].q;
499 at->atom[i].mB = at->atom[i].m;
501 nmissat += missing_atoms(&restp[resind],resind,at,i0,i);
503 gmx_residuetype_destroy(rt);
505 return nmissat;
508 static void print_top_heavy_H(FILE *out, real mHmult)
510 if (mHmult == 2.0)
511 fprintf(out,"; Using deuterium instead of hydrogen\n\n");
512 else if (mHmult == 4.0)
513 fprintf(out,"#define HEAVY_H\n\n");
514 else if (mHmult != 1.0)
515 fprintf(stderr,"WARNING: unsupported proton mass multiplier (%g) "
516 "in pdb2top\n",mHmult);
519 void print_top_comment(FILE *out,
520 const char *filename,
521 const char *generator,
522 const char *ffdir,
523 gmx_bool bITP)
525 char tmp[256];
526 char ffdir_parent[STRLEN];
527 char *p;
529 nice_header(out,filename);
530 fprintf(out,";\tThis is a %s topology file\n;\n",bITP ? "include" : "standalone");
531 fprintf(out,";\tIt was generated using program:\n;\t%s\n;\n",
532 (NULL == generator) ? "unknown" : generator);
533 fprintf(out,";\tCommand line was:\n;\t%s\n;\n",command_line());
535 if(strchr(ffdir,'/')==NULL)
537 fprintf(out,";\tForce field was read from the standard Gromacs share directory.\n;\n\n");
539 else if(ffdir[0]=='.')
541 fprintf(out,";\tForce field was read from current directory or a relative path - path added.\n;\n\n");
543 else
545 strncpy(ffdir_parent,ffdir,STRLEN-1);
546 p=strrchr(ffdir_parent,'/');
548 *p='\0';
550 fprintf(out,
551 ";\tForce field data was read from:\n"
552 ";\t%s\n"
553 ";\n"
554 ";\tNote:\n"
555 ";\tThis might be a non-standard force field location. When you use this topology, the\n"
556 ";\tforce field must either be present in the current directory, or the location\n"
557 ";\tspecified in the GMXLIB path variable or with the 'include' mdp file option.\n;\n\n",
558 ffdir_parent);
562 void print_top_header(FILE *out,const char *filename,
563 const char *title,gmx_bool bITP,const char *ffdir,real mHmult)
565 const char *p;
567 print_top_comment(out,filename,title,ffdir,bITP);
569 print_top_heavy_H(out, mHmult);
570 fprintf(out,"; Include forcefield parameters\n");
572 p=strrchr(ffdir,'/');
573 p = (ffdir[0]=='.' || p==NULL) ? ffdir : p+1;
575 fprintf(out,"#include \"%s/%s\"\n\n",p,fflib_forcefield_itp());
578 static void print_top_posre(FILE *out,const char *pr)
580 fprintf(out,"; Include Position restraint file\n");
581 fprintf(out,"#ifdef POSRES\n");
582 fprintf(out,"#include \"%s\"\n",pr);
583 fprintf(out,"#endif\n\n");
586 static void print_top_water(FILE *out,const char *ffdir,const char *water)
588 const char *p;
589 char buf[STRLEN];
591 fprintf(out,"; Include water topology\n");
593 p=strrchr(ffdir,'/');
594 p = (ffdir[0]=='.' || p==NULL) ? ffdir : p+1;
595 fprintf(out,"#include \"%s/%s.itp\"\n",p,water);
597 fprintf(out,"\n");
598 fprintf(out,"#ifdef POSRES_WATER\n");
599 fprintf(out,"; Position restraint for each water oxygen\n");
600 fprintf(out,"[ position_restraints ]\n");
601 fprintf(out,";%3s %5s %9s %10s %10s\n","i","funct","fcx","fcy","fcz");
602 fprintf(out,"%4d %4d %10g %10g %10g\n",1,1,1000.0,1000.0,1000.0);
603 fprintf(out,"#endif\n");
604 fprintf(out,"\n");
606 sprintf(buf,"%s/ions.itp",p);
608 if (fflib_fexist(buf))
610 fprintf(out,"; Include topology for ions\n");
611 fprintf(out,"#include \"%s\"\n",buf);
612 fprintf(out,"\n");
616 static void print_top_system(FILE *out, const char *title)
618 fprintf(out,"[ %s ]\n",dir2str(d_system));
619 fprintf(out,"; Name\n");
620 fprintf(out,"%s\n\n",title[0]?title:"Protein");
623 void print_top_mols(FILE *out,
624 const char *title, const char *ffdir, const char *water,
625 int nincl, char **incls, int nmol, t_mols *mols)
627 int i;
628 char *incl;
630 if (nincl>0) {
631 fprintf(out,"; Include chain topologies\n");
632 for (i=0; (i<nincl); i++) {
633 incl = strrchr(incls[i],DIR_SEPARATOR);
634 if (incl == NULL) {
635 incl = incls[i];
636 } else {
637 /* Remove the path from the include name */
638 incl = incl + 1;
640 fprintf(out,"#include \"%s\"\n",incl);
642 fprintf(out,"\n");
645 if (water)
647 print_top_water(out,ffdir,water);
649 print_top_system(out, title);
651 if (nmol) {
652 fprintf(out,"[ %s ]\n",dir2str(d_molecules));
653 fprintf(out,"; %-15s %5s\n","Compound","#mols");
654 for (i=0; (i<nmol); i++)
655 fprintf(out,"%-15s %5d\n",mols[i].name,mols[i].nr);
659 void write_top(FILE *out, char *pr,char *molname,
660 t_atoms *at,gmx_bool bRTPresname,
661 int bts[],t_params plist[],t_excls excls[],
662 gpp_atomtype_t atype,int *cgnr, int nrexcl)
663 /* NOTE: nrexcl is not the size of *excl! */
665 if (at && atype && cgnr) {
666 fprintf(out,"[ %s ]\n",dir2str(d_moleculetype));
667 fprintf(out,"; %-15s %5s\n","Name","nrexcl");
668 fprintf(out,"%-15s %5d\n\n",molname?molname:"Protein",nrexcl);
670 print_atoms(out, atype, at, cgnr, bRTPresname);
671 print_bondeds(out,at->nr,d_bonds, F_BONDS, bts[ebtsBONDS], plist);
672 print_bondeds(out,at->nr,d_constraints,F_CONSTR, 0, plist);
673 print_bondeds(out,at->nr,d_constraints,F_CONSTRNC, 0, plist);
674 print_bondeds(out,at->nr,d_pairs, F_LJ14, 0, plist);
675 print_excl(out,at->nr,excls);
676 print_bondeds(out,at->nr,d_angles, F_ANGLES, bts[ebtsANGLES],plist);
677 print_bondeds(out,at->nr,d_dihedrals, F_PDIHS, bts[ebtsPDIHS], plist);
678 print_bondeds(out,at->nr,d_dihedrals, F_IDIHS, bts[ebtsIDIHS], plist);
679 print_bondeds(out,at->nr,d_cmap, F_CMAP, bts[ebtsCMAP], plist);
680 print_bondeds(out,at->nr,d_polarization,F_POLARIZATION, 0, plist);
681 print_bondeds(out,at->nr,d_thole_polarization,F_THOLE_POL,0, plist);
682 print_bondeds(out,at->nr,d_vsites2, F_VSITE2, 0, plist);
683 print_bondeds(out,at->nr,d_vsites3, F_VSITE3, 0, plist);
684 print_bondeds(out,at->nr,d_vsites3, F_VSITE3FD, 0, plist);
685 print_bondeds(out,at->nr,d_vsites3, F_VSITE3FAD,0, plist);
686 print_bondeds(out,at->nr,d_vsites3, F_VSITE3OUT,0, plist);
687 print_bondeds(out,at->nr,d_vsites4, F_VSITE4FD, 0, plist);
688 print_bondeds(out,at->nr,d_vsites4, F_VSITE4FDN, 0, plist);
690 if (pr)
691 print_top_posre(out,pr);
695 static atom_id search_res_atom(const char *type,int resind,
696 int natom,t_atom at[],
697 char ** const *aname,
698 const char *bondtype,gmx_bool bAllowMissing)
700 int i;
702 for(i=0; (i<natom); i++)
703 if (at[i].resind == resind)
704 return search_atom(type,i,natom,at,aname,bondtype,bAllowMissing);
706 return NO_ATID;
709 static void do_ssbonds(t_params *ps,int natoms,t_atom atom[],char **aname[],
710 int nssbonds,t_ssbond *ssbonds,gmx_bool bAllowMissing)
712 int i,ri,rj;
713 atom_id ai,aj;
715 for(i=0; (i<nssbonds); i++) {
716 ri = ssbonds[i].res1;
717 rj = ssbonds[i].res2;
718 ai = search_res_atom(ssbonds[i].a1,ri,natoms,atom,aname,
719 "special bond",bAllowMissing);
720 aj = search_res_atom(ssbonds[i].a2,rj,natoms,atom,aname,
721 "special bond",bAllowMissing);
722 if ((ai == NO_ATID) || (aj == NO_ATID))
723 gmx_fatal(FARGS,"Trying to make impossible special bond (%s-%s)!",
724 ssbonds[i].a1,ssbonds[i].a2);
725 add_param(ps,ai,aj,NULL,NULL);
729 static gmx_bool inter_res_bond(const t_rbonded *b)
731 return (b->AI[0] == '-' || b->AI[0] == '+' ||
732 b->AJ[0] == '-' || b->AJ[0] == '+');
735 static void at2bonds(t_params *psb, t_hackblock *hb,
736 int natoms, t_atom atom[], char **aname[],
737 int nres, rvec x[],
738 real long_bond_dist, real short_bond_dist,
739 gmx_bool bAllowMissing)
741 int resind,i,j,k;
742 atom_id ai,aj;
743 real dist2, long_bond_dist2, short_bond_dist2;
744 const char *ptr;
746 long_bond_dist2 = sqr(long_bond_dist);
747 short_bond_dist2 = sqr(short_bond_dist);
749 if (debug)
750 ptr = "bond";
751 else
752 ptr = "check";
754 fprintf(stderr,"Making bonds...\n");
755 i=0;
756 for(resind=0; (resind < nres) && (i<natoms); resind++) {
757 /* add bonds from list of bonded interactions */
758 for(j=0; j < hb[resind].rb[ebtsBONDS].nb; j++) {
759 /* Unfortunately we can not issue errors or warnings
760 * for missing atoms in bonds, as the hydrogens and terminal atoms
761 * have not been added yet.
763 ai=search_atom(hb[resind].rb[ebtsBONDS].b[j].AI,i,natoms,atom,aname,
764 ptr,TRUE);
765 aj=search_atom(hb[resind].rb[ebtsBONDS].b[j].AJ,i,natoms,atom,aname,
766 ptr,TRUE);
767 if (ai != NO_ATID && aj != NO_ATID) {
768 dist2 = distance2(x[ai],x[aj]);
769 if (dist2 > long_bond_dist2 )
771 fprintf(stderr,"Warning: Long Bond (%d-%d = %g nm)\n",
772 ai+1,aj+1,sqrt(dist2));
774 else if (dist2 < short_bond_dist2 )
776 fprintf(stderr,"Warning: Short Bond (%d-%d = %g nm)\n",
777 ai+1,aj+1,sqrt(dist2));
779 add_param(psb,ai,aj,NULL,hb[resind].rb[ebtsBONDS].b[j].s);
782 /* add bonds from list of hacks (each added atom gets a bond) */
783 while( (i<natoms) && (atom[i].resind == resind) ) {
784 for(j=0; j < hb[resind].nhack; j++)
785 if ( ( hb[resind].hack[j].tp > 0 ||
786 hb[resind].hack[j].oname==NULL ) &&
787 strcmp(hb[resind].hack[j].AI,*(aname[i])) == 0 ) {
788 switch(hb[resind].hack[j].tp) {
789 case 9: /* COOH terminus */
790 add_param(psb,i,i+1,NULL,NULL); /* C-O */
791 add_param(psb,i,i+2,NULL,NULL); /* C-OA */
792 add_param(psb,i+2,i+3,NULL,NULL); /* OA-H */
793 break;
794 default:
795 for(k=0; (k<hb[resind].hack[j].nr); k++)
796 add_param(psb,i,i+k+1,NULL,NULL);
799 i++;
801 /* we're now at the start of the next residue */
805 static int pcompar(const void *a, const void *b)
807 t_param *pa,*pb;
808 int d;
809 pa=(t_param *)a;
810 pb=(t_param *)b;
812 d = pa->AI - pb->AI;
813 if (d == 0)
814 d = pa->AJ - pb->AJ;
815 if (d == 0)
816 return strlen(pb->s) - strlen(pa->s);
817 else
818 return d;
821 static void clean_bonds(t_params *ps)
823 int i,j;
824 atom_id a;
826 if (ps->nr > 0) {
827 /* swap atomnumbers in bond if first larger than second: */
828 for(i=0; (i<ps->nr); i++)
829 if ( ps->param[i].AJ < ps->param[i].AI ) {
830 a = ps->param[i].AI;
831 ps->param[i].AI = ps->param[i].AJ;
832 ps->param[i].AJ = a;
835 /* Sort bonds */
836 qsort(ps->param,ps->nr,(size_t)sizeof(ps->param[0]),pcompar);
838 /* remove doubles, keep the first one always. */
839 j = 1;
840 for(i=1; (i<ps->nr); i++) {
841 if ((ps->param[i].AI != ps->param[j-1].AI) ||
842 (ps->param[i].AJ != ps->param[j-1].AJ) ) {
843 if (j != i) {
844 cp_param(&(ps->param[j]),&(ps->param[i]));
846 j++;
849 fprintf(stderr,"Number of bonds was %d, now %d\n",ps->nr,j);
850 ps->nr=j;
852 else
853 fprintf(stderr,"No bonds\n");
856 void print_sums(t_atoms *atoms, gmx_bool bSystem)
858 double m,qtot;
859 int i;
860 const char *where;
862 if (bSystem)
863 where=" in system";
864 else
865 where="";
867 m=0;
868 qtot=0;
869 for(i=0; (i<atoms->nr); i++) {
870 m+=atoms->atom[i].m;
871 qtot+=atoms->atom[i].q;
874 fprintf(stderr,"Total mass%s %.3f a.m.u.\n",where,m);
875 fprintf(stderr,"Total charge%s %.3f e\n",where,qtot);
878 static void check_restp_type(const char *name,int t1,int t2)
880 if (t1 != t2)
882 gmx_fatal(FARGS,"Residues in one molecule have a different '%s' type: %d and %d",name,t1,t2);
886 static void check_restp_types(t_restp *r0,t_restp *r1)
888 int i;
890 check_restp_type("all dihedrals",r0->bAlldih,r1->bAlldih);
891 check_restp_type("nrexcl",r0->nrexcl,r1->nrexcl);
892 check_restp_type("HH14",r0->HH14,r1->HH14);
893 check_restp_type("remove dihedrals",r0->bRemoveDih,r1->bRemoveDih);
895 for(i=0; i<ebtsNR; i++)
897 check_restp_type(btsNames[i],r0->rb[i].type,r1->rb[i].type);
901 void add_atom_to_restp(t_restp *restp,int resnr,int at_start,const t_hack *hack)
903 char buf[STRLEN];
904 int k;
905 const char *Hnum="123456";
907 /*if (debug)
909 fprintf(debug,"adding atom(s) %s to atom %s in res %d%s in rtp\n",
910 hack->nname,
911 *restp->atomname[at_start], resnr, restp->resname);
913 strcpy(buf, hack->nname);
914 buf[strlen(buf)+1]='\0';
915 if ( hack->nr > 1 )
917 buf[strlen(buf)]='-';
919 /* make space */
920 restp->natom += hack->nr;
921 srenew(restp->atom, restp->natom);
922 srenew(restp->atomname, restp->natom);
923 srenew(restp->cgnr, restp->natom);
924 /* shift the rest */
925 for(k=restp->natom-1; k > at_start+hack->nr; k--)
927 restp->atom[k] =
928 restp->atom [k - hack->nr];
929 restp->atomname[k] =
930 restp->atomname[k - hack->nr];
931 restp->cgnr[k] =
932 restp->cgnr [k - hack->nr];
934 /* now add them */
935 for(k=0; k < hack->nr; k++)
937 /* set counter in atomname */
938 if ( hack->nr > 1 )
940 buf[strlen(buf)-1] = Hnum[k];
942 snew( restp->atomname[at_start+1+k], 1);
943 restp->atom [at_start+1+k] = *hack->atom;
944 *restp->atomname[at_start+1+k] = strdup(buf);
945 if ( hack->cgnr != NOTSET )
947 restp->cgnr [at_start+1+k] = hack->cgnr;
949 else
951 restp->cgnr [at_start+1+k] = restp->cgnr[at_start];
956 void get_hackblocks_rtp(t_hackblock **hb, t_restp **restp,
957 int nrtp, t_restp rtp[],
958 int nres, t_resinfo *resinfo,
959 int nterpairs,
960 t_hackblock **ntdb, t_hackblock **ctdb,
961 int *rn, int *rc)
963 int i, j, k, l;
964 char *key;
965 t_restp *res;
966 char buf[STRLEN];
967 const char *Hnum="123456";
968 int tern,terc;
969 gmx_bool bN,bC,bRM;
971 snew(*hb,nres);
972 snew(*restp,nres);
973 /* first the termini */
974 for(i=0; i<nterpairs; i++) {
975 if (rn[i] >= 0 && ntdb[i] != NULL) {
976 copy_t_hackblock(ntdb[i], &(*hb)[rn[i]]);
978 if (rc[i] >= 0 && ctdb[i] != NULL) {
979 merge_t_hackblock(ctdb[i], &(*hb)[rc[i]]);
983 /* then the whole rtp */
984 for(i=0; i < nres; i++) {
985 /* Here we allow a mismatch of one character when looking for the rtp entry.
986 * For such a mismatch there should be only one mismatching name.
987 * This is mainly useful for small molecules such as ions.
988 * Note that this will usually not work for protein, DNA and RNA,
989 * since there the residue names should be listed in residuetypes.dat
990 * and an error will have been generated earlier in the process.
992 key = *resinfo[i].rtp;
993 snew(resinfo[i].rtp,1);
994 *resinfo[i].rtp = search_rtp(key,nrtp,rtp);
995 res = get_restp(*resinfo[i].rtp,nrtp,rtp);
996 copy_t_restp(res, &(*restp)[i]);
998 /* Check that we do not have different bonded types in one molecule */
999 check_restp_types(&(*restp)[0],&(*restp)[i]);
1001 tern = -1;
1002 for(j=0; j<nterpairs && tern==-1; j++) {
1003 if (i == rn[j]) {
1004 tern = j;
1007 terc = -1;
1008 for(j=0; j<nterpairs && terc == -1; j++) {
1009 if (i == rc[j]) {
1010 terc = j;
1013 bRM = merge_t_bondeds(res->rb, (*hb)[i].rb,tern>=0,terc>=0);
1015 if (bRM && ((tern >= 0 && ntdb[tern] == NULL) ||
1016 (terc >= 0 && ctdb[terc] == NULL))) {
1017 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.");
1019 if (bRM && ((tern >= 0 && ntdb[tern]->nhack == 0) ||
1020 (terc >= 0 && ctdb[terc]->nhack == 0))) {
1021 gmx_fatal(FARGS,"There is a dangling bond at at least one of the terminal ends. Select a proper terminal entry.");
1025 /* now perform t_hack's on t_restp's,
1026 i.e. add's and deletes from termini database will be
1027 added to/removed from residue topology
1028 we'll do this on one big dirty loop, so it won't make easy reading! */
1029 for(i=0; i < nres; i++)
1031 for(j=0; j < (*hb)[i].nhack; j++)
1033 if ( (*hb)[i].hack[j].nr )
1035 /* find atom in restp */
1036 for(l=0; l < (*restp)[i].natom; l++)
1037 if ( ( (*hb)[i].hack[j].oname==NULL &&
1038 strcmp((*hb)[i].hack[j].AI, *(*restp)[i].atomname[l])==0 ) ||
1039 ( (*hb)[i].hack[j].oname!=NULL &&
1040 strcmp((*hb)[i].hack[j].oname,*(*restp)[i].atomname[l])==0 ) )
1041 break;
1042 if (l == (*restp)[i].natom)
1044 /* If we are doing an atom rename only, we don't need
1045 * to generate a fatal error if the old name is not found
1046 * in the rtp.
1048 /* Deleting can happen also only on the input atoms,
1049 * not necessarily always on the rtp entry.
1051 if (!((*hb)[i].hack[j].oname != NULL &&
1052 (*hb)[i].hack[j].nname != NULL) &&
1053 !((*hb)[i].hack[j].oname != NULL &&
1054 (*hb)[i].hack[j].nname == NULL))
1056 gmx_fatal(FARGS,
1057 "atom %s not found in buiding block %d%s "
1058 "while combining tdb and rtp",
1059 (*hb)[i].hack[j].oname!=NULL ?
1060 (*hb)[i].hack[j].oname : (*hb)[i].hack[j].AI,
1061 i+1,*resinfo[i].rtp);
1064 else
1066 if ( (*hb)[i].hack[j].oname==NULL ) {
1067 /* we're adding: */
1068 add_atom_to_restp(&(*restp)[i],resinfo[i].nr,l,
1069 &(*hb)[i].hack[j]);
1071 else
1073 /* oname != NULL */
1074 if ( (*hb)[i].hack[j].nname==NULL ) {
1075 /* we're deleting */
1076 if (debug)
1077 fprintf(debug, "deleting atom %s from res %d%s in rtp\n",
1078 *(*restp)[i].atomname[l],
1079 i+1,(*restp)[i].resname);
1080 /* shift the rest */
1081 (*restp)[i].natom--;
1082 for(k=l; k < (*restp)[i].natom; k++) {
1083 (*restp)[i].atom [k] = (*restp)[i].atom [k+1];
1084 (*restp)[i].atomname[k] = (*restp)[i].atomname[k+1];
1085 (*restp)[i].cgnr [k] = (*restp)[i].cgnr [k+1];
1087 /* give back space */
1088 srenew((*restp)[i].atom, (*restp)[i].natom);
1089 srenew((*restp)[i].atomname, (*restp)[i].natom);
1090 srenew((*restp)[i].cgnr, (*restp)[i].natom);
1091 } else { /* nname != NULL */
1092 /* we're replacing */
1093 if (debug)
1094 fprintf(debug, "replacing atom %s by %s in res %d%s in rtp\n",
1095 *(*restp)[i].atomname[l], (*hb)[i].hack[j].nname,
1096 i+1,(*restp)[i].resname);
1097 snew( (*restp)[i].atomname[l], 1);
1098 (*restp)[i].atom[l] = *(*hb)[i].hack[j].atom;
1099 *(*restp)[i].atomname[l] = strdup((*hb)[i].hack[j].nname);
1100 if ( (*hb)[i].hack[j].cgnr != NOTSET )
1101 (*restp)[i].cgnr [l] = (*hb)[i].hack[j].cgnr;
1110 static gmx_bool atomname_cmp_nr(const char *anm,t_hack *hack,int *nr)
1113 if (hack->nr == 1)
1115 *nr = 0;
1117 return (gmx_strcasecmp(anm,hack->nname) == 0);
1119 else
1121 if (isdigit(anm[strlen(anm)-1]))
1123 *nr = anm[strlen(anm)-1] - '0';
1125 else
1127 *nr = 0;
1129 if (*nr <= 0 || *nr > hack->nr)
1131 return FALSE;
1133 else
1135 return (strlen(anm) == strlen(hack->nname) + 1 &&
1136 gmx_strncasecmp(anm,hack->nname,strlen(hack->nname)) == 0);
1141 static gmx_bool match_atomnames_with_rtp_atom(t_atoms *pdba,rvec *x,int atind,
1142 t_restp *rptr,t_hackblock *hbr,
1143 gmx_bool bVerbose)
1145 int resnr;
1146 int i,j,k;
1147 char *oldnm,*newnm;
1148 int anmnr;
1149 char *start_at,buf[STRLEN];
1150 int start_nr;
1151 gmx_bool bReplaceReplace,bFoundInAdd;
1152 gmx_bool bDeleted;
1154 oldnm = *pdba->atomname[atind];
1155 resnr = pdba->resinfo[pdba->atom[atind].resind].nr;
1157 bDeleted = FALSE;
1158 for(j=0; j<hbr->nhack; j++)
1160 if (hbr->hack[j].oname != NULL && hbr->hack[j].nname != NULL &&
1161 gmx_strcasecmp(oldnm,hbr->hack[j].oname) == 0)
1163 /* This is a replace entry. */
1164 /* Check if we are not replacing a replaced atom. */
1165 bReplaceReplace = FALSE;
1166 for(k=0; k<hbr->nhack; k++) {
1167 if (k != j &&
1168 hbr->hack[k].oname != NULL && hbr->hack[k].nname != NULL &&
1169 gmx_strcasecmp(hbr->hack[k].nname,hbr->hack[j].oname) == 0)
1171 /* The replace in hack[j] replaces an atom that
1172 * was already replaced in hack[k], we do not want
1173 * second or higher level replaces at this stage.
1175 bReplaceReplace = TRUE;
1178 if (bReplaceReplace)
1180 /* Skip this replace. */
1181 continue;
1184 /* This atom still has the old name, rename it */
1185 newnm = hbr->hack[j].nname;
1186 for(k=0; k<rptr->natom; k++)
1188 if (gmx_strcasecmp(newnm,*rptr->atomname[k]) == 0)
1190 break;
1193 if (k == rptr->natom)
1195 /* The new name is not present in the rtp.
1196 * We need to apply the replace also to the rtp entry.
1199 /* We need to find the add hack that can add this atom
1200 * to find out after which atom it should be added.
1202 bFoundInAdd = FALSE;
1203 for(k=0; k<hbr->nhack; k++)
1205 if (hbr->hack[k].oname == NULL &&
1206 hbr->hack[k].nname != NULL &&
1207 atomname_cmp_nr(newnm,&hbr->hack[k],&anmnr))
1209 if (anmnr <= 1)
1211 start_at = hbr->hack[k].a[0];
1213 else
1215 sprintf(buf,"%s%d",hbr->hack[k].nname,anmnr-1);
1216 start_at = buf;
1218 for(start_nr=0; start_nr<rptr->natom; start_nr++)
1220 if (gmx_strcasecmp(start_at,(*rptr->atomname[start_nr])) == 0)
1222 break;
1225 if (start_nr == rptr->natom)
1227 gmx_fatal(FARGS,"Could not find atom '%s' in residue building block '%s' to add atom '%s' to",
1228 start_at,rptr->resname,newnm);
1230 /* We can add the atom after atom start_nr */
1231 add_atom_to_restp(rptr,resnr,start_nr,
1232 &hbr->hack[j]);
1234 bFoundInAdd = TRUE;
1238 if (!bFoundInAdd)
1240 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'",
1241 newnm,
1242 hbr->hack[j].oname,hbr->hack[j].nname,
1243 rptr->resname);
1247 if (bVerbose)
1249 printf("Renaming atom '%s' in residue '%s' %d to '%s'\n",
1250 oldnm,rptr->resname,resnr,newnm);
1252 /* Rename the atom in pdba */
1253 snew(pdba->atomname[atind],1);
1254 *pdba->atomname[atind] = strdup(newnm);
1256 else if (hbr->hack[j].oname != NULL && hbr->hack[j].nname == NULL &&
1257 gmx_strcasecmp(oldnm,hbr->hack[j].oname) == 0)
1259 /* This is a delete entry, check if this atom is present
1260 * in the rtp entry of this residue.
1262 for(k=0; k<rptr->natom; k++)
1264 if (gmx_strcasecmp(oldnm,*rptr->atomname[k]) == 0)
1266 break;
1269 if (k == rptr->natom)
1271 /* This atom is not present in the rtp entry,
1272 * delete is now from the input pdba.
1274 if (bVerbose)
1276 printf("Deleting atom '%s' in residue '%s' %d\n",
1277 oldnm,rptr->resname,resnr);
1279 /* We should free the atom name,
1280 * but it might be used multiple times in the symtab.
1281 * sfree(pdba->atomname[atind]);
1283 for(k=atind+1; k<pdba->nr; k++)
1285 pdba->atom[k-1] = pdba->atom[k];
1286 pdba->atomname[k-1] = pdba->atomname[k];
1287 copy_rvec(x[k],x[k-1]);
1289 pdba->nr--;
1290 bDeleted = TRUE;
1295 return bDeleted;
1298 void match_atomnames_with_rtp(t_restp restp[],t_hackblock hb[],
1299 t_atoms *pdba,rvec *x,
1300 gmx_bool bVerbose)
1302 int i,j,k;
1303 char *oldnm,*newnm;
1304 int resnr;
1305 t_restp *rptr;
1306 t_hackblock *hbr;
1307 int anmnr;
1308 char *start_at,buf[STRLEN];
1309 int start_nr;
1310 gmx_bool bFoundInAdd;
1312 for(i=0; i<pdba->nr; i++)
1314 oldnm = *pdba->atomname[i];
1315 resnr = pdba->resinfo[pdba->atom[i].resind].nr;
1316 rptr = &restp[pdba->atom[i].resind];
1317 for(j=0; (j<rptr->natom); j++)
1319 if (gmx_strcasecmp(oldnm,*(rptr->atomname[j])) == 0)
1321 break;
1324 if (j == rptr->natom)
1326 /* Not found yet, check if we have to rename this atom */
1327 if (match_atomnames_with_rtp_atom(pdba,x,i,
1328 rptr,&(hb[pdba->atom[i].resind]),
1329 bVerbose))
1331 /* We deleted this atom, decrease the atom counter by 1. */
1332 i--;
1338 void gen_cmap(t_params *psb, t_restp *restp, int natoms, t_atom atom[], char **aname[], int nres)
1340 int residx,i,ii,j,k;
1341 atom_id ai,aj,ak,al,am;
1342 const char *ptr;
1344 if (debug)
1345 ptr = "cmap";
1346 else
1347 ptr = "check";
1349 fprintf(stderr,"Making cmap torsions...");
1350 i=0;
1351 /* End loop at nres-1, since the very last residue does not have a +N atom, and
1352 * therefore we get a valgrind invalid 4 byte read error with atom am */
1353 for(residx=0; residx<nres-1; residx++)
1355 /* Add CMAP terms from the list of CMAP interactions */
1356 for(j=0;j<restp[residx].rb[ebtsCMAP].nb; j++)
1358 ai=search_atom(restp[residx].rb[ebtsCMAP].b[j].a[0],i,natoms,atom,aname,
1359 ptr,TRUE);
1360 aj=search_atom(restp[residx].rb[ebtsCMAP].b[j].a[1],i,natoms,atom,aname,
1361 ptr,TRUE);
1362 ak=search_atom(restp[residx].rb[ebtsCMAP].b[j].a[2],i,natoms,atom,aname,
1363 ptr,TRUE);
1364 al=search_atom(restp[residx].rb[ebtsCMAP].b[j].a[3],i,natoms,atom,aname,
1365 ptr,TRUE);
1366 am=search_atom(restp[residx].rb[ebtsCMAP].b[j].a[4],i,natoms,atom,aname,
1367 ptr,TRUE);
1369 /* The first and last residues no not have cmap torsions */
1370 if(ai!=NO_ATID && aj!=NO_ATID && ak!=NO_ATID && al!=NO_ATID && am!=NO_ATID)
1372 add_cmap_param(psb,ai,aj,ak,al,am,restp[residx].rb[ebtsCMAP].b[j].s);
1376 if(residx<nres-1)
1378 while(atom[i].resind<residx+1)
1380 i++;
1385 /* Start the next residue */
1388 static void
1389 scrub_charge_groups(int *cgnr, int natoms)
1391 int i;
1393 for(i=0;i<natoms;i++)
1395 cgnr[i]=i+1;
1400 void pdb2top(FILE *top_file, char *posre_fn, char *molname,
1401 t_atoms *atoms, rvec **x, gpp_atomtype_t atype, t_symtab *tab,
1402 int nrtp, t_restp rtp[],
1403 t_restp *restp, t_hackblock *hb,
1404 int nterpairs,t_hackblock **ntdb, t_hackblock **ctdb,
1405 gmx_bool bAllowMissing,
1406 gmx_bool bVsites, gmx_bool bVsiteAromatics,
1407 const char *ff, const char *ffdir,
1408 real mHmult,
1409 int nssbonds, t_ssbond *ssbonds,
1410 real long_bond_dist, real short_bond_dist,
1411 gmx_bool bDeuterate, gmx_bool bChargeGroups, gmx_bool bCmap,
1412 gmx_bool bRenumRes,gmx_bool bRTPresname)
1415 t_hackblock *hb;
1416 t_restp *restp;
1418 t_params plist[F_NRE];
1419 t_excls *excls;
1420 t_nextnb nnb;
1421 int *cgnr;
1422 int *vsite_type;
1423 int i,nmissat;
1424 int bts[ebtsNR];
1426 init_plist(plist);
1428 if (debug) {
1429 print_resall(debug, atoms->nres, restp, atype);
1430 dump_hb(debug, atoms->nres, hb);
1433 /* Make bonds */
1434 at2bonds(&(plist[F_BONDS]), hb,
1435 atoms->nr, atoms->atom, atoms->atomname, atoms->nres, *x,
1436 long_bond_dist, short_bond_dist, bAllowMissing);
1438 /* specbonds: disulphide bonds & heme-his */
1439 do_ssbonds(&(plist[F_BONDS]),
1440 atoms->nr, atoms->atom, atoms->atomname, nssbonds, ssbonds,
1441 bAllowMissing);
1443 nmissat = name2type(atoms, &cgnr, atype, restp);
1444 if (nmissat) {
1445 if (bAllowMissing)
1446 fprintf(stderr,"There were %d missing atoms in molecule %s\n",
1447 nmissat,molname);
1448 else
1449 gmx_fatal(FARGS,"There were %d missing atoms in molecule %s, if you want to use this incomplete topology anyhow, use the option -missing",
1450 nmissat,molname);
1453 /* Cleanup bonds (sort and rm doubles) */
1454 clean_bonds(&(plist[F_BONDS]));
1456 snew(vsite_type,atoms->nr);
1457 for(i=0; i<atoms->nr; i++)
1458 vsite_type[i]=NOTSET;
1459 if (bVsites) {
1460 /* determine which atoms will be vsites and add dummy masses
1461 also renumber atom numbers in plist[0..F_NRE]! */
1462 do_vsites(nrtp, rtp, atype, atoms, tab, x, plist,
1463 &vsite_type, &cgnr, mHmult, bVsiteAromatics, ffdir);
1466 /* Make Angles and Dihedrals */
1467 fprintf(stderr,"Generating angles, dihedrals and pairs...\n");
1468 snew(excls,atoms->nr);
1469 init_nnb(&nnb,atoms->nr,4);
1470 gen_nnb(&nnb,plist);
1471 print_nnb(&nnb,"NNB");
1472 gen_pad(&nnb,atoms,restp[0].nrexcl,restp[0].HH14,
1473 plist,excls,hb,restp[0].bAlldih,restp[0].bRemoveDih,
1474 bAllowMissing);
1475 done_nnb(&nnb);
1477 /* Make CMAP */
1478 if (TRUE == bCmap)
1480 gen_cmap(&(plist[F_CMAP]), restp, atoms->nr, atoms->atom, atoms->atomname, atoms->nres);
1481 if (plist[F_CMAP].nr > 0)
1483 fprintf(stderr, "There are %4d cmap torsion pairs\n",
1484 plist[F_CMAP].nr);
1488 /* set mass of all remaining hydrogen atoms */
1489 if (mHmult != 1.0)
1490 do_h_mass(&(plist[F_BONDS]),vsite_type,atoms,mHmult,bDeuterate);
1491 sfree(vsite_type);
1493 /* Cleanup bonds (sort and rm doubles) */
1494 /* clean_bonds(&(plist[F_BONDS]));*/
1496 fprintf(stderr,
1497 "There are %4d dihedrals, %4d impropers, %4d angles\n"
1498 " %4d pairs, %4d bonds and %4d virtual sites\n",
1499 plist[F_PDIHS].nr, plist[F_IDIHS].nr, plist[F_ANGLES].nr,
1500 plist[F_LJ14].nr, plist[F_BONDS].nr,
1501 plist[F_VSITE2].nr +
1502 plist[F_VSITE3].nr +
1503 plist[F_VSITE3FD].nr +
1504 plist[F_VSITE3FAD].nr +
1505 plist[F_VSITE3OUT].nr +
1506 plist[F_VSITE4FD].nr +
1507 plist[F_VSITE4FDN].nr );
1509 print_sums(atoms, FALSE);
1511 if (FALSE == bChargeGroups)
1513 scrub_charge_groups(cgnr, atoms->nr);
1516 if (bRenumRes)
1518 for(i=0; i<atoms->nres; i++)
1520 atoms->resinfo[i].nr = i + 1;
1521 atoms->resinfo[i].ic = ' ';
1525 if (top_file) {
1526 fprintf(stderr,"Writing topology\n");
1527 /* We can copy the bonded types from the first restp,
1528 * since the types have to be identical for all residues in one molecule.
1530 for(i=0; i<ebtsNR; i++) {
1531 bts[i] = restp[0].rb[i].type;
1533 write_top(top_file, posre_fn, molname,
1534 atoms, bRTPresname,
1535 bts, plist, excls, atype, cgnr, restp[0].nrexcl);
1538 /* cleaning up */
1539 free_t_hackblock(atoms->nres, &hb);
1540 free_t_restp(atoms->nres, &restp);
1542 /* we should clean up hb and restp here, but that is a *L*O*T* of work! */
1543 sfree(cgnr);
1544 for (i=0; i<F_NRE; i++)
1545 sfree(plist[i].param);
1546 for (i=0; i<atoms->nr; i++)
1547 sfree(excls[i].e);
1548 sfree(excls);