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.
39 #include "vsite_parm.h"
46 #include "gromacs/gmxpreprocess/add_par.h"
47 #include "gromacs/gmxpreprocess/resall.h"
48 #include "gromacs/gmxpreprocess/toputil.h"
49 #include "gromacs/legacyheaders/macros.h"
50 #include "gromacs/legacyheaders/names.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/smalloc.h"
70 vsitebondparam_t *vsbp;
78 static int vsite_bond_nrcheck(int ftype)
82 if ((interaction_function[ftype].flags & (IF_BTYPE | IF_CONSTRAINT | IF_ATYPE)) || (ftype == F_IDIHS))
84 nrcheck = NRAL(ftype);
94 static void enter_bonded(int nratoms, int *nrbonded, t_mybonded **bondeds,
99 srenew(*bondeds, *nrbonded+1);
101 /* copy atom numbers */
102 for (j = 0; j < nratoms; j++)
104 (*bondeds)[*nrbonded].a[j] = param->a[j];
107 (*bondeds)[*nrbonded].c = param->C0;
112 static void get_bondeds(int nrat, t_iatom atoms[],
113 at2vsitebond_t *at2vb,
114 int *nrbond, t_mybonded **bonds,
115 int *nrang, t_mybonded **angles,
116 int *nridih, t_mybonded **idihs )
118 int k, i, ftype, nrcheck;
121 for (k = 0; k < nrat; k++)
123 for (i = 0; i < at2vb[atoms[k]].nr; i++)
125 ftype = at2vb[atoms[k]].vsbp[i].ftype;
126 param = at2vb[atoms[k]].vsbp[i].param;
127 nrcheck = vsite_bond_nrcheck(ftype);
128 /* abuse nrcheck to see if we're adding bond, angle or idih */
131 case 2: enter_bonded(nrcheck, nrbond, bonds, param); break;
132 case 3: enter_bonded(nrcheck, nrang, angles, param); break;
133 case 4: enter_bonded(nrcheck, nridih, idihs, param); break;
139 static at2vsitebond_t *make_at2vsitebond(int natoms, t_params plist[])
142 int ftype, i, j, nrcheck, nr;
144 at2vsitebond_t *at2vb;
149 for (ftype = 0; (ftype < F_NRE); ftype++)
151 if ((interaction_function[ftype].flags & IF_VSITE) && ftype != F_VSITEN)
153 for (i = 0; (i < plist[ftype].nr); i++)
155 for (j = 0; j < NRAL(ftype); j++)
157 bVSI[plist[ftype].param[i].a[j]] = TRUE;
163 for (ftype = 0; (ftype < F_NRE); ftype++)
165 nrcheck = vsite_bond_nrcheck(ftype);
168 for (i = 0; (i < plist[ftype].nr); i++)
170 aa = plist[ftype].param[i].a;
171 for (j = 0; j < nrcheck; j++)
175 nr = at2vb[aa[j]].nr;
178 srenew(at2vb[aa[j]].vsbp, nr+10);
180 at2vb[aa[j]].vsbp[nr].ftype = ftype;
181 at2vb[aa[j]].vsbp[nr].param = &plist[ftype].param[i];
193 static void done_at2vsitebond(int natoms, at2vsitebond_t *at2vb)
197 for (i = 0; i < natoms; i++)
201 sfree(at2vb[i].vsbp);
207 static at2vsitecon_t *make_at2vsitecon(int natoms, t_params plist[])
210 int ftype, i, j, ai, aj, nr;
211 at2vsitecon_t *at2vc;
216 for (ftype = 0; (ftype < F_NRE); ftype++)
218 if ((interaction_function[ftype].flags & IF_VSITE) && ftype != F_VSITEN)
220 for (i = 0; (i < plist[ftype].nr); i++)
222 for (j = 0; j < NRAL(ftype); j++)
224 bVSI[plist[ftype].param[i].a[j]] = TRUE;
230 for (ftype = 0; (ftype < F_NRE); ftype++)
232 if (interaction_function[ftype].flags & IF_CONSTRAINT)
234 for (i = 0; (i < plist[ftype].nr); i++)
236 ai = plist[ftype].param[i].AI;
237 aj = plist[ftype].param[i].AJ;
238 if (bVSI[ai] && bVSI[aj])
240 /* Store forward direction */
244 srenew(at2vc[ai].aj, nr+10);
246 at2vc[ai].aj[nr] = aj;
248 /* Store backward direction */
252 srenew(at2vc[aj].aj, nr+10);
254 at2vc[aj].aj[nr] = ai;
265 static void done_at2vsitecon(int natoms, at2vsitecon_t *at2vc)
269 for (i = 0; i < natoms; i++)
280 static void print_bad(FILE *fp,
281 int nrbond, t_mybonded *bonds,
282 int nrang, t_mybonded *angles,
283 int nridih, t_mybonded *idihs )
289 fprintf(fp, "bonds:");
290 for (i = 0; i < nrbond; i++)
292 fprintf(fp, " %u-%u (%g)",
293 bonds[i].AI+1, bonds[i].AJ+1, bonds[i].c);
299 fprintf(fp, "angles:");
300 for (i = 0; i < nrang; i++)
302 fprintf(fp, " %u-%u-%u (%g)",
303 angles[i].AI+1, angles[i].AJ+1,
304 angles[i].AK+1, angles[i].c);
310 fprintf(fp, "idihs:");
311 for (i = 0; i < nridih; i++)
313 fprintf(fp, " %u-%u-%u-%u (%g)",
314 idihs[i].AI+1, idihs[i].AJ+1,
315 idihs[i].AK+1, idihs[i].AL+1, idihs[i].c);
321 static void print_param(FILE *fp, int ftype, int i, t_param *param)
324 static int prev_ftype = NOTSET;
325 static int prev_i = NOTSET;
328 if ( (ftype != prev_ftype) || (i != prev_i) )
334 fprintf(fp, "(%d) plist[%s].param[%d]",
335 pass, interaction_function[ftype].name, i);
336 for (j = 0; j < NRFP(ftype); j++)
338 fprintf(fp, ".c[%d]=%g ", j, param->c[j]);
344 static real get_bond_length(int nrbond, t_mybonded bonds[],
345 t_iatom ai, t_iatom aj)
351 for (i = 0; i < nrbond && (bondlen == NOTSET); i++)
353 /* check both ways */
354 if ( ( (ai == bonds[i].AI) && (aj == bonds[i].AJ) ) ||
355 ( (ai == bonds[i].AJ) && (aj == bonds[i].AI) ) )
357 bondlen = bonds[i].c; /* note: bonds[i].c might be NOTSET */
363 static real get_angle(int nrang, t_mybonded angles[],
364 t_iatom ai, t_iatom aj, t_iatom ak)
370 for (i = 0; i < nrang && (angle == NOTSET); i++)
372 /* check both ways */
373 if ( ( (ai == angles[i].AI) && (aj == angles[i].AJ) && (ak == angles[i].AK) ) ||
374 ( (ai == angles[i].AK) && (aj == angles[i].AJ) && (ak == angles[i].AI) ) )
376 angle = DEG2RAD*angles[i].c;
382 static char *get_atomtype_name_AB(t_atom *atom, gpp_atomtype_t atype)
386 name = get_atomtype_name(atom->type, atype);
388 /* When using the decoupling option, atom types are changed
389 * to decoupled for the non-bonded interactions, but the virtual
390 * sites constructions should be based on the "normal" interactions.
391 * So we return the state B atom type name if the state A atom
392 * type is the decoupled one. We should actually check for the atom
393 * type number, but that's not passed here. So we check for
394 * the decoupled atom type name. This should not cause trouble
395 * as this code is only used for topologies with v-sites without
396 * parameters generated by pdb2gmx.
398 if (strcmp(name, "decoupled") == 0)
400 name = get_atomtype_name(atom->typeB, atype);
406 static gmx_bool calc_vsite3_param(gpp_atomtype_t atype,
407 t_param *param, t_atoms *at,
408 int nrbond, t_mybonded *bonds,
409 int nrang, t_mybonded *angles )
411 /* i = virtual site | ,k
412 * j = 1st bonded heavy atom | i-j
413 * k,l = 2nd bonded atoms | `l
416 gmx_bool bXH3, bError;
417 real bjk, bjl, a = -1, b = -1;
418 /* check if this is part of a NH3 , NH2-umbrella or CH3 group,
419 * i.e. if atom k and l are dummy masses (MNH* or MCH3*) */
423 for (i = 0; i < 4; i++)
425 fprintf(debug, "atom %u type %s ",
427 get_atomtype_name_AB(&at->atom[param->a[i]], atype));
429 fprintf(debug, "\n");
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 = 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 = sqrt( sqr(bCM) - sqr(rM) );
463 /* are we dealing with the X atom? */
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 * cos(aCNH);
479 rH = bNH * 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);
496 fprintf(debug, "params for vsite3 %u: %g %g\n",
497 param->AI+1, param->C0, param->C1);
503 static gmx_bool calc_vsite3fd_param(t_param *param,
504 int nrbond, t_mybonded *bonds,
505 int nrang, t_mybonded *angles)
507 /* i = virtual site | ,k
508 * j = 1st bonded heavy atom | i-j
509 * k,l = 2nd bonded atoms | `l
513 real bij, bjk, bjl, aijk, aijl, rk, rl;
515 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
516 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
517 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
518 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
519 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
520 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) ||
521 (aijk == NOTSET) || (aijl == NOTSET);
523 rk = bjk * sin(aijk);
524 rl = bjl * sin(aijl);
525 param->C0 = rk / (rk + rl);
526 param->C1 = -bij; /* 'bond'-length for fixed distance vsite */
530 fprintf(debug, "params for vsite3fd %u: %g %g\n",
531 param->AI+1, param->C0, param->C1);
536 static gmx_bool calc_vsite3fad_param(t_param *param,
537 int nrbond, t_mybonded *bonds,
538 int nrang, t_mybonded *angles)
540 /* i = virtual site |
541 * j = 1st bonded heavy atom | i-j
542 * k = 2nd bonded heavy atom | `k-l
543 * l = 3d bonded heavy atom |
546 gmx_bool bSwapParity, bError;
549 bSwapParity = ( param->C1 == -1 );
551 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
552 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
553 bError = (bij == NOTSET) || (aijk == NOTSET);
555 param->C1 = bij; /* 'bond'-length for fixed distance vsite */
556 param->C0 = RAD2DEG*aijk; /* 'bond'-angle for fixed angle vsite */
560 param->C0 = 360 - param->C0;
565 fprintf(debug, "params for vsite3fad %u: %g %g\n",
566 param->AI+1, param->C0, param->C1);
571 static gmx_bool calc_vsite3out_param(gpp_atomtype_t atype,
572 t_param *param, t_atoms *at,
573 int nrbond, t_mybonded *bonds,
574 int nrang, t_mybonded *angles)
576 /* i = virtual site | ,k
577 * j = 1st bonded heavy atom | i-j
578 * k,l = 2nd bonded atoms | `l
579 * NOTE: i is out of the j-k-l plane!
582 gmx_bool bXH3, bError, bSwapParity;
583 real bij, bjk, bjl, aijk, aijl, akjl, pijk, pijl, a, b, c;
585 /* check if this is part of a NH2-umbrella, NH3 or CH3 group,
586 * i.e. if atom k and l are dummy masses (MNH* or MCH3*) */
590 for (i = 0; i < 4; i++)
592 fprintf(debug, "atom %u type %s ",
593 param->a[i]+1, get_atomtype_name_AB(&at->atom[param->a[i]], atype));
595 fprintf(debug, "\n");
598 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AK], atype), "MNH", 3) == 0) &&
599 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AL], atype), "MNH", 3) == 0) ) ||
600 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AK], atype), "MCH3", 4) == 0) &&
601 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AL], atype), "MCH3", 4) == 0) );
603 /* check if construction parity must be swapped */
604 bSwapParity = ( param->C1 == -1 );
606 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
607 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
608 bError = (bjk == NOTSET) || (bjl == NOTSET);
611 /* now we get some XH3 group specific construction */
612 /* note: we call the heavy atom 'C' and the X atom 'N' */
613 real bMM, bCM, bCN, bNH, aCNH, dH, rH, rHx, rHy, dM, rM;
616 /* check if bonds from heavy atom (j) to dummy masses (k,l) are equal: */
617 bError = bError || (bjk != bjl);
619 /* the X atom (C or N) in the XH3 group is the first after the masses: */
620 aN = max(param->AK, param->AL)+1;
622 /* get all bondlengths and angles: */
623 bMM = get_bond_length(nrbond, bonds, param->AK, param->AL);
625 bCN = get_bond_length(nrbond, bonds, param->AJ, aN);
626 bNH = get_bond_length(nrbond, bonds, aN, param->AI);
627 aCNH = get_angle (nrang, angles, param->AJ, aN, param->AI);
629 (bMM == NOTSET) || (bCN == NOTSET) || (bNH == NOTSET) || (aCNH == NOTSET);
632 dH = bCN - bNH * cos(aCNH);
633 rH = bNH * sin(aCNH);
634 /* we assume the H's are symmetrically distributed */
635 rHx = rH*cos(DEG2RAD*30);
636 rHy = rH*sin(DEG2RAD*30);
638 dM = sqrt( sqr(bCM) - sqr(rM) );
639 a = 0.5*( (dH/dM) - (rHy/rM) );
640 b = 0.5*( (dH/dM) + (rHy/rM) );
646 /* this is the general construction */
648 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
649 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
650 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
651 akjl = get_angle (nrang, angles, param->AK, param->AJ, param->AL);
653 (bij == NOTSET) || (aijk == NOTSET) || (aijl == NOTSET) || (akjl == NOTSET);
655 pijk = cos(aijk)*bij;
656 pijl = cos(aijl)*bij;
657 a = ( pijk + (pijk*cos(akjl)-pijl) * cos(akjl) / sqr(sin(akjl)) ) / bjk;
658 b = ( pijl + (pijl*cos(akjl)-pijk) * cos(akjl) / sqr(sin(akjl)) ) / bjl;
659 c = -sqrt( sqr(bij) -
660 ( sqr(pijk) - 2*pijk*pijl*cos(akjl) + sqr(pijl) )
662 / ( bjk*bjl*sin(akjl) );
677 fprintf(debug, "params for vsite3out %u: %g %g %g\n",
678 param->AI+1, param->C0, param->C1, param->C2);
683 static gmx_bool calc_vsite4fd_param(t_param *param,
684 int nrbond, t_mybonded *bonds,
685 int nrang, t_mybonded *angles)
687 /* i = virtual site | ,k
688 * j = 1st bonded heavy atom | i-j-m
689 * k,l,m = 2nd bonded atoms | `l
693 real bij, bjk, bjl, bjm, aijk, aijl, aijm, akjm, akjl;
694 real pk, pl, pm, cosakl, cosakm, sinakl, sinakm, cl, cm;
696 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
697 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
698 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
699 bjm = get_bond_length(nrbond, bonds, param->AJ, param->AM);
700 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
701 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
702 aijm = get_angle (nrang, angles, param->AI, param->AJ, param->AM);
703 akjm = get_angle (nrang, angles, param->AK, param->AJ, param->AM);
704 akjl = get_angle (nrang, angles, param->AK, param->AJ, param->AL);
705 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) || (bjm == NOTSET) ||
706 (aijk == NOTSET) || (aijl == NOTSET) || (aijm == NOTSET) || (akjm == NOTSET) ||
714 cosakl = (cos(akjl) - cos(aijk)*cos(aijl)) / (sin(aijk)*sin(aijl));
715 cosakm = (cos(akjm) - cos(aijk)*cos(aijm)) / (sin(aijk)*sin(aijm));
716 if (cosakl < -1 || cosakl > 1 || cosakm < -1 || cosakm > 1)
718 fprintf(stderr, "virtual site %d: angle ijk = %f, angle ijl = %f, angle ijm = %f\n",
719 param->AI+1, RAD2DEG*aijk, RAD2DEG*aijl, RAD2DEG*aijm);
720 gmx_fatal(FARGS, "invalid construction in calc_vsite4fd for atom %d: "
721 "cosakl=%f, cosakm=%f\n", param->AI+1, cosakl, cosakm);
723 sinakl = sqrt(1-sqr(cosakl));
724 sinakm = sqrt(1-sqr(cosakm));
726 /* note: there is a '+' because of the way the sines are calculated */
727 cl = -pk / ( pl*cosakl - pk + pl*sinakl*(pm*cosakm-pk)/(pm*sinakm) );
728 cm = -pk / ( pm*cosakm - pk + pm*sinakm*(pl*cosakl-pk)/(pl*sinakl) );
735 fprintf(debug, "params for vsite4fd %u: %g %g %g\n",
736 param->AI+1, param->C0, param->C1, param->C2);
745 calc_vsite4fdn_param(t_param *param,
746 int nrbond, t_mybonded *bonds,
747 int nrang, t_mybonded *angles)
749 /* i = virtual site | ,k
750 * j = 1st bonded heavy atom | i-j-m
751 * k,l,m = 2nd bonded atoms | `l
755 real bij, bjk, bjl, bjm, aijk, aijl, aijm;
756 real pk, pl, pm, a, b;
758 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
759 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
760 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
761 bjm = get_bond_length(nrbond, bonds, param->AJ, param->AM);
762 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
763 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
764 aijm = get_angle (nrang, angles, param->AI, param->AJ, param->AM);
766 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) || (bjm == NOTSET) ||
767 (aijk == NOTSET) || (aijl == NOTSET) || (aijm == NOTSET);
772 /* Calculate component of bond j-k along the direction i-j */
775 /* Calculate component of bond j-l along the direction i-j */
778 /* Calculate component of bond j-m along the direction i-j */
781 if (fabs(pl) < 1000*GMX_REAL_MIN || fabs(pm) < 1000*GMX_REAL_MIN)
783 fprintf(stderr, "virtual site %d: angle ijk = %f, angle ijl = %f, angle ijm = %f\n",
784 param->AI+1, RAD2DEG*aijk, RAD2DEG*aijl, RAD2DEG*aijm);
785 gmx_fatal(FARGS, "invalid construction in calc_vsite4fdn for atom %d: "
786 "pl=%f, pm=%f\n", param->AI+1, pl, pm);
798 fprintf(debug, "params for vsite4fdn %u: %g %g %g\n",
799 param->AI+1, param->C0, param->C1, param->C2);
808 int set_vsites(gmx_bool bVerbose, t_atoms *atoms, gpp_atomtype_t atype,
812 int nvsite, nrbond, nrang, nridih, nrset;
813 gmx_bool bFirst, bSet, bERROR;
814 at2vsitebond_t *at2vb;
824 fprintf(debug, "\nCalculating parameters for virtual sites\n");
827 /* Make a reverse list to avoid ninteractions^2 operations */
828 at2vb = make_at2vsitebond(atoms->nr, plist);
830 for (ftype = 0; (ftype < F_NRE); ftype++)
832 if ((interaction_function[ftype].flags & IF_VSITE) && ftype != F_VSITEN)
835 nvsite += plist[ftype].nr;
836 for (i = 0; (i < plist[ftype].nr); i++)
838 /* check if all parameters are set */
840 for (j = 0; j < NRFP(ftype) && bSet; j++)
842 bSet = plist[ftype].param[i].c[j] != NOTSET;
847 fprintf(debug, "bSet=%s ", bool_names[bSet]);
848 print_param(debug, ftype, i, &plist[ftype].param[i]);
852 if (bVerbose && bFirst)
854 fprintf(stderr, "Calculating parameters for virtual sites\n");
858 nrbond = nrang = nridih = 0;
863 /* now set the vsite parameters: */
864 get_bondeds(NRAL(ftype), plist[ftype].param[i].a, at2vb,
865 &nrbond, &bonds, &nrang, &angles, &nridih, &idihs);
868 fprintf(debug, "Found %d bonds, %d angles and %d idihs "
869 "for virtual site %u (%s)\n", nrbond, nrang, nridih,
870 plist[ftype].param[i].AI+1,
871 interaction_function[ftype].longname);
872 print_bad(debug, nrbond, bonds, nrang, angles, nridih, idihs);
878 calc_vsite3_param(atype, &(plist[ftype].param[i]), atoms,
879 nrbond, bonds, nrang, angles);
883 calc_vsite3fd_param(&(plist[ftype].param[i]),
884 nrbond, bonds, nrang, angles);
888 calc_vsite3fad_param(&(plist[ftype].param[i]),
889 nrbond, bonds, nrang, angles);
893 calc_vsite3out_param(atype, &(plist[ftype].param[i]), atoms,
894 nrbond, bonds, nrang, angles);
898 calc_vsite4fd_param(&(plist[ftype].param[i]),
899 nrbond, bonds, nrang, angles);
903 calc_vsite4fdn_param(&(plist[ftype].param[i]),
904 nrbond, bonds, nrang, angles);
907 gmx_fatal(FARGS, "Automatic parameter generation not supported "
909 interaction_function[ftype].longname,
910 plist[ftype].param[i].AI+1);
914 gmx_fatal(FARGS, "Automatic parameter generation not supported "
915 "for %s atom %d for this bonding configuration",
916 interaction_function[ftype].longname,
917 plist[ftype].param[i].AI+1);
924 if (debug && plist[ftype].nr)
926 fprintf(stderr, "Calculated parameters for %d out of %d %s atoms\n",
927 nrset, plist[ftype].nr, interaction_function[ftype].longname);
932 done_at2vsitebond(atoms->nr, at2vb);
937 void set_vsites_ptype(gmx_bool bVerbose, gmx_moltype_t *molt)
946 fprintf(stderr, "Setting particle type to V for virtual sites\n");
950 fprintf(stderr, "checking %d functypes\n", F_NRE);
952 for (ftype = 0; ftype < F_NRE; ftype++)
954 il = &molt->ilist[ftype];
955 if (interaction_function[ftype].flags & IF_VSITE)
957 nra = interaction_function[ftype].nratoms;
963 fprintf(stderr, "doing %d %s virtual sites\n",
964 (nrd / (nra+1)), interaction_function[ftype].longname);
967 for (i = 0; (i < nrd); )
969 /* The virtual site */
971 molt->atoms.atom[avsite].ptype = eptVSite;
985 static void check_vsite_constraints(t_params *plist,
986 int cftype, int vsite_type[])
993 ps = &(plist[cftype]);
994 for (i = 0; (i < ps->nr); i++)
996 for (k = 0; k < 2; k++)
998 atom = ps->param[i].a[k];
999 if (vsite_type[atom] != NOTSET)
1001 fprintf(stderr, "ERROR: Cannot have constraint (%u-%u) with virtual site (%u)\n",
1002 ps->param[i].AI+1, ps->param[i].AJ+1, atom+1);
1009 gmx_fatal(FARGS, "There were %d virtual sites involved in constraints", n);
1013 static void clean_vsite_bonds(t_params *plist, t_pindex pindex[],
1014 int cftype, int vsite_type[])
1016 int ftype, i, j, parnr, k, l, m, n, nvsite, nOut, kept_i, vsnral, vsitetype;
1017 int nconverted, nremoved;
1018 atom_id atom, oatom, constr, at1, at2;
1019 atom_id vsiteatoms[MAXATOMLIST];
1020 gmx_bool bKeep, bRemove, bUsed, bPresent, bThisFD, bThisOUT, bAllFD, bFirstTwo;
1023 if (cftype == F_CONNBONDS)
1028 ps = &(plist[cftype]);
1034 for (i = 0; (i < ps->nr); i++) /* for all bonds in the plist */
1039 /* check if all virtual sites are constructed from the same atoms */
1043 fprintf(debug, "constr %u %u:", ps->param[i].AI+1, ps->param[i].AJ+1);
1045 for (k = 0; (k < 2) && !bKeep && !bRemove; k++)
1047 /* for all atoms in the bond */
1048 atom = ps->param[i].a[k];
1049 if (vsite_type[atom] != NOTSET)
1053 fprintf(debug, " d%d[%d: %d %d %d]", k, atom+1,
1054 plist[pindex[atom].ftype].param[pindex[atom].parnr].AJ+1,
1055 plist[pindex[atom].ftype].param[pindex[atom].parnr].AK+1,
1056 plist[pindex[atom].ftype].param[pindex[atom].parnr].AL+1);
1059 bThisFD = ( (pindex[atom].ftype == F_VSITE3FD ) ||
1060 (pindex[atom].ftype == F_VSITE3FAD) ||
1061 (pindex[atom].ftype == F_VSITE4FD ) ||
1062 (pindex[atom].ftype == F_VSITE4FDN ) );
1063 bThisOUT = ( (pindex[atom].ftype == F_VSITE3OUT) &&
1064 (interaction_function[cftype].flags & IF_CONSTRAINT) );
1065 bAllFD = bAllFD && bThisFD;
1066 if (bThisFD || bThisOUT)
1070 fprintf(debug, " %s", bThisOUT ? "out" : "fd");
1072 oatom = ps->param[i].a[1-k]; /* the other atom */
1073 if (vsite_type[oatom] == NOTSET &&
1074 oatom == plist[pindex[atom].ftype].param[pindex[atom].parnr].AJ)
1076 /* if the other atom isn't a vsite, and it is AI */
1084 fprintf(debug, " D-AI");
1092 /* if this is the first vsite we encounter then
1093 store construction atoms */
1094 vsnral = NRAL(pindex[atom].ftype)-1;
1095 for (m = 0; (m < vsnral); m++)
1098 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1103 /* if it is not the first then
1104 check if this vsite is constructed from the same atoms */
1105 if (vsnral == NRAL(pindex[atom].ftype)-1)
1107 for (m = 0; (m < vsnral) && !bKeep; m++)
1111 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1112 for (n = 0; (n < vsnral) && !bPresent; n++)
1114 if (constr == vsiteatoms[n])
1124 fprintf(debug, " !present");
1134 fprintf(debug, " !same#at");
1148 /* if we have no virtual sites in this bond, keep it */
1153 fprintf(debug, " no vsite");
1158 /* check if all non-vsite atoms are used in construction: */
1160 for (k = 0; (k < 2) && !bKeep; k++) /* for all atoms in the bond */
1162 atom = ps->param[i].a[k];
1163 if (vsite_type[atom] == NOTSET)
1166 for (m = 0; (m < vsnral) && !bUsed; m++)
1168 if (atom == vsiteatoms[m])
1171 bFirstTwo = bFirstTwo && m < 2;
1179 fprintf(debug, " !used");
1185 if (!( bAllFD && bFirstTwo ) )
1187 /* check if all constructing atoms are constrained together */
1188 for (m = 0; m < vsnral && !bKeep; m++) /* all constr. atoms */
1190 at1 = vsiteatoms[m];
1191 at2 = vsiteatoms[(m+1) % vsnral];
1193 for (ftype = 0; ftype < F_NRE; ftype++)
1195 if (interaction_function[ftype].flags & IF_CONSTRAINT)
1197 for (j = 0; (j < plist[ftype].nr) && !bPresent; j++)
1199 /* all constraints until one matches */
1200 bPresent = ( ( (plist[ftype].param[j].AI == at1) &&
1201 (plist[ftype].param[j].AJ == at2) ) ||
1202 ( (plist[ftype].param[j].AI == at2) &&
1203 (plist[ftype].param[j].AJ == at1) ) );
1212 fprintf(debug, " !bonded");
1223 fprintf(debug, " keeping");
1225 /* now copy the bond to the new array */
1226 ps->param[kept_i] = ps->param[i];
1229 else if (IS_CHEMBOND(cftype))
1231 srenew(plist[F_CONNBONDS].param, plist[F_CONNBONDS].nr+1);
1232 plist[F_CONNBONDS].param[plist[F_CONNBONDS].nr] = ps->param[i];
1233 plist[F_CONNBONDS].nr++;
1242 fprintf(debug, "\n");
1248 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1249 nremoved, interaction_function[cftype].longname, kept_i);
1253 fprintf(stderr, "Converted %4d %15ss with virtual sites to connections, %5d left\n",
1254 nconverted, interaction_function[cftype].longname, kept_i);
1258 fprintf(stderr, "Warning: removed %d %ss with vsite with %s construction\n"
1259 " This vsite construction does not guarantee constant "
1261 " If the constructions were generated by pdb2gmx ignore "
1263 nOut, interaction_function[cftype].longname,
1264 interaction_function[F_VSITE3OUT].longname );
1269 static void clean_vsite_angles(t_params *plist, t_pindex pindex[],
1270 int cftype, int vsite_type[],
1271 at2vsitecon_t *at2vc)
1273 int i, j, parnr, k, l, m, n, nvsite, kept_i, vsnral, vsitetype;
1274 atom_id atom, constr, at1, at2;
1275 atom_id vsiteatoms[MAXATOMLIST];
1276 gmx_bool bKeep, bUsed, bPresent, bAll3FAD, bFirstTwo;
1279 ps = &(plist[cftype]);
1282 for (i = 0; (i < ps->nr); i++) /* for all angles in the plist */
1286 /* check if all virtual sites are constructed from the same atoms */
1288 for (k = 0; (k < 3) && !bKeep; k++) /* for all atoms in the angle */
1290 atom = ps->param[i].a[k];
1291 if (vsite_type[atom] != NOTSET)
1294 bAll3FAD = bAll3FAD && (pindex[atom].ftype == F_VSITE3FAD);
1297 /* store construction atoms of first vsite */
1298 vsnral = NRAL(pindex[atom].ftype)-1;
1299 for (m = 0; (m < vsnral); m++)
1302 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1306 /* check if this vsite is constructed from the same atoms */
1307 if (vsnral == NRAL(pindex[atom].ftype)-1)
1309 for (m = 0; (m < vsnral) && !bKeep; m++)
1313 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1314 for (n = 0; (n < vsnral) && !bPresent; n++)
1316 if (constr == vsiteatoms[n])
1334 /* keep all angles with no virtual sites in them or
1335 with virtual sites with more than 3 constr. atoms */
1336 if (nvsite == 0 && vsnral > 3)
1341 /* check if all non-vsite atoms are used in construction: */
1343 for (k = 0; (k < 3) && !bKeep; k++) /* for all atoms in the angle */
1345 atom = ps->param[i].a[k];
1346 if (vsite_type[atom] == NOTSET)
1349 for (m = 0; (m < vsnral) && !bUsed; m++)
1351 if (atom == vsiteatoms[m])
1354 bFirstTwo = bFirstTwo && m < 2;
1364 if (!( bAll3FAD && bFirstTwo ) )
1366 /* check if all constructing atoms are constrained together */
1367 for (m = 0; m < vsnral && !bKeep; m++) /* all constr. atoms */
1369 at1 = vsiteatoms[m];
1370 at2 = vsiteatoms[(m+1) % vsnral];
1372 for (j = 0; j < at2vc[at1].nr; j++)
1374 if (at2vc[at1].aj[j] == at2)
1388 /* now copy the angle to the new array */
1389 ps->param[kept_i] = ps->param[i];
1394 if (ps->nr != kept_i)
1396 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1397 ps->nr-kept_i, interaction_function[cftype].longname, kept_i);
1402 static void clean_vsite_dihs(t_params *plist, t_pindex pindex[],
1403 int cftype, int vsite_type[])
1408 ps = &(plist[cftype]);
1411 for (i = 0; (i < ps->nr); i++) /* for all dihedrals in the plist */
1413 int ftype, parnr, k, l, m, n, nvsite;
1414 int vsnral = 0; /* keep the compiler happy */
1415 atom_id atom, constr;
1416 atom_id vsiteatoms[4] = { 0 }; /* init to zero to make gcc4.8 happy */
1417 gmx_bool bKeep, bUsed, bPresent;
1421 /* check if all virtual sites are constructed from the same atoms */
1423 for (k = 0; (k < 4) && !bKeep; k++) /* for all atoms in the dihedral */
1425 atom = ps->param[i].a[k];
1426 if (vsite_type[atom] != NOTSET)
1430 /* store construction atoms of first vsite */
1431 vsnral = NRAL(pindex[atom].ftype)-1;
1432 assert(vsnral <= 4);
1433 for (m = 0; (m < vsnral); m++)
1436 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1440 fprintf(debug, "dih w. vsite: %u %u %u %u\n",
1441 ps->param[i].AI+1, ps->param[i].AJ+1,
1442 ps->param[i].AK+1, ps->param[i].AL+1);
1443 fprintf(debug, "vsite %u from: %u %u %u\n",
1444 atom+1, vsiteatoms[0]+1, vsiteatoms[1]+1, vsiteatoms[2]+1);
1449 /* check if this vsite is constructed from the same atoms */
1450 if (vsnral == NRAL(pindex[atom].ftype)-1)
1452 for (m = 0; (m < vsnral) && !bKeep; m++)
1456 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1457 for (n = 0; (n < vsnral) && !bPresent; n++)
1459 if (constr == vsiteatoms[n])
1475 /* keep all dihedrals with no virtual sites in them */
1481 /* check if all atoms in dihedral are either virtual sites, or used in
1482 construction of virtual sites. If so, keep it, if not throw away: */
1483 for (k = 0; (k < 4) && !bKeep; k++) /* for all atoms in the dihedral */
1485 atom = ps->param[i].a[k];
1486 if (vsite_type[atom] == NOTSET)
1488 /* vsnral will be set here, we don't get here with nvsite==0 */
1490 for (m = 0; (m < vsnral) && !bUsed; m++)
1492 if (atom == vsiteatoms[m])
1502 fprintf(debug, "unused atom in dih: %u\n", atom+1);
1510 ps->param[kept_i] = ps->param[i];
1515 if (ps->nr != kept_i)
1517 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1518 ps->nr-kept_i, interaction_function[cftype].longname, kept_i);
1523 void clean_vsite_bondeds(t_params *plist, int natoms, gmx_bool bRmVSiteBds)
1525 int i, k, nvsite, ftype, vsite, parnr;
1528 at2vsitecon_t *at2vc;
1530 pindex = 0; /* avoid warnings */
1531 /* make vsite_type array */
1532 snew(vsite_type, natoms);
1533 for (i = 0; i < natoms; i++)
1535 vsite_type[i] = NOTSET;
1538 for (ftype = 0; ftype < F_NRE; ftype++)
1540 if (interaction_function[ftype].flags & IF_VSITE)
1542 nvsite += plist[ftype].nr;
1544 while (i < plist[ftype].nr)
1546 vsite = plist[ftype].param[i].AI;
1547 if (vsite_type[vsite] == NOTSET)
1549 vsite_type[vsite] = ftype;
1553 gmx_fatal(FARGS, "multiple vsite constructions for atom %d", vsite+1);
1555 if (ftype == F_VSITEN)
1557 while (i < plist[ftype].nr && plist[ftype].param[i].AI == vsite)
1570 /* the rest only if we have virtual sites: */
1573 fprintf(stderr, "Cleaning up constraints %swith virtual sites\n",
1574 bRmVSiteBds ? "and constant bonded interactions " : "");
1576 /* Make a reverse list to avoid ninteractions^2 operations */
1577 at2vc = make_at2vsitecon(natoms, plist);
1579 snew(pindex, natoms);
1580 for (ftype = 0; ftype < F_NRE; ftype++)
1582 if ((interaction_function[ftype].flags & IF_VSITE) &&
1585 for (parnr = 0; (parnr < plist[ftype].nr); parnr++)
1587 k = plist[ftype].param[parnr].AI;
1588 pindex[k].ftype = ftype;
1589 pindex[k].parnr = parnr;
1596 for (i = 0; i < natoms; i++)
1598 fprintf(debug, "atom %d vsite_type %s\n", i,
1599 vsite_type[i] == NOTSET ? "NOTSET" :
1600 interaction_function[vsite_type[i]].name);
1604 /* remove things with vsite atoms */
1605 for (ftype = 0; ftype < F_NRE; ftype++)
1607 if ( ( ( interaction_function[ftype].flags & IF_BOND ) && bRmVSiteBds ) ||
1608 ( interaction_function[ftype].flags & IF_CONSTRAINT ) )
1610 if (interaction_function[ftype].flags & (IF_BTYPE | IF_CONSTRAINT) )
1612 clean_vsite_bonds (plist, pindex, ftype, vsite_type);
1614 else if (interaction_function[ftype].flags & IF_ATYPE)
1616 clean_vsite_angles(plist, pindex, ftype, vsite_type, at2vc);
1618 else if ( (ftype == F_PDIHS) || (ftype == F_IDIHS) )
1620 clean_vsite_dihs (plist, pindex, ftype, vsite_type);
1624 /* check if we have constraints left with virtual sites in them */
1625 for (ftype = 0; ftype < F_NRE; ftype++)
1627 if (interaction_function[ftype].flags & IF_CONSTRAINT)
1629 check_vsite_constraints(plist, ftype, vsite_type);
1633 done_at2vsitecon(natoms, at2vc);