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