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,2015,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.
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.
39 #include "vsite_parm.h"
48 #include "gromacs/gmxpreprocess/add_par.h"
49 #include "gromacs/gmxpreprocess/notset.h"
50 #include "gromacs/gmxpreprocess/resall.h"
51 #include "gromacs/gmxpreprocess/toputil.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdtypes/md_enums.h"
56 #include "gromacs/topology/ifunc.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/gmxassert.h"
61 #include "gromacs/utility/smalloc.h"
62 #include "gromacs/utility/strconvert.h"
67 t_iatom &ai() { return a[0]; }
68 t_iatom &aj() { return a[1]; }
69 t_iatom &ak() { return a[2]; }
70 t_iatom &al() { return a[3]; }
81 vsitebondparam_t *vsbp;
89 static int vsite_bond_nrcheck(int ftype)
93 if ((interaction_function[ftype].flags & (IF_BTYPE | IF_CONSTRAINT | IF_ATYPE)) || (ftype == F_IDIHS))
95 nrcheck = NRAL(ftype);
105 static void enter_bonded(int nratoms, int *nrbonded, t_mybonded **bondeds,
110 srenew(*bondeds, *nrbonded+1);
112 /* copy atom numbers */
113 for (j = 0; j < nratoms; j++)
115 (*bondeds)[*nrbonded].a[j] = param->a[j];
118 (*bondeds)[*nrbonded].c = param->c0();
123 static void get_bondeds(int nrat, const t_iatom atoms[],
124 at2vsitebond_t *at2vb,
125 int *nrbond, t_mybonded **bonds,
126 int *nrang, t_mybonded **angles,
127 int *nridih, t_mybonded **idihs )
129 int k, i, ftype, nrcheck;
132 for (k = 0; k < nrat; k++)
134 for (i = 0; i < at2vb[atoms[k]].nr; i++)
136 ftype = at2vb[atoms[k]].vsbp[i].ftype;
137 param = at2vb[atoms[k]].vsbp[i].param;
138 nrcheck = vsite_bond_nrcheck(ftype);
139 /* abuse nrcheck to see if we're adding bond, angle or idih */
142 case 2: enter_bonded(nrcheck, nrbond, bonds, param); break;
143 case 3: enter_bonded(nrcheck, nrang, angles, param); break;
144 case 4: enter_bonded(nrcheck, nridih, idihs, param); break;
150 static at2vsitebond_t *make_at2vsitebond(int natoms, t_params plist[])
153 int ftype, i, j, nrcheck, nr;
155 at2vsitebond_t *at2vb;
160 for (ftype = 0; (ftype < F_NRE); ftype++)
162 if ((interaction_function[ftype].flags & IF_VSITE) && ftype != F_VSITEN)
164 for (i = 0; (i < plist[ftype].nr); i++)
166 for (j = 0; j < NRAL(ftype); j++)
168 bVSI[plist[ftype].param[i].a[j]] = TRUE;
174 for (ftype = 0; (ftype < F_NRE); ftype++)
176 nrcheck = vsite_bond_nrcheck(ftype);
179 for (i = 0; (i < plist[ftype].nr); i++)
181 aa = plist[ftype].param[i].a;
182 for (j = 0; j < nrcheck; j++)
186 nr = at2vb[aa[j]].nr;
189 srenew(at2vb[aa[j]].vsbp, nr+10);
191 at2vb[aa[j]].vsbp[nr].ftype = ftype;
192 at2vb[aa[j]].vsbp[nr].param = &plist[ftype].param[i];
204 static void done_at2vsitebond(int natoms, at2vsitebond_t *at2vb)
208 for (i = 0; i < natoms; i++)
212 sfree(at2vb[i].vsbp);
218 static at2vsitecon_t *make_at2vsitecon(int natoms, t_params plist[])
221 int ftype, i, j, ai, aj, nr;
222 at2vsitecon_t *at2vc;
227 for (ftype = 0; (ftype < F_NRE); ftype++)
229 if ((interaction_function[ftype].flags & IF_VSITE) && ftype != F_VSITEN)
231 for (i = 0; (i < plist[ftype].nr); i++)
233 for (j = 0; j < NRAL(ftype); j++)
235 bVSI[plist[ftype].param[i].a[j]] = TRUE;
241 for (ftype = 0; (ftype < F_NRE); ftype++)
243 if (interaction_function[ftype].flags & IF_CONSTRAINT)
245 for (i = 0; (i < plist[ftype].nr); i++)
247 ai = plist[ftype].param[i].ai();
248 aj = plist[ftype].param[i].aj();
249 if (bVSI[ai] && bVSI[aj])
251 /* Store forward direction */
255 srenew(at2vc[ai].aj, nr+10);
257 at2vc[ai].aj[nr] = aj;
259 /* Store backward direction */
263 srenew(at2vc[aj].aj, nr+10);
265 at2vc[aj].aj[nr] = ai;
276 static void done_at2vsitecon(int natoms, at2vsitecon_t *at2vc)
280 for (i = 0; i < natoms; i++)
291 static void print_bad(FILE *fp,
292 int nrbond, t_mybonded *bonds,
293 int nrang, t_mybonded *angles,
294 int nridih, t_mybonded *idihs )
300 fprintf(fp, "bonds:");
301 for (i = 0; i < nrbond; i++)
303 fprintf(fp, " %d-%d (%g)",
304 bonds[i].ai()+1, bonds[i].aj()+1, bonds[i].c);
310 fprintf(fp, "angles:");
311 for (i = 0; i < nrang; i++)
313 fprintf(fp, " %d-%d-%d (%g)",
314 angles[i].ai()+1, angles[i].aj()+1,
315 angles[i].ak()+1, angles[i].c);
321 fprintf(fp, "idihs:");
322 for (i = 0; i < nridih; i++)
324 fprintf(fp, " %d-%d-%d-%d (%g)",
325 idihs[i].ai()+1, idihs[i].aj()+1,
326 idihs[i].ak()+1, idihs[i].al()+1, idihs[i].c);
332 static void print_param(FILE *fp, int ftype, int i, t_param *param)
335 static int prev_ftype = NOTSET;
336 static int prev_i = NOTSET;
339 if ( (ftype != prev_ftype) || (i != prev_i) )
345 fprintf(fp, "(%d) plist[%s].param[%d]",
346 pass, interaction_function[ftype].name, i);
347 for (j = 0; j < NRFP(ftype); j++)
349 fprintf(fp, ".c[%d]=%g ", j, param->c[j]);
355 static real get_bond_length(int nrbond, t_mybonded bonds[],
356 t_iatom ai, t_iatom aj)
362 for (i = 0; i < nrbond && (bondlen == NOTSET); i++)
364 /* check both ways */
365 if ( ( (ai == bonds[i].ai()) && (aj == bonds[i].aj()) ) ||
366 ( (ai == bonds[i].aj()) && (aj == bonds[i].ai()) ) )
368 bondlen = bonds[i].c; /* note: bonds[i].c might be NOTSET */
374 static real get_angle(int nrang, t_mybonded angles[],
375 t_iatom ai, t_iatom aj, t_iatom ak)
381 for (i = 0; i < nrang && (angle == NOTSET); i++)
383 /* check both ways */
384 if ( ( (ai == angles[i].ai()) && (aj == angles[i].aj()) && (ak == angles[i].ak()) ) ||
385 ( (ai == angles[i].ak()) && (aj == angles[i].aj()) && (ak == angles[i].ai()) ) )
387 angle = DEG2RAD*angles[i].c;
393 static char *get_atomtype_name_AB(t_atom *atom, gpp_atomtype_t atype)
397 name = get_atomtype_name(atom->type, atype);
399 /* When using the decoupling option, atom types are changed
400 * to decoupled for the non-bonded interactions, but the virtual
401 * sites constructions should be based on the "normal" interactions.
402 * So we return the state B atom type name if the state A atom
403 * type is the decoupled one. We should actually check for the atom
404 * type number, but that's not passed here. So we check for
405 * the decoupled atom type name. This should not cause trouble
406 * as this code is only used for topologies with v-sites without
407 * parameters generated by pdb2gmx.
409 if (strcmp(name, "decoupled") == 0)
411 name = get_atomtype_name(atom->typeB, atype);
417 static bool calc_vsite3_param(gpp_atomtype_t atype,
418 t_param *param, t_atoms *at,
419 int nrbond, t_mybonded *bonds,
420 int nrang, t_mybonded *angles )
422 /* i = virtual site | ,k
423 * j = 1st bonded heavy atom | i-j
424 * k,l = 2nd bonded atoms | `l
428 real bjk, bjl, a = -1, b = -1;
429 /* check if this is part of a NH3 , NH2-umbrella or CH3 group,
430 * i.e. if atom k and l are dummy masses (MNH* or MCH3*) */
432 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->ak()], atype), "MNH", 3) == 0) &&
433 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->al()], atype), "MNH", 3) == 0) ) ||
434 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->ak()], atype), "MCH3", 4) == 0) &&
435 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->al()], atype), "MCH3", 4) == 0) );
437 bjk = get_bond_length(nrbond, bonds, param->aj(), param->ak());
438 bjl = get_bond_length(nrbond, bonds, param->aj(), param->al());
439 bError = (bjk == NOTSET) || (bjl == NOTSET);
442 /* now we get some XH2/XH3 group specific construction */
443 /* note: we call the heavy atom 'C' and the X atom 'N' */
444 real bMM, bCM, bCN, bNH, aCNH, dH, rH, dM, rM;
447 /* check if bonds from heavy atom (j) to dummy masses (k,l) are equal: */
448 bError = bError || (bjk != bjl);
450 /* the X atom (C or N) in the XH2/XH3 group is the first after the masses: */
451 aN = std::max(param->ak(), param->al())+1;
453 /* get common bonds */
454 bMM = get_bond_length(nrbond, bonds, param->ak(), param->al());
456 bCN = get_bond_length(nrbond, bonds, param->aj(), aN);
457 bError = bError || (bMM == NOTSET) || (bCN == NOTSET);
459 /* calculate common things */
461 dM = std::sqrt( gmx::square(bCM) - gmx::square(rM) );
463 /* are we dealing with the X atom? */
464 if (param->ai() == aN)
466 /* this is trivial */
467 a = b = 0.5 * bCN/dM;
472 /* get other bondlengths and angles: */
473 bNH = get_bond_length(nrbond, bonds, aN, param->ai());
474 aCNH = get_angle (nrang, angles, param->aj(), aN, param->ai());
475 bError = bError || (bNH == NOTSET) || (aCNH == NOTSET);
478 dH = bCN - bNH * std::cos(aCNH);
479 rH = bNH * std::sin(aCNH);
481 a = 0.5 * ( dH/dM + rH/rM );
482 b = 0.5 * ( dH/dM - rH/rM );
487 gmx_fatal(FARGS, "calc_vsite3_param not implemented for the general case "
488 "(atom %d)", param->ai()+1);
497 static bool calc_vsite3fd_param(t_param *param,
498 int nrbond, t_mybonded *bonds,
499 int nrang, t_mybonded *angles)
501 /* i = virtual site | ,k
502 * j = 1st bonded heavy atom | i-j
503 * k,l = 2nd bonded atoms | `l
507 real bij, bjk, bjl, aijk, aijl, rk, rl;
509 bij = get_bond_length(nrbond, bonds, param->ai(), param->aj());
510 bjk = get_bond_length(nrbond, bonds, param->aj(), param->ak());
511 bjl = get_bond_length(nrbond, bonds, param->aj(), param->al());
512 aijk = get_angle (nrang, angles, param->ai(), param->aj(), param->ak());
513 aijl = get_angle (nrang, angles, param->ai(), param->aj(), param->al());
514 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) ||
515 (aijk == NOTSET) || (aijl == NOTSET);
517 rk = bjk * std::sin(aijk);
518 rl = bjl * std::sin(aijl);
519 param->c0() = rk / (rk + rl);
520 param->c1() = -bij; /* 'bond'-length for fixed distance vsite */
525 static bool calc_vsite3fad_param(t_param *param,
526 int nrbond, t_mybonded *bonds,
527 int nrang, t_mybonded *angles)
529 /* i = virtual site |
530 * j = 1st bonded heavy atom | i-j
531 * k = 2nd bonded heavy atom | `k-l
532 * l = 3d bonded heavy atom |
535 bool bSwapParity, bError;
538 bSwapParity = ( param->c1() == -1 );
540 bij = get_bond_length(nrbond, bonds, param->ai(), param->aj());
541 aijk = get_angle (nrang, angles, param->ai(), param->aj(), param->ak());
542 bError = (bij == NOTSET) || (aijk == NOTSET);
544 param->c1() = bij; /* 'bond'-length for fixed distance vsite */
545 param->c0() = RAD2DEG*aijk; /* 'bond'-angle for fixed angle vsite */
549 param->c0() = 360 - param->c0();
555 static bool calc_vsite3out_param(gpp_atomtype_t atype,
556 t_param *param, t_atoms *at,
557 int nrbond, t_mybonded *bonds,
558 int nrang, t_mybonded *angles)
560 /* i = virtual site | ,k
561 * j = 1st bonded heavy atom | i-j
562 * k,l = 2nd bonded atoms | `l
563 * NOTE: i is out of the j-k-l plane!
566 bool bXH3, bError, bSwapParity;
567 real bij, bjk, bjl, aijk, aijl, akjl, pijk, pijl, a, b, c;
569 /* check if this is part of a NH2-umbrella, NH3 or CH3 group,
570 * i.e. if atom k and l are dummy masses (MNH* or MCH3*) */
572 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->ak()], atype), "MNH", 3) == 0) &&
573 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->al()], atype), "MNH", 3) == 0) ) ||
574 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->ak()], atype), "MCH3", 4) == 0) &&
575 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->al()], atype), "MCH3", 4) == 0) );
577 /* check if construction parity must be swapped */
578 bSwapParity = ( param->c1() == -1 );
580 bjk = get_bond_length(nrbond, bonds, param->aj(), param->ak());
581 bjl = get_bond_length(nrbond, bonds, param->aj(), param->al());
582 bError = (bjk == NOTSET) || (bjl == NOTSET);
585 /* now we get some XH3 group specific construction */
586 /* note: we call the heavy atom 'C' and the X atom 'N' */
587 real bMM, bCM, bCN, bNH, aCNH, dH, rH, rHx, rHy, dM, rM;
590 /* check if bonds from heavy atom (j) to dummy masses (k,l) are equal: */
591 bError = bError || (bjk != bjl);
593 /* the X atom (C or N) in the XH3 group is the first after the masses: */
594 aN = std::max(param->ak(), param->al())+1;
596 /* get all bondlengths and angles: */
597 bMM = get_bond_length(nrbond, bonds, param->ak(), param->al());
599 bCN = get_bond_length(nrbond, bonds, param->aj(), aN);
600 bNH = get_bond_length(nrbond, bonds, aN, param->ai());
601 aCNH = get_angle (nrang, angles, param->aj(), aN, param->ai());
603 (bMM == NOTSET) || (bCN == NOTSET) || (bNH == NOTSET) || (aCNH == NOTSET);
606 dH = bCN - bNH * std::cos(aCNH);
607 rH = bNH * std::sin(aCNH);
608 /* we assume the H's are symmetrically distributed */
609 rHx = rH*std::cos(DEG2RAD*30);
610 rHy = rH*std::sin(DEG2RAD*30);
612 dM = std::sqrt( gmx::square(bCM) - gmx::square(rM) );
613 a = 0.5*( (dH/dM) - (rHy/rM) );
614 b = 0.5*( (dH/dM) + (rHy/rM) );
620 /* this is the general construction */
622 bij = get_bond_length(nrbond, bonds, param->ai(), param->aj());
623 aijk = get_angle (nrang, angles, param->ai(), param->aj(), param->ak());
624 aijl = get_angle (nrang, angles, param->ai(), param->aj(), param->al());
625 akjl = get_angle (nrang, angles, param->ak(), param->aj(), param->al());
627 (bij == NOTSET) || (aijk == NOTSET) || (aijl == NOTSET) || (akjl == NOTSET);
629 pijk = std::cos(aijk)*bij;
630 pijl = std::cos(aijl)*bij;
631 a = ( pijk + (pijk*std::cos(akjl)-pijl) * std::cos(akjl) / gmx::square(std::sin(akjl)) ) / bjk;
632 b = ( pijl + (pijl*std::cos(akjl)-pijk) * std::cos(akjl) / gmx::square(std::sin(akjl)) ) / bjl;
633 c = -std::sqrt( gmx::square(bij) -
634 ( gmx::square(pijk) - 2*pijk*pijl*std::cos(akjl) + gmx::square(pijl) )
635 / gmx::square(std::sin(akjl)) )
636 / ( bjk*bjl*std::sin(akjl) );
652 static bool calc_vsite4fd_param(t_param *param,
653 int nrbond, t_mybonded *bonds,
654 int nrang, t_mybonded *angles)
656 /* i = virtual site | ,k
657 * j = 1st bonded heavy atom | i-j-m
658 * k,l,m = 2nd bonded atoms | `l
662 real bij, bjk, bjl, bjm, aijk, aijl, aijm, akjm, akjl;
663 real pk, pl, pm, cosakl, cosakm, sinakl, sinakm, cl, cm;
665 bij = get_bond_length(nrbond, bonds, param->ai(), param->aj());
666 bjk = get_bond_length(nrbond, bonds, param->aj(), param->ak());
667 bjl = get_bond_length(nrbond, bonds, param->aj(), param->al());
668 bjm = get_bond_length(nrbond, bonds, param->aj(), param->am());
669 aijk = get_angle (nrang, angles, param->ai(), param->aj(), param->ak());
670 aijl = get_angle (nrang, angles, param->ai(), param->aj(), param->al());
671 aijm = get_angle (nrang, angles, param->ai(), param->aj(), param->am());
672 akjm = get_angle (nrang, angles, param->ak(), param->aj(), param->am());
673 akjl = get_angle (nrang, angles, param->ak(), param->aj(), param->al());
674 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) || (bjm == NOTSET) ||
675 (aijk == NOTSET) || (aijl == NOTSET) || (aijm == NOTSET) || (akjm == NOTSET) ||
680 pk = bjk*std::sin(aijk);
681 pl = bjl*std::sin(aijl);
682 pm = bjm*std::sin(aijm);
683 cosakl = (std::cos(akjl) - std::cos(aijk)*std::cos(aijl)) / (std::sin(aijk)*std::sin(aijl));
684 cosakm = (std::cos(akjm) - std::cos(aijk)*std::cos(aijm)) / (std::sin(aijk)*std::sin(aijm));
685 if (cosakl < -1 || cosakl > 1 || cosakm < -1 || cosakm > 1)
687 fprintf(stderr, "virtual site %d: angle ijk = %f, angle ijl = %f, angle ijm = %f\n",
688 param->ai()+1, RAD2DEG*aijk, RAD2DEG*aijl, RAD2DEG*aijm);
689 gmx_fatal(FARGS, "invalid construction in calc_vsite4fd for atom %d: "
690 "cosakl=%f, cosakm=%f\n", param->ai()+1, cosakl, cosakm);
692 sinakl = std::sqrt(1-gmx::square(cosakl));
693 sinakm = std::sqrt(1-gmx::square(cosakm));
695 /* note: there is a '+' because of the way the sines are calculated */
696 cl = -pk / ( pl*cosakl - pk + pl*sinakl*(pm*cosakm-pk)/(pm*sinakm) );
697 cm = -pk / ( pm*cosakm - pk + pm*sinakm*(pl*cosakl-pk)/(pl*sinakl) );
709 calc_vsite4fdn_param(t_param *param,
710 int nrbond, t_mybonded *bonds,
711 int nrang, t_mybonded *angles)
713 /* i = virtual site | ,k
714 * j = 1st bonded heavy atom | i-j-m
715 * k,l,m = 2nd bonded atoms | `l
719 real bij, bjk, bjl, bjm, aijk, aijl, aijm;
720 real pk, pl, pm, a, b;
722 bij = get_bond_length(nrbond, bonds, param->ai(), param->aj());
723 bjk = get_bond_length(nrbond, bonds, param->aj(), param->ak());
724 bjl = get_bond_length(nrbond, bonds, param->aj(), param->al());
725 bjm = get_bond_length(nrbond, bonds, param->aj(), param->am());
726 aijk = get_angle (nrang, angles, param->ai(), param->aj(), param->ak());
727 aijl = get_angle (nrang, angles, param->ai(), param->aj(), param->al());
728 aijm = get_angle (nrang, angles, param->ai(), param->aj(), param->am());
730 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) || (bjm == NOTSET) ||
731 (aijk == NOTSET) || (aijl == NOTSET) || (aijm == NOTSET);
736 /* Calculate component of bond j-k along the direction i-j */
737 pk = -bjk*std::cos(aijk);
739 /* Calculate component of bond j-l along the direction i-j */
740 pl = -bjl*std::cos(aijl);
742 /* Calculate component of bond j-m along the direction i-j */
743 pm = -bjm*std::cos(aijm);
745 if (fabs(pl) < 1000*GMX_REAL_MIN || fabs(pm) < 1000*GMX_REAL_MIN)
747 fprintf(stderr, "virtual site %d: angle ijk = %f, angle ijl = %f, angle ijm = %f\n",
748 param->ai()+1, RAD2DEG*aijk, RAD2DEG*aijl, RAD2DEG*aijm);
749 gmx_fatal(FARGS, "invalid construction in calc_vsite4fdn for atom %d: "
750 "pl=%f, pm=%f\n", param->ai()+1, pl, pm);
767 int set_vsites(bool bVerbose, t_atoms *atoms, gpp_atomtype_t atype,
771 int nvsite, nrbond, nrang, nridih, nrset;
772 bool bFirst, bSet, bERROR;
773 at2vsitebond_t *at2vb;
781 /* Make a reverse list to avoid ninteractions^2 operations */
782 at2vb = make_at2vsitebond(atoms->nr, plist);
784 for (ftype = 0; (ftype < F_NRE); ftype++)
786 if (interaction_function[ftype].flags & IF_VSITE)
788 nvsite += plist[ftype].nr;
790 if (ftype == F_VSITEN)
792 /* We don't calculate parameters for VSITEN */
797 for (i = 0; (i < plist[ftype].nr); i++)
799 /* check if all parameters are set */
801 for (j = 0; j < NRFP(ftype) && bSet; j++)
803 bSet = plist[ftype].param[i].c[j] != NOTSET;
808 fprintf(debug, "bSet=%s ", gmx::boolToString(bSet));
809 print_param(debug, ftype, i, &plist[ftype].param[i]);
813 if (bVerbose && bFirst)
815 fprintf(stderr, "Calculating parameters for virtual sites\n");
819 nrbond = nrang = nridih = 0;
824 /* now set the vsite parameters: */
825 get_bondeds(NRAL(ftype), plist[ftype].param[i].a, at2vb,
826 &nrbond, &bonds, &nrang, &angles, &nridih, &idihs);
829 fprintf(debug, "Found %d bonds, %d angles and %d idihs "
830 "for virtual site %d (%s)\n", nrbond, nrang, nridih,
831 plist[ftype].param[i].ai()+1,
832 interaction_function[ftype].longname);
833 print_bad(debug, nrbond, bonds, nrang, angles, nridih, idihs);
839 calc_vsite3_param(atype, &(plist[ftype].param[i]), atoms,
840 nrbond, bonds, nrang, angles);
844 calc_vsite3fd_param(&(plist[ftype].param[i]),
845 nrbond, bonds, nrang, angles);
849 calc_vsite3fad_param(&(plist[ftype].param[i]),
850 nrbond, bonds, nrang, angles);
854 calc_vsite3out_param(atype, &(plist[ftype].param[i]), atoms,
855 nrbond, bonds, nrang, angles);
859 calc_vsite4fd_param(&(plist[ftype].param[i]),
860 nrbond, bonds, nrang, angles);
864 calc_vsite4fdn_param(&(plist[ftype].param[i]),
865 nrbond, bonds, nrang, angles);
868 gmx_fatal(FARGS, "Automatic parameter generation not supported "
870 interaction_function[ftype].longname,
871 plist[ftype].param[i].ai()+1);
876 gmx_fatal(FARGS, "Automatic parameter generation not supported "
877 "for %s atom %d for this bonding configuration",
878 interaction_function[ftype].longname,
879 plist[ftype].param[i].ai()+1);
889 done_at2vsitebond(atoms->nr, at2vb);
894 void set_vsites_ptype(bool bVerbose, gmx_moltype_t *molt)
903 fprintf(stderr, "Setting particle type to V for virtual sites\n");
905 for (ftype = 0; ftype < F_NRE; ftype++)
907 il = &molt->ilist[ftype];
908 if (interaction_function[ftype].flags & IF_VSITE)
910 nra = interaction_function[ftype].nratoms;
916 fprintf(stderr, "doing %d %s virtual sites\n",
917 (nrd / (nra+1)), interaction_function[ftype].longname);
920 for (i = 0; (i < nrd); )
922 /* The virtual site */
924 molt->atoms.atom[avsite].ptype = eptVSite;
938 static void check_vsite_constraints(t_params *plist,
939 int cftype, const int vsite_type[])
946 ps = &(plist[cftype]);
947 for (i = 0; (i < ps->nr); i++)
949 for (k = 0; k < 2; k++)
951 atom = ps->param[i].a[k];
952 if (vsite_type[atom] != NOTSET)
954 fprintf(stderr, "ERROR: Cannot have constraint (%d-%d) with virtual site (%d)\n",
955 ps->param[i].ai()+1, ps->param[i].aj()+1, atom+1);
962 gmx_fatal(FARGS, "There were %d virtual sites involved in constraints", n);
966 static void clean_vsite_bonds(t_params *plist, t_pindex pindex[],
967 int cftype, const int vsite_type[])
969 int ftype, i, j, k, m, n, nvsite, nOut, kept_i;
970 int nconverted, nremoved;
971 int atom, oatom, at1, at2;
972 bool bKeep, bRemove, bUsed, bPresent, bThisFD, bThisOUT, bAllFD, bFirstTwo;
975 if (cftype == F_CONNBONDS)
980 ps = &(plist[cftype]);
985 for (i = 0; (i < ps->nr); i++) /* for all bonds in the plist */
988 const int *first_atoms = nullptr;
993 /* check if all virtual sites are constructed from the same atoms */
995 for (k = 0; (k < 2) && !bKeep && !bRemove; k++)
997 /* for all atoms in the bond */
998 atom = ps->param[i].a[k];
999 if (vsite_type[atom] != NOTSET && vsite_type[atom] != F_VSITEN)
1002 bThisFD = ( (pindex[atom].ftype == F_VSITE3FD ) ||
1003 (pindex[atom].ftype == F_VSITE3FAD) ||
1004 (pindex[atom].ftype == F_VSITE4FD ) ||
1005 (pindex[atom].ftype == F_VSITE4FDN ) );
1006 bThisOUT = ( (pindex[atom].ftype == F_VSITE3OUT) &&
1007 (interaction_function[cftype].flags & IF_CONSTRAINT) );
1008 bAllFD = bAllFD && bThisFD;
1009 if (bThisFD || bThisOUT)
1011 oatom = ps->param[i].a[1-k]; /* the other atom */
1012 if (vsite_type[oatom] == NOTSET &&
1013 oatom == plist[pindex[atom].ftype].param[pindex[atom].parnr].aj())
1015 /* if the other atom isn't a vsite, and it is AI */
1025 /* TODO This fragment, and corresponding logic in
1026 clean_vsite_angles and clean_vsite_dihs should
1027 be refactored into a common function */
1030 /* if this is the first vsite we encounter then
1031 store construction atoms */
1032 /* TODO This would be nicer to implement with
1033 a C++ "vector view" class" with an
1034 STL-container-like interface. */
1035 vsnral = NRAL(pindex[atom].ftype) - 1;
1036 first_atoms = plist[pindex[atom].ftype].param[pindex[atom].parnr].a + 1;
1040 GMX_ASSERT(vsnral != 0, "nvsite > 1 must have vsnral != 0");
1041 GMX_ASSERT(first_atoms != nullptr, "nvsite > 1 must have first_atoms != NULL");
1042 /* if it is not the first then
1043 check if this vsite is constructed from the same atoms */
1044 if (vsnral == NRAL(pindex[atom].ftype)-1)
1046 for (m = 0; (m < vsnral) && !bKeep; m++)
1051 atoms = plist[pindex[atom].ftype].param[pindex[atom].parnr].a + 1;
1052 for (n = 0; (n < vsnral) && !bPresent; n++)
1054 if (atoms[m] == first_atoms[n])
1080 /* if we have no virtual sites in this bond, keep it */
1086 /* TODO This loop and the corresponding loop in
1087 check_vsite_angles should be refactored into a common
1089 /* check if all non-vsite atoms are used in construction: */
1091 for (k = 0; (k < 2) && !bKeep; k++) /* for all atoms in the bond */
1093 atom = ps->param[i].a[k];
1094 if (vsite_type[atom] == NOTSET)
1097 for (m = 0; (m < vsnral) && !bUsed; m++)
1099 GMX_ASSERT(first_atoms != nullptr, "If we've seen a vsite before, we know what its first atom index was");
1101 if (atom == first_atoms[m])
1104 bFirstTwo = bFirstTwo && m < 2;
1114 if (!( bAllFD && bFirstTwo ) )
1116 /* Two atom bonded interactions include constraints.
1117 * We need to remove constraints between vsite pairs that have
1118 * a fixed distance due to being constructed from the same
1119 * atoms, since this can be numerically unstable.
1121 for (m = 0; m < vsnral && !bKeep; m++) /* all constr. atoms */
1123 at1 = first_atoms[m];
1124 at2 = first_atoms[(m+1) % vsnral];
1126 for (ftype = 0; ftype < F_NRE; ftype++)
1128 if (interaction_function[ftype].flags & IF_CONSTRAINT)
1130 for (j = 0; (j < plist[ftype].nr) && !bPresent; j++)
1132 /* all constraints until one matches */
1133 bPresent = ( ( (plist[ftype].param[j].ai() == at1) &&
1134 (plist[ftype].param[j].aj() == at2) ) ||
1135 ( (plist[ftype].param[j].ai() == at2) &&
1136 (plist[ftype].param[j].aj() == at1) ) );
1150 /* now copy the bond to the new array */
1151 ps->param[kept_i] = ps->param[i];
1154 else if (IS_CHEMBOND(cftype))
1156 srenew(plist[F_CONNBONDS].param, plist[F_CONNBONDS].nr+1);
1157 plist[F_CONNBONDS].param[plist[F_CONNBONDS].nr] = ps->param[i];
1158 plist[F_CONNBONDS].nr++;
1169 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1170 nremoved, interaction_function[cftype].longname, kept_i);
1174 fprintf(stderr, "Converted %4d %15ss with virtual sites to connections, %5d left\n",
1175 nconverted, interaction_function[cftype].longname, kept_i);
1179 fprintf(stderr, "Warning: removed %d %ss with vsite with %s construction\n"
1180 " This vsite construction does not guarantee constant "
1182 " If the constructions were generated by pdb2gmx ignore "
1184 nOut, interaction_function[cftype].longname,
1185 interaction_function[F_VSITE3OUT].longname );
1190 static void clean_vsite_angles(t_params *plist, t_pindex pindex[],
1191 int cftype, const int vsite_type[],
1192 at2vsitecon_t *at2vc)
1194 int i, j, k, m, n, nvsite, kept_i;
1196 bool bKeep, bUsed, bPresent, bAll3FAD, bFirstTwo;
1199 ps = &(plist[cftype]);
1201 for (i = 0; (i < ps->nr); i++) /* for all angles in the plist */
1204 const int *first_atoms = nullptr;
1208 /* check if all virtual sites are constructed from the same atoms */
1210 for (k = 0; (k < 3) && !bKeep; k++) /* for all atoms in the angle */
1212 atom = ps->param[i].a[k];
1213 if (vsite_type[atom] != NOTSET && vsite_type[atom] != F_VSITEN)
1216 bAll3FAD = bAll3FAD && (pindex[atom].ftype == F_VSITE3FAD);
1219 /* store construction atoms of first vsite */
1220 vsnral = NRAL(pindex[atom].ftype) - 1;
1221 first_atoms = plist[pindex[atom].ftype].param[pindex[atom].parnr].a + 1;
1225 GMX_ASSERT(vsnral != 0, "If we've seen a vsite before, we know how many constructing atoms it had");
1226 GMX_ASSERT(first_atoms != nullptr, "If we've seen a vsite before, we know what its first atom index was");
1227 /* check if this vsite is constructed from the same atoms */
1228 if (vsnral == NRAL(pindex[atom].ftype)-1)
1230 for (m = 0; (m < vsnral) && !bKeep; m++)
1235 atoms = plist[pindex[atom].ftype].param[pindex[atom].parnr].a + 1;
1236 for (n = 0; (n < vsnral) && !bPresent; n++)
1238 if (atoms[m] == first_atoms[n])
1257 /* keep all angles with no virtual sites in them or
1258 with virtual sites with more than 3 constr. atoms */
1259 if (nvsite == 0 && vsnral > 3)
1264 /* check if all non-vsite atoms are used in construction: */
1266 for (k = 0; (k < 3) && !bKeep; k++) /* for all atoms in the angle */
1268 atom = ps->param[i].a[k];
1269 if (vsite_type[atom] == NOTSET)
1272 for (m = 0; (m < vsnral) && !bUsed; m++)
1274 GMX_ASSERT(first_atoms != nullptr, "If we've seen a vsite before, we know what its first atom index was");
1276 if (atom == first_atoms[m])
1279 bFirstTwo = bFirstTwo && m < 2;
1289 if (!( bAll3FAD && bFirstTwo ) )
1291 /* check if all constructing atoms are constrained together */
1292 for (m = 0; m < vsnral && !bKeep; m++) /* all constr. atoms */
1294 at1 = first_atoms[m];
1295 at2 = first_atoms[(m+1) % vsnral];
1297 for (j = 0; j < at2vc[at1].nr; j++)
1299 if (at2vc[at1].aj[j] == at2)
1313 /* now copy the angle to the new array */
1314 ps->param[kept_i] = ps->param[i];
1319 if (ps->nr != kept_i)
1321 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1322 ps->nr-kept_i, interaction_function[cftype].longname, kept_i);
1327 static void clean_vsite_dihs(t_params *plist, t_pindex pindex[],
1328 int cftype, const int vsite_type[])
1333 ps = &(plist[cftype]);
1336 for (i = 0; (i < ps->nr); i++) /* for all dihedrals in the plist */
1338 int k, m, n, nvsite;
1340 const int *first_atoms = nullptr;
1342 bool bKeep, bUsed, bPresent;
1346 /* check if all virtual sites are constructed from the same atoms */
1348 for (k = 0; (k < 4) && !bKeep; k++) /* for all atoms in the dihedral */
1350 atom = ps->param[i].a[k];
1351 if (vsite_type[atom] != NOTSET && vsite_type[atom] != F_VSITEN)
1355 /* store construction atoms of first vsite */
1356 vsnral = NRAL(pindex[atom].ftype) - 1;
1357 first_atoms = plist[pindex[atom].ftype].param[pindex[atom].parnr].a + 1;
1361 GMX_ASSERT(vsnral != 0, "If we've seen a vsite before, we know how many constructing atoms it had");
1362 GMX_ASSERT(first_atoms != nullptr, "If we've seen a vsite before, we know what its first atom index was");
1363 /* check if this vsite is constructed from the same atoms */
1364 if (vsnral == NRAL(pindex[atom].ftype)-1)
1366 for (m = 0; (m < vsnral) && !bKeep; m++)
1371 atoms = plist[pindex[atom].ftype].param[pindex[atom].parnr].a + 1;
1372 for (n = 0; (n < vsnral) && !bPresent; n++)
1374 if (atoms[m] == first_atoms[n])
1386 /* TODO clean_site_bonds and _angles do this increment
1387 at the top of the loop. Refactor this for
1393 /* keep all dihedrals with no virtual sites in them */
1399 /* check if all atoms in dihedral are either virtual sites, or used in
1400 construction of virtual sites. If so, keep it, if not throw away: */
1401 for (k = 0; (k < 4) && !bKeep; k++) /* for all atoms in the dihedral */
1403 GMX_ASSERT(vsnral != 0, "If we've seen a vsite before, we know how many constructing atoms it had");
1404 GMX_ASSERT(first_atoms != nullptr, "If we've seen a vsite before, we know what its first atom index was");
1405 atom = ps->param[i].a[k];
1406 if (vsite_type[atom] == NOTSET)
1408 /* vsnral will be set here, we don't get here with nvsite==0 */
1410 for (m = 0; (m < vsnral) && !bUsed; m++)
1412 if (atom == first_atoms[m])
1426 ps->param[kept_i] = ps->param[i];
1431 if (ps->nr != kept_i)
1433 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1434 ps->nr-kept_i, interaction_function[cftype].longname, kept_i);
1439 void clean_vsite_bondeds(t_params *plist, int natoms, bool bRmVSiteBds)
1441 int i, k, nvsite, ftype, vsite, parnr;
1444 at2vsitecon_t *at2vc;
1446 pindex = nullptr; /* avoid warnings */
1447 /* make vsite_type array */
1448 snew(vsite_type, natoms);
1449 for (i = 0; i < natoms; i++)
1451 vsite_type[i] = NOTSET;
1454 for (ftype = 0; ftype < F_NRE; ftype++)
1456 if (interaction_function[ftype].flags & IF_VSITE)
1458 nvsite += plist[ftype].nr;
1460 while (i < plist[ftype].nr)
1462 vsite = plist[ftype].param[i].ai();
1463 if (vsite_type[vsite] == NOTSET)
1465 vsite_type[vsite] = ftype;
1469 gmx_fatal(FARGS, "multiple vsite constructions for atom %d", vsite+1);
1471 if (ftype == F_VSITEN)
1473 while (i < plist[ftype].nr && plist[ftype].param[i].ai() == vsite)
1486 /* the rest only if we have virtual sites: */
1489 fprintf(stderr, "Cleaning up constraints %swith virtual sites\n",
1490 bRmVSiteBds ? "and constant bonded interactions " : "");
1492 /* Make a reverse list to avoid ninteractions^2 operations */
1493 at2vc = make_at2vsitecon(natoms, plist);
1495 snew(pindex, natoms);
1496 for (ftype = 0; ftype < F_NRE; ftype++)
1498 /* Here we skip VSITEN. In neary all practical use cases this
1499 * is not an issue, since VSITEN is intended for constructing
1500 * additional interaction sites, not for replacing normal atoms
1501 * with bonded interactions. Thus we do not expect constant
1502 * bonded interactions. If VSITEN does get used with constant
1503 * bonded interactions, these are not removed which only leads
1504 * to very minor extra computation and constant energy.
1505 * The only problematic case is onstraints between VSITEN
1506 * constructions with fixed distance (which is anyhow useless).
1507 * This will generate a fatal error in check_vsite_constraints.
1509 if ((interaction_function[ftype].flags & IF_VSITE) &&
1512 for (parnr = 0; (parnr < plist[ftype].nr); parnr++)
1514 k = plist[ftype].param[parnr].ai();
1515 pindex[k].ftype = ftype;
1516 pindex[k].parnr = parnr;
1521 /* remove interactions that include virtual sites */
1522 for (ftype = 0; ftype < F_NRE; ftype++)
1524 if ( ( ( interaction_function[ftype].flags & IF_BOND ) && bRmVSiteBds ) ||
1525 ( interaction_function[ftype].flags & IF_CONSTRAINT ) )
1527 if (interaction_function[ftype].flags & (IF_BTYPE | IF_CONSTRAINT) )
1529 clean_vsite_bonds (plist, pindex, ftype, vsite_type);
1531 else if (interaction_function[ftype].flags & IF_ATYPE)
1533 clean_vsite_angles(plist, pindex, ftype, vsite_type, at2vc);
1535 else if ( (ftype == F_PDIHS) || (ftype == F_IDIHS) )
1537 clean_vsite_dihs (plist, pindex, ftype, vsite_type);
1541 /* check that no remaining constraints include virtual sites */
1542 for (ftype = 0; ftype < F_NRE; ftype++)
1544 if (interaction_function[ftype].flags & IF_CONSTRAINT)
1546 check_vsite_constraints(plist, ftype, vsite_type);
1550 done_at2vsitecon(natoms, at2vc);