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.
44 #include "vsite_parm.h"
45 #include "gromacs/utility/smalloc.h"
48 #include "gromacs/math/vec.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/legacyheaders/macros.h"
69 vsitebondparam_t *vsbp;
77 static int vsite_bond_nrcheck(int ftype)
81 if ((interaction_function[ftype].flags & (IF_BTYPE | IF_CONSTRAINT | IF_ATYPE)) || (ftype == F_IDIHS))
83 nrcheck = NRAL(ftype);
93 static void enter_bonded(int nratoms, int *nrbonded, t_mybonded **bondeds,
98 srenew(*bondeds, *nrbonded+1);
100 /* copy atom numbers */
101 for (j = 0; j < nratoms; j++)
103 (*bondeds)[*nrbonded].a[j] = param->a[j];
106 (*bondeds)[*nrbonded].c = param->C0;
111 static void get_bondeds(int nrat, t_iatom atoms[],
112 at2vsitebond_t *at2vb,
113 int *nrbond, t_mybonded **bonds,
114 int *nrang, t_mybonded **angles,
115 int *nridih, t_mybonded **idihs )
117 int k, i, ftype, nrcheck;
120 for (k = 0; k < nrat; k++)
122 for (i = 0; i < at2vb[atoms[k]].nr; i++)
124 ftype = at2vb[atoms[k]].vsbp[i].ftype;
125 param = at2vb[atoms[k]].vsbp[i].param;
126 nrcheck = vsite_bond_nrcheck(ftype);
127 /* abuse nrcheck to see if we're adding bond, angle or idih */
130 case 2: enter_bonded(nrcheck, nrbond, bonds, param); break;
131 case 3: enter_bonded(nrcheck, nrang, angles, param); break;
132 case 4: enter_bonded(nrcheck, nridih, idihs, param); break;
138 static at2vsitebond_t *make_at2vsitebond(int natoms, t_params plist[])
141 int ftype, i, j, nrcheck, nr;
143 at2vsitebond_t *at2vb;
148 for (ftype = 0; (ftype < F_NRE); ftype++)
150 if ((interaction_function[ftype].flags & IF_VSITE) && ftype != F_VSITEN)
152 for (i = 0; (i < plist[ftype].nr); i++)
154 for (j = 0; j < NRAL(ftype); j++)
156 bVSI[plist[ftype].param[i].a[j]] = TRUE;
162 for (ftype = 0; (ftype < F_NRE); ftype++)
164 nrcheck = vsite_bond_nrcheck(ftype);
167 for (i = 0; (i < plist[ftype].nr); i++)
169 aa = plist[ftype].param[i].a;
170 for (j = 0; j < nrcheck; j++)
174 nr = at2vb[aa[j]].nr;
177 srenew(at2vb[aa[j]].vsbp, nr+10);
179 at2vb[aa[j]].vsbp[nr].ftype = ftype;
180 at2vb[aa[j]].vsbp[nr].param = &plist[ftype].param[i];
192 static void done_at2vsitebond(int natoms, at2vsitebond_t *at2vb)
196 for (i = 0; i < natoms; i++)
200 sfree(at2vb[i].vsbp);
206 static at2vsitecon_t *make_at2vsitecon(int natoms, t_params plist[])
209 int ftype, i, j, ai, aj, nr;
210 at2vsitecon_t *at2vc;
215 for (ftype = 0; (ftype < F_NRE); ftype++)
217 if ((interaction_function[ftype].flags & IF_VSITE) && ftype != F_VSITEN)
219 for (i = 0; (i < plist[ftype].nr); i++)
221 for (j = 0; j < NRAL(ftype); j++)
223 bVSI[plist[ftype].param[i].a[j]] = TRUE;
229 for (ftype = 0; (ftype < F_NRE); ftype++)
231 if (interaction_function[ftype].flags & IF_CONSTRAINT)
233 for (i = 0; (i < plist[ftype].nr); i++)
235 ai = plist[ftype].param[i].AI;
236 aj = plist[ftype].param[i].AJ;
237 if (bVSI[ai] && bVSI[aj])
239 /* Store forward direction */
243 srenew(at2vc[ai].aj, nr+10);
245 at2vc[ai].aj[nr] = aj;
247 /* Store backward direction */
251 srenew(at2vc[aj].aj, nr+10);
253 at2vc[aj].aj[nr] = ai;
264 static void done_at2vsitecon(int natoms, at2vsitecon_t *at2vc)
268 for (i = 0; i < natoms; i++)
279 static void print_bad(FILE *fp,
280 int nrbond, t_mybonded *bonds,
281 int nrang, t_mybonded *angles,
282 int nridih, t_mybonded *idihs )
288 fprintf(fp, "bonds:");
289 for (i = 0; i < nrbond; i++)
291 fprintf(fp, " %u-%u (%g)",
292 bonds[i].AI+1, bonds[i].AJ+1, bonds[i].c);
298 fprintf(fp, "angles:");
299 for (i = 0; i < nrang; i++)
301 fprintf(fp, " %u-%u-%u (%g)",
302 angles[i].AI+1, angles[i].AJ+1,
303 angles[i].AK+1, angles[i].c);
309 fprintf(fp, "idihs:");
310 for (i = 0; i < nridih; i++)
312 fprintf(fp, " %u-%u-%u-%u (%g)",
313 idihs[i].AI+1, idihs[i].AJ+1,
314 idihs[i].AK+1, idihs[i].AL+1, idihs[i].c);
320 static void print_param(FILE *fp, int ftype, int i, t_param *param)
323 static int prev_ftype = NOTSET;
324 static int prev_i = NOTSET;
327 if ( (ftype != prev_ftype) || (i != prev_i) )
333 fprintf(fp, "(%d) plist[%s].param[%d]",
334 pass, interaction_function[ftype].name, i);
335 for (j = 0; j < NRFP(ftype); j++)
337 fprintf(fp, ".c[%d]=%g ", j, param->c[j]);
343 static real get_bond_length(int nrbond, t_mybonded bonds[],
344 t_iatom ai, t_iatom aj)
350 for (i = 0; i < nrbond && (bondlen == NOTSET); i++)
352 /* check both ways */
353 if ( ( (ai == bonds[i].AI) && (aj == bonds[i].AJ) ) ||
354 ( (ai == bonds[i].AJ) && (aj == bonds[i].AI) ) )
356 bondlen = bonds[i].c; /* note: bonds[i].c might be NOTSET */
362 static real get_angle(int nrang, t_mybonded angles[],
363 t_iatom ai, t_iatom aj, t_iatom ak)
369 for (i = 0; i < nrang && (angle == NOTSET); i++)
371 /* check both ways */
372 if ( ( (ai == angles[i].AI) && (aj == angles[i].AJ) && (ak == angles[i].AK) ) ||
373 ( (ai == angles[i].AK) && (aj == angles[i].AJ) && (ak == angles[i].AI) ) )
375 angle = DEG2RAD*angles[i].c;
381 static char *get_atomtype_name_AB(t_atom *atom, gpp_atomtype_t atype)
385 name = get_atomtype_name(atom->type, atype);
387 /* When using the decoupling option, atom types are changed
388 * to decoupled for the non-bonded interactions, but the virtual
389 * sites constructions should be based on the "normal" interactions.
390 * So we return the state B atom type name if the state A atom
391 * type is the decoupled one. We should actually check for the atom
392 * type number, but that's not passed here. So we check for
393 * the decoupled atom type name. This should not cause trouble
394 * as this code is only used for topologies with v-sites without
395 * parameters generated by pdb2gmx.
397 if (strcmp(name, "decoupled") == 0)
399 name = get_atomtype_name(atom->typeB, atype);
405 static gmx_bool calc_vsite3_param(gpp_atomtype_t atype,
406 t_param *param, t_atoms *at,
407 int nrbond, t_mybonded *bonds,
408 int nrang, t_mybonded *angles )
410 /* i = virtual site | ,k
411 * j = 1st bonded heavy atom | i-j
412 * k,l = 2nd bonded atoms | `l
415 gmx_bool bXH3, bError;
416 real bjk, bjl, a = -1, b = -1;
417 /* check if this is part of a NH3 , NH2-umbrella or CH3 group,
418 * i.e. if atom k and l are dummy masses (MNH* or MCH3*) */
422 for (i = 0; i < 4; i++)
424 fprintf(debug, "atom %u type %s ",
426 get_atomtype_name_AB(&at->atom[param->a[i]], atype));
428 fprintf(debug, "\n");
431 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AK], atype), "MNH", 3) == 0) &&
432 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AL], atype), "MNH", 3) == 0) ) ||
433 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AK], atype), "MCH3", 4) == 0) &&
434 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AL], atype), "MCH3", 4) == 0) );
436 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
437 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
438 bError = (bjk == NOTSET) || (bjl == NOTSET);
441 /* now we get some XH2/XH3 group specific construction */
442 /* note: we call the heavy atom 'C' and the X atom 'N' */
443 real bMM, bCM, bCN, bNH, aCNH, dH, rH, dM, rM;
446 /* check if bonds from heavy atom (j) to dummy masses (k,l) are equal: */
447 bError = bError || (bjk != bjl);
449 /* the X atom (C or N) in the XH2/XH3 group is the first after the masses: */
450 aN = max(param->AK, param->AL)+1;
452 /* get common bonds */
453 bMM = get_bond_length(nrbond, bonds, param->AK, param->AL);
455 bCN = get_bond_length(nrbond, bonds, param->AJ, aN);
456 bError = bError || (bMM == NOTSET) || (bCN == NOTSET);
458 /* calculate common things */
460 dM = sqrt( sqr(bCM) - sqr(rM) );
462 /* are we dealing with the X atom? */
465 /* this is trivial */
466 a = b = 0.5 * bCN/dM;
471 /* get other bondlengths and angles: */
472 bNH = get_bond_length(nrbond, bonds, aN, param->AI);
473 aCNH = get_angle (nrang, angles, param->AJ, aN, param->AI);
474 bError = bError || (bNH == NOTSET) || (aCNH == NOTSET);
477 dH = bCN - bNH * cos(aCNH);
478 rH = bNH * sin(aCNH);
480 a = 0.5 * ( dH/dM + rH/rM );
481 b = 0.5 * ( dH/dM - rH/rM );
486 gmx_fatal(FARGS, "calc_vsite3_param not implemented for the general case "
487 "(atom %d)", param->AI+1);
495 fprintf(debug, "params for vsite3 %u: %g %g\n",
496 param->AI+1, param->C0, param->C1);
502 static gmx_bool calc_vsite3fd_param(t_param *param,
503 int nrbond, t_mybonded *bonds,
504 int nrang, t_mybonded *angles)
506 /* i = virtual site | ,k
507 * j = 1st bonded heavy atom | i-j
508 * k,l = 2nd bonded atoms | `l
512 real bij, bjk, bjl, aijk, aijl, rk, rl;
514 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
515 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
516 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
517 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
518 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
519 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) ||
520 (aijk == NOTSET) || (aijl == NOTSET);
522 rk = bjk * sin(aijk);
523 rl = bjl * sin(aijl);
524 param->C0 = rk / (rk + rl);
525 param->C1 = -bij; /* 'bond'-length for fixed distance vsite */
529 fprintf(debug, "params for vsite3fd %u: %g %g\n",
530 param->AI+1, param->C0, param->C1);
535 static gmx_bool calc_vsite3fad_param(t_param *param,
536 int nrbond, t_mybonded *bonds,
537 int nrang, t_mybonded *angles)
539 /* i = virtual site |
540 * j = 1st bonded heavy atom | i-j
541 * k = 2nd bonded heavy atom | `k-l
542 * l = 3d bonded heavy atom |
545 gmx_bool bSwapParity, bError;
548 bSwapParity = ( param->C1 == -1 );
550 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
551 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
552 bError = (bij == NOTSET) || (aijk == NOTSET);
554 param->C1 = bij; /* 'bond'-length for fixed distance vsite */
555 param->C0 = RAD2DEG*aijk; /* 'bond'-angle for fixed angle vsite */
559 param->C0 = 360 - param->C0;
564 fprintf(debug, "params for vsite3fad %u: %g %g\n",
565 param->AI+1, param->C0, param->C1);
570 static gmx_bool calc_vsite3out_param(gpp_atomtype_t atype,
571 t_param *param, t_atoms *at,
572 int nrbond, t_mybonded *bonds,
573 int nrang, t_mybonded *angles)
575 /* i = virtual site | ,k
576 * j = 1st bonded heavy atom | i-j
577 * k,l = 2nd bonded atoms | `l
578 * NOTE: i is out of the j-k-l plane!
581 gmx_bool bXH3, bError, bSwapParity;
582 real bij, bjk, bjl, aijk, aijl, akjl, pijk, pijl, a, b, c;
584 /* check if this is part of a NH2-umbrella, NH3 or CH3 group,
585 * i.e. if atom k and l are dummy masses (MNH* or MCH3*) */
589 for (i = 0; i < 4; i++)
591 fprintf(debug, "atom %u type %s ",
592 param->a[i]+1, get_atomtype_name_AB(&at->atom[param->a[i]], atype));
594 fprintf(debug, "\n");
597 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AK], atype), "MNH", 3) == 0) &&
598 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AL], atype), "MNH", 3) == 0) ) ||
599 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AK], atype), "MCH3", 4) == 0) &&
600 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AL], atype), "MCH3", 4) == 0) );
602 /* check if construction parity must be swapped */
603 bSwapParity = ( param->C1 == -1 );
605 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
606 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
607 bError = (bjk == NOTSET) || (bjl == NOTSET);
610 /* now we get some XH3 group specific construction */
611 /* note: we call the heavy atom 'C' and the X atom 'N' */
612 real bMM, bCM, bCN, bNH, aCNH, dH, rH, rHx, rHy, dM, rM;
615 /* check if bonds from heavy atom (j) to dummy masses (k,l) are equal: */
616 bError = bError || (bjk != bjl);
618 /* the X atom (C or N) in the XH3 group is the first after the masses: */
619 aN = max(param->AK, param->AL)+1;
621 /* get all bondlengths and angles: */
622 bMM = get_bond_length(nrbond, bonds, param->AK, param->AL);
624 bCN = get_bond_length(nrbond, bonds, param->AJ, aN);
625 bNH = get_bond_length(nrbond, bonds, aN, param->AI);
626 aCNH = get_angle (nrang, angles, param->AJ, aN, param->AI);
628 (bMM == NOTSET) || (bCN == NOTSET) || (bNH == NOTSET) || (aCNH == NOTSET);
631 dH = bCN - bNH * cos(aCNH);
632 rH = bNH * sin(aCNH);
633 /* we assume the H's are symmetrically distributed */
634 rHx = rH*cos(DEG2RAD*30);
635 rHy = rH*sin(DEG2RAD*30);
637 dM = sqrt( sqr(bCM) - sqr(rM) );
638 a = 0.5*( (dH/dM) - (rHy/rM) );
639 b = 0.5*( (dH/dM) + (rHy/rM) );
645 /* this is the general construction */
647 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
648 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
649 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
650 akjl = get_angle (nrang, angles, param->AK, param->AJ, param->AL);
652 (bij == NOTSET) || (aijk == NOTSET) || (aijl == NOTSET) || (akjl == NOTSET);
654 pijk = cos(aijk)*bij;
655 pijl = cos(aijl)*bij;
656 a = ( pijk + (pijk*cos(akjl)-pijl) * cos(akjl) / sqr(sin(akjl)) ) / bjk;
657 b = ( pijl + (pijl*cos(akjl)-pijk) * cos(akjl) / sqr(sin(akjl)) ) / bjl;
658 c = -sqrt( sqr(bij) -
659 ( sqr(pijk) - 2*pijk*pijl*cos(akjl) + sqr(pijl) )
661 / ( bjk*bjl*sin(akjl) );
676 fprintf(debug, "params for vsite3out %u: %g %g %g\n",
677 param->AI+1, param->C0, param->C1, param->C2);
682 static gmx_bool calc_vsite4fd_param(t_param *param,
683 int nrbond, t_mybonded *bonds,
684 int nrang, t_mybonded *angles)
686 /* i = virtual site | ,k
687 * j = 1st bonded heavy atom | i-j-m
688 * k,l,m = 2nd bonded atoms | `l
692 real bij, bjk, bjl, bjm, aijk, aijl, aijm, akjm, akjl;
693 real pk, pl, pm, cosakl, cosakm, sinakl, sinakm, cl, cm;
695 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
696 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
697 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
698 bjm = get_bond_length(nrbond, bonds, param->AJ, param->AM);
699 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
700 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
701 aijm = get_angle (nrang, angles, param->AI, param->AJ, param->AM);
702 akjm = get_angle (nrang, angles, param->AK, param->AJ, param->AM);
703 akjl = get_angle (nrang, angles, param->AK, param->AJ, param->AL);
704 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) || (bjm == NOTSET) ||
705 (aijk == NOTSET) || (aijl == NOTSET) || (aijm == NOTSET) || (akjm == NOTSET) ||
713 cosakl = (cos(akjl) - cos(aijk)*cos(aijl)) / (sin(aijk)*sin(aijl));
714 cosakm = (cos(akjm) - cos(aijk)*cos(aijm)) / (sin(aijk)*sin(aijm));
715 if (cosakl < -1 || cosakl > 1 || cosakm < -1 || cosakm > 1)
717 fprintf(stderr, "virtual site %d: angle ijk = %f, angle ijl = %f, angle ijm = %f\n",
718 param->AI+1, RAD2DEG*aijk, RAD2DEG*aijl, RAD2DEG*aijm);
719 gmx_fatal(FARGS, "invalid construction in calc_vsite4fd for atom %d: "
720 "cosakl=%f, cosakm=%f\n", param->AI+1, cosakl, cosakm);
722 sinakl = sqrt(1-sqr(cosakl));
723 sinakm = sqrt(1-sqr(cosakm));
725 /* note: there is a '+' because of the way the sines are calculated */
726 cl = -pk / ( pl*cosakl - pk + pl*sinakl*(pm*cosakm-pk)/(pm*sinakm) );
727 cm = -pk / ( pm*cosakm - pk + pm*sinakm*(pl*cosakl-pk)/(pl*sinakl) );
734 fprintf(debug, "params for vsite4fd %u: %g %g %g\n",
735 param->AI+1, param->C0, param->C1, param->C2);
744 calc_vsite4fdn_param(t_param *param,
745 int nrbond, t_mybonded *bonds,
746 int nrang, t_mybonded *angles)
748 /* i = virtual site | ,k
749 * j = 1st bonded heavy atom | i-j-m
750 * k,l,m = 2nd bonded atoms | `l
754 real bij, bjk, bjl, bjm, aijk, aijl, aijm;
755 real pk, pl, pm, a, b;
757 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
758 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
759 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
760 bjm = get_bond_length(nrbond, bonds, param->AJ, param->AM);
761 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
762 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
763 aijm = get_angle (nrang, angles, param->AI, param->AJ, param->AM);
765 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) || (bjm == NOTSET) ||
766 (aijk == NOTSET) || (aijl == NOTSET) || (aijm == NOTSET);
771 /* Calculate component of bond j-k along the direction i-j */
774 /* Calculate component of bond j-l along the direction i-j */
777 /* Calculate component of bond j-m along the direction i-j */
780 if (fabs(pl) < 1000*GMX_REAL_MIN || fabs(pm) < 1000*GMX_REAL_MIN)
782 fprintf(stderr, "virtual site %d: angle ijk = %f, angle ijl = %f, angle ijm = %f\n",
783 param->AI+1, RAD2DEG*aijk, RAD2DEG*aijl, RAD2DEG*aijm);
784 gmx_fatal(FARGS, "invalid construction in calc_vsite4fdn for atom %d: "
785 "pl=%f, pm=%f\n", param->AI+1, pl, pm);
797 fprintf(debug, "params for vsite4fdn %u: %g %g %g\n",
798 param->AI+1, param->C0, param->C1, param->C2);
807 int set_vsites(gmx_bool bVerbose, t_atoms *atoms, gpp_atomtype_t atype,
811 int nvsite, nrbond, nrang, nridih, nrset;
812 gmx_bool bFirst, bSet, bERROR;
813 at2vsitebond_t *at2vb;
823 fprintf(debug, "\nCalculating parameters for virtual sites\n");
826 /* Make a reverse list to avoid ninteractions^2 operations */
827 at2vb = make_at2vsitebond(atoms->nr, plist);
829 for (ftype = 0; (ftype < F_NRE); ftype++)
831 if ((interaction_function[ftype].flags & IF_VSITE) && ftype != F_VSITEN)
834 nvsite += plist[ftype].nr;
835 for (i = 0; (i < plist[ftype].nr); i++)
837 /* check if all parameters are set */
839 for (j = 0; j < NRFP(ftype) && bSet; j++)
841 bSet = plist[ftype].param[i].c[j] != NOTSET;
846 fprintf(debug, "bSet=%s ", bool_names[bSet]);
847 print_param(debug, ftype, i, &plist[ftype].param[i]);
851 if (bVerbose && bFirst)
853 fprintf(stderr, "Calculating parameters for virtual sites\n");
857 nrbond = nrang = nridih = 0;
862 /* now set the vsite parameters: */
863 get_bondeds(NRAL(ftype), plist[ftype].param[i].a, at2vb,
864 &nrbond, &bonds, &nrang, &angles, &nridih, &idihs);
867 fprintf(debug, "Found %d bonds, %d angles and %d idihs "
868 "for virtual site %u (%s)\n", nrbond, nrang, nridih,
869 plist[ftype].param[i].AI+1,
870 interaction_function[ftype].longname);
871 print_bad(debug, nrbond, bonds, nrang, angles, nridih, idihs);
877 calc_vsite3_param(atype, &(plist[ftype].param[i]), atoms,
878 nrbond, bonds, nrang, angles);
882 calc_vsite3fd_param(&(plist[ftype].param[i]),
883 nrbond, bonds, nrang, angles);
887 calc_vsite3fad_param(&(plist[ftype].param[i]),
888 nrbond, bonds, nrang, angles);
892 calc_vsite3out_param(atype, &(plist[ftype].param[i]), atoms,
893 nrbond, bonds, nrang, angles);
897 calc_vsite4fd_param(&(plist[ftype].param[i]),
898 nrbond, bonds, nrang, angles);
902 calc_vsite4fdn_param(&(plist[ftype].param[i]),
903 nrbond, bonds, nrang, angles);
906 gmx_fatal(FARGS, "Automatic parameter generation not supported "
908 interaction_function[ftype].longname,
909 plist[ftype].param[i].AI+1);
913 gmx_fatal(FARGS, "Automatic parameter generation not supported "
914 "for %s atom %d for this bonding configuration",
915 interaction_function[ftype].longname,
916 plist[ftype].param[i].AI+1);
923 if (debug && plist[ftype].nr)
925 fprintf(stderr, "Calculated parameters for %d out of %d %s atoms\n",
926 nrset, plist[ftype].nr, interaction_function[ftype].longname);
931 done_at2vsitebond(atoms->nr, at2vb);
936 void set_vsites_ptype(gmx_bool bVerbose, gmx_moltype_t *molt)
945 fprintf(stderr, "Setting particle type to V for virtual sites\n");
949 fprintf(stderr, "checking %d functypes\n", F_NRE);
951 for (ftype = 0; ftype < F_NRE; ftype++)
953 il = &molt->ilist[ftype];
954 if (interaction_function[ftype].flags & IF_VSITE)
956 nra = interaction_function[ftype].nratoms;
962 fprintf(stderr, "doing %d %s virtual sites\n",
963 (nrd / (nra+1)), interaction_function[ftype].longname);
966 for (i = 0; (i < nrd); )
968 /* The virtual site */
970 molt->atoms.atom[avsite].ptype = eptVSite;
984 static void check_vsite_constraints(t_params *plist,
985 int cftype, int vsite_type[])
992 ps = &(plist[cftype]);
993 for (i = 0; (i < ps->nr); i++)
995 for (k = 0; k < 2; k++)
997 atom = ps->param[i].a[k];
998 if (vsite_type[atom] != NOTSET)
1000 fprintf(stderr, "ERROR: Cannot have constraint (%u-%u) with virtual site (%u)\n",
1001 ps->param[i].AI+1, ps->param[i].AJ+1, atom+1);
1008 gmx_fatal(FARGS, "There were %d virtual sites involved in constraints", n);
1012 static void clean_vsite_bonds(t_params *plist, t_pindex pindex[],
1013 int cftype, int vsite_type[])
1015 int ftype, i, j, parnr, k, l, m, n, nvsite, nOut, kept_i, vsnral, vsitetype;
1016 int nconverted, nremoved;
1017 atom_id atom, oatom, constr, at1, at2;
1018 atom_id vsiteatoms[MAXATOMLIST];
1019 gmx_bool bKeep, bRemove, bUsed, bPresent, bThisFD, bThisOUT, bAllFD, bFirstTwo;
1022 if (cftype == F_CONNBONDS)
1027 ps = &(plist[cftype]);
1033 for (i = 0; (i < ps->nr); i++) /* for all bonds in the plist */
1038 /* check if all virtual sites are constructed from the same atoms */
1042 fprintf(debug, "constr %u %u:", ps->param[i].AI+1, ps->param[i].AJ+1);
1044 for (k = 0; (k < 2) && !bKeep && !bRemove; k++)
1046 /* for all atoms in the bond */
1047 atom = ps->param[i].a[k];
1048 if (vsite_type[atom] != NOTSET)
1052 fprintf(debug, " d%d[%d: %d %d %d]", k, atom+1,
1053 plist[pindex[atom].ftype].param[pindex[atom].parnr].AJ+1,
1054 plist[pindex[atom].ftype].param[pindex[atom].parnr].AK+1,
1055 plist[pindex[atom].ftype].param[pindex[atom].parnr].AL+1);
1058 bThisFD = ( (pindex[atom].ftype == F_VSITE3FD ) ||
1059 (pindex[atom].ftype == F_VSITE3FAD) ||
1060 (pindex[atom].ftype == F_VSITE4FD ) ||
1061 (pindex[atom].ftype == F_VSITE4FDN ) );
1062 bThisOUT = ( (pindex[atom].ftype == F_VSITE3OUT) &&
1063 (interaction_function[cftype].flags & IF_CONSTRAINT) );
1064 bAllFD = bAllFD && bThisFD;
1065 if (bThisFD || bThisOUT)
1069 fprintf(debug, " %s", bThisOUT ? "out" : "fd");
1071 oatom = ps->param[i].a[1-k]; /* the other atom */
1072 if (vsite_type[oatom] == NOTSET &&
1073 oatom == plist[pindex[atom].ftype].param[pindex[atom].parnr].AJ)
1075 /* if the other atom isn't a vsite, and it is AI */
1083 fprintf(debug, " D-AI");
1091 /* if this is the first vsite we encounter then
1092 store construction atoms */
1093 vsnral = NRAL(pindex[atom].ftype)-1;
1094 for (m = 0; (m < vsnral); m++)
1097 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1102 /* if it is not the first then
1103 check if this vsite is constructed from the same atoms */
1104 if (vsnral == NRAL(pindex[atom].ftype)-1)
1106 for (m = 0; (m < vsnral) && !bKeep; m++)
1110 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1111 for (n = 0; (n < vsnral) && !bPresent; n++)
1113 if (constr == vsiteatoms[n])
1123 fprintf(debug, " !present");
1133 fprintf(debug, " !same#at");
1147 /* if we have no virtual sites in this bond, keep it */
1152 fprintf(debug, " no vsite");
1157 /* check if all non-vsite atoms are used in construction: */
1159 for (k = 0; (k < 2) && !bKeep; k++) /* for all atoms in the bond */
1161 atom = ps->param[i].a[k];
1162 if (vsite_type[atom] == NOTSET)
1165 for (m = 0; (m < vsnral) && !bUsed; m++)
1167 if (atom == vsiteatoms[m])
1170 bFirstTwo = bFirstTwo && m < 2;
1178 fprintf(debug, " !used");
1184 if (!( bAllFD && bFirstTwo ) )
1186 /* check if all constructing atoms are constrained together */
1187 for (m = 0; m < vsnral && !bKeep; m++) /* all constr. atoms */
1189 at1 = vsiteatoms[m];
1190 at2 = vsiteatoms[(m+1) % vsnral];
1192 for (ftype = 0; ftype < F_NRE; ftype++)
1194 if (interaction_function[ftype].flags & IF_CONSTRAINT)
1196 for (j = 0; (j < plist[ftype].nr) && !bPresent; j++)
1198 /* all constraints until one matches */
1199 bPresent = ( ( (plist[ftype].param[j].AI == at1) &&
1200 (plist[ftype].param[j].AJ == at2) ) ||
1201 ( (plist[ftype].param[j].AI == at2) &&
1202 (plist[ftype].param[j].AJ == at1) ) );
1211 fprintf(debug, " !bonded");
1222 fprintf(debug, " keeping");
1224 /* now copy the bond to the new array */
1225 ps->param[kept_i] = ps->param[i];
1228 else if (IS_CHEMBOND(cftype))
1230 srenew(plist[F_CONNBONDS].param, plist[F_CONNBONDS].nr+1);
1231 plist[F_CONNBONDS].param[plist[F_CONNBONDS].nr] = ps->param[i];
1232 plist[F_CONNBONDS].nr++;
1241 fprintf(debug, "\n");
1247 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1248 nremoved, interaction_function[cftype].longname, kept_i);
1252 fprintf(stderr, "Converted %4d %15ss with virtual sites to connections, %5d left\n",
1253 nconverted, interaction_function[cftype].longname, kept_i);
1257 fprintf(stderr, "Warning: removed %d %ss with vsite with %s construction\n"
1258 " This vsite construction does not guarantee constant "
1260 " If the constructions were generated by pdb2gmx ignore "
1262 nOut, interaction_function[cftype].longname,
1263 interaction_function[F_VSITE3OUT].longname );
1268 static void clean_vsite_angles(t_params *plist, t_pindex pindex[],
1269 int cftype, int vsite_type[],
1270 at2vsitecon_t *at2vc)
1272 int i, j, parnr, k, l, m, n, nvsite, kept_i, vsnral, vsitetype;
1273 atom_id atom, constr, at1, at2;
1274 atom_id vsiteatoms[MAXATOMLIST];
1275 gmx_bool bKeep, bUsed, bPresent, bAll3FAD, bFirstTwo;
1278 ps = &(plist[cftype]);
1281 for (i = 0; (i < ps->nr); i++) /* for all angles in the plist */
1285 /* check if all virtual sites are constructed from the same atoms */
1287 for (k = 0; (k < 3) && !bKeep; k++) /* for all atoms in the angle */
1289 atom = ps->param[i].a[k];
1290 if (vsite_type[atom] != NOTSET)
1293 bAll3FAD = bAll3FAD && (pindex[atom].ftype == F_VSITE3FAD);
1296 /* store construction atoms of first vsite */
1297 vsnral = NRAL(pindex[atom].ftype)-1;
1298 for (m = 0; (m < vsnral); m++)
1301 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1305 /* check if this vsite is constructed from the same atoms */
1306 if (vsnral == NRAL(pindex[atom].ftype)-1)
1308 for (m = 0; (m < vsnral) && !bKeep; m++)
1312 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1313 for (n = 0; (n < vsnral) && !bPresent; n++)
1315 if (constr == vsiteatoms[n])
1333 /* keep all angles with no virtual sites in them or
1334 with virtual sites with more than 3 constr. atoms */
1335 if (nvsite == 0 && vsnral > 3)
1340 /* check if all non-vsite atoms are used in construction: */
1342 for (k = 0; (k < 3) && !bKeep; k++) /* for all atoms in the angle */
1344 atom = ps->param[i].a[k];
1345 if (vsite_type[atom] == NOTSET)
1348 for (m = 0; (m < vsnral) && !bUsed; m++)
1350 if (atom == vsiteatoms[m])
1353 bFirstTwo = bFirstTwo && m < 2;
1363 if (!( bAll3FAD && bFirstTwo ) )
1365 /* check if all constructing atoms are constrained together */
1366 for (m = 0; m < vsnral && !bKeep; m++) /* all constr. atoms */
1368 at1 = vsiteatoms[m];
1369 at2 = vsiteatoms[(m+1) % vsnral];
1371 for (j = 0; j < at2vc[at1].nr; j++)
1373 if (at2vc[at1].aj[j] == at2)
1387 /* now copy the angle to the new array */
1388 ps->param[kept_i] = ps->param[i];
1393 if (ps->nr != kept_i)
1395 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1396 ps->nr-kept_i, interaction_function[cftype].longname, kept_i);
1401 static void clean_vsite_dihs(t_params *plist, t_pindex pindex[],
1402 int cftype, int vsite_type[])
1407 ps = &(plist[cftype]);
1410 for (i = 0; (i < ps->nr); i++) /* for all dihedrals in the plist */
1412 int ftype, parnr, k, l, m, n, nvsite;
1413 int vsnral = 0; /* keep the compiler happy */
1414 atom_id atom, constr;
1415 atom_id vsiteatoms[4] = { 0 }; /* init to zero to make gcc4.8 happy */
1416 gmx_bool bKeep, bUsed, bPresent;
1420 /* check if all virtual sites are constructed from the same atoms */
1422 for (k = 0; (k < 4) && !bKeep; k++) /* for all atoms in the dihedral */
1424 atom = ps->param[i].a[k];
1425 if (vsite_type[atom] != NOTSET)
1429 /* store construction atoms of first vsite */
1430 vsnral = NRAL(pindex[atom].ftype)-1;
1431 assert(vsnral <= 4);
1432 for (m = 0; (m < vsnral); m++)
1435 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1439 fprintf(debug, "dih w. vsite: %u %u %u %u\n",
1440 ps->param[i].AI+1, ps->param[i].AJ+1,
1441 ps->param[i].AK+1, ps->param[i].AL+1);
1442 fprintf(debug, "vsite %u from: %u %u %u\n",
1443 atom+1, vsiteatoms[0]+1, vsiteatoms[1]+1, vsiteatoms[2]+1);
1448 /* check if this vsite is constructed from the same atoms */
1449 if (vsnral == NRAL(pindex[atom].ftype)-1)
1451 for (m = 0; (m < vsnral) && !bKeep; m++)
1455 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1456 for (n = 0; (n < vsnral) && !bPresent; n++)
1458 if (constr == vsiteatoms[n])
1474 /* keep all dihedrals with no virtual sites in them */
1480 /* check if all atoms in dihedral are either virtual sites, or used in
1481 construction of virtual sites. If so, keep it, if not throw away: */
1482 for (k = 0; (k < 4) && !bKeep; k++) /* for all atoms in the dihedral */
1484 atom = ps->param[i].a[k];
1485 if (vsite_type[atom] == NOTSET)
1487 /* vsnral will be set here, we don't get here with nvsite==0 */
1489 for (m = 0; (m < vsnral) && !bUsed; m++)
1491 if (atom == vsiteatoms[m])
1501 fprintf(debug, "unused atom in dih: %u\n", atom+1);
1509 ps->param[kept_i] = ps->param[i];
1514 if (ps->nr != kept_i)
1516 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1517 ps->nr-kept_i, interaction_function[cftype].longname, kept_i);
1522 void clean_vsite_bondeds(t_params *plist, int natoms, gmx_bool bRmVSiteBds)
1524 int i, k, nvsite, ftype, vsite, parnr;
1527 at2vsitecon_t *at2vc;
1529 pindex = 0; /* avoid warnings */
1530 /* make vsite_type array */
1531 snew(vsite_type, natoms);
1532 for (i = 0; i < natoms; i++)
1534 vsite_type[i] = NOTSET;
1537 for (ftype = 0; ftype < F_NRE; ftype++)
1539 if (interaction_function[ftype].flags & IF_VSITE)
1541 nvsite += plist[ftype].nr;
1543 while (i < plist[ftype].nr)
1545 vsite = plist[ftype].param[i].AI;
1546 if (vsite_type[vsite] == NOTSET)
1548 vsite_type[vsite] = ftype;
1552 gmx_fatal(FARGS, "multiple vsite constructions for atom %d", vsite+1);
1554 if (ftype == F_VSITEN)
1556 while (i < plist[ftype].nr && plist[ftype].param[i].AI == vsite)
1569 /* the rest only if we have virtual sites: */
1572 fprintf(stderr, "Cleaning up constraints %swith virtual sites\n",
1573 bRmVSiteBds ? "and constant bonded interactions " : "");
1575 /* Make a reverse list to avoid ninteractions^2 operations */
1576 at2vc = make_at2vsitecon(natoms, plist);
1578 snew(pindex, natoms);
1579 for (ftype = 0; ftype < F_NRE; ftype++)
1581 if ((interaction_function[ftype].flags & IF_VSITE) &&
1584 for (parnr = 0; (parnr < plist[ftype].nr); parnr++)
1586 k = plist[ftype].param[parnr].AI;
1587 pindex[k].ftype = ftype;
1588 pindex[k].parnr = parnr;
1595 for (i = 0; i < natoms; i++)
1597 fprintf(debug, "atom %d vsite_type %s\n", i,
1598 vsite_type[i] == NOTSET ? "NOTSET" :
1599 interaction_function[vsite_type[i]].name);
1603 /* remove things with vsite atoms */
1604 for (ftype = 0; ftype < F_NRE; ftype++)
1606 if ( ( ( interaction_function[ftype].flags & IF_BOND ) && bRmVSiteBds ) ||
1607 ( interaction_function[ftype].flags & IF_CONSTRAINT ) )
1609 if (interaction_function[ftype].flags & (IF_BTYPE | IF_CONSTRAINT) )
1611 clean_vsite_bonds (plist, pindex, ftype, vsite_type);
1613 else if (interaction_function[ftype].flags & IF_ATYPE)
1615 clean_vsite_angles(plist, pindex, ftype, vsite_type, at2vc);
1617 else if ( (ftype == F_PDIHS) || (ftype == F_IDIHS) )
1619 clean_vsite_dihs (plist, pindex, ftype, vsite_type);
1623 /* check if we have constraints left with virtual sites in them */
1624 for (ftype = 0; ftype < F_NRE; ftype++)
1626 if (interaction_function[ftype].flags & IF_CONSTRAINT)
1628 check_vsite_constraints(plist, ftype, vsite_type);
1632 done_at2vsitecon(natoms, at2vc);