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