2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, 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.
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.
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.
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.
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.
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.
45 #include "gen_vsite.h"
46 #include "gromacs/utility/smalloc.h"
54 #include "gromacs/fileio/futil.h"
55 #include "gpp_atomtype.h"
56 #include "fflibutil.h"
60 #define OPENDIR '[' /* starting sign for directive */
61 #define CLOSEDIR ']' /* ending sign for directive */
64 char atomtype[MAXNAME]; /* Type for the XH3/XH2 atom */
65 gmx_bool isplanar; /* If true, the atomtype above and the three connected
66 * ones are in a planar geometry. The two next entries
67 * are undefined in that case
69 int nhydrogens; /* number of connected hydrogens */
70 char nextheavytype[MAXNAME]; /* Type for the heavy atom bonded to XH2/XH3 */
71 char dummymass[MAXNAME]; /* The type of MNH* or MCH3* dummy mass to use */
75 /* Structure to represent average bond and angles values in vsite aromatic
76 * residues. Note that these are NOT necessarily the bonds and angles from the
77 * forcefield; many forcefields (like Amber, OPLS) have some inherent strain in
78 * 5-rings (i.e. the sum of angles is !=540, but impropers keep it planar)
81 char resname[MAXNAME];
84 struct vsitetop_bond {
88 } *bond; /* list of bonds */
89 struct vsitetop_angle {
94 } *angle; /* list of angles */
99 DDB_CH3, DDB_NH3, DDB_NH2, DDB_PHE, DDB_TYR,
100 DDB_TRP, DDB_HISA, DDB_HISB, DDB_HISH, DDB_DIR_NR
103 typedef char t_dirname[STRLEN];
105 static const t_dirname ddb_dirnames[DDB_DIR_NR] = {
117 static int ddb_name2dir(char *name)
119 /* Translate a directive name to the number of the directive.
120 * HID/HIE/HIP names are translated to the ones we use in Gromacs.
127 for (i = 0; i < DDB_DIR_NR && index < 0; i++)
129 if (!gmx_strcasecmp(name, ddb_dirnames[i]))
139 static void read_vsite_database(const char *ddbname,
140 t_vsiteconf **pvsiteconflist, int *nvsiteconf,
141 t_vsitetop **pvsitetoplist, int *nvsitetop)
143 /* This routine is a quick hack to fix the problem with hardcoded atomtypes
144 * and aromatic vsite parameters by reading them from a ff???.vsd file.
146 * The file can contain sections [ NH3 ], [ CH3 ], [ NH2 ], and ring residue names.
147 * For the NH3 and CH3 section each line has three fields. The first is the atomtype
148 * (nb: not bonded type) of the N/C atom to be replaced, the second field is
149 * the type of the next heavy atom it is bonded to, and the third field the type
150 * of dummy mass that will be used for this group.
152 * If the NH2 group planar (sp2 N) a different vsite construct is used, so in this
153 * case the second field should just be the word planar.
159 int i, j, n, k, nvsite, ntop, curdir, prevdir;
160 t_vsiteconf *vsiteconflist;
161 t_vsitetop *vsitetoplist;
163 char s1[MAXNAME], s2[MAXNAME], s3[MAXNAME], s4[MAXNAME];
165 ddb = libopen(ddbname);
167 nvsite = *nvsiteconf;
168 vsiteconflist = *pvsiteconflist;
170 vsitetoplist = *pvsitetoplist;
174 snew(vsiteconflist, 1);
175 snew(vsitetoplist, 1);
177 while (fgets2(pline, STRLEN-2, ddb) != NULL)
179 strip_comment(pline);
181 if (strlen(pline) > 0)
183 if (pline[0] == OPENDIR)
185 strncpy(dirstr, pline+1, STRLEN-2);
186 if ((ch = strchr (dirstr, CLOSEDIR)) != NULL)
192 if (!gmx_strcasecmp(dirstr, "HID") ||
193 !gmx_strcasecmp(dirstr, "HISD"))
195 sprintf(dirstr, "HISA");
197 else if (!gmx_strcasecmp(dirstr, "HIE") ||
198 !gmx_strcasecmp(dirstr, "HISE"))
200 sprintf(dirstr, "HISB");
202 else if (!gmx_strcasecmp(dirstr, "HIP"))
204 sprintf(dirstr, "HISH");
207 curdir = ddb_name2dir(dirstr);
210 gmx_fatal(FARGS, "Invalid directive %s in vsite database %s",
219 gmx_fatal(FARGS, "First entry in vsite database must be a directive.\n");
224 n = sscanf(pline, "%s%s%s", s1, s2, s3);
225 if (n < 3 && !gmx_strcasecmp(s2, "planar"))
227 srenew(vsiteconflist, nvsite+1);
228 strncpy(vsiteconflist[nvsite].atomtype, s1, MAXNAME-1);
229 vsiteconflist[nvsite].isplanar = TRUE;
230 vsiteconflist[nvsite].nextheavytype[0] = 0;
231 vsiteconflist[nvsite].dummymass[0] = 0;
232 vsiteconflist[nvsite].nhydrogens = 2;
237 srenew(vsiteconflist, (nvsite+1));
238 strncpy(vsiteconflist[nvsite].atomtype, s1, MAXNAME-1);
239 vsiteconflist[nvsite].isplanar = FALSE;
240 strncpy(vsiteconflist[nvsite].nextheavytype, s2, MAXNAME-1);
241 strncpy(vsiteconflist[nvsite].dummymass, s3, MAXNAME-1);
242 if (curdir == DDB_NH2)
244 vsiteconflist[nvsite].nhydrogens = 2;
248 vsiteconflist[nvsite].nhydrogens = 3;
254 gmx_fatal(FARGS, "Not enough directives in vsite database line: %s\n", pline);
264 while ((i < ntop) && gmx_strcasecmp(dirstr, vsitetoplist[i].resname))
268 /* Allocate a new topology entry if this is a new residue */
271 srenew(vsitetoplist, ntop+1);
272 ntop++; /* i still points to current vsite topology entry */
273 strncpy(vsitetoplist[i].resname, dirstr, MAXNAME-1);
274 vsitetoplist[i].nbonds = vsitetoplist[i].nangles = 0;
275 snew(vsitetoplist[i].bond, 1);
276 snew(vsitetoplist[i].angle, 1);
278 n = sscanf(pline, "%s%s%s%s", s1, s2, s3, s4);
282 k = vsitetoplist[i].nbonds++;
283 srenew(vsitetoplist[i].bond, k+1);
284 strncpy(vsitetoplist[i].bond[k].atom1, s1, MAXNAME-1);
285 strncpy(vsitetoplist[i].bond[k].atom2, s2, MAXNAME-1);
286 vsitetoplist[i].bond[k].value = strtod(s3, NULL);
291 k = vsitetoplist[i].nangles++;
292 srenew(vsitetoplist[i].angle, k+1);
293 strncpy(vsitetoplist[i].angle[k].atom1, s1, MAXNAME-1);
294 strncpy(vsitetoplist[i].angle[k].atom2, s2, MAXNAME-1);
295 strncpy(vsitetoplist[i].angle[k].atom3, s3, MAXNAME-1);
296 vsitetoplist[i].angle[k].value = strtod(s4, NULL);
300 gmx_fatal(FARGS, "Need 3 or 4 values to specify bond/angle values in %s: %s\n", ddbname, pline);
304 gmx_fatal(FARGS, "Didnt find a case for directive %s in read_vsite_database\n", dirstr);
310 *pvsiteconflist = vsiteconflist;
311 *pvsitetoplist = vsitetoplist;
312 *nvsiteconf = nvsite;
318 static int nitrogen_is_planar(t_vsiteconf vsiteconflist[], int nvsiteconf, char atomtype[])
320 /* Return 1 if atomtype exists in database list and is planar, 0 if not,
321 * and -1 if not found.
324 gmx_bool found = FALSE;
325 for (i = 0; i < nvsiteconf && !found; i++)
327 found = (!gmx_strcasecmp(vsiteconflist[i].atomtype, atomtype) && (vsiteconflist[i].nhydrogens == 2));
331 res = (vsiteconflist[i-1].isplanar == TRUE);
341 static char *get_dummymass_name(t_vsiteconf vsiteconflist[], int nvsiteconf, char atom[], char nextheavy[])
343 /* Return the dummy mass name if found, or NULL if not set in ddb database */
345 gmx_bool found = FALSE;
346 for (i = 0; i < nvsiteconf && !found; i++)
348 found = (!gmx_strcasecmp(vsiteconflist[i].atomtype, atom) &&
349 !gmx_strcasecmp(vsiteconflist[i].nextheavytype, nextheavy));
353 return vsiteconflist[i-1].dummymass;
363 static real get_ddb_bond(t_vsitetop *vsitetop, int nvsitetop,
365 const char atom1[], const char atom2[])
370 while (i < nvsitetop && gmx_strcasecmp(res, vsitetop[i].resname))
376 gmx_fatal(FARGS, "No vsite information for residue %s found in vsite database.\n", res);
379 while (j < vsitetop[i].nbonds &&
380 ( strcmp(atom1, vsitetop[i].bond[j].atom1) || strcmp(atom2, vsitetop[i].bond[j].atom2)) &&
381 ( strcmp(atom2, vsitetop[i].bond[j].atom1) || strcmp(atom1, vsitetop[i].bond[j].atom2)))
385 if (j == vsitetop[i].nbonds)
387 gmx_fatal(FARGS, "Couldnt find bond %s-%s for residue %s in vsite database.\n", atom1, atom2, res);
390 return vsitetop[i].bond[j].value;
394 static real get_ddb_angle(t_vsitetop *vsitetop, int nvsitetop,
395 const char res[], const char atom1[],
396 const char atom2[], const char atom3[])
401 while (i < nvsitetop && gmx_strcasecmp(res, vsitetop[i].resname))
407 gmx_fatal(FARGS, "No vsite information for residue %s found in vsite database.\n", res);
410 while (j < vsitetop[i].nangles &&
411 ( strcmp(atom1, vsitetop[i].angle[j].atom1) ||
412 strcmp(atom2, vsitetop[i].angle[j].atom2) ||
413 strcmp(atom3, vsitetop[i].angle[j].atom3)) &&
414 ( strcmp(atom3, vsitetop[i].angle[j].atom1) ||
415 strcmp(atom2, vsitetop[i].angle[j].atom2) ||
416 strcmp(atom1, vsitetop[i].angle[j].atom3)))
420 if (j == vsitetop[i].nangles)
422 gmx_fatal(FARGS, "Couldnt find angle %s-%s-%s for residue %s in vsite database.\n", atom1, atom2, atom3, res);
425 return vsitetop[i].angle[j].value;
429 static void count_bonds(int atom, t_params *psb, char ***atomname,
430 int *nrbonds, int *nrHatoms, int Hatoms[], int *Heavy,
431 int *nrheavies, int heavies[])
433 int i, heavy, other, nrb, nrH, nrhv;
435 /* find heavy atom bound to this hydrogen */
437 for (i = 0; (i < psb->nr) && (heavy == NOTSET); i++)
439 if (psb->param[i].AI == atom)
441 heavy = psb->param[i].AJ;
443 else if (psb->param[i].AJ == atom)
445 heavy = psb->param[i].AI;
450 gmx_fatal(FARGS, "unbound hydrogen atom %d", atom+1);
452 /* find all atoms bound to heavy atom */
457 for (i = 0; i < psb->nr; i++)
459 if (psb->param[i].AI == heavy)
461 other = psb->param[i].AJ;
463 else if (psb->param[i].AJ == heavy)
465 other = psb->param[i].AI;
470 if (is_hydrogen(*(atomname[other])))
477 heavies[nrhv] = other;
489 static void print_bonds(FILE *fp, int o2n[],
490 int nrHatoms, int Hatoms[], int Heavy,
491 int nrheavies, int heavies[])
495 fprintf(fp, "Found: %d Hatoms: ", nrHatoms);
496 for (i = 0; i < nrHatoms; i++)
498 fprintf(fp, " %d", o2n[Hatoms[i]]+1);
500 fprintf(fp, "; %d Heavy atoms: %d", nrheavies+1, o2n[Heavy]+1);
501 for (i = 0; i < nrheavies; i++)
503 fprintf(fp, " %d", o2n[heavies[i]]+1);
508 static int get_atype(int atom, t_atoms *at, int nrtp, t_restp rtp[],
509 gmx_residuetype_t rt)
516 if (at->atom[atom].m)
518 type = at->atom[atom].type;
522 /* get type from rtp */
523 rtpp = get_restp(*(at->resinfo[at->atom[atom].resind].name), nrtp, rtp);
524 bNterm = gmx_residuetype_is_protein(rt, *(at->resinfo[at->atom[atom].resind].name)) &&
525 (at->atom[atom].resind == 0);
526 j = search_jtype(rtpp, *(at->atomname[atom]), bNterm);
527 type = rtpp->atom[j].type;
532 static int vsite_nm2type(const char *name, gpp_atomtype_t atype)
536 tp = get_atomtype_type(name, atype);
539 gmx_fatal(FARGS, "Dummy mass type (%s) not found in atom type database",
546 static real get_amass(int atom, t_atoms *at, int nrtp, t_restp rtp[],
547 gmx_residuetype_t rt)
554 if (at->atom[atom].m)
556 mass = at->atom[atom].m;
560 /* get mass from rtp */
561 rtpp = get_restp(*(at->resinfo[at->atom[atom].resind].name), nrtp, rtp);
562 bNterm = gmx_residuetype_is_protein(rt, *(at->resinfo[at->atom[atom].resind].name)) &&
563 (at->atom[atom].resind == 0);
564 j = search_jtype(rtpp, *(at->atomname[atom]), bNterm);
565 mass = rtpp->atom[j].m;
570 static void my_add_param(t_params *plist, int ai, int aj, real b)
572 static real c[MAXFORCEPARAM] =
573 { NOTSET, NOTSET, NOTSET, NOTSET, NOTSET, NOTSET };
576 add_param(plist, ai, aj, c, NULL);
579 static void add_vsites(t_params plist[], int vsite_type[],
580 int Heavy, int nrHatoms, int Hatoms[],
581 int nrheavies, int heavies[])
583 int i, j, ftype, other, moreheavy, bb;
584 gmx_bool bSwapParity;
586 for (i = 0; i < nrHatoms; i++)
588 ftype = vsite_type[Hatoms[i]];
589 /* Errors in setting the vsite_type should really be caugth earlier,
590 * because here it's not possible to print any useful error message.
591 * But it's still better to print a message than to segfault.
595 gmx_incons("Undetected error in setting up virtual sites");
597 bSwapParity = (ftype < 0);
598 vsite_type[Hatoms[i]] = ftype = abs(ftype);
599 if (ftype == F_BONDS)
601 if ( (nrheavies != 1) && (nrHatoms != 1) )
603 gmx_fatal(FARGS, "cannot make constraint in add_vsites for %d heavy "
604 "atoms and %d hydrogen atoms", nrheavies, nrHatoms);
606 my_add_param(&(plist[F_CONSTRNC]), Hatoms[i], heavies[0], NOTSET);
617 gmx_fatal(FARGS, "Not enough heavy atoms (%d) for %s (min 3)",
619 interaction_function[vsite_type[Hatoms[i]]].name);
621 add_vsite3_atoms(&plist[ftype], Hatoms[i], Heavy, heavies[0], heavies[1],
628 moreheavy = heavies[1];
632 /* find more heavy atoms */
633 other = moreheavy = NOTSET;
634 for (j = 0; (j < plist[F_BONDS].nr) && (moreheavy == NOTSET); j++)
636 if (plist[F_BONDS].param[j].AI == heavies[0])
638 other = plist[F_BONDS].param[j].AJ;
640 else if (plist[F_BONDS].param[j].AJ == heavies[0])
642 other = plist[F_BONDS].param[j].AI;
644 if ( (other != NOTSET) && (other != Heavy) )
649 if (moreheavy == NOTSET)
651 gmx_fatal(FARGS, "Unbound molecule part %d-%d", Heavy+1, Hatoms[0]+1);
654 add_vsite3_atoms(&plist[ftype], Hatoms[i], Heavy, heavies[0], moreheavy,
662 gmx_fatal(FARGS, "Not enough heavy atoms (%d) for %s (min 4)",
664 interaction_function[vsite_type[Hatoms[i]]].name);
666 add_vsite4_atoms(&plist[ftype],
667 Hatoms[0], Heavy, heavies[0], heavies[1], heavies[2]);
671 gmx_fatal(FARGS, "can't use add_vsites for interaction function %s",
672 interaction_function[vsite_type[Hatoms[i]]].name);
678 #define ANGLE_6RING (DEG2RAD*120)
680 /* cosine rule: a^2 = b^2 + c^2 - 2 b c cos(alpha) */
681 /* get a^2 when a, b and alpha are given: */
682 #define cosrule(b, c, alpha) ( sqr(b) + sqr(c) - 2*b*c*cos(alpha) )
683 /* get cos(alpha) when a, b and c are given: */
684 #define acosrule(a, b, c) ( (sqr(b)+sqr(c)-sqr(a))/(2*b*c) )
686 static int gen_vsites_6ring(t_atoms *at, int *vsite_type[], t_params plist[],
687 int nrfound, int *ats, real bond_cc, real bond_ch,
688 real xcom, gmx_bool bDoZ)
690 /* these MUST correspond to the atnms array in do_vsite_aromatics! */
692 atCG, atCD1, atHD1, atCD2, atHD2, atCE1, atHE1, atCE2, atHE2,
697 real a, b, dCGCE, tmp1, tmp2, mtot, mG, mrest;
698 real xCG, yCG, xCE1, yCE1, xCE2, yCE2;
699 /* CG, CE1 and CE2 stay and each get a part of the total mass,
700 * so the c-o-m stays the same.
707 gmx_incons("Generating vsites on 6-rings");
711 /* constraints between CG, CE1 and CE2: */
712 dCGCE = sqrt( cosrule(bond_cc, bond_cc, ANGLE_6RING) );
713 my_add_param(&(plist[F_CONSTRNC]), ats[atCG], ats[atCE1], dCGCE);
714 my_add_param(&(plist[F_CONSTRNC]), ats[atCG], ats[atCE2], dCGCE);
715 my_add_param(&(plist[F_CONSTRNC]), ats[atCE1], ats[atCE2], dCGCE);
717 /* rest will be vsite3 */
720 for (i = 0; i < (bDoZ ? atNR : atHZ); i++)
722 mtot += at->atom[ats[i]].m;
723 if (i != atCG && i != atCE1 && i != atCE2 && (bDoZ || (i != atHZ && i != atCZ) ) )
725 at->atom[ats[i]].m = at->atom[ats[i]].mB = 0;
726 (*vsite_type)[ats[i]] = F_VSITE3;
730 /* Distribute mass so center-of-mass stays the same.
731 * The center-of-mass in the call is defined with x=0 at
732 * the CE1-CE2 bond and y=0 at the line from CG to the middle of CE1-CE2 bond.
734 xCG = -bond_cc+bond_cc*cos(ANGLE_6RING);
737 yCE1 = bond_cc*sin(0.5*ANGLE_6RING);
739 yCE2 = -bond_cc*sin(0.5*ANGLE_6RING);
741 mG = at->atom[ats[atCG]].m = at->atom[ats[atCG]].mB = xcom*mtot/xCG;
743 at->atom[ats[atCE1]].m = at->atom[ats[atCE1]].mB =
744 at->atom[ats[atCE2]].m = at->atom[ats[atCE2]].mB = mrest / 2;
746 /* vsite3 construction: r_d = r_i + a r_ij + b r_ik */
747 tmp1 = dCGCE*sin(ANGLE_6RING*0.5);
748 tmp2 = bond_cc*cos(0.5*ANGLE_6RING) + tmp1;
750 a = b = -bond_ch / tmp1;
752 add_vsite3_param(&plist[F_VSITE3],
753 ats[atHE1], ats[atCE1], ats[atCE2], ats[atCG], a, b);
754 add_vsite3_param(&plist[F_VSITE3],
755 ats[atHE2], ats[atCE2], ats[atCE1], ats[atCG], a, b);
756 /* CD1, CD2 and CZ: */
758 add_vsite3_param(&plist[F_VSITE3],
759 ats[atCD1], ats[atCE2], ats[atCE1], ats[atCG], a, b);
760 add_vsite3_param(&plist[F_VSITE3],
761 ats[atCD2], ats[atCE1], ats[atCE2], ats[atCG], a, b);
764 add_vsite3_param(&plist[F_VSITE3],
765 ats[atCZ], ats[atCG], ats[atCE1], ats[atCE2], a, b);
767 /* HD1, HD2 and HZ: */
768 a = b = ( bond_ch + tmp2 ) / tmp1;
769 add_vsite3_param(&plist[F_VSITE3],
770 ats[atHD1], ats[atCE2], ats[atCE1], ats[atCG], a, b);
771 add_vsite3_param(&plist[F_VSITE3],
772 ats[atHD2], ats[atCE1], ats[atCE2], ats[atCG], a, b);
775 add_vsite3_param(&plist[F_VSITE3],
776 ats[atHZ], ats[atCG], ats[atCE1], ats[atCE2], a, b);
782 static int gen_vsites_phe(t_atoms *at, int *vsite_type[], t_params plist[],
783 int nrfound, int *ats, t_vsitetop *vsitetop, int nvsitetop)
785 real bond_cc, bond_ch;
788 /* these MUST correspond to the atnms array in do_vsite_aromatics! */
790 atCG, atCD1, atHD1, atCD2, atHD2, atCE1, atHE1, atCE2, atHE2,
793 real x[atNR], y[atNR];
794 /* Aromatic rings have 6-fold symmetry, so we only need one bond length.
795 * (angle is always 120 degrees).
797 bond_cc = get_ddb_bond(vsitetop, nvsitetop, "PHE", "CD1", "CE1");
798 bond_ch = get_ddb_bond(vsitetop, nvsitetop, "PHE", "CD1", "HD1");
800 x[atCG] = -bond_cc+bond_cc*cos(ANGLE_6RING);
803 y[atCD1] = bond_cc*sin(0.5*ANGLE_6RING);
804 x[atHD1] = x[atCD1]+bond_ch*cos(ANGLE_6RING);
805 y[atHD1] = y[atCD1]+bond_ch*sin(ANGLE_6RING);
808 x[atHE1] = x[atCE1]-bond_ch*cos(ANGLE_6RING);
809 y[atHE1] = y[atCE1]+bond_ch*sin(ANGLE_6RING);
811 y[atCD2] = -y[atCD1];
813 y[atHD2] = -y[atHD1];
815 y[atCE2] = -y[atCE1];
817 y[atHE2] = -y[atHE1];
818 x[atCZ] = bond_cc*cos(0.5*ANGLE_6RING);
820 x[atHZ] = x[atCZ]+bond_ch;
824 for (i = 0; i < atNR; i++)
826 xcom += x[i]*at->atom[ats[i]].m;
827 mtot += at->atom[ats[i]].m;
831 return gen_vsites_6ring(at, vsite_type, plist, nrfound, ats, bond_cc, bond_ch, xcom, TRUE);
834 static void calc_vsite3_param(real xd, real yd, real xi, real yi, real xj, real yj,
835 real xk, real yk, real *a, real *b)
837 /* determine parameters by solving the equation system, since we know the
838 * virtual site coordinates here.
840 real dx_ij, dx_ik, dy_ij, dy_ik;
847 b_ij = sqrt(dx_ij*dx_ij+dy_ij*dy_ij);
848 b_ik = sqrt(dx_ik*dx_ik+dy_ik*dy_ik);
850 *a = ( (xd-xi)*dy_ik - dx_ik*(yd-yi) ) / (dx_ij*dy_ik - dx_ik*dy_ij);
851 *b = ( yd - yi - (*a)*dy_ij ) / dy_ik;
855 static int gen_vsites_trp(gpp_atomtype_t atype, rvec *newx[],
856 t_atom *newatom[], char ***newatomname[],
857 int *o2n[], int *newvsite_type[], int *newcgnr[],
858 t_symtab *symtab, int *nadd, rvec x[], int *cgnr[],
859 t_atoms *at, int *vsite_type[], t_params plist[],
860 int nrfound, int *ats, int add_shift,
861 t_vsitetop *vsitetop, int nvsitetop)
864 /* these MUST correspond to the atnms array in do_vsite_aromatics! */
866 atCB, atCG, atCD1, atHD1, atCD2, atNE1, atHE1, atCE2, atCE3, atHE3,
867 atCZ2, atHZ2, atCZ3, atHZ3, atCH2, atHH2, atNR
869 /* weights for determining the COM's of both rings (M1 and M2): */
870 real mw[NMASS][atNR] = {
871 { 0, 1, 1, 1, 0.5, 1, 1, 0.5, 0, 0,
873 { 0, 0, 0, 0, 0.5, 0, 0, 0.5, 1, 1,
877 real xi[atNR], yi[atNR];
878 real xcom[NMASS], ycom[NMASS], I, alpha;
879 real lineA, lineB, dist;
880 real b_CD2_CE2, b_NE1_CE2, b_CG_CD2, b_CH2_HH2, b_CE2_CZ2;
881 real b_NE1_HE1, b_CD2_CE3, b_CE3_CZ3, b_CB_CG;
882 real b_CZ2_CH2, b_CZ2_HZ2, b_CD1_HD1, b_CE3_HE3;
883 real b_CG_CD1, b_CZ3_HZ3;
884 real a_NE1_CE2_CD2, a_CE2_CD2_CG, a_CB_CG_CD2, a_CE2_CD2_CE3;
885 real a_CB_CG_CD1, a_CD2_CG_CD1, a_CE2_CZ2_HZ2, a_CZ2_CH2_HH2;
886 real a_CD2_CE2_CZ2, a_CD2_CE3_CZ3, a_CE3_CZ3_HZ3, a_CG_CD1_HD1;
887 real a_CE2_CZ2_CH2, a_HE1_NE1_CE2, a_CD2_CE3_HE3;
889 int atM[NMASS], tpM, i, i0, j, nvsite;
890 real mwtot, mtot, mM[NMASS], dCBM1, dCBM2, dM1M2;
891 real a, b, c[MAXFORCEPARAM];
892 rvec r_ij, r_ik, t1, t2;
897 gmx_incons("atom types in gen_vsites_trp");
899 /* Get geometry from database */
900 b_CD2_CE2 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CD2", "CE2");
901 b_NE1_CE2 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "NE1", "CE2");
902 b_CG_CD1 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CG", "CD1");
903 b_CG_CD2 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CG", "CD2");
904 b_CB_CG = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CB", "CG");
905 b_CE2_CZ2 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CE2", "CZ2");
906 b_CD2_CE3 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CD2", "CE3");
907 b_CE3_CZ3 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CE3", "CZ3");
908 b_CZ2_CH2 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CZ2", "CH2");
910 b_CD1_HD1 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CD1", "HD1");
911 b_CZ2_HZ2 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CZ2", "HZ2");
912 b_NE1_HE1 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "NE1", "HE1");
913 b_CH2_HH2 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CH2", "HH2");
914 b_CE3_HE3 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CE3", "HE3");
915 b_CZ3_HZ3 = get_ddb_bond(vsitetop, nvsitetop, "TRP", "CZ3", "HZ3");
917 a_NE1_CE2_CD2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "NE1", "CE2", "CD2");
918 a_CE2_CD2_CG = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CE2", "CD2", "CG");
919 a_CB_CG_CD2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CB", "CG", "CD2");
920 a_CD2_CG_CD1 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CD2", "CG", "CD1");
921 a_CB_CG_CD1 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CB", "CG", "CD1");
923 a_CE2_CD2_CE3 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CE2", "CD2", "CE3");
924 a_CD2_CE2_CZ2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CD2", "CE2", "CZ2");
925 a_CD2_CE3_CZ3 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CD2", "CE3", "CZ3");
926 a_CE3_CZ3_HZ3 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CE3", "CZ3", "HZ3");
927 a_CZ2_CH2_HH2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CZ2", "CH2", "HH2");
928 a_CE2_CZ2_HZ2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CE2", "CZ2", "HZ2");
929 a_CE2_CZ2_CH2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CE2", "CZ2", "CH2");
930 a_CG_CD1_HD1 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CG", "CD1", "HD1");
931 a_HE1_NE1_CE2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "HE1", "NE1", "CE2");
932 a_CD2_CE3_HE3 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TRP", "CD2", "CE3", "HE3");
934 /* Calculate local coordinates.
935 * y-axis (x=0) is the bond CD2-CE2.
936 * x-axis (y=0) is perpendicular to the bond CD2-CE2 and
937 * intersects the middle of the bond.
940 yi[atCD2] = -0.5*b_CD2_CE2;
943 yi[atCE2] = 0.5*b_CD2_CE2;
945 xi[atNE1] = -b_NE1_CE2*sin(a_NE1_CE2_CD2);
946 yi[atNE1] = yi[atCE2]-b_NE1_CE2*cos(a_NE1_CE2_CD2);
948 xi[atCG] = -b_CG_CD2*sin(a_CE2_CD2_CG);
949 yi[atCG] = yi[atCD2]+b_CG_CD2*cos(a_CE2_CD2_CG);
951 alpha = a_CE2_CD2_CG + M_PI - a_CB_CG_CD2;
952 xi[atCB] = xi[atCG]-b_CB_CG*sin(alpha);
953 yi[atCB] = yi[atCG]+b_CB_CG*cos(alpha);
955 alpha = a_CE2_CD2_CG + a_CD2_CG_CD1 - M_PI;
956 xi[atCD1] = xi[atCG]-b_CG_CD1*sin(alpha);
957 yi[atCD1] = yi[atCG]+b_CG_CD1*cos(alpha);
959 xi[atCE3] = b_CD2_CE3*sin(a_CE2_CD2_CE3);
960 yi[atCE3] = yi[atCD2]+b_CD2_CE3*cos(a_CE2_CD2_CE3);
962 xi[atCZ2] = b_CE2_CZ2*sin(a_CD2_CE2_CZ2);
963 yi[atCZ2] = yi[atCE2]-b_CE2_CZ2*cos(a_CD2_CE2_CZ2);
965 alpha = a_CE2_CD2_CE3 + a_CD2_CE3_CZ3 - M_PI;
966 xi[atCZ3] = xi[atCE3]+b_CE3_CZ3*sin(alpha);
967 yi[atCZ3] = yi[atCE3]+b_CE3_CZ3*cos(alpha);
969 alpha = a_CD2_CE2_CZ2 + a_CE2_CZ2_CH2 - M_PI;
970 xi[atCH2] = xi[atCZ2]+b_CZ2_CH2*sin(alpha);
971 yi[atCH2] = yi[atCZ2]-b_CZ2_CH2*cos(alpha);
974 alpha = a_CE2_CD2_CG + a_CD2_CG_CD1 - a_CG_CD1_HD1;
975 xi[atHD1] = xi[atCD1]-b_CD1_HD1*sin(alpha);
976 yi[atHD1] = yi[atCD1]+b_CD1_HD1*cos(alpha);
978 alpha = a_NE1_CE2_CD2 + M_PI - a_HE1_NE1_CE2;
979 xi[atHE1] = xi[atNE1]-b_NE1_HE1*sin(alpha);
980 yi[atHE1] = yi[atNE1]-b_NE1_HE1*cos(alpha);
982 alpha = a_CE2_CD2_CE3 + M_PI - a_CD2_CE3_HE3;
983 xi[atHE3] = xi[atCE3]+b_CE3_HE3*sin(alpha);
984 yi[atHE3] = yi[atCE3]+b_CE3_HE3*cos(alpha);
986 alpha = a_CD2_CE2_CZ2 + M_PI - a_CE2_CZ2_HZ2;
987 xi[atHZ2] = xi[atCZ2]+b_CZ2_HZ2*sin(alpha);
988 yi[atHZ2] = yi[atCZ2]-b_CZ2_HZ2*cos(alpha);
990 alpha = a_CD2_CE2_CZ2 + a_CE2_CZ2_CH2 - a_CZ2_CH2_HH2;
991 xi[atHZ3] = xi[atCZ3]+b_CZ3_HZ3*sin(alpha);
992 yi[atHZ3] = yi[atCZ3]+b_CZ3_HZ3*cos(alpha);
994 alpha = a_CE2_CD2_CE3 + a_CD2_CE3_CZ3 - a_CE3_CZ3_HZ3;
995 xi[atHH2] = xi[atCH2]+b_CH2_HH2*sin(alpha);
996 yi[atHH2] = yi[atCH2]-b_CH2_HH2*cos(alpha);
998 /* Determine coeff. for the line CB-CG */
999 lineA = (yi[atCB]-yi[atCG])/(xi[atCB]-xi[atCG]);
1000 lineB = yi[atCG]-lineA*xi[atCG];
1002 /* Calculate masses for each ring and put it on the dummy masses */
1003 for (j = 0; j < NMASS; j++)
1005 mM[j] = xcom[j] = ycom[j] = 0;
1007 for (i = 0; i < atNR; i++)
1011 for (j = 0; j < NMASS; j++)
1013 mM[j] += mw[j][i] * at->atom[ats[i]].m;
1014 xcom[j] += xi[i] * mw[j][i] * at->atom[ats[i]].m;
1015 ycom[j] += yi[i] * mw[j][i] * at->atom[ats[i]].m;
1019 for (j = 0; j < NMASS; j++)
1025 /* get dummy mass type */
1026 tpM = vsite_nm2type("MW", atype);
1027 /* make space for 2 masses: shift all atoms starting with CB */
1029 for (j = 0; j < NMASS; j++)
1031 atM[j] = i0+*nadd+j;
1035 fprintf(stderr, "Inserting %d dummy masses at %d\n", NMASS, (*o2n)[i0]+1);
1038 for (j = i0; j < at->nr; j++)
1040 (*o2n)[j] = j+*nadd;
1042 srenew(*newx, at->nr+*nadd);
1043 srenew(*newatom, at->nr+*nadd);
1044 srenew(*newatomname, at->nr+*nadd);
1045 srenew(*newvsite_type, at->nr+*nadd);
1046 srenew(*newcgnr, at->nr+*nadd);
1047 for (j = 0; j < NMASS; j++)
1049 (*newatomname)[at->nr+*nadd-1-j] = NULL;
1052 /* Dummy masses will be placed at the center-of-mass in each ring. */
1054 /* calc initial position for dummy masses in real (non-local) coordinates.
1055 * Cheat by using the routine to calculate virtual site parameters. It is
1056 * much easier when we have the coordinates expressed in terms of
1059 rvec_sub(x[ats[atCB]], x[ats[atCG]], r_ij);
1060 rvec_sub(x[ats[atCD2]], x[ats[atCG]], r_ik);
1061 calc_vsite3_param(xcom[0], ycom[0], xi[atCG], yi[atCG], xi[atCB], yi[atCB],
1062 xi[atCD2], yi[atCD2], &a, &b);
1065 rvec_add(t1, t2, t1);
1066 rvec_add(t1, x[ats[atCG]], (*newx)[atM[0]]);
1068 calc_vsite3_param(xcom[1], ycom[1], xi[atCG], yi[atCG], xi[atCB], yi[atCB],
1069 xi[atCD2], yi[atCD2], &a, &b);
1072 rvec_add(t1, t2, t1);
1073 rvec_add(t1, x[ats[atCG]], (*newx)[atM[1]]);
1075 /* set parameters for the masses */
1076 for (j = 0; j < NMASS; j++)
1078 sprintf(name, "MW%d", j+1);
1079 (*newatomname) [atM[j]] = put_symtab(symtab, name);
1080 (*newatom) [atM[j]].m = (*newatom)[atM[j]].mB = mM[j];
1081 (*newatom) [atM[j]].q = (*newatom)[atM[j]].qB = 0.0;
1082 (*newatom) [atM[j]].type = (*newatom)[atM[j]].typeB = tpM;
1083 (*newatom) [atM[j]].ptype = eptAtom;
1084 (*newatom) [atM[j]].resind = at->atom[i0].resind;
1085 (*newatom) [atM[j]].elem[0] = 'M';
1086 (*newatom) [atM[j]].elem[1] = '\0';
1087 (*newvsite_type)[atM[j]] = NOTSET;
1088 (*newcgnr) [atM[j]] = (*cgnr)[i0];
1090 /* renumber cgnr: */
1091 for (i = i0; i < at->nr; i++)
1096 /* constraints between CB, M1 and M2 */
1097 /* 'add_shift' says which atoms won't be renumbered afterwards */
1098 dCBM1 = sqrt( sqr(xcom[0]-xi[atCB]) + sqr(ycom[0]-yi[atCB]) );
1099 dM1M2 = sqrt( sqr(xcom[0]-xcom[1]) + sqr(ycom[0]-ycom[1]) );
1100 dCBM2 = sqrt( sqr(xcom[1]-xi[atCB]) + sqr(ycom[1]-yi[atCB]) );
1101 my_add_param(&(plist[F_CONSTRNC]), ats[atCB], add_shift+atM[0], dCBM1);
1102 my_add_param(&(plist[F_CONSTRNC]), ats[atCB], add_shift+atM[1], dCBM2);
1103 my_add_param(&(plist[F_CONSTRNC]), add_shift+atM[0], add_shift+atM[1], dM1M2);
1105 /* rest will be vsite3 */
1107 for (i = 0; i < atNR; i++)
1111 at->atom[ats[i]].m = at->atom[ats[i]].mB = 0;
1112 (*vsite_type)[ats[i]] = F_VSITE3;
1117 /* now define all vsites from M1, M2, CB, ie:
1118 r_d = r_M1 + a r_M1_M2 + b r_M1_CB */
1119 for (i = 0; i < atNR; i++)
1121 if ( (*vsite_type)[ats[i]] == F_VSITE3)
1123 calc_vsite3_param(xi[i], yi[i], xcom[0], ycom[0], xcom[1], ycom[1], xi[atCB], yi[atCB], &a, &b);
1124 add_vsite3_param(&plist[F_VSITE3],
1125 ats[i], add_shift+atM[0], add_shift+atM[1], ats[atCB], a, b);
1133 static int gen_vsites_tyr(gpp_atomtype_t atype, rvec *newx[],
1134 t_atom *newatom[], char ***newatomname[],
1135 int *o2n[], int *newvsite_type[], int *newcgnr[],
1136 t_symtab *symtab, int *nadd, rvec x[], int *cgnr[],
1137 t_atoms *at, int *vsite_type[], t_params plist[],
1138 int nrfound, int *ats, int add_shift,
1139 t_vsitetop *vsitetop, int nvsitetop)
1141 int nvsite, i, i0, j, atM, tpM;
1142 real dCGCE, dCEOH, dCGM, tmp1, a, b;
1143 real bond_cc, bond_ch, bond_co, bond_oh, angle_coh;
1145 real vmass, vdist, mM;
1149 /* these MUST correspond to the atnms array in do_vsite_aromatics! */
1151 atCG, atCD1, atHD1, atCD2, atHD2, atCE1, atHE1, atCE2, atHE2,
1152 atCZ, atOH, atHH, atNR
1154 real xi[atNR], yi[atNR];
1155 /* CG, CE1, CE2 (as in general 6-ring) and OH and HH stay,
1156 rest gets virtualized.
1157 Now we have two linked triangles with one improper keeping them flat */
1158 if (atNR != nrfound)
1160 gmx_incons("Number of atom types in gen_vsites_tyr");
1163 /* Aromatic rings have 6-fold symmetry, so we only need one bond length
1164 * for the ring part (angle is always 120 degrees).
1166 bond_cc = get_ddb_bond(vsitetop, nvsitetop, "TYR", "CD1", "CE1");
1167 bond_ch = get_ddb_bond(vsitetop, nvsitetop, "TYR", "CD1", "HD1");
1168 bond_co = get_ddb_bond(vsitetop, nvsitetop, "TYR", "CZ", "OH");
1169 bond_oh = get_ddb_bond(vsitetop, nvsitetop, "TYR", "OH", "HH");
1170 angle_coh = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, "TYR", "CZ", "OH", "HH");
1172 xi[atCG] = -bond_cc+bond_cc*cos(ANGLE_6RING);
1174 xi[atCD1] = -bond_cc;
1175 yi[atCD1] = bond_cc*sin(0.5*ANGLE_6RING);
1176 xi[atHD1] = xi[atCD1]+bond_ch*cos(ANGLE_6RING);
1177 yi[atHD1] = yi[atCD1]+bond_ch*sin(ANGLE_6RING);
1179 yi[atCE1] = yi[atCD1];
1180 xi[atHE1] = xi[atCE1]-bond_ch*cos(ANGLE_6RING);
1181 yi[atHE1] = yi[atCE1]+bond_ch*sin(ANGLE_6RING);
1182 xi[atCD2] = xi[atCD1];
1183 yi[atCD2] = -yi[atCD1];
1184 xi[atHD2] = xi[atHD1];
1185 yi[atHD2] = -yi[atHD1];
1186 xi[atCE2] = xi[atCE1];
1187 yi[atCE2] = -yi[atCE1];
1188 xi[atHE2] = xi[atHE1];
1189 yi[atHE2] = -yi[atHE1];
1190 xi[atCZ] = bond_cc*cos(0.5*ANGLE_6RING);
1192 xi[atOH] = xi[atCZ]+bond_co;
1196 for (i = 0; i < atOH; i++)
1198 xcom += xi[i]*at->atom[ats[i]].m;
1199 mtot += at->atom[ats[i]].m;
1203 /* first do 6 ring as default,
1204 except CZ (we'll do that different) and HZ (we don't have that): */
1205 nvsite = gen_vsites_6ring(at, vsite_type, plist, nrfound, ats, bond_cc, bond_ch, xcom, FALSE);
1207 /* then construct CZ from the 2nd triangle */
1208 /* vsite3 construction: r_d = r_i + a r_ij + b r_ik */
1209 a = b = 0.5 * bond_co / ( bond_co - bond_cc*cos(ANGLE_6RING) );
1210 add_vsite3_param(&plist[F_VSITE3],
1211 ats[atCZ], ats[atOH], ats[atCE1], ats[atCE2], a, b);
1212 at->atom[ats[atCZ]].m = at->atom[ats[atCZ]].mB = 0;
1214 /* constraints between CE1, CE2 and OH */
1215 dCGCE = sqrt( cosrule(bond_cc, bond_cc, ANGLE_6RING) );
1216 dCEOH = sqrt( cosrule(bond_cc, bond_co, ANGLE_6RING) );
1217 my_add_param(&(plist[F_CONSTRNC]), ats[atCE1], ats[atOH], dCEOH);
1218 my_add_param(&(plist[F_CONSTRNC]), ats[atCE2], ats[atOH], dCEOH);
1220 /* We also want to constrain the angle C-O-H, but since CZ is constructed
1221 * we need to introduce a constraint to CG.
1222 * CG is much further away, so that will lead to instabilities in LINCS
1223 * when we constrain both CG-HH and OH-HH distances. Instead of requiring
1224 * the use of lincs_order=8 we introduce a dummy mass three times further
1225 * away from OH than HH. The mass is accordingly a third, with the remaining
1226 * 2/3 moved to OH. This shouldnt cause any problems since the forces will
1227 * apply to the HH constructed atom and not directly on the virtual mass.
1230 vdist = 2.0*bond_oh;
1231 mM = at->atom[ats[atHH]].m/2.0;
1232 at->atom[ats[atOH]].m += mM; /* add 1/2 of original H mass */
1233 at->atom[ats[atOH]].mB += mM; /* add 1/2 of original H mass */
1234 at->atom[ats[atHH]].m = at->atom[ats[atHH]].mB = 0;
1236 /* get dummy mass type */
1237 tpM = vsite_nm2type("MW", atype);
1238 /* make space for 1 mass: shift HH only */
1243 fprintf(stderr, "Inserting 1 dummy mass at %d\n", (*o2n)[i0]+1);
1246 for (j = i0; j < at->nr; j++)
1248 (*o2n)[j] = j+*nadd;
1250 srenew(*newx, at->nr+*nadd);
1251 srenew(*newatom, at->nr+*nadd);
1252 srenew(*newatomname, at->nr+*nadd);
1253 srenew(*newvsite_type, at->nr+*nadd);
1254 srenew(*newcgnr, at->nr+*nadd);
1255 (*newatomname)[at->nr+*nadd-1] = NULL;
1257 /* Calc the dummy mass initial position */
1258 rvec_sub(x[ats[atHH]], x[ats[atOH]], r1);
1260 rvec_add(r1, x[ats[atHH]], (*newx)[atM]);
1262 strcpy(name, "MW1");
1263 (*newatomname) [atM] = put_symtab(symtab, name);
1264 (*newatom) [atM].m = (*newatom)[atM].mB = mM;
1265 (*newatom) [atM].q = (*newatom)[atM].qB = 0.0;
1266 (*newatom) [atM].type = (*newatom)[atM].typeB = tpM;
1267 (*newatom) [atM].ptype = eptAtom;
1268 (*newatom) [atM].resind = at->atom[i0].resind;
1269 (*newatom) [atM].elem[0] = 'M';
1270 (*newatom) [atM].elem[1] = '\0';
1271 (*newvsite_type)[atM] = NOTSET;
1272 (*newcgnr) [atM] = (*cgnr)[i0];
1273 /* renumber cgnr: */
1274 for (i = i0; i < at->nr; i++)
1279 (*vsite_type)[ats[atHH]] = F_VSITE2;
1281 /* assume we also want the COH angle constrained: */
1282 tmp1 = bond_cc*cos(0.5*ANGLE_6RING) + dCGCE*sin(ANGLE_6RING*0.5) + bond_co;
1283 dCGM = sqrt( cosrule(tmp1, vdist, angle_coh) );
1284 my_add_param(&(plist[F_CONSTRNC]), ats[atCG], add_shift+atM, dCGM);
1285 my_add_param(&(plist[F_CONSTRNC]), ats[atOH], add_shift+atM, vdist);
1287 add_vsite2_param(&plist[F_VSITE2],
1288 ats[atHH], ats[atOH], add_shift+atM, 1.0/2.0);
1292 static int gen_vsites_his(t_atoms *at, int *vsite_type[], t_params plist[],
1293 int nrfound, int *ats, t_vsitetop *vsitetop, int nvsitetop)
1296 real a, b, alpha, dCGCE1, dCGNE2;
1297 real sinalpha, cosalpha;
1298 real xcom, ycom, mtot;
1299 real mG, mrest, mCE1, mNE2;
1300 real b_CG_ND1, b_ND1_CE1, b_CE1_NE2, b_CG_CD2, b_CD2_NE2;
1301 real b_ND1_HD1, b_NE2_HE2, b_CE1_HE1, b_CD2_HD2;
1302 real a_CG_ND1_CE1, a_CG_CD2_NE2, a_ND1_CE1_NE2, a_CE1_NE2_CD2;
1303 real a_NE2_CE1_HE1, a_NE2_CD2_HD2, a_CE1_ND1_HD1, a_CE1_NE2_HE2;
1306 /* these MUST correspond to the atnms array in do_vsite_aromatics! */
1308 atCG, atND1, atHD1, atCD2, atHD2, atCE1, atHE1, atNE2, atHE2, atNR
1310 real x[atNR], y[atNR];
1312 /* CG, CE1 and NE2 stay, each gets part of the total mass,
1313 rest gets virtualized */
1314 /* check number of atoms, 3 hydrogens may be missing: */
1315 /* assert( nrfound >= atNR-3 || nrfound <= atNR );
1316 * Don't understand the above logic. Shouldn't it be && rather than || ???
1318 if ((nrfound < atNR-3) || (nrfound > atNR))
1320 gmx_incons("Generating vsites for HIS");
1323 /* avoid warnings about uninitialized variables */
1324 b_ND1_HD1 = b_NE2_HE2 = b_CE1_HE1 = b_CD2_HD2 = a_NE2_CE1_HE1 =
1325 a_NE2_CD2_HD2 = a_CE1_ND1_HD1 = a_CE1_NE2_HE2 = 0;
1327 if (ats[atHD1] != NOTSET)
1329 if (ats[atHE2] != NOTSET)
1331 sprintf(resname, "HISH");
1335 sprintf(resname, "HISA");
1340 sprintf(resname, "HISB");
1343 /* Get geometry from database */
1344 b_CG_ND1 = get_ddb_bond(vsitetop, nvsitetop, resname, "CG", "ND1");
1345 b_ND1_CE1 = get_ddb_bond(vsitetop, nvsitetop, resname, "ND1", "CE1");
1346 b_CE1_NE2 = get_ddb_bond(vsitetop, nvsitetop, resname, "CE1", "NE2");
1347 b_CG_CD2 = get_ddb_bond(vsitetop, nvsitetop, resname, "CG", "CD2");
1348 b_CD2_NE2 = get_ddb_bond(vsitetop, nvsitetop, resname, "CD2", "NE2");
1349 a_CG_ND1_CE1 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, resname, "CG", "ND1", "CE1");
1350 a_CG_CD2_NE2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, resname, "CG", "CD2", "NE2");
1351 a_ND1_CE1_NE2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, resname, "ND1", "CE1", "NE2");
1352 a_CE1_NE2_CD2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, resname, "CE1", "NE2", "CD2");
1354 if (ats[atHD1] != NOTSET)
1356 b_ND1_HD1 = get_ddb_bond(vsitetop, nvsitetop, resname, "ND1", "HD1");
1357 a_CE1_ND1_HD1 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, resname, "CE1", "ND1", "HD1");
1359 if (ats[atHE2] != NOTSET)
1361 b_NE2_HE2 = get_ddb_bond(vsitetop, nvsitetop, resname, "NE2", "HE2");
1362 a_CE1_NE2_HE2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, resname, "CE1", "NE2", "HE2");
1364 if (ats[atHD2] != NOTSET)
1366 b_CD2_HD2 = get_ddb_bond(vsitetop, nvsitetop, resname, "CD2", "HD2");
1367 a_NE2_CD2_HD2 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, resname, "NE2", "CD2", "HD2");
1369 if (ats[atHE1] != NOTSET)
1371 b_CE1_HE1 = get_ddb_bond(vsitetop, nvsitetop, resname, "CE1", "HE1");
1372 a_NE2_CE1_HE1 = DEG2RAD*get_ddb_angle(vsitetop, nvsitetop, resname, "NE2", "CE1", "HE1");
1375 /* constraints between CG, CE1 and NE1 */
1376 dCGCE1 = sqrt( cosrule(b_CG_ND1, b_ND1_CE1, a_CG_ND1_CE1) );
1377 dCGNE2 = sqrt( cosrule(b_CG_CD2, b_CD2_NE2, a_CG_CD2_NE2) );
1379 my_add_param(&(plist[F_CONSTRNC]), ats[atCG], ats[atCE1], dCGCE1);
1380 my_add_param(&(plist[F_CONSTRNC]), ats[atCG], ats[atNE2], dCGNE2);
1381 /* we already have a constraint CE1-NE2, so we don't add it again */
1383 /* calculate the positions in a local frame of reference.
1384 * The x-axis is the line from CG that makes a right angle
1385 * with the bond CE1-NE2, and the y-axis the bond CE1-NE2.
1387 /* First calculate the x-axis intersection with y-axis (=yCE1).
1388 * Get cos(angle CG-CE1-NE2) :
1390 cosalpha = acosrule(dCGNE2, dCGCE1, b_CE1_NE2);
1392 y[atCE1] = cosalpha*dCGCE1;
1394 y[atNE2] = y[atCE1]-b_CE1_NE2;
1395 sinalpha = sqrt(1-cosalpha*cosalpha);
1396 x[atCG] = -sinalpha*dCGCE1;
1398 x[atHE1] = x[atHE2] = x[atHD1] = x[atHD2] = 0;
1399 y[atHE1] = y[atHE2] = y[atHD1] = y[atHD2] = 0;
1401 /* calculate ND1 and CD2 positions from CE1 and NE2 */
1403 x[atND1] = -b_ND1_CE1*sin(a_ND1_CE1_NE2);
1404 y[atND1] = y[atCE1]-b_ND1_CE1*cos(a_ND1_CE1_NE2);
1406 x[atCD2] = -b_CD2_NE2*sin(a_CE1_NE2_CD2);
1407 y[atCD2] = y[atNE2]+b_CD2_NE2*cos(a_CE1_NE2_CD2);
1409 /* And finally the hydrogen positions */
1410 if (ats[atHE1] != NOTSET)
1412 x[atHE1] = x[atCE1] + b_CE1_HE1*sin(a_NE2_CE1_HE1);
1413 y[atHE1] = y[atCE1] - b_CE1_HE1*cos(a_NE2_CE1_HE1);
1415 /* HD2 - first get (ccw) angle from (positive) y-axis */
1416 if (ats[atHD2] != NOTSET)
1418 alpha = a_CE1_NE2_CD2 + M_PI - a_NE2_CD2_HD2;
1419 x[atHD2] = x[atCD2] - b_CD2_HD2*sin(alpha);
1420 y[atHD2] = y[atCD2] + b_CD2_HD2*cos(alpha);
1422 if (ats[atHD1] != NOTSET)
1424 /* HD1 - first get (cw) angle from (positive) y-axis */
1425 alpha = a_ND1_CE1_NE2 + M_PI - a_CE1_ND1_HD1;
1426 x[atHD1] = x[atND1] - b_ND1_HD1*sin(alpha);
1427 y[atHD1] = y[atND1] - b_ND1_HD1*cos(alpha);
1429 if (ats[atHE2] != NOTSET)
1431 x[atHE2] = x[atNE2] + b_NE2_HE2*sin(a_CE1_NE2_HE2);
1432 y[atHE2] = y[atNE2] + b_NE2_HE2*cos(a_CE1_NE2_HE2);
1434 /* Have all coordinates now */
1436 /* calc center-of-mass; keep atoms CG, CE1, NE2 and
1437 * set the rest to vsite3
1439 mtot = xcom = ycom = 0;
1441 for (i = 0; i < atNR; i++)
1443 if (ats[i] != NOTSET)
1445 mtot += at->atom[ats[i]].m;
1446 xcom += x[i]*at->atom[ats[i]].m;
1447 ycom += y[i]*at->atom[ats[i]].m;
1448 if (i != atCG && i != atCE1 && i != atNE2)
1450 at->atom[ats[i]].m = at->atom[ats[i]].mB = 0;
1451 (*vsite_type)[ats[i]] = F_VSITE3;
1456 if (nvsite+3 != nrfound)
1458 gmx_incons("Generating vsites for HIS");
1464 /* distribute mass so that com stays the same */
1465 mG = xcom*mtot/x[atCG];
1467 mCE1 = (ycom-y[atNE2])*mrest/(y[atCE1]-y[atNE2]);
1470 at->atom[ats[atCG]].m = at->atom[ats[atCG]].mB = mG;
1471 at->atom[ats[atCE1]].m = at->atom[ats[atCE1]].mB = mCE1;
1472 at->atom[ats[atNE2]].m = at->atom[ats[atNE2]].mB = mNE2;
1475 if (ats[atHE1] != NOTSET)
1477 calc_vsite3_param(x[atHE1], y[atHE1], x[atCE1], y[atCE1], x[atNE2], y[atNE2],
1478 x[atCG], y[atCG], &a, &b);
1479 add_vsite3_param(&plist[F_VSITE3],
1480 ats[atHE1], ats[atCE1], ats[atNE2], ats[atCG], a, b);
1483 if (ats[atHE2] != NOTSET)
1485 calc_vsite3_param(x[atHE2], y[atHE2], x[atNE2], y[atNE2], x[atCE1], y[atCE1],
1486 x[atCG], y[atCG], &a, &b);
1487 add_vsite3_param(&plist[F_VSITE3],
1488 ats[atHE2], ats[atNE2], ats[atCE1], ats[atCG], a, b);
1492 calc_vsite3_param(x[atND1], y[atND1], x[atNE2], y[atNE2], x[atCE1], y[atCE1],
1493 x[atCG], y[atCG], &a, &b);
1494 add_vsite3_param(&plist[F_VSITE3],
1495 ats[atND1], ats[atNE2], ats[atCE1], ats[atCG], a, b);
1498 calc_vsite3_param(x[atCD2], y[atCD2], x[atCE1], y[atCE1], x[atNE2], y[atNE2],
1499 x[atCG], y[atCG], &a, &b);
1500 add_vsite3_param(&plist[F_VSITE3],
1501 ats[atCD2], ats[atCE1], ats[atNE2], ats[atCG], a, b);
1504 if (ats[atHD1] != NOTSET)
1506 calc_vsite3_param(x[atHD1], y[atHD1], x[atNE2], y[atNE2], x[atCE1], y[atCE1],
1507 x[atCG], y[atCG], &a, &b);
1508 add_vsite3_param(&plist[F_VSITE3],
1509 ats[atHD1], ats[atNE2], ats[atCE1], ats[atCG], a, b);
1512 if (ats[atHD2] != NOTSET)
1514 calc_vsite3_param(x[atHD2], y[atHD2], x[atCE1], y[atCE1], x[atNE2], y[atNE2],
1515 x[atCG], y[atCG], &a, &b);
1516 add_vsite3_param(&plist[F_VSITE3],
1517 ats[atHD2], ats[atCE1], ats[atNE2], ats[atCG], a, b);
1522 static gmx_bool is_vsite(int vsite_type)
1524 if (vsite_type == NOTSET)
1528 switch (abs(vsite_type) )
1542 static char atomnamesuffix[] = "1234";
1544 void do_vsites(int nrtp, t_restp rtp[], gpp_atomtype_t atype,
1545 t_atoms *at, t_symtab *symtab, rvec *x[],
1546 t_params plist[], int *vsite_type[], int *cgnr[],
1547 real mHmult, gmx_bool bVsiteAromatics,
1550 #define MAXATOMSPERRESIDUE 16
1551 int i, j, k, m, i0, ni0, whatres, resind, add_shift, ftype, nvsite, nadd;
1553 int nrfound = 0, needed, nrbonds, nrHatoms, Heavy, nrheavies, tpM, tpHeavy;
1554 int Hatoms[4], heavies[4], bb;
1555 gmx_bool bWARNING, bAddVsiteParam, bFirstWater;
1557 gmx_bool *bResProcessed;
1558 real mHtot, mtot, fact, fact2;
1559 rvec rpar, rperp, temp;
1560 char name[10], tpname[32], nexttpname[32], *ch;
1562 int *o2n, *newvsite_type, *newcgnr, ats[MAXATOMSPERRESIDUE];
1565 char ***newatomname;
1569 int nvsiteconf, nvsitetop, cmplength;
1570 gmx_bool isN, planarN, bFound;
1571 gmx_residuetype_t rt;
1573 t_vsiteconf *vsiteconflist;
1574 /* pointer to a list of CH3/NH3/NH2 configuration entries.
1575 * See comments in read_vsite_database. It isnt beautiful,
1576 * but it had to be fixed, and I dont even want to try to
1577 * maintain this part of the code...
1579 t_vsitetop *vsitetop;
1580 /* Pointer to a list of geometry (bond/angle) entries for
1581 * residues like PHE, TRP, TYR, HIS, etc., where we need
1582 * to know the geometry to construct vsite aromatics.
1583 * Note that equilibrium geometry isnt necessarily the same
1584 * as the individual bond and angle values given in the
1585 * force field (rings can be strained).
1588 /* if bVsiteAromatics=TRUE do_vsites will specifically convert atoms in
1589 PHE, TRP, TYR and HIS to a construction of virtual sites */
1591 resPHE, resTRP, resTYR, resHIS, resNR
1593 const char *resnms[resNR] = { "PHE", "TRP", "TYR", "HIS" };
1594 /* Amber03 alternative names for termini */
1595 const char *resnmsN[resNR] = { "NPHE", "NTRP", "NTYR", "NHIS" };
1596 const char *resnmsC[resNR] = { "CPHE", "CTRP", "CTYR", "CHIS" };
1597 /* HIS can be known as HISH, HIS1, HISA, HID, HIE, HIP, etc. too */
1598 gmx_bool bPartial[resNR] = { FALSE, FALSE, FALSE, TRUE };
1599 /* the atnms for every residue MUST correspond to the enums in the
1600 gen_vsites_* (one for each residue) routines! */
1601 /* also the atom names in atnms MUST be in the same order as in the .rtp! */
1602 const char *atnms[resNR][MAXATOMSPERRESIDUE+1] = {
1604 "CD1", "HD1", "CD2", "HD2",
1605 "CE1", "HE1", "CE2", "HE2",
1609 "CD1", "HD1", "CD2",
1610 "NE1", "HE1", "CE2", "CE3", "HE3",
1611 "CZ2", "HZ2", "CZ3", "HZ3",
1612 "CH2", "HH2", NULL },
1614 "CD1", "HD1", "CD2", "HD2",
1615 "CE1", "HE1", "CE2", "HE2",
1616 "CZ", "OH", "HH", NULL },
1618 "ND1", "HD1", "CD2", "HD2",
1619 "CE1", "HE1", "NE2", "HE2", NULL }
1624 printf("Searching for atoms to make virtual sites ...\n");
1625 fprintf(debug, "# # # VSITES # # #\n");
1628 ndb = fflib_search_file_end(ffdir, ".vsd", FALSE, &db);
1630 vsiteconflist = NULL;
1633 for (f = 0; f < ndb; f++)
1635 read_vsite_database(db[f], &vsiteconflist, &nvsiteconf, &vsitetop, &nvsitetop);
1643 /* we need a marker for which atoms should *not* be renumbered afterwards */
1644 add_shift = 10*at->nr;
1645 /* make arrays where masses can be inserted into */
1647 snew(newatom, at->nr);
1648 snew(newatomname, at->nr);
1649 snew(newvsite_type, at->nr);
1650 snew(newcgnr, at->nr);
1651 /* make index array to tell where the atoms go to when masses are inserted */
1653 for (i = 0; i < at->nr; i++)
1657 /* make index to tell which residues were already processed */
1658 snew(bResProcessed, at->nres);
1660 gmx_residuetype_init(&rt);
1662 /* generate vsite constructions */
1663 /* loop over all atoms */
1665 for (i = 0; (i < at->nr); i++)
1667 if (at->atom[i].resind != resind)
1669 resind = at->atom[i].resind;
1670 resnm = *(at->resinfo[resind].name);
1672 /* first check for aromatics to virtualize */
1673 /* don't waste our effort on DNA, water etc. */
1674 /* Only do the vsite aromatic stuff when we reach the
1675 * CA atom, since there might be an X2/X3 group on the
1676 * N-terminus that must be treated first.
1678 if (bVsiteAromatics &&
1679 !strcmp(*(at->atomname[i]), "CA") &&
1680 !bResProcessed[resind] &&
1681 gmx_residuetype_is_protein(rt, *(at->resinfo[resind].name)) )
1683 /* mark this residue */
1684 bResProcessed[resind] = TRUE;
1685 /* find out if this residue needs converting */
1687 for (j = 0; j < resNR && whatres == NOTSET; j++)
1690 cmplength = bPartial[j] ? strlen(resnm)-1 : strlen(resnm);
1692 bFound = ((gmx_strncasecmp(resnm, resnms[j], cmplength) == 0) ||
1693 (gmx_strncasecmp(resnm, resnmsN[j], cmplength) == 0) ||
1694 (gmx_strncasecmp(resnm, resnmsC[j], cmplength) == 0));
1699 /* get atoms we will be needing for the conversion */
1701 for (k = 0; atnms[j][k]; k++)
1704 for (m = i; m < at->nr && at->atom[m].resind == resind && ats[k] == NOTSET; m++)
1706 if (gmx_strcasecmp(*(at->atomname[m]), atnms[j][k]) == 0)
1714 /* now k is number of atom names in atnms[j] */
1723 if (nrfound < needed)
1725 gmx_fatal(FARGS, "not enough atoms found (%d, need %d) in "
1726 "residue %s %d while\n "
1727 "generating aromatics virtual site construction",
1728 nrfound, needed, resnm, at->resinfo[resind].nr);
1730 /* Advance overall atom counter */
1734 /* the enums for every residue MUST correspond to atnms[residue] */
1740 fprintf(stderr, "PHE at %d\n", o2n[ats[0]]+1);
1742 nvsite += gen_vsites_phe(at, vsite_type, plist, nrfound, ats, vsitetop, nvsitetop);
1747 fprintf(stderr, "TRP at %d\n", o2n[ats[0]]+1);
1749 nvsite += gen_vsites_trp(atype, &newx, &newatom, &newatomname, &o2n,
1750 &newvsite_type, &newcgnr, symtab, &nadd, *x, cgnr,
1751 at, vsite_type, plist, nrfound, ats, add_shift, vsitetop, nvsitetop);
1756 fprintf(stderr, "TYR at %d\n", o2n[ats[0]]+1);
1758 nvsite += gen_vsites_tyr(atype, &newx, &newatom, &newatomname, &o2n,
1759 &newvsite_type, &newcgnr, symtab, &nadd, *x, cgnr,
1760 at, vsite_type, plist, nrfound, ats, add_shift, vsitetop, nvsitetop);
1765 fprintf(stderr, "HIS at %d\n", o2n[ats[0]]+1);
1767 nvsite += gen_vsites_his(at, vsite_type, plist, nrfound, ats, vsitetop, nvsitetop);
1770 /* this means this residue won't be processed */
1773 gmx_fatal(FARGS, "DEATH HORROR in do_vsites (%s:%d)",
1774 __FILE__, __LINE__);
1775 } /* switch whatres */
1776 /* skip back to beginning of residue */
1777 while (i > 0 && at->atom[i-1].resind == resind)
1781 } /* if bVsiteAromatics & is protein */
1783 /* now process the rest of the hydrogens */
1784 /* only process hydrogen atoms which are not already set */
1785 if ( ((*vsite_type)[i] == NOTSET) && is_hydrogen(*(at->atomname[i])))
1787 /* find heavy atom, count #bonds from it and #H atoms bound to it
1788 and return H atom numbers (Hatoms) and heavy atom numbers (heavies) */
1789 count_bonds(i, &plist[F_BONDS], at->atomname,
1790 &nrbonds, &nrHatoms, Hatoms, &Heavy, &nrheavies, heavies);
1791 /* get Heavy atom type */
1792 tpHeavy = get_atype(Heavy, at, nrtp, rtp, rt);
1793 strcpy(tpname, get_atomtype_name(tpHeavy, atype));
1796 bAddVsiteParam = TRUE;
1797 /* nested if's which check nrHatoms, nrbonds and atomname */
1803 (*vsite_type)[i] = F_BONDS;
1805 case 3: /* =CH-, -NH- or =NH+- */
1806 (*vsite_type)[i] = F_VSITE3FD;
1808 case 4: /* --CH- (tert) */
1809 /* The old type 4FD had stability issues, so
1810 * all new constructs should use 4FDN
1812 (*vsite_type)[i] = F_VSITE4FDN;
1814 /* Check parity of heavy atoms from coordinates */
1819 rvec_sub((*x)[aj], (*x)[ai], tmpmat[0]);
1820 rvec_sub((*x)[ak], (*x)[ai], tmpmat[1]);
1821 rvec_sub((*x)[al], (*x)[ai], tmpmat[2]);
1823 if (det(tmpmat) > 0)
1831 default: /* nrbonds != 2, 3 or 4 */
1836 else if ( /*(nrHatoms == 2) && (nrbonds == 2) && REMOVED this test
1838 (gmx_strncasecmp(*at->atomname[Heavy], "OW", 2) == 0) )
1840 bAddVsiteParam = FALSE; /* this is water: skip these hydrogens */
1843 bFirstWater = FALSE;
1847 "Not converting hydrogens in water to virtual sites\n");
1851 else if ( (nrHatoms == 2) && (nrbonds == 4) )
1853 /* -CH2- , -NH2+- */
1854 (*vsite_type)[Hatoms[0]] = F_VSITE3OUT;
1855 (*vsite_type)[Hatoms[1]] = -F_VSITE3OUT;
1859 /* 2 or 3 hydrogen atom, with 3 or 4 bonds in total to the heavy atom.
1860 * If it is a nitrogen, first check if it is planar.
1862 isN = planarN = FALSE;
1863 if ((nrHatoms == 2) && ((*at->atomname[Heavy])[0] == 'N'))
1866 j = nitrogen_is_planar(vsiteconflist, nvsiteconf, tpname);
1869 gmx_fatal(FARGS, "No vsite database NH2 entry for type %s\n", tpname);
1873 if ( (nrHatoms == 2) && (nrbonds == 3) && ( !isN || planarN ) )
1875 /* =CH2 or, if it is a nitrogen NH2, it is a planar one */
1876 (*vsite_type)[Hatoms[0]] = F_VSITE3FAD;
1877 (*vsite_type)[Hatoms[1]] = -F_VSITE3FAD;
1879 else if ( ( (nrHatoms == 2) && (nrbonds == 3) &&
1880 ( isN && !planarN ) ) ||
1881 ( (nrHatoms == 3) && (nrbonds == 4) ) )
1883 /* CH3, NH3 or non-planar NH2 group */
1884 int Hat_vsite_type[3] = { F_VSITE3, F_VSITE3OUT, F_VSITE3OUT };
1885 gmx_bool Hat_SwapParity[3] = { FALSE, TRUE, FALSE };
1889 fprintf(stderr, "-XH3 or nonplanar NH2 group at %d\n", i+1);
1891 bAddVsiteParam = FALSE; /* we'll do this ourselves! */
1892 /* -NH2 (umbrella), -NH3+ or -CH3 */
1893 (*vsite_type)[Heavy] = F_VSITE3;
1894 for (j = 0; j < nrHatoms; j++)
1896 (*vsite_type)[Hatoms[j]] = Hat_vsite_type[j];
1898 /* get dummy mass type from first char of heavy atom type (N or C) */
1900 strcpy(nexttpname, get_atomtype_name(get_atype(heavies[0], at, nrtp, rtp, rt), atype));
1901 ch = get_dummymass_name(vsiteconflist, nvsiteconf, tpname, nexttpname);
1907 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);
1911 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);
1919 tpM = vsite_nm2type(name, atype);
1920 /* make space for 2 masses: shift all atoms starting with 'Heavy' */
1926 fprintf(stderr, "Inserting %d dummy masses at %d\n", NMASS, o2n[i0]+1);
1929 for (j = i0; j < at->nr; j++)
1934 srenew(newx, at->nr+nadd);
1935 srenew(newatom, at->nr+nadd);
1936 srenew(newatomname, at->nr+nadd);
1937 srenew(newvsite_type, at->nr+nadd);
1938 srenew(newcgnr, at->nr+nadd);
1940 for (j = 0; j < NMASS; j++)
1942 newatomname[at->nr+nadd-1-j] = NULL;
1945 /* calculate starting position for the masses */
1947 /* get atom masses, and set Heavy and Hatoms mass to zero */
1948 for (j = 0; j < nrHatoms; j++)
1950 mHtot += get_amass(Hatoms[j], at, nrtp, rtp, rt);
1951 at->atom[Hatoms[j]].m = at->atom[Hatoms[j]].mB = 0;
1953 mtot = mHtot + get_amass(Heavy, at, nrtp, rtp, rt);
1954 at->atom[Heavy].m = at->atom[Heavy].mB = 0;
1961 /* generate vectors parallel and perpendicular to rotational axis:
1962 * rpar = Heavy -> Hcom
1963 * rperp = Hcom -> H1 */
1965 for (j = 0; j < nrHatoms; j++)
1967 rvec_inc(rpar, (*x)[Hatoms[j]]);
1969 svmul(1.0/nrHatoms, rpar, rpar); /* rpar = ( H1+H2+H3 ) / 3 */
1970 rvec_dec(rpar, (*x)[Heavy]); /* - Heavy */
1971 rvec_sub((*x)[Hatoms[0]], (*x)[Heavy], rperp);
1972 rvec_dec(rperp, rpar); /* rperp = H1 - Heavy - rpar */
1973 /* calc mass positions */
1974 svmul(fact2, rpar, temp);
1975 for (j = 0; (j < NMASS); j++) /* xM = xN + fact2 * rpar +/- fact * rperp */
1977 rvec_add((*x)[Heavy], temp, newx[ni0+j]);
1979 svmul(fact, rperp, temp);
1980 rvec_inc(newx[ni0 ], temp);
1981 rvec_dec(newx[ni0+1], temp);
1982 /* set atom parameters for the masses */
1983 for (j = 0; (j < NMASS); j++)
1985 /* make name: "M??#" or "M?#" (? is atomname, # is number) */
1987 for (k = 0; (*at->atomname[Heavy])[k] && ( k < NMASS ); k++)
1989 name[k+1] = (*at->atomname[Heavy])[k];
1991 name[k+1] = atomnamesuffix[j];
1993 newatomname[ni0+j] = put_symtab(symtab, name);
1994 newatom[ni0+j].m = newatom[ni0+j].mB = mtot/NMASS;
1995 newatom[ni0+j].q = newatom[ni0+j].qB = 0.0;
1996 newatom[ni0+j].type = newatom[ni0+j].typeB = tpM;
1997 newatom[ni0+j].ptype = eptAtom;
1998 newatom[ni0+j].resind = at->atom[i0].resind;
1999 newatom[ni0+j].elem[0] = 'M';
2000 newatom[ni0+j].elem[1] = '\0';
2001 newvsite_type[ni0+j] = NOTSET;
2002 newcgnr[ni0+j] = (*cgnr)[i0];
2004 /* add constraints between dummy masses and to heavies[0] */
2005 /* 'add_shift' says which atoms won't be renumbered afterwards */
2006 my_add_param(&(plist[F_CONSTRNC]), heavies[0], add_shift+ni0, NOTSET);
2007 my_add_param(&(plist[F_CONSTRNC]), heavies[0], add_shift+ni0+1, NOTSET);
2008 my_add_param(&(plist[F_CONSTRNC]), add_shift+ni0, add_shift+ni0+1, NOTSET);
2010 /* generate Heavy, H1, H2 and H3 from M1, M2 and heavies[0] */
2011 /* note that vsite_type cannot be NOTSET, because we just set it */
2012 add_vsite3_atoms (&plist[(*vsite_type)[Heavy]],
2013 Heavy, heavies[0], add_shift+ni0, add_shift+ni0+1,
2015 for (j = 0; j < nrHatoms; j++)
2017 add_vsite3_atoms(&plist[(*vsite_type)[Hatoms[j]]],
2018 Hatoms[j], heavies[0], add_shift+ni0, add_shift+ni0+1,
2032 "Warning: cannot convert atom %d %s (bound to a heavy atom "
2034 " %d bonds and %d bound hydrogens atoms) to virtual site\n",
2035 i+1, *(at->atomname[i]), tpname, nrbonds, nrHatoms);
2039 /* add vsite parameters to topology,
2040 also get rid of negative vsite_types */
2041 add_vsites(plist, (*vsite_type), Heavy, nrHatoms, Hatoms,
2042 nrheavies, heavies);
2043 /* transfer mass of virtual site to Heavy atom */
2044 for (j = 0; j < nrHatoms; j++)
2046 if (is_vsite((*vsite_type)[Hatoms[j]]))
2048 at->atom[Heavy].m += at->atom[Hatoms[j]].m;
2049 at->atom[Heavy].mB = at->atom[Heavy].m;
2050 at->atom[Hatoms[j]].m = at->atom[Hatoms[j]].mB = 0;
2057 fprintf(debug, "atom %d: ", o2n[i]+1);
2058 print_bonds(debug, o2n, nrHatoms, Hatoms, Heavy, nrheavies, heavies);
2060 } /* if vsite NOTSET & is hydrogen */
2062 } /* for i < at->nr */
2064 gmx_residuetype_destroy(rt);
2068 fprintf(debug, "Before inserting new atoms:\n");
2069 for (i = 0; i < at->nr; i++)
2071 fprintf(debug, "%4d %4d %4s %4d %4s %6d %-10s\n", i+1, o2n[i]+1,
2072 at->atomname[i] ? *(at->atomname[i]) : "(NULL)",
2073 at->resinfo[at->atom[i].resind].nr,
2074 at->resinfo[at->atom[i].resind].name ?
2075 *(at->resinfo[at->atom[i].resind].name) : "(NULL)",
2077 ((*vsite_type)[i] == NOTSET) ?
2078 "NOTSET" : interaction_function[(*vsite_type)[i]].name);
2080 fprintf(debug, "new atoms to be inserted:\n");
2081 for (i = 0; i < at->nr+nadd; i++)
2085 fprintf(debug, "%4d %4s %4d %6d %-10s\n", i+1,
2086 newatomname[i] ? *(newatomname[i]) : "(NULL)",
2087 newatom[i].resind, newcgnr[i],
2088 (newvsite_type[i] == NOTSET) ?
2089 "NOTSET" : interaction_function[newvsite_type[i]].name);
2094 /* add all original atoms to the new arrays, using o2n index array */
2095 for (i = 0; i < at->nr; i++)
2097 newatomname [o2n[i]] = at->atomname [i];
2098 newatom [o2n[i]] = at->atom [i];
2099 newvsite_type[o2n[i]] = (*vsite_type)[i];
2100 newcgnr [o2n[i]] = (*cgnr) [i];
2101 copy_rvec((*x)[i], newx[o2n[i]]);
2103 /* throw away old atoms */
2105 sfree(at->atomname);
2109 /* put in the new ones */
2112 at->atomname = newatomname;
2113 *vsite_type = newvsite_type;
2116 if (at->nr > add_shift)
2118 gmx_fatal(FARGS, "Added impossible amount of dummy masses "
2119 "(%d on a total of %d atoms)\n", nadd, at->nr-nadd);
2124 fprintf(debug, "After inserting new atoms:\n");
2125 for (i = 0; i < at->nr; i++)
2127 fprintf(debug, "%4d %4s %4d %4s %6d %-10s\n", i+1,
2128 at->atomname[i] ? *(at->atomname[i]) : "(NULL)",
2129 at->resinfo[at->atom[i].resind].nr,
2130 at->resinfo[at->atom[i].resind].name ?
2131 *(at->resinfo[at->atom[i].resind].name) : "(NULL)",
2133 ((*vsite_type)[i] == NOTSET) ?
2134 "NOTSET" : interaction_function[(*vsite_type)[i]].name);
2138 /* now renumber all the interactions because of the added atoms */
2139 for (ftype = 0; ftype < F_NRE; ftype++)
2141 params = &(plist[ftype]);
2144 fprintf(debug, "Renumbering %d %s\n", params->nr,
2145 interaction_function[ftype].longname);
2147 for (i = 0; i < params->nr; i++)
2149 for (j = 0; j < NRAL(ftype); j++)
2151 if (params->param[i].a[j] >= add_shift)
2155 fprintf(debug, " [%u -> %u]", params->param[i].a[j],
2156 params->param[i].a[j]-add_shift);
2158 params->param[i].a[j] = params->param[i].a[j]-add_shift;
2164 fprintf(debug, " [%u -> %d]", params->param[i].a[j],
2165 o2n[params->param[i].a[j]]);
2167 params->param[i].a[j] = o2n[params->param[i].a[j]];
2172 fprintf(debug, "\n");
2176 /* now check if atoms in the added constraints are in increasing order */
2177 params = &(plist[F_CONSTRNC]);
2178 for (i = 0; i < params->nr; i++)
2180 if (params->param[i].AI > params->param[i].AJ)
2182 j = params->param[i].AJ;
2183 params->param[i].AJ = params->param[i].AI;
2184 params->param[i].AI = j;
2191 /* tell the user what we did */
2192 fprintf(stderr, "Marked %d virtual sites\n", nvsite);
2193 fprintf(stderr, "Added %d dummy masses\n", nadd);
2194 fprintf(stderr, "Added %d new constraints\n", plist[F_CONSTRNC].nr);
2197 void do_h_mass(t_params *psb, int vsite_type[], t_atoms *at, real mHmult,
2198 gmx_bool bDeuterate)
2202 /* loop over all atoms */
2203 for (i = 0; i < at->nr; i++)
2205 /* adjust masses if i is hydrogen and not a virtual site */
2206 if (!is_vsite(vsite_type[i]) && is_hydrogen(*(at->atomname[i])) )
2208 /* find bonded heavy atom */
2210 for (j = 0; (j < psb->nr) && (a == NOTSET); j++)
2212 /* if other atom is not a virtual site, it is the one we want */
2213 if ( (psb->param[j].AI == i) &&
2214 !is_vsite(vsite_type[psb->param[j].AJ]) )
2216 a = psb->param[j].AJ;
2218 else if ( (psb->param[j].AJ == i) &&
2219 !is_vsite(vsite_type[psb->param[j].AI]) )
2221 a = psb->param[j].AI;
2226 gmx_fatal(FARGS, "Unbound hydrogen atom (%d) found while adjusting mass",
2230 /* adjust mass of i (hydrogen) with mHmult
2231 and correct mass of a (bonded atom) with same amount */
2234 at->atom[a].m -= (mHmult-1.0)*at->atom[i].m;
2235 at->atom[a].mB -= (mHmult-1.0)*at->atom[i].m;
2237 at->atom[i].m *= mHmult;
2238 at->atom[i].mB *= mHmult;