f3fffdd685cd0c02d78710bd03c277d6cdebca02
[alexxy/gromacs.git] / src / kernel / gen_vsite.c
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  * 
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  * 
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "string2.h"
40 #include <stdio.h>
41 #include <math.h>
42 #include <string.h>
43 #include "gen_vsite.h"
44 #include "smalloc.h"
45 #include "resall.h"
46 #include "add_par.h"
47 #include "vec.h"
48 #include "toputil.h"
49 #include "physics.h"
50 #include "index.h"
51 #include "names.h"
52 #include "futil.h"
53 #include "gpp_atomtype.h"
54 #include "fflibutil.h"
55
56 #define MAXNAME 32
57 #define OPENDIR         '['     /* starting sign for directive          */
58 #define CLOSEDIR        ']'     /* ending sign for directive            */
59
60 typedef struct {
61   char   atomtype[MAXNAME];      /* Type for the XH3/XH2 atom */
62   gmx_bool   isplanar;               /* If true, the atomtype above and the three connected
63                                   * ones are in a planar geometry. The two next entries
64                                   * are undefined in that case 
65                                   */
66   int    nhydrogens;             /* number of connected hydrogens */
67   char   nextheavytype[MAXNAME]; /* Type for the heavy atom bonded to XH2/XH3 */
68   char   dummymass[MAXNAME];     /* The type of MNH* or MCH3* dummy mass to use */
69 } t_vsiteconf;
70
71
72 /* Structure to represent average bond and angles values in vsite aromatic
73  * residues. Note that these are NOT necessarily the bonds and angles from the
74  * forcefield; many forcefields (like Amber, OPLS) have some inherent strain in
75  * 5-rings (i.e. the sum of angles is !=540, but impropers keep it planar)
76  */
77 typedef struct {
78   char resname[MAXNAME];
79   int nbonds;
80   int nangles;
81   struct vsitetop_bond {
82     char   atom1[MAXNAME];
83     char   atom2[MAXNAME];
84     float  value;
85   } *bond; /* list of bonds */
86   struct vsitetop_angle {
87     char   atom1[MAXNAME];
88     char   atom2[MAXNAME];
89     char   atom3[MAXNAME];
90     float  value;
91   } *angle; /* list of angles */
92 } t_vsitetop;
93
94
95 enum { DDB_CH3, DDB_NH3, DDB_NH2, DDB_PHE, DDB_TYR, 
96        DDB_TRP, DDB_HISA, DDB_HISB, DDB_HISH , DDB_DIR_NR };
97
98 typedef char t_dirname[STRLEN];
99
100 static const t_dirname ddb_dirnames[DDB_DIR_NR] = {
101   "CH3",
102   "NH3",
103   "NH2",
104   "PHE",
105   "TYR",
106   "TRP",
107   "HISA",
108   "HISB",
109   "HISH"
110 };
111
112 static int ddb_name2dir(char *name) 
113 {
114   /* Translate a directive name to the number of the directive.
115    * HID/HIE/HIP names are translated to the ones we use in Gromacs.
116    */
117
118   int i,index;
119
120   index=-1;
121
122   for(i=0;i<DDB_DIR_NR && index<0;i++)
123     if(!gmx_strcasecmp(name,ddb_dirnames[i]))
124       index=i;
125   
126   return index;
127 }
128
129
130 static void read_vsite_database(const char *ddbname,
131                                 t_vsiteconf **pvsiteconflist, int *nvsiteconf,
132                                 t_vsitetop **pvsitetoplist, int *nvsitetop)
133 {
134   /* This routine is a quick hack to fix the problem with hardcoded atomtypes
135    * and aromatic vsite parameters by reading them from a ff???.vsd file.
136    *
137    * The file can contain sections [ NH3 ], [ CH3 ], [ NH2 ], and ring residue names.
138    * For the NH3 and CH3 section each line has three fields. The first is the atomtype 
139    * (nb: not bonded type) of the N/C atom to be replaced, the second field is
140    * the type of the next heavy atom it is bonded to, and the third field the type
141    * of dummy mass that will be used for this group. 
142    * 
143    * If the NH2 group planar (sp2 N) a different vsite construct is used, so in this
144    * case the second field should just be the word planar.
145    */
146
147   FILE *ddb;
148   char dirstr[STRLEN];
149   char pline[STRLEN];
150   int i,j,n,k,nvsite,ntop,curdir,prevdir;
151   t_vsiteconf *vsiteconflist;
152   t_vsitetop *vsitetoplist;
153   char *ch;
154   char s1[MAXNAME],s2[MAXNAME],s3[MAXNAME],s4[MAXNAME];
155
156   ddb = libopen(ddbname);
157
158   nvsite        = *nvsiteconf;
159   vsiteconflist = *pvsiteconflist;
160   ntop          = *nvsitetop;
161   vsitetoplist  = *pvsitetoplist;
162
163   curdir=-1;
164   
165   snew(vsiteconflist,1);
166   snew(vsitetoplist,1);
167
168   while(fgets2(pline,STRLEN-2,ddb) != NULL) {
169     strip_comment(pline);
170     trim(pline);
171     if(strlen(pline)>0) {
172       if(pline[0] == OPENDIR ) {
173         strncpy(dirstr,pline+1,STRLEN-2);
174         if ((ch = strchr (dirstr,CLOSEDIR)) != NULL)
175           (*ch) = 0;
176         trim (dirstr);
177
178         if(!gmx_strcasecmp(dirstr,"HID") ||
179            !gmx_strcasecmp(dirstr,"HISD"))
180           sprintf(dirstr,"HISA");
181         else if(!gmx_strcasecmp(dirstr,"HIE") ||
182                 !gmx_strcasecmp(dirstr,"HISE"))
183           sprintf(dirstr,"HISB");
184         else if(!gmx_strcasecmp(dirstr,"HIP"))
185           sprintf(dirstr,"HISH");
186
187         curdir=ddb_name2dir(dirstr);
188         if (curdir < 0) {
189           gmx_fatal(FARGS,"Invalid directive %s in vsite database %s",
190                     dirstr,ddbname);
191         }
192       } else {
193         switch(curdir) {
194         case -1:
195           gmx_fatal(FARGS,"First entry in vsite database must be a directive.\n");
196           break;
197         case DDB_CH3:
198         case DDB_NH3:
199         case DDB_NH2:
200           n = sscanf(pline,"%s%s%s",s1,s2,s3);
201           if(n<3 && !gmx_strcasecmp(s2,"planar")) {
202             srenew(vsiteconflist,nvsite+1);
203             strncpy(vsiteconflist[nvsite].atomtype,s1,MAXNAME-1);
204             vsiteconflist[nvsite].isplanar=TRUE;
205             vsiteconflist[nvsite].nextheavytype[0]=0;
206             vsiteconflist[nvsite].dummymass[0]=0;
207             vsiteconflist[nvsite].nhydrogens=2;
208             nvsite++;
209           } else if(n==3) {
210             srenew(vsiteconflist,(nvsite+1));
211             strncpy(vsiteconflist[nvsite].atomtype,s1,MAXNAME-1);
212             vsiteconflist[nvsite].isplanar=FALSE;
213             strncpy(vsiteconflist[nvsite].nextheavytype,s2,MAXNAME-1);
214             strncpy(vsiteconflist[nvsite].dummymass,s3,MAXNAME-1);
215             if(curdir==DDB_NH2)
216               vsiteconflist[nvsite].nhydrogens=2;
217             else
218               vsiteconflist[nvsite].nhydrogens=3;
219             nvsite++;
220           } else {
221             gmx_fatal(FARGS,"Not enough directives in vsite database line: %s\n",pline);
222           }
223           break;
224         case DDB_PHE:
225         case DDB_TYR:
226         case DDB_TRP:
227         case DDB_HISA:
228         case DDB_HISB:
229         case DDB_HISH:
230           i=0;
231           while((i<ntop) && gmx_strcasecmp(dirstr,vsitetoplist[i].resname))
232             i++;
233           /* Allocate a new topology entry if this is a new residue */
234           if(i==ntop) {
235             srenew(vsitetoplist,ntop+1);            
236             ntop++; /* i still points to current vsite topology entry */
237             strncpy(vsitetoplist[i].resname,dirstr,MAXNAME-1);
238             vsitetoplist[i].nbonds=vsitetoplist[i].nangles=0;
239             snew(vsitetoplist[i].bond,1);
240             snew(vsitetoplist[i].angle,1);
241           }
242           n = sscanf(pline,"%s%s%s%s",s1,s2,s3,s4);
243           if(n==3) {
244             /* bond */
245             k=vsitetoplist[i].nbonds++;
246             srenew(vsitetoplist[i].bond,k+1);
247             strncpy(vsitetoplist[i].bond[k].atom1,s1,MAXNAME-1);
248             strncpy(vsitetoplist[i].bond[k].atom2,s2,MAXNAME-1);
249             vsitetoplist[i].bond[k].value=strtod(s3,NULL);
250           } else if (n==4) {
251             /* angle */
252             k=vsitetoplist[i].nangles++;
253             srenew(vsitetoplist[i].angle,k+1);
254             strncpy(vsitetoplist[i].angle[k].atom1,s1,MAXNAME-1);
255             strncpy(vsitetoplist[i].angle[k].atom2,s2,MAXNAME-1);
256             strncpy(vsitetoplist[i].angle[k].atom3,s3,MAXNAME-1);
257             vsitetoplist[i].angle[k].value=strtod(s4,NULL);
258           } else {
259             gmx_fatal(FARGS,"Need 3 or 4 values to specify bond/angle values in %s: %s\n",ddbname,pline);
260           }
261           break;
262         default:
263           gmx_fatal(FARGS,"Didnt find a case for directive %s in read_vsite_database\n",dirstr);
264         }
265       }
266     }
267   }
268
269   *pvsiteconflist=vsiteconflist;
270   *pvsitetoplist=vsitetoplist;
271   *nvsiteconf=nvsite;
272   *nvsitetop=ntop;
273   
274   ffclose(ddb);
275 }
276
277 static int nitrogen_is_planar(t_vsiteconf vsiteconflist[],int nvsiteconf,char atomtype[])
278 {
279   /* Return 1 if atomtype exists in database list and is planar, 0 if not,
280    * and -1 if not found.
281    */
282   int i,res;
283   gmx_bool found=FALSE;
284   for(i=0;i<nvsiteconf && !found;i++) {
285     found=(!gmx_strcasecmp(vsiteconflist[i].atomtype,atomtype) && (vsiteconflist[i].nhydrogens==2));
286   }
287   if(found)
288     res=(vsiteconflist[i-1].isplanar==TRUE);
289   else
290     res=-1;
291
292   return res;
293 }
294   
295 static char *get_dummymass_name(t_vsiteconf vsiteconflist[],int nvsiteconf,char atom[], char nextheavy[])
296 {
297   /* Return the dummy mass name if found, or NULL if not set in ddb database */
298   int i;
299   gmx_bool found=FALSE;
300   for(i=0;i<nvsiteconf && !found;i++) {
301     found=(!gmx_strcasecmp(vsiteconflist[i].atomtype,atom) &&
302            !gmx_strcasecmp(vsiteconflist[i].nextheavytype,nextheavy));
303   }
304   if(found)
305     return vsiteconflist[i-1].dummymass;
306   else
307     return NULL;
308 }
309   
310
311
312 static real get_ddb_bond(t_vsitetop *vsitetop, int nvsitetop,
313                          const char res[],
314                          const char atom1[], const char atom2[])
315 {
316   int i,j;
317   
318   i=0;
319   while(i<nvsitetop && gmx_strcasecmp(res,vsitetop[i].resname))
320     i++;
321   if(i==nvsitetop)
322     gmx_fatal(FARGS,"No vsite information for residue %s found in vsite database.\n",res);
323   j=0;
324   while(j<vsitetop[i].nbonds && 
325         ( strcmp(atom1,vsitetop[i].bond[j].atom1) || strcmp(atom2,vsitetop[i].bond[j].atom2)) &&
326         ( strcmp(atom2,vsitetop[i].bond[j].atom1) || strcmp(atom1,vsitetop[i].bond[j].atom2)))
327     j++;
328   if(j==vsitetop[i].nbonds)
329     gmx_fatal(FARGS,"Couldnt find bond %s-%s for residue %s in vsite database.\n",atom1,atom2,res);
330   
331   return vsitetop[i].bond[j].value;
332 }
333       
334
335 static real get_ddb_angle(t_vsitetop *vsitetop, int nvsitetop,
336                           const char res[], const char atom1[],
337                           const char atom2[], const char atom3[])
338 {
339   int i,j;
340   
341   i=0;
342   while(i<nvsitetop && gmx_strcasecmp(res,vsitetop[i].resname))
343     i++;
344   if(i==nvsitetop)
345     gmx_fatal(FARGS,"No vsite information for residue %s found in vsite database.\n",res);
346   j=0;
347   while(j<vsitetop[i].nangles && 
348         ( strcmp(atom1,vsitetop[i].angle[j].atom1) || 
349           strcmp(atom2,vsitetop[i].angle[j].atom2) || 
350           strcmp(atom3,vsitetop[i].angle[j].atom3)) &&
351         ( strcmp(atom3,vsitetop[i].angle[j].atom1) || 
352           strcmp(atom2,vsitetop[i].angle[j].atom2) || 
353           strcmp(atom1,vsitetop[i].angle[j].atom3)))
354     j++;
355   if(j==vsitetop[i].nangles)
356     gmx_fatal(FARGS,"Couldnt find angle %s-%s-%s for residue %s in vsite database.\n",atom1,atom2,atom3,res);
357   
358   return vsitetop[i].angle[j].value;
359 }
360       
361
362 static void count_bonds(int atom, t_params *psb, char ***atomname, 
363                         int *nrbonds, int *nrHatoms, int Hatoms[], int *Heavy,
364                         int *nrheavies, int heavies[])
365 {
366   int i,heavy,other,nrb,nrH,nrhv;
367   
368   /* find heavy atom bound to this hydrogen */
369   heavy=NOTSET;
370   for(i=0; (i<psb->nr) && (heavy==NOTSET); i++)
371     if (psb->param[i].AI==atom)
372       heavy=psb->param[i].AJ;
373     else if (psb->param[i].AJ==atom)
374       heavy=psb->param[i].AI;
375   if (heavy==NOTSET)
376     gmx_fatal(FARGS,"unbound hydrogen atom %d",atom+1);
377   /* find all atoms bound to heavy atom */
378   other=NOTSET;
379   nrb=0;
380   nrH=0;
381   nrhv=0;
382   for(i=0; i<psb->nr; i++) {
383     if (psb->param[i].AI==heavy)
384       other=psb->param[i].AJ;
385     else if (psb->param[i].AJ==heavy)
386       other=psb->param[i].AI;
387     if (other!=NOTSET) {
388       nrb++;
389       if (is_hydrogen(*(atomname[other]))) {
390         Hatoms[nrH]=other;
391         nrH++;
392       } else {
393         heavies[nrhv]=other;
394         nrhv++;
395       }
396       other=NOTSET;
397     }
398   }
399   *Heavy   =heavy;
400   *nrbonds =nrb;
401   *nrHatoms=nrH;
402   *nrheavies=nrhv;
403 }
404
405 static void print_bonds(FILE *fp, int o2n[],
406                         int nrHatoms, int Hatoms[], int Heavy,
407                         int nrheavies, int heavies[])
408 {
409   int i;
410   
411   fprintf(fp,"Found: %d Hatoms: ",nrHatoms);
412   for(i=0; i<nrHatoms; i++)
413     fprintf(fp," %d",o2n[Hatoms[i]]+1);
414   fprintf(fp,"; %d Heavy atoms: %d",nrheavies+1,o2n[Heavy]+1);
415   for(i=0; i<nrheavies; i++)
416     fprintf(fp," %d",o2n[heavies[i]]+1);
417   fprintf(fp,"\n");
418 }
419
420 static int get_atype(int atom, t_atoms *at, int nrtp, t_restp rtp[],
421                      gmx_residuetype_t rt)
422 {
423   int type;
424   gmx_bool bNterm;
425   int  j;
426   t_restp *rtpp;
427   
428   if (at->atom[atom].m)
429     type=at->atom[atom].type;
430   else {
431     /* get type from rtp */
432     rtpp = get_restp(*(at->resinfo[at->atom[atom].resind].name),nrtp,rtp);
433     bNterm = gmx_residuetype_is_protein(rt,*(at->resinfo[at->atom[atom].resind].name)) && 
434       (at->atom[atom].resind == 0);
435     j=search_jtype(rtpp,*(at->atomname[atom]),bNterm);
436     type=rtpp->atom[j].type;
437   }
438   return type;
439 }
440
441 static int vsite_nm2type(const char *name, gpp_atomtype_t atype)
442 {
443   int tp;
444   
445   tp = get_atomtype_type(name,atype);
446   if (tp == NOTSET)
447     gmx_fatal(FARGS,"Dummy mass type (%s) not found in atom type database",
448               name);
449               
450   return tp;
451 }
452
453 static real get_amass(int atom, t_atoms *at, int nrtp, t_restp rtp[],
454                       gmx_residuetype_t rt)
455 {
456   real mass;
457   gmx_bool bNterm;
458   int  j;
459   t_restp *rtpp;
460   
461   if (at->atom[atom].m)
462     mass=at->atom[atom].m;
463   else {
464     /* get mass from rtp */
465     rtpp = get_restp(*(at->resinfo[at->atom[atom].resind].name),nrtp,rtp);
466     bNterm = gmx_residuetype_is_protein(rt,*(at->resinfo[at->atom[atom].resind].name)) && 
467       (at->atom[atom].resind == 0);
468     j=search_jtype(rtpp,*(at->atomname[atom]),bNterm);
469     mass=rtpp->atom[j].m;
470   }
471   return mass;
472 }
473
474 static void my_add_param(t_params *plist, int ai, int aj, real b)
475 {
476   static real c[MAXFORCEPARAM] = 
477   { NOTSET, NOTSET, NOTSET, NOTSET, NOTSET, NOTSET };
478   
479   c[0]=b;
480   add_param(plist,ai,aj,c,NULL);
481 }
482
483 static void add_vsites(t_params plist[], int vsite_type[], 
484                        int Heavy, int nrHatoms, int Hatoms[], 
485                        int nrheavies, int heavies[])
486 {
487   int i,j,ftype,other,moreheavy,bb;
488   gmx_bool bSwapParity;
489   
490   for(i=0; i<nrHatoms; i++) {
491     ftype=vsite_type[Hatoms[i]];
492     /* Errors in setting the vsite_type should really be caugth earlier,
493      * because here it's not possible to print any useful error message.
494      * But it's still better to print a message than to segfault.
495      */
496     if (ftype == NOTSET) {
497       gmx_incons("Undetected error in setting up virtual sites");
498     }
499     bSwapParity = (ftype<0);
500     vsite_type[Hatoms[i]] = ftype = abs(ftype);
501     if (ftype == F_BONDS) {
502       if ( (nrheavies != 1) && (nrHatoms != 1) )
503         gmx_fatal(FARGS,"cannot make constraint in add_vsites for %d heavy "
504                     "atoms and %d hydrogen atoms",nrheavies,nrHatoms);
505       my_add_param(&(plist[F_CONSTRNC]),Hatoms[i],heavies[0],NOTSET);
506     } else {
507       switch (ftype) {
508       case F_VSITE3:
509       case F_VSITE3FD:
510       case F_VSITE3OUT:
511         if (nrheavies < 2) 
512           gmx_fatal(FARGS,"Not enough heavy atoms (%d) for %s (min 3)",
513                       nrheavies+1,
514                       interaction_function[vsite_type[Hatoms[i]]].name);
515         add_vsite3_atoms(&plist[ftype],Hatoms[i],Heavy,heavies[0],heavies[1],
516                        bSwapParity);
517         break;
518       case F_VSITE3FAD: {
519         if (nrheavies > 1)
520           moreheavy=heavies[1];
521         else {
522           /* find more heavy atoms */
523           other=moreheavy=NOTSET;
524           for(j=0; (j<plist[F_BONDS].nr) && (moreheavy==NOTSET); j++) {
525             if (plist[F_BONDS].param[j].AI==heavies[0])
526               other=plist[F_BONDS].param[j].AJ;
527             else if (plist[F_BONDS].param[j].AJ==heavies[0])
528               other=plist[F_BONDS].param[j].AI;
529             if ( (other != NOTSET) && (other != Heavy) ) 
530               moreheavy=other;
531           }
532           if (moreheavy==NOTSET)
533             gmx_fatal(FARGS,"Unbound molecule part %d-%d",Heavy+1,Hatoms[0]+1);
534         }
535         add_vsite3_atoms(&plist[ftype],Hatoms[i],Heavy,heavies[0],moreheavy,
536                        bSwapParity);
537         break;
538       }
539       case F_VSITE4FD: 
540       case F_VSITE4FDN:
541           if (nrheavies < 3) 
542           {
543               gmx_fatal(FARGS,"Not enough heavy atoms (%d) for %s (min 4)",
544                         nrheavies+1,
545                         interaction_function[vsite_type[Hatoms[i]]].name);
546           }
547           add_vsite4_atoms(&plist[ftype],  
548                            Hatoms[0], Heavy, heavies[0], heavies[1], heavies[2]);
549           break;
550       
551       default:
552         gmx_fatal(FARGS,"can't use add_vsites for interaction function %s",
553                     interaction_function[vsite_type[Hatoms[i]]].name);
554       } /* switch ftype */
555     } /* else */
556   } /* for i */
557 }
558
559 #define ANGLE_6RING (DEG2RAD*120)
560
561 /* cosine rule: a^2 = b^2 + c^2 - 2 b c cos(alpha) */
562 /* get a^2 when a, b and alpha are given: */
563 #define cosrule(b,c,alpha) ( sqr(b) + sqr(c) - 2*b*c*cos(alpha) )
564 /* get cos(alpha) when a, b and c are given: */
565 #define acosrule(a,b,c) ( (sqr(b)+sqr(c)-sqr(a))/(2*b*c) )
566
567 static int gen_vsites_6ring(t_atoms *at, int *vsite_type[], t_params plist[], 
568                           int nrfound, int *ats, real bond_cc, real bond_ch,
569                           real xcom, real ycom, gmx_bool bDoZ)
570 {
571   /* these MUST correspond to the atnms array in do_vsite_aromatics! */
572   enum { atCG, atCD1, atHD1, atCD2, atHD2, atCE1, atHE1, atCE2, atHE2, 
573          atCZ, atHZ, atNR };
574   
575   int i,nvsite;
576    real a,b,dCGCE,tmp1,tmp2,mtot,mG,mrest;
577    real xCG,yCG,xCE1,yCE1,xCE2,yCE2;
578   /* CG, CE1 and CE2 stay and each get a part of the total mass, 
579    * so the c-o-m stays the same.
580    */
581
582   if (bDoZ) {
583     if (atNR != nrfound)
584       gmx_incons("Generating vsites on 6-rings");
585   }
586
587   /* constraints between CG, CE1 and CE2: */
588   dCGCE = sqrt( cosrule(bond_cc,bond_cc,ANGLE_6RING) );
589   my_add_param(&(plist[F_CONSTRNC]),ats[atCG] ,ats[atCE1],dCGCE);
590   my_add_param(&(plist[F_CONSTRNC]),ats[atCG] ,ats[atCE2],dCGCE);
591   my_add_param(&(plist[F_CONSTRNC]),ats[atCE1],ats[atCE2],dCGCE);
592   
593   /* rest will be vsite3 */
594   mtot=0;
595   nvsite=0;
596   for(i=0; i<  (bDoZ ? atNR : atHZ) ; i++) {
597     mtot+=at->atom[ats[i]].m;
598     if ( i!=atCG && i!=atCE1 && i!=atCE2 && (bDoZ || (i!=atHZ && i!=atCZ) ) ) {
599       at->atom[ats[i]].m = at->atom[ats[i]].mB = 0;
600       (*vsite_type)[ats[i]]=F_VSITE3;
601       nvsite++;
602     }
603   }
604   /* Distribute mass so center-of-mass stays the same. 
605    * The center-of-mass in the call is defined with x=0 at
606    * the CE1-CE2 bond and y=0 at the line from CG to the middle of CE1-CE2 bond.
607    */
608   xCG=-bond_cc+bond_cc*cos(ANGLE_6RING);
609   yCG=0;
610   xCE1=0;
611   yCE1=bond_cc*sin(0.5*ANGLE_6RING);
612   xCE2=0;
613   yCE2=-bond_cc*sin(0.5*ANGLE_6RING);  
614   
615   mG = at->atom[ats[atCG]].m = at->atom[ats[atCG]].mB = xcom*mtot/xCG;
616   mrest=mtot-mG;
617   at->atom[ats[atCE1]].m = at->atom[ats[atCE1]].mB = 
618     at->atom[ats[atCE2]].m = at->atom[ats[atCE2]].mB = mrest / 2;
619   
620   /* vsite3 construction: r_d = r_i + a r_ij + b r_ik */
621   tmp1 = dCGCE*sin(ANGLE_6RING*0.5);
622   tmp2 = bond_cc*cos(0.5*ANGLE_6RING) + tmp1;
623   tmp1 *= 2;
624   a = b = - bond_ch / tmp1;
625   /* HE1 and HE2: */
626   add_vsite3_param(&plist[F_VSITE3],
627                  ats[atHE1],ats[atCE1],ats[atCE2],ats[atCG], a,b);
628   add_vsite3_param(&plist[F_VSITE3],
629                  ats[atHE2],ats[atCE2],ats[atCE1],ats[atCG], a,b);
630   /* CD1, CD2 and CZ: */
631   a = b = tmp2 / tmp1;
632   add_vsite3_param(&plist[F_VSITE3],
633                  ats[atCD1],ats[atCE2],ats[atCE1],ats[atCG], a,b);
634   add_vsite3_param(&plist[F_VSITE3],
635                  ats[atCD2],ats[atCE1],ats[atCE2],ats[atCG], a,b);
636   if (bDoZ)
637     add_vsite3_param(&plist[F_VSITE3],
638                    ats[atCZ], ats[atCG], ats[atCE1],ats[atCE2],a,b);
639   /* HD1, HD2 and HZ: */
640   a = b = ( bond_ch + tmp2 ) / tmp1;
641   add_vsite3_param(&plist[F_VSITE3],
642                  ats[atHD1],ats[atCE2],ats[atCE1],ats[atCG], a,b);
643   add_vsite3_param(&plist[F_VSITE3],
644                  ats[atHD2],ats[atCE1],ats[atCE2],ats[atCG], a,b);
645   if (bDoZ)
646     add_vsite3_param(&plist[F_VSITE3],
647                    ats[atHZ], ats[atCG], ats[atCE1],ats[atCE2],a,b);
648   
649   return nvsite;
650 }
651
652 static int gen_vsites_phe(t_atoms *at, int *vsite_type[], t_params plist[], 
653                          int nrfound, int *ats, t_vsitetop *vsitetop, int nvsitetop)
654 {
655   real bond_cc,bond_ch;
656   real xcom,ycom,mtot;
657   int i;
658   /* these MUST correspond to the atnms array in do_vsite_aromatics! */
659   enum { atCG, atCD1, atHD1, atCD2, atHD2, atCE1, atHE1, atCE2, atHE2, 
660          atCZ, atHZ, atNR };
661   real x[atNR],y[atNR];
662   /* Aromatic rings have 6-fold symmetry, so we only need one bond length.
663    * (angle is always 120 degrees).
664    */   
665   bond_cc=get_ddb_bond(vsitetop,nvsitetop,"PHE","CD1","CE1");
666   bond_ch=get_ddb_bond(vsitetop,nvsitetop,"PHE","CD1","HD1");
667   
668   x[atCG]=-bond_cc+bond_cc*cos(ANGLE_6RING);
669   y[atCG]=0;
670   x[atCD1]=-bond_cc;
671   y[atCD1]=bond_cc*sin(0.5*ANGLE_6RING);
672   x[atHD1]=x[atCD1]+bond_ch*cos(ANGLE_6RING);
673   y[atHD1]=y[atCD1]+bond_ch*sin(ANGLE_6RING);
674   x[atCE1]=0;
675   y[atCE1]=y[atCD1];
676   x[atHE1]=x[atCE1]-bond_ch*cos(ANGLE_6RING);
677   y[atHE1]=y[atCE1]+bond_ch*sin(ANGLE_6RING);
678   x[atCD2]=x[atCD1];
679   y[atCD2]=-y[atCD1];
680   x[atHD2]=x[atHD1];
681   y[atHD2]=-y[atHD1];
682   x[atCE2]=x[atCE1];
683   y[atCE2]=-y[atCE1];
684   x[atHE2]=x[atHE1];
685   y[atHE2]=-y[atHE1];
686   x[atCZ]=bond_cc*cos(0.5*ANGLE_6RING);
687   y[atCZ]=0;
688   x[atHZ]=x[atCZ]+bond_ch;
689   y[atHZ]=0;
690
691   xcom=ycom=mtot=0;
692   for(i=0;i<atNR;i++) {
693     xcom+=x[i]*at->atom[ats[i]].m;
694     ycom+=y[i]*at->atom[ats[i]].m;
695     mtot+=at->atom[ats[i]].m;
696   }
697   xcom/=mtot;
698   ycom/=mtot;
699
700   return gen_vsites_6ring(at, vsite_type, plist, nrfound, ats, bond_cc, bond_ch, xcom, ycom, TRUE);
701 }
702
703 static void calc_vsite3_param(real xd,real yd,real xi,real yi,real xj,real yj,
704                               real xk, real yk, real *a, real *b)
705 {
706   /* determine parameters by solving the equation system, since we know the
707    * virtual site coordinates here.
708    */
709   real dx_ij,dx_ik,dy_ij,dy_ik;
710   real b_ij,b_ik;
711
712   dx_ij=xj-xi;
713   dy_ij=yj-yi;
714   dx_ik=xk-xi;
715   dy_ik=yk-yi;
716   b_ij=sqrt(dx_ij*dx_ij+dy_ij*dy_ij);
717   b_ik=sqrt(dx_ik*dx_ik+dy_ik*dy_ik);
718
719   *a = ( (xd-xi)*dy_ik - dx_ik*(yd-yi) ) / (dx_ij*dy_ik - dx_ik*dy_ij);
720   *b = ( yd - yi - (*a)*dy_ij ) / dy_ik;
721 }
722
723
724 static int gen_vsites_trp(gpp_atomtype_t atype, rvec *newx[],
725                           t_atom *newatom[], char ***newatomname[], 
726                           int *o2n[], int *newvsite_type[], int *newcgnr[],
727                           t_symtab *symtab, int *nadd, rvec x[], int *cgnr[],
728                           t_atoms *at, int *vsite_type[], t_params plist[], 
729                           int nrfound, int *ats, int add_shift,
730                           t_vsitetop *vsitetop, int nvsitetop)
731 {
732 #define NMASS 2
733   /* these MUST correspond to the atnms array in do_vsite_aromatics! */
734   enum { 
735     atCB,  atCG,  atCD1, atHD1, atCD2, atNE1, atHE1, atCE2, atCE3, atHE3, 
736     atCZ2, atHZ2, atCZ3, atHZ3, atCH2, atHH2, atNR };
737   /* weights for determining the COM's of both rings (M1 and M2): */
738   real mw[NMASS][atNR] = {
739     {   0,     1,     1,     1,   0.5,     1,     1,   0.5,     0,     0,
740         0,     0,     0,     0,     0,     0 },
741     {   0,     0,     0,     0,   0.5,     0,     0,   0.5,     1,     1,
742         1,     1,     1,     1,     1,     1 }
743   };
744  
745   real xi[atNR],yi[atNR];
746   real xcom[NMASS],ycom[NMASS],I,alpha;
747   real lineA,lineB,dist;
748   real b_CD2_CE2,b_NE1_CE2,b_CG_CD2,b_CH2_HH2,b_CE2_CZ2;
749   real b_NE1_HE1,b_CD2_CE3,b_CE3_CZ3,b_CB_CG;
750   real b_CZ2_CH2,b_CZ2_HZ2,b_CD1_HD1,b_CE3_HE3;
751   real b_CG_CD1,b_CZ3_HZ3;
752   real a_NE1_CE2_CD2,a_CE2_CD2_CG,a_CB_CG_CD2,a_CE2_CD2_CE3;
753   real a_CB_CG_CD1,a_CD2_CG_CD1,a_CE2_CZ2_HZ2,a_CZ2_CH2_HH2;
754   real a_CD2_CE2_CZ2,a_CD2_CE3_CZ3,a_CE3_CZ3_HZ3,a_CG_CD1_HD1;
755   real a_CE2_CZ2_CH2,a_HE1_NE1_CE2,a_CD2_CE3_HE3;
756   real xM[NMASS];
757   int  atM[NMASS],tpM,i,i0,j,nvsite;
758   real mwtot,mtot,mM[NMASS],dCBM1,dCBM2,dM1M2;
759   real a,b,c[MAXFORCEPARAM];
760   rvec r_ij,r_ik,t1,t2;
761   char name[10];
762
763   if (atNR != nrfound)
764     gmx_incons("atom types in gen_vsites_trp");
765   /* Get geometry from database */
766   b_CD2_CE2=get_ddb_bond(vsitetop,nvsitetop,"TRP","CD2","CE2");
767   b_NE1_CE2=get_ddb_bond(vsitetop,nvsitetop,"TRP","NE1","CE2");
768   b_CG_CD1=get_ddb_bond(vsitetop,nvsitetop,"TRP","CG","CD1");
769   b_CG_CD2=get_ddb_bond(vsitetop,nvsitetop,"TRP","CG","CD2");
770   b_CB_CG=get_ddb_bond(vsitetop,nvsitetop,"TRP","CB","CG");
771   b_CE2_CZ2=get_ddb_bond(vsitetop,nvsitetop,"TRP","CE2","CZ2");
772   b_CD2_CE3=get_ddb_bond(vsitetop,nvsitetop,"TRP","CD2","CE3");
773   b_CE3_CZ3=get_ddb_bond(vsitetop,nvsitetop,"TRP","CE3","CZ3");
774   b_CZ2_CH2=get_ddb_bond(vsitetop,nvsitetop,"TRP","CZ2","CH2");
775
776   b_CD1_HD1=get_ddb_bond(vsitetop,nvsitetop,"TRP","CD1","HD1");
777   b_CZ2_HZ2=get_ddb_bond(vsitetop,nvsitetop,"TRP","CZ2","HZ2");
778   b_NE1_HE1=get_ddb_bond(vsitetop,nvsitetop,"TRP","NE1","HE1");
779   b_CH2_HH2=get_ddb_bond(vsitetop,nvsitetop,"TRP","CH2","HH2");
780   b_CE3_HE3=get_ddb_bond(vsitetop,nvsitetop,"TRP","CE3","HE3");
781   b_CZ3_HZ3=get_ddb_bond(vsitetop,nvsitetop,"TRP","CZ3","HZ3");
782
783   a_NE1_CE2_CD2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","NE1","CE2","CD2");
784   a_CE2_CD2_CG=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CE2","CD2","CG");
785   a_CB_CG_CD2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CB","CG","CD2");
786   a_CD2_CG_CD1=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CD2","CG","CD1");
787   a_CB_CG_CD1=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CB","CG","CD1");
788  
789   a_CE2_CD2_CE3=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CE2","CD2","CE3");
790   a_CD2_CE2_CZ2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CD2","CE2","CZ2");
791   a_CD2_CE3_CZ3=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CD2","CE3","CZ3");
792   a_CE3_CZ3_HZ3=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CE3","CZ3","HZ3");
793   a_CZ2_CH2_HH2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CZ2","CH2","HH2");
794   a_CE2_CZ2_HZ2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CE2","CZ2","HZ2");
795   a_CE2_CZ2_CH2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CE2","CZ2","CH2");
796   a_CG_CD1_HD1=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CG","CD1","HD1");
797   a_HE1_NE1_CE2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","HE1","NE1","CE2");
798   a_CD2_CE3_HE3=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TRP","CD2","CE3","HE3");
799  
800   /* Calculate local coordinates.
801    * y-axis (x=0) is the bond CD2-CE2.
802    * x-axis (y=0) is perpendicular to the bond CD2-CE2 and
803    * intersects the middle of the bond.
804    */
805   xi[atCD2]=0;
806   yi[atCD2]=-0.5*b_CD2_CE2;
807   
808   xi[atCE2]=0;
809   yi[atCE2]=0.5*b_CD2_CE2;
810
811   xi[atNE1]=-b_NE1_CE2*sin(a_NE1_CE2_CD2);
812   yi[atNE1]=yi[atCE2]-b_NE1_CE2*cos(a_NE1_CE2_CD2);
813
814   xi[atCG]=-b_CG_CD2*sin(a_CE2_CD2_CG);
815   yi[atCG]=yi[atCD2]+b_CG_CD2*cos(a_CE2_CD2_CG);
816
817   alpha = a_CE2_CD2_CG + M_PI - a_CB_CG_CD2;
818   xi[atCB]=xi[atCG]-b_CB_CG*sin(alpha);
819   yi[atCB]=yi[atCG]+b_CB_CG*cos(alpha);
820
821   alpha = a_CE2_CD2_CG + a_CD2_CG_CD1 - M_PI;
822   xi[atCD1]=xi[atCG]-b_CG_CD1*sin(alpha);
823   yi[atCD1]=yi[atCG]+b_CG_CD1*cos(alpha);
824
825   xi[atCE3]=b_CD2_CE3*sin(a_CE2_CD2_CE3);
826   yi[atCE3]=yi[atCD2]+b_CD2_CE3*cos(a_CE2_CD2_CE3);
827
828   xi[atCZ2]=b_CE2_CZ2*sin(a_CD2_CE2_CZ2);
829   yi[atCZ2]=yi[atCE2]-b_CE2_CZ2*cos(a_CD2_CE2_CZ2);
830
831   alpha = a_CE2_CD2_CE3 + a_CD2_CE3_CZ3 - M_PI;
832   xi[atCZ3]=xi[atCE3]+b_CE3_CZ3*sin(alpha);
833   yi[atCZ3]=yi[atCE3]+b_CE3_CZ3*cos(alpha);
834
835   alpha = a_CD2_CE2_CZ2 + a_CE2_CZ2_CH2 - M_PI;
836   xi[atCH2]=xi[atCZ2]+b_CZ2_CH2*sin(alpha);
837   yi[atCH2]=yi[atCZ2]-b_CZ2_CH2*cos(alpha);
838
839   /* hydrogens */
840   alpha = a_CE2_CD2_CG + a_CD2_CG_CD1 - a_CG_CD1_HD1;
841   xi[atHD1]=xi[atCD1]-b_CD1_HD1*sin(alpha);
842   yi[atHD1]=yi[atCD1]+b_CD1_HD1*cos(alpha);
843
844   alpha = a_NE1_CE2_CD2 + M_PI - a_HE1_NE1_CE2;
845   xi[atHE1]=xi[atNE1]-b_NE1_HE1*sin(alpha);
846   yi[atHE1]=yi[atNE1]-b_NE1_HE1*cos(alpha);
847   
848   alpha = a_CE2_CD2_CE3 + M_PI - a_CD2_CE3_HE3;
849   xi[atHE3]=xi[atCE3]+b_CE3_HE3*sin(alpha);
850   yi[atHE3]=yi[atCE3]+b_CE3_HE3*cos(alpha);
851
852   alpha = a_CD2_CE2_CZ2 + M_PI - a_CE2_CZ2_HZ2;
853   xi[atHZ2]=xi[atCZ2]+b_CZ2_HZ2*sin(alpha);
854   yi[atHZ2]=yi[atCZ2]-b_CZ2_HZ2*cos(alpha);
855
856   alpha = a_CD2_CE2_CZ2 + a_CE2_CZ2_CH2 - a_CZ2_CH2_HH2;
857   xi[atHZ3]=xi[atCZ3]+b_CZ3_HZ3*sin(alpha);
858   yi[atHZ3]=yi[atCZ3]+b_CZ3_HZ3*cos(alpha);
859
860   alpha = a_CE2_CD2_CE3 + a_CD2_CE3_CZ3 - a_CE3_CZ3_HZ3;
861   xi[atHH2]=xi[atCH2]+b_CH2_HH2*sin(alpha);
862   yi[atHH2]=yi[atCH2]-b_CH2_HH2*cos(alpha);
863
864   /* Determine coeff. for the line CB-CG */
865   lineA=(yi[atCB]-yi[atCG])/(xi[atCB]-xi[atCG]);
866   lineB=yi[atCG]-lineA*xi[atCG];
867   
868   /* Calculate masses for each ring and put it on the dummy masses */
869   for(j=0; j<NMASS; j++)
870     mM[j]=xcom[j]=ycom[j]=0;
871   for(i=0; i<atNR; i++) {
872     if (i!=atCB) {
873       for(j=0; j<NMASS; j++) {
874          mM[j] += mw[j][i] * at->atom[ats[i]].m;
875          xcom[j] += xi[i] * mw[j][i] * at->atom[ats[i]].m;
876          ycom[j] += yi[i] * mw[j][i] * at->atom[ats[i]].m;
877       }  
878     }
879   }
880   for(j=0; j<NMASS; j++) {
881     xcom[j]/=mM[j];
882     ycom[j]/=mM[j];
883   }
884   
885   /* get dummy mass type */
886   tpM=vsite_nm2type("MW",atype);
887   /* make space for 2 masses: shift all atoms starting with CB */
888   i0=ats[atCB];
889   for(j=0; j<NMASS; j++)
890     atM[j]=i0+*nadd+j;
891   if (debug)
892      fprintf(stderr,"Inserting %d dummy masses at %d\n",NMASS,(*o2n)[i0]+1);
893   *nadd+=NMASS;
894   for(j=i0; j<at->nr; j++)
895     (*o2n)[j]=j+*nadd;
896   srenew(*newx,at->nr+*nadd);
897   srenew(*newatom,at->nr+*nadd);
898   srenew(*newatomname,at->nr+*nadd);
899   srenew(*newvsite_type,at->nr+*nadd);
900   srenew(*newcgnr,at->nr+*nadd);
901   for(j=0; j<NMASS; j++)
902     (*newatomname)[at->nr+*nadd-1-j]=NULL;
903   
904   /* Dummy masses will be placed at the center-of-mass in each ring. */
905
906   /* calc initial position for dummy masses in real (non-local) coordinates.
907    * Cheat by using the routine to calculate virtual site parameters. It is
908    * much easier when we have the coordinates expressed in terms of
909    * CB, CG, CD2.
910    */
911   rvec_sub(x[ats[atCB]],x[ats[atCG]],r_ij);
912   rvec_sub(x[ats[atCD2]],x[ats[atCG]],r_ik);
913   calc_vsite3_param(xcom[0],ycom[0],xi[atCG],yi[atCG],xi[atCB],yi[atCB],
914                     xi[atCD2],yi[atCD2],&a,&b);
915   svmul(a,r_ij,t1);
916   svmul(b,r_ik,t2);
917   rvec_add(t1,t2,t1);
918   rvec_add(t1,x[ats[atCG]],(*newx)[atM[0]]);
919   
920   calc_vsite3_param(xcom[1],ycom[1],xi[atCG],yi[atCG],xi[atCB],yi[atCB],
921                     xi[atCD2],yi[atCD2],&a,&b);
922   svmul(a,r_ij,t1);
923   svmul(b,r_ik,t2);
924   rvec_add(t1,t2,t1);
925   rvec_add(t1,x[ats[atCG]],(*newx)[atM[1]]);
926
927   /* set parameters for the masses */
928   for(j=0; j<NMASS; j++) {
929     sprintf(name,"MW%d",j+1);
930     (*newatomname)  [atM[j]]       = put_symtab(symtab,name);
931     (*newatom)      [atM[j]].m     = (*newatom)[atM[j]].mB    = mM[j];
932     (*newatom)      [atM[j]].q     = (*newatom)[atM[j]].qB    = 0.0;
933     (*newatom)      [atM[j]].type  = (*newatom)[atM[j]].typeB = tpM;
934     (*newatom)      [atM[j]].ptype = eptAtom;
935     (*newatom)      [atM[j]].resind= at->atom[i0].resind;
936     (*newvsite_type)[atM[j]]       = NOTSET;
937     (*newcgnr)      [atM[j]]       = (*cgnr)[i0];
938   }
939   /* renumber cgnr: */
940   for(i=i0; i<at->nr; i++)
941     (*cgnr)[i]++;
942
943   /* constraints between CB, M1 and M2 */
944   /* 'add_shift' says which atoms won't be renumbered afterwards */
945   dCBM1 = sqrt( sqr(xcom[0]-xi[atCB]) + sqr(ycom[0]-yi[atCB]) );
946   dM1M2 = sqrt( sqr(xcom[0]-xcom[1]) + sqr(ycom[0]-ycom[1]) );
947   dCBM2 = sqrt( sqr(xcom[1]-xi[atCB]) + sqr(ycom[1]-yi[atCB]) );
948   my_add_param(&(plist[F_CONSTRNC]),ats[atCB],       add_shift+atM[0],dCBM1);
949   my_add_param(&(plist[F_CONSTRNC]),ats[atCB],       add_shift+atM[1],dCBM2);
950   my_add_param(&(plist[F_CONSTRNC]),add_shift+atM[0],add_shift+atM[1],dM1M2);
951   
952   /* rest will be vsite3 */
953   nvsite=0;
954   for(i=0; i<atNR; i++)
955     if (i!=atCB) {
956       at->atom[ats[i]].m = at->atom[ats[i]].mB = 0;
957       (*vsite_type)[ats[i]] = F_VSITE3;
958       nvsite++;
959     }
960   
961   /* now define all vsites from M1, M2, CB, ie:
962      r_d = r_M1 + a r_M1_M2 + b r_M1_CB */
963   for(i=0; i<atNR; i++)
964     if ( (*vsite_type)[ats[i]] == F_VSITE3 ) {
965       calc_vsite3_param(xi[i],yi[i],xcom[0],ycom[0],xcom[1],ycom[1],xi[atCB],yi[atCB],&a,&b);
966       add_vsite3_param(&plist[F_VSITE3],
967                      ats[i],add_shift+atM[0],add_shift+atM[1],ats[atCB],a,b);
968     }
969   return nvsite;
970 #undef NMASS
971 }
972
973
974 static int gen_vsites_tyr(gpp_atomtype_t atype, rvec *newx[],
975                           t_atom *newatom[], char ***newatomname[], 
976                         int *o2n[], int *newvsite_type[], int *newcgnr[],
977                           t_symtab *symtab, int *nadd, rvec x[], int *cgnr[],
978                           t_atoms *at, int *vsite_type[], t_params plist[], 
979                           int nrfound, int *ats, int add_shift,
980                           t_vsitetop *vsitetop, int nvsitetop)
981 {
982   int nvsite,i,i0,j,atM,tpM;
983   real dCGCE,dCEOH,dCGM,tmp1,a,b;
984   real bond_cc,bond_ch,bond_co,bond_oh,angle_coh;
985   real xcom,ycom,mtot;
986   real vmass,vdist,mM;
987   rvec r1;
988   char name[10];
989
990   /* these MUST correspond to the atnms array in do_vsite_aromatics! */
991   enum { atCG, atCD1, atHD1, atCD2, atHD2, atCE1, atHE1, atCE2, atHE2, 
992          atCZ, atOH, atHH, atNR };
993   real xi[atNR],yi[atNR];
994   /* CG, CE1, CE2 (as in general 6-ring) and OH and HH stay, 
995      rest gets virtualized.
996      Now we have two linked triangles with one improper keeping them flat */
997   if (atNR != nrfound)
998     gmx_incons("Number of atom types in gen_vsites_tyr");
999
1000   /* Aromatic rings have 6-fold symmetry, so we only need one bond length
1001    * for the ring part (angle is always 120 degrees).
1002    */
1003   bond_cc=get_ddb_bond(vsitetop,nvsitetop,"TYR","CD1","CE1");
1004   bond_ch=get_ddb_bond(vsitetop,nvsitetop,"TYR","CD1","HD1");
1005   bond_co=get_ddb_bond(vsitetop,nvsitetop,"TYR","CZ","OH");
1006   bond_oh=get_ddb_bond(vsitetop,nvsitetop,"TYR","OH","HH");
1007   angle_coh=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,"TYR","CZ","OH","HH");
1008   
1009   xi[atCG]=-bond_cc+bond_cc*cos(ANGLE_6RING);
1010   yi[atCG]=0;
1011   xi[atCD1]=-bond_cc;
1012   yi[atCD1]=bond_cc*sin(0.5*ANGLE_6RING);
1013   xi[atHD1]=xi[atCD1]+bond_ch*cos(ANGLE_6RING);
1014   yi[atHD1]=yi[atCD1]+bond_ch*sin(ANGLE_6RING);
1015   xi[atCE1]=0;
1016   yi[atCE1]=yi[atCD1];
1017   xi[atHE1]=xi[atCE1]-bond_ch*cos(ANGLE_6RING);
1018   yi[atHE1]=yi[atCE1]+bond_ch*sin(ANGLE_6RING);
1019   xi[atCD2]=xi[atCD1];
1020   yi[atCD2]=-yi[atCD1];
1021   xi[atHD2]=xi[atHD1];
1022   yi[atHD2]=-yi[atHD1];
1023   xi[atCE2]=xi[atCE1];
1024   yi[atCE2]=-yi[atCE1];
1025   xi[atHE2]=xi[atHE1];
1026   yi[atHE2]=-yi[atHE1];
1027   xi[atCZ]=bond_cc*cos(0.5*ANGLE_6RING);
1028   yi[atCZ]=0;
1029   xi[atOH]=xi[atCZ]+bond_co;
1030   yi[atOH]=0;
1031
1032   xcom=ycom=mtot=0;
1033   for(i=0;i<atOH;i++) {
1034     xcom+=xi[i]*at->atom[ats[i]].m;
1035     ycom+=yi[i]*at->atom[ats[i]].m;
1036     mtot+=at->atom[ats[i]].m;
1037   }
1038   xcom/=mtot;
1039   ycom/=mtot;
1040
1041   /* first do 6 ring as default, 
1042      except CZ (we'll do that different) and HZ (we don't have that): */
1043   nvsite = gen_vsites_6ring(at, vsite_type, plist, nrfound, ats, bond_cc, bond_ch, xcom, ycom, FALSE);
1044   
1045   /* then construct CZ from the 2nd triangle */
1046   /* vsite3 construction: r_d = r_i + a r_ij + b r_ik */
1047   a = b = 0.5 * bond_co / ( bond_co - bond_cc*cos(ANGLE_6RING) );
1048   add_vsite3_param(&plist[F_VSITE3],
1049                  ats[atCZ],ats[atOH],ats[atCE1],ats[atCE2],a,b);
1050   at->atom[ats[atCZ]].m = at->atom[ats[atCZ]].mB = 0;
1051   
1052   /* constraints between CE1, CE2 and OH */
1053   dCGCE = sqrt( cosrule(bond_cc,bond_cc,ANGLE_6RING) );
1054   dCEOH = sqrt( cosrule(bond_cc,bond_co,ANGLE_6RING) );
1055   my_add_param(&(plist[F_CONSTRNC]),ats[atCE1],ats[atOH],dCEOH);
1056   my_add_param(&(plist[F_CONSTRNC]),ats[atCE2],ats[atOH],dCEOH);
1057
1058   /* We also want to constrain the angle C-O-H, but since CZ is constructed
1059    * we need to introduce a constraint to CG.
1060    * CG is much further away, so that will lead to instabilities in LINCS
1061    * when we constrain both CG-HH and OH-HH distances. Instead of requiring
1062    * the use of lincs_order=8 we introduce a dummy mass three times further
1063    * away from OH than HH. The mass is accordingly a third, with the remaining
1064    * 2/3 moved to OH. This shouldnt cause any problems since the forces will
1065    * apply to the HH constructed atom and not directly on the virtual mass.
1066    */
1067
1068   vdist=2.0*bond_oh;
1069   mM=at->atom[ats[atHH]].m/2.0;
1070   at->atom[ats[atOH]].m+=mM;  /* add 1/2 of original H mass */
1071   at->atom[ats[atOH]].mB+=mM; /* add 1/2 of original H mass */
1072   at->atom[ats[atHH]].m=at->atom[ats[atHH]].mB=0;
1073   
1074   /* get dummy mass type */
1075   tpM=vsite_nm2type("MW",atype);
1076   /* make space for 1 mass: shift HH only */
1077   i0=ats[atHH];
1078   atM=i0+*nadd;
1079   if (debug)
1080      fprintf(stderr,"Inserting 1 dummy mass at %d\n",(*o2n)[i0]+1);
1081   (*nadd)++;
1082   for(j=i0; j<at->nr; j++)
1083     (*o2n)[j]=j+*nadd;
1084   srenew(*newx,at->nr+*nadd);
1085   srenew(*newatom,at->nr+*nadd);
1086   srenew(*newatomname,at->nr+*nadd);
1087   srenew(*newvsite_type,at->nr+*nadd);
1088   srenew(*newcgnr,at->nr+*nadd);
1089   (*newatomname)[at->nr+*nadd-1]=NULL; 
1090   
1091   /* Calc the dummy mass initial position */
1092   rvec_sub(x[ats[atHH]],x[ats[atOH]],r1);
1093   svmul(2.0,r1,r1);
1094   rvec_add(r1,x[ats[atHH]],(*newx)[atM]);
1095   
1096   strcpy(name,"MW1");
1097   (*newatomname)  [atM]       = put_symtab(symtab,name);
1098   (*newatom)      [atM].m     = (*newatom)[atM].mB    = mM;
1099   (*newatom)      [atM].q     = (*newatom)[atM].qB    = 0.0;
1100   (*newatom)      [atM].type  = (*newatom)[atM].typeB = tpM;
1101   (*newatom)      [atM].ptype = eptAtom;
1102   (*newatom)      [atM].resind= at->atom[i0].resind;
1103   (*newvsite_type)[atM]       = NOTSET;
1104   (*newcgnr)      [atM]       = (*cgnr)[i0]; 
1105   /* renumber cgnr: */
1106   for(i=i0; i<at->nr; i++)
1107     (*cgnr)[i]++;
1108
1109   (*vsite_type)[ats[atHH]] = F_VSITE2;
1110   nvsite++;
1111   /* assume we also want the COH angle constrained: */
1112   tmp1 = bond_cc*cos(0.5*ANGLE_6RING) + dCGCE*sin(ANGLE_6RING*0.5) + bond_co;
1113   dCGM = sqrt( cosrule(tmp1,vdist,angle_coh) );
1114   my_add_param(&(plist[F_CONSTRNC]),ats[atCG],add_shift+atM,dCGM);
1115   my_add_param(&(plist[F_CONSTRNC]),ats[atOH],add_shift+atM,vdist);
1116
1117   add_vsite2_param(&plist[F_VSITE2],
1118                  ats[atHH],ats[atOH],add_shift+atM,1.0/2.0);
1119   return nvsite;
1120 }
1121       
1122 static int gen_vsites_his(t_atoms *at, int *vsite_type[], t_params plist[], 
1123                         int nrfound, int *ats, t_vsitetop *vsitetop, int nvsitetop)
1124 {
1125   int nvsite,i;
1126   real a,b,alpha,dCGCE1,dCGNE2;
1127   real sinalpha,cosalpha;
1128   real xcom,ycom,mtot;
1129   real mG,mrest,mCE1,mNE2;
1130   real b_CG_ND1,b_ND1_CE1,b_CE1_NE2,b_CG_CD2,b_CD2_NE2;
1131   real b_ND1_HD1,b_NE2_HE2,b_CE1_HE1,b_CD2_HD2;
1132   real a_CG_ND1_CE1,a_CG_CD2_NE2,a_ND1_CE1_NE2,a_CE1_NE2_CD2;
1133   real a_NE2_CE1_HE1,a_NE2_CD2_HD2,a_CE1_ND1_HD1,a_CE1_NE2_HE2;
1134   char resname[10];
1135
1136   /* these MUST correspond to the atnms array in do_vsite_aromatics! */
1137   enum { atCG, atND1, atHD1, atCD2, atHD2, atCE1, atHE1, atNE2, atHE2, atNR };
1138   real x[atNR],y[atNR];
1139
1140   /* CG, CE1 and NE2 stay, each gets part of the total mass,
1141      rest gets virtualized */
1142   /* check number of atoms, 3 hydrogens may be missing: */
1143   /* assert( nrfound >= atNR-3 || nrfound <= atNR );
1144    * Don't understand the above logic. Shouldn't it be && rather than || ???
1145    */
1146   if ((nrfound < atNR-3) || (nrfound > atNR))
1147     gmx_incons("Generating vsites for HIS");
1148   
1149   /* avoid warnings about uninitialized variables */
1150   b_ND1_HD1=b_NE2_HE2=b_CE1_HE1=b_CD2_HD2=a_NE2_CE1_HE1=
1151     a_NE2_CD2_HD2=a_CE1_ND1_HD1=a_CE1_NE2_HE2=0;
1152
1153   if(ats[atHD1]!=NOTSET) {
1154     if(ats[atHE2]!=NOTSET) 
1155       sprintf(resname,"HISH");
1156     else 
1157       sprintf(resname,"HISA");
1158   } else {
1159     sprintf(resname,"HISB");
1160   }
1161   
1162   /* Get geometry from database */
1163   b_CG_ND1=get_ddb_bond(vsitetop,nvsitetop,resname,"CG","ND1");
1164   b_ND1_CE1=get_ddb_bond(vsitetop,nvsitetop,resname,"ND1","CE1");
1165   b_CE1_NE2=get_ddb_bond(vsitetop,nvsitetop,resname,"CE1","NE2");
1166   b_CG_CD2 =get_ddb_bond(vsitetop,nvsitetop,resname,"CG","CD2");
1167   b_CD2_NE2=get_ddb_bond(vsitetop,nvsitetop,resname,"CD2","NE2");
1168   a_CG_ND1_CE1=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,resname,"CG","ND1","CE1");
1169   a_CG_CD2_NE2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,resname,"CG","CD2","NE2");
1170   a_ND1_CE1_NE2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,resname,"ND1","CE1","NE2");
1171   a_CE1_NE2_CD2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,resname,"CE1","NE2","CD2");
1172
1173   if(ats[atHD1]!=NOTSET) {
1174     b_ND1_HD1=get_ddb_bond(vsitetop,nvsitetop,resname,"ND1","HD1");
1175     a_CE1_ND1_HD1=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,resname,"CE1","ND1","HD1");
1176   }
1177   if(ats[atHE2]!=NOTSET) {
1178     b_NE2_HE2=get_ddb_bond(vsitetop,nvsitetop,resname,"NE2","HE2");
1179     a_CE1_NE2_HE2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,resname,"CE1","NE2","HE2");
1180   }
1181   if(ats[atHD2]!=NOTSET) {
1182     b_CD2_HD2=get_ddb_bond(vsitetop,nvsitetop,resname,"CD2","HD2");
1183     a_NE2_CD2_HD2=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,resname,"NE2","CD2","HD2");
1184   }
1185   if(ats[atHE1]!=NOTSET) {
1186     b_CE1_HE1=get_ddb_bond(vsitetop,nvsitetop,resname,"CE1","HE1");
1187     a_NE2_CE1_HE1=DEG2RAD*get_ddb_angle(vsitetop,nvsitetop,resname,"NE2","CE1","HE1");
1188   }
1189
1190   /* constraints between CG, CE1 and NE1 */
1191   dCGCE1   = sqrt( cosrule(b_CG_ND1,b_ND1_CE1,a_CG_ND1_CE1) );
1192   dCGNE2   = sqrt( cosrule(b_CG_CD2,b_CD2_NE2,a_CG_CD2_NE2) );
1193
1194   my_add_param(&(plist[F_CONSTRNC]),ats[atCG] ,ats[atCE1],dCGCE1);
1195   my_add_param(&(plist[F_CONSTRNC]),ats[atCG] ,ats[atNE2],dCGNE2);
1196   /* we already have a constraint CE1-NE2, so we don't add it again */
1197   
1198   /* calculate the positions in a local frame of reference. 
1199    * The x-axis is the line from CG that makes a right angle
1200    * with the bond CE1-NE2, and the y-axis the bond CE1-NE2.
1201    */
1202   /* First calculate the x-axis intersection with y-axis (=yCE1).
1203    * Get cos(angle CG-CE1-NE2) :
1204    */
1205   cosalpha=acosrule(dCGNE2,dCGCE1,b_CE1_NE2);
1206   x[atCE1] = 0;
1207   y[atCE1] = cosalpha*dCGCE1;
1208   x[atNE2] = 0;
1209   y[atNE2] = y[atCE1]-b_CE1_NE2;
1210   sinalpha=sqrt(1-cosalpha*cosalpha);
1211   x[atCG]  = - sinalpha*dCGCE1;
1212   y[atCG]  = 0;
1213   x[atHE1] = x[atHE2] = x[atHD1] = x[atHD2] = 0;
1214   y[atHE1] = y[atHE2] = y[atHD1] = y[atHD2] = 0;
1215   
1216   /* calculate ND1 and CD2 positions from CE1 and NE2 */
1217
1218   x[atND1] = -b_ND1_CE1*sin(a_ND1_CE1_NE2);
1219   y[atND1] = y[atCE1]-b_ND1_CE1*cos(a_ND1_CE1_NE2);
1220
1221   x[atCD2] = -b_CD2_NE2*sin(a_CE1_NE2_CD2);
1222   y[atCD2] = y[atNE2]+b_CD2_NE2*cos(a_CE1_NE2_CD2);
1223
1224   /* And finally the hydrogen positions */
1225   if(ats[atHE1]!=NOTSET) {  
1226     x[atHE1] = x[atCE1] + b_CE1_HE1*sin(a_NE2_CE1_HE1);
1227     y[atHE1] = y[atCE1] - b_CE1_HE1*cos(a_NE2_CE1_HE1);
1228   }
1229   /* HD2 - first get (ccw) angle from (positive) y-axis */
1230   if(ats[atHD2]!=NOTSET) {  
1231     alpha = a_CE1_NE2_CD2 + M_PI - a_NE2_CD2_HD2;
1232     x[atHD2] = x[atCD2] - b_CD2_HD2*sin(alpha);
1233     y[atHD2] = y[atCD2] + b_CD2_HD2*cos(alpha);
1234   }
1235   if(ats[atHD1]!=NOTSET) {
1236     /* HD1 - first get (cw) angle from (positive) y-axis */
1237     alpha = a_ND1_CE1_NE2 + M_PI - a_CE1_ND1_HD1;
1238     x[atHD1] = x[atND1] - b_ND1_HD1*sin(alpha);
1239     y[atHD1] = y[atND1] - b_ND1_HD1*cos(alpha);
1240   }
1241   if(ats[atHE2]!=NOTSET) {
1242     x[atHE2] = x[atNE2] + b_NE2_HE2*sin(a_CE1_NE2_HE2);
1243     y[atHE2] = y[atNE2] + b_NE2_HE2*cos(a_CE1_NE2_HE2);
1244   }
1245   /* Have all coordinates now */
1246   
1247   /* calc center-of-mass; keep atoms CG, CE1, NE2 and
1248    * set the rest to vsite3
1249    */
1250   mtot=xcom=ycom=0;
1251   nvsite=0;
1252   for(i=0; i<atNR; i++) 
1253     if (ats[i]!=NOTSET) {
1254       mtot+=at->atom[ats[i]].m;
1255       xcom+=x[i]*at->atom[ats[i]].m;
1256       ycom+=y[i]*at->atom[ats[i]].m;
1257       if (i!=atCG && i!=atCE1 && i!=atNE2) {
1258         at->atom[ats[i]].m = at->atom[ats[i]].mB = 0;
1259         (*vsite_type)[ats[i]]=F_VSITE3;
1260         nvsite++;
1261       }
1262     } 
1263   if (nvsite+3 != nrfound)
1264     gmx_incons("Generating vsites for HIS");
1265
1266   xcom/=mtot;
1267   ycom/=mtot;
1268   
1269   /* distribute mass so that com stays the same */
1270   mG = xcom*mtot/x[atCG];
1271   mrest = mtot-mG;
1272   mCE1 = (ycom-y[atNE2])*mrest/(y[atCE1]-y[atNE2]);
1273   mNE2 = mrest-mCE1;
1274   
1275   at->atom[ats[atCG]].m = at->atom[ats[atCG]].mB = mG;
1276   at->atom[ats[atCE1]].m = at->atom[ats[atCE1]].mB = mCE1;
1277   at->atom[ats[atNE2]].m = at->atom[ats[atNE2]].mB = mNE2;
1278   
1279   /* HE1 */
1280   if (ats[atHE1]!=NOTSET) {
1281     calc_vsite3_param(x[atHE1],y[atHE1],x[atCE1],y[atCE1],x[atNE2],y[atNE2],
1282                       x[atCG],y[atCG],&a,&b);
1283     add_vsite3_param(&plist[F_VSITE3],
1284                    ats[atHE1],ats[atCE1],ats[atNE2],ats[atCG],a,b);
1285   }
1286   /* HE2 */
1287   if (ats[atHE2]!=NOTSET) {
1288     calc_vsite3_param(x[atHE2],y[atHE2],x[atNE2],y[atNE2],x[atCE1],y[atCE1],
1289                       x[atCG],y[atCG],&a,&b);
1290     add_vsite3_param(&plist[F_VSITE3],
1291                    ats[atHE2],ats[atNE2],ats[atCE1],ats[atCG],a,b);
1292   }
1293   
1294   /* ND1 */
1295   calc_vsite3_param(x[atND1],y[atND1],x[atNE2],y[atNE2],x[atCE1],y[atCE1],
1296                     x[atCG],y[atCG],&a,&b);
1297   add_vsite3_param(&plist[F_VSITE3],
1298                  ats[atND1],ats[atNE2],ats[atCE1],ats[atCG],a,b);
1299   
1300   /* CD2 */
1301   calc_vsite3_param(x[atCD2],y[atCD2],x[atCE1],y[atCE1],x[atNE2],y[atNE2],
1302                     x[atCG],y[atCG],&a,&b);
1303   add_vsite3_param(&plist[F_VSITE3],
1304                  ats[atCD2],ats[atCE1],ats[atNE2],ats[atCG],a,b);
1305   
1306   /* HD1 */
1307   if (ats[atHD1]!=NOTSET) {
1308     calc_vsite3_param(x[atHD1],y[atHD1],x[atNE2],y[atNE2],x[atCE1],y[atCE1],
1309                       x[atCG],y[atCG],&a,&b);
1310     add_vsite3_param(&plist[F_VSITE3],
1311                    ats[atHD1],ats[atNE2],ats[atCE1],ats[atCG],a,b);
1312   }
1313   /* HD2 */
1314   if (ats[atHD2]!=NOTSET) {
1315     calc_vsite3_param(x[atHD2],y[atHD2],x[atCE1],y[atCE1],x[atNE2],y[atNE2],
1316                       x[atCG],y[atCG],&a,&b);
1317     add_vsite3_param(&plist[F_VSITE3],
1318                    ats[atHD2],ats[atCE1],ats[atNE2],ats[atCG],a,b);
1319   }  
1320   return nvsite;
1321 }
1322
1323 static gmx_bool is_vsite(int vsite_type)
1324 {
1325   if (vsite_type == NOTSET)
1326     return FALSE;
1327   switch ( abs(vsite_type) ) {
1328   case F_VSITE3:
1329   case F_VSITE3FD:
1330   case F_VSITE3OUT:
1331   case F_VSITE3FAD:
1332   case F_VSITE4FD:
1333   case F_VSITE4FDN:
1334     return TRUE;
1335   default:
1336     return FALSE;
1337   }
1338 }
1339
1340 static char atomnamesuffix[] = "1234";
1341
1342 void do_vsites(int nrtp, t_restp rtp[], gpp_atomtype_t atype, 
1343                t_atoms *at, t_symtab *symtab, rvec *x[], 
1344                t_params plist[], int *vsite_type[], int *cgnr[], 
1345                real mHmult, gmx_bool bVsiteAromatics,
1346                const char *ffdir)
1347 {
1348 #define MAXATOMSPERRESIDUE 16
1349   int  i,j,k,m,i0,ni0,whatres,resind,add_shift,ftype,nvsite,nadd;
1350   int  ai,aj,ak,al;
1351   int  nrfound=0,needed,nrbonds,nrHatoms,Heavy,nrheavies,tpM,tpHeavy;
1352   int  Hatoms[4],heavies[4],bb;
1353   gmx_bool bWARNING,bAddVsiteParam,bFirstWater;
1354   matrix tmpmat;
1355   gmx_bool *bResProcessed;
1356   real mHtot,mtot,fact,fact2;
1357   rvec rpar,rperp,temp;
1358   char name[10],tpname[32],nexttpname[32],*ch;
1359   rvec *newx;
1360   int  *o2n,*newvsite_type,*newcgnr,ats[MAXATOMSPERRESIDUE];
1361   t_atom *newatom;
1362   t_params *params;
1363   char ***newatomname;
1364   char *resnm=NULL;
1365   int  ndb,f;
1366   char **db;
1367   int  nvsiteconf,nvsitetop,cmplength;
1368   gmx_bool isN,planarN,bFound;
1369   gmx_residuetype_t rt;
1370     
1371   t_vsiteconf *vsiteconflist; 
1372   /* pointer to a list of CH3/NH3/NH2 configuration entries.
1373    * See comments in read_vsite_database. It isnt beautiful,
1374    * but it had to be fixed, and I dont even want to try to 
1375    * maintain this part of the code...
1376    */
1377   t_vsitetop *vsitetop;       
1378   /* Pointer to a list of geometry (bond/angle) entries for
1379    * residues like PHE, TRP, TYR, HIS, etc., where we need
1380    * to know the geometry to construct vsite aromatics.
1381    * Note that equilibrium geometry isnt necessarily the same
1382    * as the individual bond and angle values given in the
1383    * force field (rings can be strained).
1384    */
1385
1386   /* if bVsiteAromatics=TRUE do_vsites will specifically convert atoms in 
1387      PHE, TRP, TYR and HIS to a construction of virtual sites */
1388   enum                    { resPHE, resTRP, resTYR, resHIS, resNR };
1389   const char *resnms[resNR]   = {   "PHE",  "TRP",  "TYR",  "HIS" };
1390         /* Amber03 alternative names for termini */
1391   const char *resnmsN[resNR]  = {  "NPHE", "NTRP", "NTYR", "NHIS" };
1392   const char *resnmsC[resNR]  = {  "CPHE", "CTRP", "CTYR", "CHIS" };
1393   /* HIS can be known as HISH, HIS1, HISA, HID, HIE, HIP, etc. too */
1394   gmx_bool bPartial[resNR]  = {  FALSE,  FALSE,  FALSE,   TRUE  };
1395   /* the atnms for every residue MUST correspond to the enums in the 
1396      gen_vsites_* (one for each residue) routines! */
1397   /* also the atom names in atnms MUST be in the same order as in the .rtp! */
1398   const char *atnms[resNR][MAXATOMSPERRESIDUE+1] = { 
1399     { "CG", /* PHE */
1400       "CD1", "HD1", "CD2", "HD2", 
1401       "CE1", "HE1", "CE2", "HE2", 
1402       "CZ", "HZ", NULL },
1403     { "CB", /* TRP */
1404       "CG",
1405       "CD1", "HD1", "CD2", 
1406       "NE1", "HE1", "CE2", "CE3", "HE3", 
1407       "CZ2", "HZ2", "CZ3", "HZ3", 
1408       "CH2", "HH2", NULL },
1409     { "CG", /* TYR */
1410       "CD1", "HD1", "CD2", "HD2", 
1411       "CE1", "HE1", "CE2", "HE2", 
1412       "CZ", "OH", "HH", NULL },
1413     { "CG", /* HIS */
1414       "ND1", "HD1", "CD2", "HD2",
1415       "CE1", "HE1", "NE2", "HE2", NULL }
1416   };
1417     
1418   if (debug) {
1419     printf("Searching for atoms to make virtual sites ...\n");
1420     fprintf(debug,"# # # VSITES # # #\n");
1421   }
1422
1423   ndb = fflib_search_file_end(ffdir,".vsd",FALSE,&db);
1424   nvsiteconf    = 0;
1425   vsiteconflist = NULL;
1426   nvsitetop     = 0;
1427   vsitetop      = NULL;
1428   for(f=0; f<ndb; f++) {
1429     read_vsite_database(db[f],&vsiteconflist,&nvsiteconf,&vsitetop,&nvsitetop);
1430     sfree(db[f]);
1431   }
1432   sfree(db);
1433   
1434   bFirstWater=TRUE;
1435   nvsite=0;
1436   nadd=0;
1437   /* we need a marker for which atoms should *not* be renumbered afterwards */
1438   add_shift = 10*at->nr;
1439   /* make arrays where masses can be inserted into */ 
1440   snew(newx,at->nr); 
1441   snew(newatom,at->nr);
1442   snew(newatomname,at->nr);
1443   snew(newvsite_type,at->nr);
1444   snew(newcgnr,at->nr);
1445   /* make index array to tell where the atoms go to when masses are inserted */
1446   snew(o2n,at->nr);
1447   for(i=0; i<at->nr; i++)
1448     o2n[i]=i;
1449   /* make index to tell which residues were already processed */
1450   snew(bResProcessed,at->nres);
1451     
1452   gmx_residuetype_init(&rt);
1453     
1454   /* generate vsite constructions */
1455   /* loop over all atoms */
1456   resind = -1;
1457   for(i=0; (i<at->nr); i++) {
1458     if (at->atom[i].resind != resind) {
1459       resind = at->atom[i].resind;
1460       resnm  = *(at->resinfo[resind].name);
1461     }
1462     /* first check for aromatics to virtualize */
1463     /* don't waste our effort on DNA, water etc. */
1464     /* Only do the vsite aromatic stuff when we reach the 
1465      * CA atom, since there might be an X2/X3 group on the
1466      * N-terminus that must be treated first.
1467      */
1468     if ( bVsiteAromatics &&
1469          !strcmp(*(at->atomname[i]),"CA") &&
1470          !bResProcessed[resind] &&
1471          gmx_residuetype_is_protein(rt,*(at->resinfo[resind].name)) ) {
1472       /* mark this residue */
1473       bResProcessed[resind] = TRUE;
1474       /* find out if this residue needs converting */
1475       whatres=NOTSET;
1476       for(j=0; j<resNR && whatres==NOTSET; j++) {
1477                   
1478         cmplength = bPartial[j] ? strlen(resnm)-1 : strlen(resnm);
1479         
1480         bFound = ((gmx_strncasecmp(resnm,resnms[j], cmplength)==0) ||
1481                   (gmx_strncasecmp(resnm,resnmsN[j],cmplength)==0) ||
1482                   (gmx_strncasecmp(resnm,resnmsC[j],cmplength)==0));
1483         
1484         if ( bFound ) {
1485           whatres=j;
1486           /* get atoms we will be needing for the conversion */
1487           nrfound=0;
1488           for (k=0; atnms[j][k]; k++) 
1489       {
1490             ats[k]=NOTSET;
1491           for(m=i; m<at->nr && at->atom[m].resind==resind && ats[k]==NOTSET;m++) 
1492         {
1493               if (gmx_strcasecmp(*(at->atomname[m]),atnms[j][k])==0) 
1494           {
1495               ats[k]=m;
1496               nrfound++;
1497               }
1498             }
1499           }
1500         
1501       /* now k is number of atom names in atnms[j] */
1502       if (j==resHIS)
1503       {
1504           needed = k-3;
1505       }
1506           else
1507           {
1508           needed = k;
1509       }
1510       if (nrfound<needed)
1511           {
1512           gmx_fatal(FARGS,"not enough atoms found (%d, need %d) in "
1513                     "residue %s %d while\n             "
1514                     "generating aromatics virtual site construction",
1515                     nrfound,needed,resnm,at->resinfo[resind].nr);
1516       }
1517       /* Advance overall atom counter */
1518       i++;
1519     }
1520       }
1521       /* the enums for every residue MUST correspond to atnms[residue] */
1522       switch (whatres) {
1523       case resPHE: 
1524         if (debug) fprintf(stderr,"PHE at %d\n",o2n[ats[0]]+1);
1525         nvsite+=gen_vsites_phe(at, vsite_type, plist, nrfound, ats, vsitetop, nvsitetop);
1526         break;
1527       case resTRP: 
1528         if (debug) fprintf(stderr,"TRP at %d\n",o2n[ats[0]]+1);
1529         nvsite+=gen_vsites_trp(atype, &newx, &newatom, &newatomname, &o2n, 
1530                            &newvsite_type, &newcgnr, symtab, &nadd, *x, cgnr,
1531                            at, vsite_type, plist, nrfound, ats, add_shift, vsitetop, nvsitetop);
1532         break;
1533       case resTYR: 
1534         if (debug) fprintf(stderr,"TYR at %d\n",o2n[ats[0]]+1);
1535         nvsite+=gen_vsites_tyr(atype, &newx, &newatom, &newatomname, &o2n, 
1536                            &newvsite_type, &newcgnr, symtab, &nadd, *x, cgnr,
1537                            at, vsite_type, plist, nrfound, ats, add_shift, vsitetop, nvsitetop);
1538         break;
1539       case resHIS: 
1540         if (debug) fprintf(stderr,"HIS at %d\n",o2n[ats[0]]+1);
1541         nvsite+=gen_vsites_his(at, vsite_type, plist, nrfound, ats, vsitetop, nvsitetop);
1542         break;
1543       case NOTSET:
1544         /* this means this residue won't be processed */
1545         break;
1546       default:
1547         gmx_fatal(FARGS,"DEATH HORROR in do_vsites (%s:%d)",
1548                     __FILE__,__LINE__);
1549       } /* switch whatres */
1550       /* skip back to beginning of residue */
1551       while(i>0 && at->atom[i-1].resind == resind)
1552         i--;
1553     } /* if bVsiteAromatics & is protein */
1554     
1555     /* now process the rest of the hydrogens */
1556     /* only process hydrogen atoms which are not already set */
1557     if ( ((*vsite_type)[i]==NOTSET) && is_hydrogen(*(at->atomname[i]))) {
1558       /* find heavy atom, count #bonds from it and #H atoms bound to it
1559          and return H atom numbers (Hatoms) and heavy atom numbers (heavies) */
1560       count_bonds(i, &plist[F_BONDS], at->atomname, 
1561                   &nrbonds, &nrHatoms, Hatoms, &Heavy, &nrheavies, heavies);
1562       /* get Heavy atom type */
1563       tpHeavy=get_atype(Heavy,at,nrtp,rtp,rt);
1564       strcpy(tpname,get_atomtype_name(tpHeavy,atype));
1565
1566       bWARNING=FALSE;
1567       bAddVsiteParam=TRUE;
1568       /* nested if's which check nrHatoms, nrbonds and atomname */
1569       if (nrHatoms == 1) {
1570         switch(nrbonds) {
1571         case 2: /* -O-H */
1572           (*vsite_type)[i]=F_BONDS;
1573           break;
1574         case 3: /* =CH-, -NH- or =NH+- */
1575           (*vsite_type)[i]=F_VSITE3FD;
1576           break;
1577         case 4: /* --CH- (tert) */
1578         /* The old type 4FD had stability issues, so 
1579          * all new constructs should use 4FDN
1580          */
1581           (*vsite_type)[i]=F_VSITE4FDN;
1582         
1583         /* Check parity of heavy atoms from coordinates */
1584         ai = Heavy;
1585         aj = heavies[0];
1586         ak = heavies[1];
1587         al = heavies[2];
1588         rvec_sub((*x)[aj],(*x)[ai],tmpmat[0]);
1589         rvec_sub((*x)[ak],(*x)[ai],tmpmat[1]);
1590         rvec_sub((*x)[al],(*x)[ai],tmpmat[2]);
1591         
1592         if(det(tmpmat)>0)
1593         {
1594             /* swap parity */
1595             heavies[1] = aj;
1596             heavies[0] = ak;
1597         }
1598         
1599           break;
1600         default: /* nrbonds != 2, 3 or 4 */
1601           bWARNING=TRUE;
1602         }
1603         
1604       } else if ( /*(nrHatoms == 2) && (nrbonds == 2) && REMOVED this test
1605                    DvdS 19-01-04 */
1606                   (gmx_strncasecmp(*at->atomname[Heavy],"OW",2)==0) ) {
1607         bAddVsiteParam=FALSE; /* this is water: skip these hydrogens */
1608         if (bFirstWater) {
1609           bFirstWater=FALSE;
1610           if (debug)
1611             fprintf(debug,
1612                     "Not converting hydrogens in water to virtual sites\n");
1613         }
1614       } else if ( (nrHatoms == 2) && (nrbonds == 4) ) {
1615         /* -CH2- , -NH2+- */
1616         (*vsite_type)[Hatoms[0]] = F_VSITE3OUT;
1617         (*vsite_type)[Hatoms[1]] =-F_VSITE3OUT;
1618       } else {
1619         /* 2 or 3 hydrogen atom, with 3 or 4 bonds in total to the heavy atom.
1620          * If it is a nitrogen, first check if it is planar.
1621          */
1622         isN=planarN=FALSE;
1623         if((nrHatoms == 2) && ((*at->atomname[Heavy])[0]=='N')) {
1624           isN=TRUE;
1625           j=nitrogen_is_planar(vsiteconflist,nvsiteconf,tpname);
1626           if(j<0) 
1627             gmx_fatal(FARGS,"No vsite database NH2 entry for type %s\n",tpname);
1628           planarN=(j==1);
1629         }
1630         if ( (nrHatoms == 2) && (nrbonds == 3) && ( !isN || planarN ) ) {
1631           /* =CH2 or, if it is a nitrogen NH2, it is a planar one */
1632           (*vsite_type)[Hatoms[0]] = F_VSITE3FAD;
1633           (*vsite_type)[Hatoms[1]] =-F_VSITE3FAD;
1634         } else if ( ( (nrHatoms == 2) && (nrbonds == 3) && 
1635                       ( isN && !planarN ) ) || 
1636                     ( (nrHatoms == 3) && (nrbonds == 4) ) ) {
1637           /* CH3, NH3 or non-planar NH2 group */
1638         int  Hat_vsite_type[3] = { F_VSITE3, F_VSITE3OUT, F_VSITE3OUT };
1639         gmx_bool Hat_SwapParity[3] = { FALSE,    TRUE,        FALSE };
1640         
1641         if (debug) fprintf(stderr,"-XH3 or nonplanar NH2 group at %d\n",i+1);
1642         bAddVsiteParam=FALSE; /* we'll do this ourselves! */
1643         /* -NH2 (umbrella), -NH3+ or -CH3 */
1644         (*vsite_type)[Heavy]       = F_VSITE3;
1645         for (j=0; j<nrHatoms; j++)
1646           (*vsite_type)[Hatoms[j]] = Hat_vsite_type[j];
1647         /* get dummy mass type from first char of heavy atom type (N or C) */
1648         
1649         strcpy(nexttpname,get_atomtype_name(get_atype(heavies[0],at,nrtp,rtp,rt),atype));
1650         ch=get_dummymass_name(vsiteconflist,nvsiteconf,tpname,nexttpname);
1651
1652         if (ch == NULL) {
1653           if (ndb > 0) {
1654             gmx_fatal(FARGS,"Can't find dummy mass for type %s bonded to type %s in the virtual site database (.vsd files). Add it to the database!\n",tpname,nexttpname);
1655           } else {
1656              gmx_fatal(FARGS,"A dummy mass for type %s bonded to type %s is required, but no virtual site database (.vsd) files where found.\n",tpname,nexttpname);
1657           }
1658         } else {
1659           strcpy(name,ch);
1660         }
1661
1662         tpM=vsite_nm2type(name,atype);
1663         /* make space for 2 masses: shift all atoms starting with 'Heavy' */
1664 #define NMASS 2
1665         i0=Heavy;
1666         ni0=i0+nadd;
1667         if (debug) 
1668           fprintf(stderr,"Inserting %d dummy masses at %d\n",NMASS,o2n[i0]+1);
1669         nadd+=NMASS;
1670         for(j=i0; j<at->nr; j++)
1671           o2n[j]=j+nadd; 
1672
1673         srenew(newx,at->nr+nadd);
1674         srenew(newatom,at->nr+nadd);
1675         srenew(newatomname,at->nr+nadd);
1676         srenew(newvsite_type,at->nr+nadd);
1677         srenew(newcgnr,at->nr+nadd);
1678
1679         for(j=0; j<NMASS; j++)
1680           newatomname[at->nr+nadd-1-j]=NULL;
1681
1682         /* calculate starting position for the masses */
1683         mHtot=0;
1684         /* get atom masses, and set Heavy and Hatoms mass to zero */
1685         for(j=0; j<nrHatoms; j++) {
1686           mHtot += get_amass(Hatoms[j],at,nrtp,rtp,rt);
1687           at->atom[Hatoms[j]].m = at->atom[Hatoms[j]].mB = 0;
1688         }
1689         mtot = mHtot + get_amass(Heavy,at,nrtp,rtp,rt);
1690         at->atom[Heavy].m = at->atom[Heavy].mB = 0;
1691         if (mHmult != 1.0)
1692           mHtot *= mHmult;
1693         fact2=mHtot/mtot;
1694         fact=sqrt(fact2);
1695         /* generate vectors parallel and perpendicular to rotational axis:
1696          * rpar  = Heavy -> Hcom
1697          * rperp = Hcom  -> H1   */
1698         clear_rvec(rpar);
1699         for(j=0; j<nrHatoms; j++)
1700           rvec_inc(rpar,(*x)[Hatoms[j]]);
1701         svmul(1.0/nrHatoms,rpar,rpar); /* rpar = ( H1+H2+H3 ) / 3 */
1702         rvec_dec(rpar,(*x)[Heavy]);    /*        - Heavy          */
1703         rvec_sub((*x)[Hatoms[0]],(*x)[Heavy],rperp);
1704         rvec_dec(rperp,rpar);          /* rperp = H1 - Heavy - rpar */
1705         /* calc mass positions */
1706         svmul(fact2,rpar,temp);
1707         for (j=0; (j<NMASS); j++) /* xM = xN + fact2 * rpar +/- fact * rperp */
1708           rvec_add((*x)[Heavy],temp,newx[ni0+j]);
1709         svmul(fact,rperp,temp);
1710         rvec_inc(newx[ni0  ],temp);
1711         rvec_dec(newx[ni0+1],temp);
1712         /* set atom parameters for the masses */
1713         for(j=0; (j<NMASS); j++) {
1714           /* make name: "M??#" or "M?#" (? is atomname, # is number) */
1715           name[0]='M';
1716           for(k=0; (*at->atomname[Heavy])[k] && ( k < NMASS ); k++ )
1717             name[k+1]=(*at->atomname[Heavy])[k];
1718           name[k+1]=atomnamesuffix[j];
1719           name[k+2]='\0';
1720           newatomname[ni0+j]   = put_symtab(symtab,name);
1721           newatom[ni0+j].m     = newatom[ni0+j].mB    = mtot/NMASS;
1722           newatom[ni0+j].q     = newatom[ni0+j].qB    = 0.0;
1723           newatom[ni0+j].type  = newatom[ni0+j].typeB = tpM;
1724           newatom[ni0+j].ptype = eptAtom;
1725           newatom[ni0+j].resind= at->atom[i0].resind;
1726           newvsite_type[ni0+j] = NOTSET;
1727           newcgnr[ni0+j]       = (*cgnr)[i0];
1728         }
1729         /* add constraints between dummy masses and to heavies[0] */
1730         /* 'add_shift' says which atoms won't be renumbered afterwards */
1731         my_add_param(&(plist[F_CONSTRNC]),heavies[0],  add_shift+ni0,  NOTSET);
1732         my_add_param(&(plist[F_CONSTRNC]),heavies[0],  add_shift+ni0+1,NOTSET);
1733         my_add_param(&(plist[F_CONSTRNC]),add_shift+ni0,add_shift+ni0+1,NOTSET);
1734         
1735         /* generate Heavy, H1, H2 and H3 from M1, M2 and heavies[0] */
1736         /* note that vsite_type cannot be NOTSET, because we just set it */
1737         add_vsite3_atoms  (&plist[(*vsite_type)[Heavy]],
1738                          Heavy,     heavies[0], add_shift+ni0, add_shift+ni0+1,
1739                          FALSE);
1740         for(j=0; j<nrHatoms; j++)
1741           add_vsite3_atoms(&plist[(*vsite_type)[Hatoms[j]]], 
1742                          Hatoms[j], heavies[0], add_shift+ni0, add_shift+ni0+1,
1743                          Hat_SwapParity[j]);
1744 #undef NMASS
1745         } else
1746           bWARNING=TRUE;
1747     
1748       }
1749       if (bWARNING)
1750         fprintf(stderr,
1751                 "Warning: cannot convert atom %d %s (bound to a heavy atom "
1752                 "%s with \n"
1753                 "         %d bonds and %d bound hydrogens atoms) to virtual site\n",
1754                 i+1,*(at->atomname[i]),tpname,nrbonds,nrHatoms);
1755       if (bAddVsiteParam) {
1756         /* add vsite parameters to topology, 
1757            also get rid of negative vsite_types */
1758         add_vsites(plist, (*vsite_type), Heavy, nrHatoms, Hatoms,
1759                    nrheavies, heavies);
1760         /* transfer mass of virtual site to Heavy atom */
1761         for(j=0; j<nrHatoms; j++) 
1762           if (is_vsite((*vsite_type)[Hatoms[j]])) {
1763             at->atom[Heavy].m += at->atom[Hatoms[j]].m;
1764             at->atom[Heavy].mB = at->atom[Heavy].m;
1765             at->atom[Hatoms[j]].m = at->atom[Hatoms[j]].mB = 0;
1766           }
1767       }
1768       nvsite+=nrHatoms;
1769       if (debug) {
1770         fprintf(debug,"atom %d: ",o2n[i]+1);
1771         print_bonds(debug,o2n,nrHatoms,Hatoms,Heavy,nrheavies,heavies);
1772       }
1773     } /* if vsite NOTSET & is hydrogen */
1774     
1775   } /* for i < at->nr */
1776
1777   gmx_residuetype_destroy(rt);
1778     
1779   if (debug) {
1780     fprintf(debug,"Before inserting new atoms:\n");
1781     for(i=0; i<at->nr; i++)
1782       fprintf(debug,"%4d %4d %4s %4d %4s %6d %-10s\n",i+1,o2n[i]+1,
1783               at->atomname[i]?*(at->atomname[i]):"(NULL)",
1784               at->resinfo[at->atom[i].resind].nr,
1785               at->resinfo[at->atom[i].resind].name ?
1786               *(at->resinfo[at->atom[i].resind].name):"(NULL)",
1787               (*cgnr)[i],
1788               ((*vsite_type)[i]==NOTSET) ? 
1789               "NOTSET" : interaction_function[(*vsite_type)[i]].name);
1790     fprintf(debug,"new atoms to be inserted:\n");
1791     for(i=0; i<at->nr+nadd; i++)
1792       if (newatomname[i])
1793         fprintf(debug,"%4d %4s %4d %6d %-10s\n",i+1,
1794                 newatomname[i]?*(newatomname[i]):"(NULL)",
1795                 newatom[i].resind,newcgnr[i],
1796                 (newvsite_type[i]==NOTSET) ? 
1797                 "NOTSET" : interaction_function[newvsite_type[i]].name);
1798   }
1799   
1800   /* add all original atoms to the new arrays, using o2n index array */
1801   for(i=0; i<at->nr; i++) {
1802     newatomname  [o2n[i]] = at->atomname [i];
1803     newatom      [o2n[i]] = at->atom     [i];
1804     newvsite_type[o2n[i]] = (*vsite_type)[i];
1805     newcgnr      [o2n[i]] = (*cgnr)      [i];
1806     copy_rvec((*x)[i],newx[o2n[i]]);
1807   }
1808   /* throw away old atoms */
1809   sfree(at->atom);
1810   sfree(at->atomname);
1811   sfree(*vsite_type);
1812   sfree(*cgnr);
1813   sfree(*x);
1814   /* put in the new ones */
1815   at->nr      += nadd;
1816   at->atom     = newatom;
1817   at->atomname = newatomname;
1818   *vsite_type  = newvsite_type;
1819   *cgnr        = newcgnr;
1820   *x           = newx;
1821   if (at->nr > add_shift)
1822     gmx_fatal(FARGS,"Added impossible amount of dummy masses "
1823                 "(%d on a total of %d atoms)\n",nadd,at->nr-nadd);
1824   
1825   if (debug) {
1826     fprintf(debug,"After inserting new atoms:\n");
1827     for(i=0; i<at->nr; i++)
1828       fprintf(debug,"%4d %4s %4d %4s %6d %-10s\n",i+1,
1829               at->atomname[i]?*(at->atomname[i]):"(NULL)",
1830               at->resinfo[at->atom[i].resind].nr, 
1831               at->resinfo[at->atom[i].resind].name ?
1832               *(at->resinfo[at->atom[i].resind].name):"(NULL)",
1833               (*cgnr)[i],
1834               ((*vsite_type)[i]==NOTSET) ? 
1835               "NOTSET" : interaction_function[(*vsite_type)[i]].name);
1836   }
1837   
1838   /* now renumber all the interactions because of the added atoms */
1839   for (ftype=0; ftype<F_NRE; ftype++) {
1840     params=&(plist[ftype]);
1841     if (debug)
1842       fprintf(debug,"Renumbering %d %s\n", params->nr,
1843               interaction_function[ftype].longname);
1844     for (i=0; i<params->nr; i++) {
1845       for (j=0; j<NRAL(ftype); j++)
1846         if (params->param[i].a[j]>=add_shift) {
1847           if (debug) fprintf(debug," [%u -> %u]",params->param[i].a[j],
1848                              params->param[i].a[j]-add_shift);
1849           params->param[i].a[j]=params->param[i].a[j]-add_shift;
1850         } else {
1851           if (debug) fprintf(debug," [%u -> %d]",params->param[i].a[j],
1852                              o2n[params->param[i].a[j]]);
1853           params->param[i].a[j]=o2n[params->param[i].a[j]];
1854         }
1855       if (debug) fprintf(debug,"\n");
1856     }
1857   }
1858   /* now check if atoms in the added constraints are in increasing order */
1859   params=&(plist[F_CONSTRNC]);
1860   for(i=0; i<params->nr; i++)
1861     if ( params->param[i].AI > params->param[i].AJ ) {
1862       j                   = params->param[i].AJ;
1863       params->param[i].AJ = params->param[i].AI;
1864       params->param[i].AI = j;
1865     }
1866   
1867   /* clean up */
1868   sfree(o2n);
1869   
1870   /* tell the user what we did */
1871   fprintf(stderr,"Marked %d virtual sites\n",nvsite);
1872   fprintf(stderr,"Added %d dummy masses\n",nadd);
1873   fprintf(stderr,"Added %d new constraints\n",plist[F_CONSTRNC].nr);
1874 }
1875   
1876 void do_h_mass(t_params *psb, int vsite_type[], t_atoms *at, real mHmult,
1877                gmx_bool bDeuterate)
1878 {
1879   int i,j,a;
1880   
1881   /* loop over all atoms */
1882   for (i=0; i<at->nr; i++)
1883     /* adjust masses if i is hydrogen and not a virtual site */
1884     if ( !is_vsite(vsite_type[i]) && is_hydrogen(*(at->atomname[i])) ) {
1885       /* find bonded heavy atom */
1886       a=NOTSET;
1887       for(j=0; (j<psb->nr) && (a==NOTSET); j++) {
1888         /* if other atom is not a virtual site, it is the one we want */
1889         if ( (psb->param[j].AI==i) && 
1890              !is_vsite(vsite_type[psb->param[j].AJ]) )
1891           a=psb->param[j].AJ;
1892         else if ( (psb->param[j].AJ==i) && 
1893                   !is_vsite(vsite_type[psb->param[j].AI]) )
1894           a=psb->param[j].AI;
1895       }
1896       if (a==NOTSET)
1897         gmx_fatal(FARGS,"Unbound hydrogen atom (%d) found while adjusting mass",
1898                     i+1);
1899       
1900       /* adjust mass of i (hydrogen) with mHmult
1901          and correct mass of a (bonded atom) with same amount */
1902       if (!bDeuterate) {
1903         at->atom[a].m -= (mHmult-1.0)*at->atom[i].m;
1904         at->atom[a].mB-= (mHmult-1.0)*at->atom[i].m;
1905       }
1906       at->atom[i].m *= mHmult;
1907       at->atom[i].mB*= mHmult;
1908     }
1909 }
1910