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, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "vsite_parm.h"
54 #include "gmx_fatal.h"
72 vsitebondparam_t *vsbp;
80 static int vsite_bond_nrcheck(int ftype)
84 if ((interaction_function[ftype].flags & (IF_BTYPE | IF_CONSTRAINT | IF_ATYPE)) || (ftype == F_IDIHS))
86 nrcheck = NRAL(ftype);
96 static void enter_bonded(int nratoms, int *nrbonded, t_mybonded **bondeds,
101 srenew(*bondeds, *nrbonded+1);
103 /* copy atom numbers */
104 for (j = 0; j < nratoms; j++)
106 (*bondeds)[*nrbonded].a[j] = param->a[j];
109 (*bondeds)[*nrbonded].c = param->C0;
114 static void get_bondeds(int nrat, t_iatom atoms[],
115 at2vsitebond_t *at2vb,
116 int *nrbond, t_mybonded **bonds,
117 int *nrang, t_mybonded **angles,
118 int *nridih, t_mybonded **idihs )
120 int k, i, ftype, nrcheck;
123 for (k = 0; k < nrat; k++)
125 for (i = 0; i < at2vb[atoms[k]].nr; i++)
127 ftype = at2vb[atoms[k]].vsbp[i].ftype;
128 param = at2vb[atoms[k]].vsbp[i].param;
129 nrcheck = vsite_bond_nrcheck(ftype);
130 /* abuse nrcheck to see if we're adding bond, angle or idih */
133 case 2: enter_bonded(nrcheck, nrbond, bonds, param); break;
134 case 3: enter_bonded(nrcheck, nrang, angles, param); break;
135 case 4: enter_bonded(nrcheck, nridih, idihs, param); break;
141 static at2vsitebond_t *make_at2vsitebond(int natoms, t_params plist[])
144 int ftype, i, j, nrcheck, nr;
146 at2vsitebond_t *at2vb;
151 for (ftype = 0; (ftype < F_NRE); ftype++)
153 if ((interaction_function[ftype].flags & IF_VSITE) && ftype != F_VSITEN)
155 for (i = 0; (i < plist[ftype].nr); i++)
157 for (j = 0; j < NRAL(ftype); j++)
159 bVSI[plist[ftype].param[i].a[j]] = TRUE;
165 for (ftype = 0; (ftype < F_NRE); ftype++)
167 nrcheck = vsite_bond_nrcheck(ftype);
170 for (i = 0; (i < plist[ftype].nr); i++)
172 aa = plist[ftype].param[i].a;
173 for (j = 0; j < nrcheck; j++)
177 nr = at2vb[aa[j]].nr;
180 srenew(at2vb[aa[j]].vsbp, nr+10);
182 at2vb[aa[j]].vsbp[nr].ftype = ftype;
183 at2vb[aa[j]].vsbp[nr].param = &plist[ftype].param[i];
195 static void done_at2vsitebond(int natoms, at2vsitebond_t *at2vb)
199 for (i = 0; i < natoms; i++)
203 sfree(at2vb[i].vsbp);
209 static at2vsitecon_t *make_at2vsitecon(int natoms, t_params plist[])
212 int ftype, i, j, ai, aj, nr;
213 at2vsitecon_t *at2vc;
218 for (ftype = 0; (ftype < F_NRE); ftype++)
220 if ((interaction_function[ftype].flags & IF_VSITE) && ftype != F_VSITEN)
222 for (i = 0; (i < plist[ftype].nr); i++)
224 for (j = 0; j < NRAL(ftype); j++)
226 bVSI[plist[ftype].param[i].a[j]] = TRUE;
232 for (ftype = 0; (ftype < F_NRE); ftype++)
234 if (interaction_function[ftype].flags & IF_CONSTRAINT)
236 for (i = 0; (i < plist[ftype].nr); i++)
238 ai = plist[ftype].param[i].AI;
239 aj = plist[ftype].param[i].AJ;
240 if (bVSI[ai] && bVSI[aj])
242 /* Store forward direction */
246 srenew(at2vc[ai].aj, nr+10);
248 at2vc[ai].aj[nr] = aj;
250 /* Store backward direction */
254 srenew(at2vc[aj].aj, nr+10);
256 at2vc[aj].aj[nr] = ai;
267 static void done_at2vsitecon(int natoms, at2vsitecon_t *at2vc)
271 for (i = 0; i < natoms; i++)
282 static void print_bad(FILE *fp,
283 int nrbond, t_mybonded *bonds,
284 int nrang, t_mybonded *angles,
285 int nridih, t_mybonded *idihs )
291 fprintf(fp, "bonds:");
292 for (i = 0; i < nrbond; i++)
294 fprintf(fp, " %u-%u (%g)",
295 bonds[i].AI+1, bonds[i].AJ+1, bonds[i].c);
301 fprintf(fp, "angles:");
302 for (i = 0; i < nrang; i++)
304 fprintf(fp, " %u-%u-%u (%g)",
305 angles[i].AI+1, angles[i].AJ+1,
306 angles[i].AK+1, angles[i].c);
312 fprintf(fp, "idihs:");
313 for (i = 0; i < nridih; i++)
315 fprintf(fp, " %u-%u-%u-%u (%g)",
316 idihs[i].AI+1, idihs[i].AJ+1,
317 idihs[i].AK+1, idihs[i].AL+1, idihs[i].c);
323 static void print_param(FILE *fp, int ftype, int i, t_param *param)
326 static int prev_ftype = NOTSET;
327 static int prev_i = NOTSET;
330 if ( (ftype != prev_ftype) || (i != prev_i) )
336 fprintf(fp, "(%d) plist[%s].param[%d]",
337 pass, interaction_function[ftype].name, i);
338 for (j = 0; j < NRFP(ftype); j++)
340 fprintf(fp, ".c[%d]=%g ", j, param->c[j]);
346 static real get_bond_length(int nrbond, t_mybonded bonds[],
347 t_iatom ai, t_iatom aj)
353 for (i = 0; i < nrbond && (bondlen == NOTSET); i++)
355 /* check both ways */
356 if ( ( (ai == bonds[i].AI) && (aj == bonds[i].AJ) ) ||
357 ( (ai == bonds[i].AJ) && (aj == bonds[i].AI) ) )
359 bondlen = bonds[i].c; /* note: bonds[i].c might be NOTSET */
365 static real get_angle(int nrang, t_mybonded angles[],
366 t_iatom ai, t_iatom aj, t_iatom ak)
372 for (i = 0; i < nrang && (angle == NOTSET); i++)
374 /* check both ways */
375 if ( ( (ai == angles[i].AI) && (aj == angles[i].AJ) && (ak == angles[i].AK) ) ||
376 ( (ai == angles[i].AK) && (aj == angles[i].AJ) && (ak == angles[i].AI) ) )
378 angle = DEG2RAD*angles[i].c;
384 static char *get_atomtype_name_AB(t_atom *atom, gpp_atomtype_t atype)
388 name = get_atomtype_name(atom->type, atype);
390 /* When using the decoupling option, atom types are changed
391 * to decoupled for the non-bonded interactions, but the virtual
392 * sites constructions should be based on the "normal" interactions.
393 * So we return the state B atom type name if the state A atom
394 * type is the decoupled one. We should actually check for the atom
395 * type number, but that's not passed here. So we check for
396 * the decoupled atom type name. This should not cause trouble
397 * as this code is only used for topologies with v-sites without
398 * parameters generated by pdb2gmx.
400 if (strcmp(name, "decoupled") == 0)
402 name = get_atomtype_name(atom->typeB, atype);
408 static gmx_bool calc_vsite3_param(gpp_atomtype_t atype,
409 t_param *param, t_atoms *at,
410 int nrbond, t_mybonded *bonds,
411 int nrang, t_mybonded *angles )
413 /* i = virtual site | ,k
414 * j = 1st bonded heavy atom | i-j
415 * k,l = 2nd bonded atoms | `l
418 gmx_bool bXH3, bError;
419 real bjk, bjl, a = -1, b = -1;
420 /* check if this is part of a NH3 , NH2-umbrella or CH3 group,
421 * i.e. if atom k and l are dummy masses (MNH* or MCH3*) */
425 for (i = 0; i < 4; i++)
427 fprintf(debug, "atom %u type %s ",
429 get_atomtype_name_AB(&at->atom[param->a[i]], atype));
431 fprintf(debug, "\n");
434 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AK], atype), "MNH", 3) == 0) &&
435 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AL], atype), "MNH", 3) == 0) ) ||
436 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AK], atype), "MCH3", 4) == 0) &&
437 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AL], atype), "MCH3", 4) == 0) );
439 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
440 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
441 bError = (bjk == NOTSET) || (bjl == NOTSET);
444 /* now we get some XH2/XH3 group specific construction */
445 /* note: we call the heavy atom 'C' and the X atom 'N' */
446 real bMM, bCM, bCN, bNH, aCNH, dH, rH, dM, rM;
449 /* check if bonds from heavy atom (j) to dummy masses (k,l) are equal: */
450 bError = bError || (bjk != bjl);
452 /* the X atom (C or N) in the XH2/XH3 group is the first after the masses: */
453 aN = max(param->AK, param->AL)+1;
455 /* get common bonds */
456 bMM = get_bond_length(nrbond, bonds, param->AK, param->AL);
458 bCN = get_bond_length(nrbond, bonds, param->AJ, aN);
459 bError = bError || (bMM == NOTSET) || (bCN == NOTSET);
461 /* calculate common things */
463 dM = sqrt( sqr(bCM) - sqr(rM) );
465 /* are we dealing with the X atom? */
468 /* this is trivial */
469 a = b = 0.5 * bCN/dM;
474 /* get other bondlengths and angles: */
475 bNH = get_bond_length(nrbond, bonds, aN, param->AI);
476 aCNH = get_angle (nrang, angles, param->AJ, aN, param->AI);
477 bError = bError || (bNH == NOTSET) || (aCNH == NOTSET);
480 dH = bCN - bNH * cos(aCNH);
481 rH = bNH * sin(aCNH);
483 a = 0.5 * ( dH/dM + rH/rM );
484 b = 0.5 * ( dH/dM - rH/rM );
489 gmx_fatal(FARGS, "calc_vsite3_param not implemented for the general case "
490 "(atom %d)", param->AI+1);
498 fprintf(debug, "params for vsite3 %u: %g %g\n",
499 param->AI+1, param->C0, param->C1);
505 static gmx_bool calc_vsite3fd_param(t_param *param,
506 int nrbond, t_mybonded *bonds,
507 int nrang, t_mybonded *angles)
509 /* i = virtual site | ,k
510 * j = 1st bonded heavy atom | i-j
511 * k,l = 2nd bonded atoms | `l
515 real bij, bjk, bjl, aijk, aijl, rk, rl;
517 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
518 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
519 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
520 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
521 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
522 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) ||
523 (aijk == NOTSET) || (aijl == NOTSET);
525 rk = bjk * sin(aijk);
526 rl = bjl * sin(aijl);
527 param->C0 = rk / (rk + rl);
528 param->C1 = -bij; /* 'bond'-length for fixed distance vsite */
532 fprintf(debug, "params for vsite3fd %u: %g %g\n",
533 param->AI+1, param->C0, param->C1);
538 static gmx_bool calc_vsite3fad_param(t_param *param,
539 int nrbond, t_mybonded *bonds,
540 int nrang, t_mybonded *angles)
542 /* i = virtual site |
543 * j = 1st bonded heavy atom | i-j
544 * k = 2nd bonded heavy atom | `k-l
545 * l = 3d bonded heavy atom |
548 gmx_bool bSwapParity, bError;
551 bSwapParity = ( param->C1 == -1 );
553 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
554 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
555 bError = (bij == NOTSET) || (aijk == NOTSET);
557 param->C1 = bij; /* 'bond'-length for fixed distance vsite */
558 param->C0 = RAD2DEG*aijk; /* 'bond'-angle for fixed angle vsite */
562 param->C0 = 360 - param->C0;
567 fprintf(debug, "params for vsite3fad %u: %g %g\n",
568 param->AI+1, param->C0, param->C1);
573 static gmx_bool calc_vsite3out_param(gpp_atomtype_t atype,
574 t_param *param, t_atoms *at,
575 int nrbond, t_mybonded *bonds,
576 int nrang, t_mybonded *angles)
578 /* i = virtual site | ,k
579 * j = 1st bonded heavy atom | i-j
580 * k,l = 2nd bonded atoms | `l
581 * NOTE: i is out of the j-k-l plane!
584 gmx_bool bXH3, bError, bSwapParity;
585 real bij, bjk, bjl, aijk, aijl, akjl, pijk, pijl, a, b, c;
587 /* check if this is part of a NH2-umbrella, NH3 or CH3 group,
588 * i.e. if atom k and l are dummy masses (MNH* or MCH3*) */
592 for (i = 0; i < 4; i++)
594 fprintf(debug, "atom %u type %s ",
595 param->a[i]+1, get_atomtype_name_AB(&at->atom[param->a[i]], atype));
597 fprintf(debug, "\n");
600 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AK], atype), "MNH", 3) == 0) &&
601 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AL], atype), "MNH", 3) == 0) ) ||
602 ( (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AK], atype), "MCH3", 4) == 0) &&
603 (gmx_strncasecmp(get_atomtype_name_AB(&at->atom[param->AL], atype), "MCH3", 4) == 0) );
605 /* check if construction parity must be swapped */
606 bSwapParity = ( param->C1 == -1 );
608 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
609 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
610 bError = (bjk == NOTSET) || (bjl == NOTSET);
613 /* now we get some XH3 group specific construction */
614 /* note: we call the heavy atom 'C' and the X atom 'N' */
615 real bMM, bCM, bCN, bNH, aCNH, dH, rH, rHx, rHy, dM, rM;
618 /* check if bonds from heavy atom (j) to dummy masses (k,l) are equal: */
619 bError = bError || (bjk != bjl);
621 /* the X atom (C or N) in the XH3 group is the first after the masses: */
622 aN = max(param->AK, param->AL)+1;
624 /* get all bondlengths and angles: */
625 bMM = get_bond_length(nrbond, bonds, param->AK, param->AL);
627 bCN = get_bond_length(nrbond, bonds, param->AJ, aN);
628 bNH = get_bond_length(nrbond, bonds, aN, param->AI);
629 aCNH = get_angle (nrang, angles, param->AJ, aN, param->AI);
631 (bMM == NOTSET) || (bCN == NOTSET) || (bNH == NOTSET) || (aCNH == NOTSET);
634 dH = bCN - bNH * cos(aCNH);
635 rH = bNH * sin(aCNH);
636 /* we assume the H's are symmetrically distributed */
637 rHx = rH*cos(DEG2RAD*30);
638 rHy = rH*sin(DEG2RAD*30);
640 dM = sqrt( sqr(bCM) - sqr(rM) );
641 a = 0.5*( (dH/dM) - (rHy/rM) );
642 b = 0.5*( (dH/dM) + (rHy/rM) );
648 /* this is the general construction */
650 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
651 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
652 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
653 akjl = get_angle (nrang, angles, param->AK, param->AJ, param->AL);
655 (bij == NOTSET) || (aijk == NOTSET) || (aijl == NOTSET) || (akjl == NOTSET);
657 pijk = cos(aijk)*bij;
658 pijl = cos(aijl)*bij;
659 a = ( pijk + (pijk*cos(akjl)-pijl) * cos(akjl) / sqr(sin(akjl)) ) / bjk;
660 b = ( pijl + (pijl*cos(akjl)-pijk) * cos(akjl) / sqr(sin(akjl)) ) / bjl;
661 c = -sqrt( sqr(bij) -
662 ( sqr(pijk) - 2*pijk*pijl*cos(akjl) + sqr(pijl) )
664 / ( bjk*bjl*sin(akjl) );
679 fprintf(debug, "params for vsite3out %u: %g %g %g\n",
680 param->AI+1, param->C0, param->C1, param->C2);
685 static gmx_bool calc_vsite4fd_param(t_param *param,
686 int nrbond, t_mybonded *bonds,
687 int nrang, t_mybonded *angles)
689 /* i = virtual site | ,k
690 * j = 1st bonded heavy atom | i-j-m
691 * k,l,m = 2nd bonded atoms | `l
695 real bij, bjk, bjl, bjm, aijk, aijl, aijm, akjm, akjl;
696 real pk, pl, pm, cosakl, cosakm, sinakl, sinakm, cl, cm;
698 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
699 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
700 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
701 bjm = get_bond_length(nrbond, bonds, param->AJ, param->AM);
702 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
703 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
704 aijm = get_angle (nrang, angles, param->AI, param->AJ, param->AM);
705 akjm = get_angle (nrang, angles, param->AK, param->AJ, param->AM);
706 akjl = get_angle (nrang, angles, param->AK, param->AJ, param->AL);
707 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) || (bjm == NOTSET) ||
708 (aijk == NOTSET) || (aijl == NOTSET) || (aijm == NOTSET) || (akjm == NOTSET) ||
716 cosakl = (cos(akjl) - cos(aijk)*cos(aijl)) / (sin(aijk)*sin(aijl));
717 cosakm = (cos(akjm) - cos(aijk)*cos(aijm)) / (sin(aijk)*sin(aijm));
718 if (cosakl < -1 || cosakl > 1 || cosakm < -1 || cosakm > 1)
720 fprintf(stderr, "virtual site %d: angle ijk = %f, angle ijl = %f, angle ijm = %f\n",
721 param->AI+1, RAD2DEG*aijk, RAD2DEG*aijl, RAD2DEG*aijm);
722 gmx_fatal(FARGS, "invalid construction in calc_vsite4fd for atom %d: "
723 "cosakl=%f, cosakm=%f\n", param->AI+1, cosakl, cosakm);
725 sinakl = sqrt(1-sqr(cosakl));
726 sinakm = sqrt(1-sqr(cosakm));
728 /* note: there is a '+' because of the way the sines are calculated */
729 cl = -pk / ( pl*cosakl - pk + pl*sinakl*(pm*cosakm-pk)/(pm*sinakm) );
730 cm = -pk / ( pm*cosakm - pk + pm*sinakm*(pl*cosakl-pk)/(pl*sinakl) );
737 fprintf(debug, "params for vsite4fd %u: %g %g %g\n",
738 param->AI+1, param->C0, param->C1, param->C2);
747 calc_vsite4fdn_param(t_param *param,
748 int nrbond, t_mybonded *bonds,
749 int nrang, t_mybonded *angles)
751 /* i = virtual site | ,k
752 * j = 1st bonded heavy atom | i-j-m
753 * k,l,m = 2nd bonded atoms | `l
757 real bij, bjk, bjl, bjm, aijk, aijl, aijm;
758 real pk, pl, pm, a, b;
760 bij = get_bond_length(nrbond, bonds, param->AI, param->AJ);
761 bjk = get_bond_length(nrbond, bonds, param->AJ, param->AK);
762 bjl = get_bond_length(nrbond, bonds, param->AJ, param->AL);
763 bjm = get_bond_length(nrbond, bonds, param->AJ, param->AM);
764 aijk = get_angle (nrang, angles, param->AI, param->AJ, param->AK);
765 aijl = get_angle (nrang, angles, param->AI, param->AJ, param->AL);
766 aijm = get_angle (nrang, angles, param->AI, param->AJ, param->AM);
768 bError = (bij == NOTSET) || (bjk == NOTSET) || (bjl == NOTSET) || (bjm == NOTSET) ||
769 (aijk == NOTSET) || (aijl == NOTSET) || (aijm == NOTSET);
774 /* Calculate component of bond j-k along the direction i-j */
777 /* Calculate component of bond j-l along the direction i-j */
780 /* Calculate component of bond j-m along the direction i-j */
783 if (fabs(pl) < 1000*GMX_REAL_MIN || fabs(pm) < 1000*GMX_REAL_MIN)
785 fprintf(stderr, "virtual site %d: angle ijk = %f, angle ijl = %f, angle ijm = %f\n",
786 param->AI+1, RAD2DEG*aijk, RAD2DEG*aijl, RAD2DEG*aijm);
787 gmx_fatal(FARGS, "invalid construction in calc_vsite4fdn for atom %d: "
788 "pl=%f, pm=%f\n", param->AI+1, pl, pm);
800 fprintf(debug, "params for vsite4fdn %u: %g %g %g\n",
801 param->AI+1, param->C0, param->C1, param->C2);
810 int set_vsites(gmx_bool bVerbose, t_atoms *atoms, gpp_atomtype_t atype,
814 int nvsite, nrbond, nrang, nridih, nrset;
815 gmx_bool bFirst, bSet, bERROR;
816 at2vsitebond_t *at2vb;
826 fprintf(debug, "\nCalculating parameters for virtual sites\n");
829 /* Make a reverse list to avoid ninteractions^2 operations */
830 at2vb = make_at2vsitebond(atoms->nr, plist);
832 for (ftype = 0; (ftype < F_NRE); ftype++)
834 if ((interaction_function[ftype].flags & IF_VSITE) && ftype != F_VSITEN)
837 nvsite += plist[ftype].nr;
838 for (i = 0; (i < plist[ftype].nr); i++)
840 /* check if all parameters are set */
842 for (j = 0; j < NRFP(ftype) && bSet; j++)
844 bSet = plist[ftype].param[i].c[j] != NOTSET;
849 fprintf(debug, "bSet=%s ", bool_names[bSet]);
850 print_param(debug, ftype, i, &plist[ftype].param[i]);
854 if (bVerbose && bFirst)
856 fprintf(stderr, "Calculating parameters for virtual sites\n");
860 nrbond = nrang = nridih = 0;
865 /* now set the vsite parameters: */
866 get_bondeds(NRAL(ftype), plist[ftype].param[i].a, at2vb,
867 &nrbond, &bonds, &nrang, &angles, &nridih, &idihs);
870 fprintf(debug, "Found %d bonds, %d angles and %d idihs "
871 "for virtual site %u (%s)\n", nrbond, nrang, nridih,
872 plist[ftype].param[i].AI+1,
873 interaction_function[ftype].longname);
874 print_bad(debug, nrbond, bonds, nrang, angles, nridih, idihs);
880 calc_vsite3_param(atype, &(plist[ftype].param[i]), atoms,
881 nrbond, bonds, nrang, angles);
885 calc_vsite3fd_param(&(plist[ftype].param[i]),
886 nrbond, bonds, nrang, angles);
890 calc_vsite3fad_param(&(plist[ftype].param[i]),
891 nrbond, bonds, nrang, angles);
895 calc_vsite3out_param(atype, &(plist[ftype].param[i]), atoms,
896 nrbond, bonds, nrang, angles);
900 calc_vsite4fd_param(&(plist[ftype].param[i]),
901 nrbond, bonds, nrang, angles);
905 calc_vsite4fdn_param(&(plist[ftype].param[i]),
906 nrbond, bonds, nrang, angles);
909 gmx_fatal(FARGS, "Automatic parameter generation not supported "
911 interaction_function[ftype].longname,
912 plist[ftype].param[i].AI+1);
916 gmx_fatal(FARGS, "Automatic parameter generation not supported "
917 "for %s atom %d for this bonding configuration",
918 interaction_function[ftype].longname,
919 plist[ftype].param[i].AI+1);
926 if (debug && plist[ftype].nr)
928 fprintf(stderr, "Calculated parameters for %d out of %d %s atoms\n",
929 nrset, plist[ftype].nr, interaction_function[ftype].longname);
934 done_at2vsitebond(atoms->nr, at2vb);
939 void set_vsites_ptype(gmx_bool bVerbose, gmx_moltype_t *molt)
948 fprintf(stderr, "Setting particle type to V for virtual sites\n");
952 fprintf(stderr, "checking %d functypes\n", F_NRE);
954 for (ftype = 0; ftype < F_NRE; ftype++)
956 il = &molt->ilist[ftype];
957 if (interaction_function[ftype].flags & IF_VSITE)
959 nra = interaction_function[ftype].nratoms;
965 fprintf(stderr, "doing %d %s virtual sites\n",
966 (nrd / (nra+1)), interaction_function[ftype].longname);
969 for (i = 0; (i < nrd); )
971 /* The virtual site */
973 molt->atoms.atom[avsite].ptype = eptVSite;
987 static void check_vsite_constraints(t_params *plist,
988 int cftype, int vsite_type[])
995 ps = &(plist[cftype]);
996 for (i = 0; (i < ps->nr); i++)
998 for (k = 0; k < 2; k++)
1000 atom = ps->param[i].a[k];
1001 if (vsite_type[atom] != NOTSET)
1003 fprintf(stderr, "ERROR: Cannot have constraint (%u-%u) with virtual site (%u)\n",
1004 ps->param[i].AI+1, ps->param[i].AJ+1, atom+1);
1011 gmx_fatal(FARGS, "There were %d virtual sites involved in constraints", n);
1015 static void clean_vsite_bonds(t_params *plist, t_pindex pindex[],
1016 int cftype, int vsite_type[])
1018 int ftype, i, j, parnr, k, l, m, n, nvsite, nOut, kept_i, vsnral, vsitetype;
1019 int nconverted, nremoved;
1020 atom_id atom, oatom, constr, at1, at2;
1021 atom_id vsiteatoms[MAXATOMLIST];
1022 gmx_bool bKeep, bRemove, bUsed, bPresent, bThisFD, bThisOUT, bAllFD, bFirstTwo;
1025 if (cftype == F_CONNBONDS)
1030 ps = &(plist[cftype]);
1036 for (i = 0; (i < ps->nr); i++) /* for all bonds in the plist */
1041 /* check if all virtual sites are constructed from the same atoms */
1045 fprintf(debug, "constr %u %u:", ps->param[i].AI+1, ps->param[i].AJ+1);
1047 for (k = 0; (k < 2) && !bKeep && !bRemove; k++)
1049 /* for all atoms in the bond */
1050 atom = ps->param[i].a[k];
1051 if (vsite_type[atom] != NOTSET)
1055 fprintf(debug, " d%d[%d: %d %d %d]", k, atom+1,
1056 plist[pindex[atom].ftype].param[pindex[atom].parnr].AJ+1,
1057 plist[pindex[atom].ftype].param[pindex[atom].parnr].AK+1,
1058 plist[pindex[atom].ftype].param[pindex[atom].parnr].AL+1);
1061 bThisFD = ( (pindex[atom].ftype == F_VSITE3FD ) ||
1062 (pindex[atom].ftype == F_VSITE3FAD) ||
1063 (pindex[atom].ftype == F_VSITE4FD ) ||
1064 (pindex[atom].ftype == F_VSITE4FDN ) );
1065 bThisOUT = ( (pindex[atom].ftype == F_VSITE3OUT) &&
1066 (interaction_function[cftype].flags & IF_CONSTRAINT) );
1067 bAllFD = bAllFD && bThisFD;
1068 if (bThisFD || bThisOUT)
1072 fprintf(debug, " %s", bThisOUT ? "out" : "fd");
1074 oatom = ps->param[i].a[1-k]; /* the other atom */
1075 if (vsite_type[oatom] == NOTSET &&
1076 oatom == plist[pindex[atom].ftype].param[pindex[atom].parnr].AJ)
1078 /* if the other atom isn't a vsite, and it is AI */
1086 fprintf(debug, " D-AI");
1094 /* if this is the first vsite we encounter then
1095 store construction atoms */
1096 vsnral = NRAL(pindex[atom].ftype)-1;
1097 for (m = 0; (m < vsnral); m++)
1100 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1105 /* if it is not the first then
1106 check if this vsite is constructed from the same atoms */
1107 if (vsnral == NRAL(pindex[atom].ftype)-1)
1109 for (m = 0; (m < vsnral) && !bKeep; m++)
1113 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1114 for (n = 0; (n < vsnral) && !bPresent; n++)
1116 if (constr == vsiteatoms[n])
1126 fprintf(debug, " !present");
1136 fprintf(debug, " !same#at");
1150 /* if we have no virtual sites in this bond, keep it */
1155 fprintf(debug, " no vsite");
1160 /* check if all non-vsite atoms are used in construction: */
1162 for (k = 0; (k < 2) && !bKeep; k++) /* for all atoms in the bond */
1164 atom = ps->param[i].a[k];
1165 if (vsite_type[atom] == NOTSET)
1168 for (m = 0; (m < vsnral) && !bUsed; m++)
1170 if (atom == vsiteatoms[m])
1173 bFirstTwo = bFirstTwo && m < 2;
1181 fprintf(debug, " !used");
1187 if (!( bAllFD && bFirstTwo ) )
1189 /* check if all constructing atoms are constrained together */
1190 for (m = 0; m < vsnral && !bKeep; m++) /* all constr. atoms */
1192 at1 = vsiteatoms[m];
1193 at2 = vsiteatoms[(m+1) % vsnral];
1195 for (ftype = 0; ftype < F_NRE; ftype++)
1197 if (interaction_function[ftype].flags & IF_CONSTRAINT)
1199 for (j = 0; (j < plist[ftype].nr) && !bPresent; j++)
1201 /* all constraints until one matches */
1202 bPresent = ( ( (plist[ftype].param[j].AI == at1) &&
1203 (plist[ftype].param[j].AJ == at2) ) ||
1204 ( (plist[ftype].param[j].AI == at2) &&
1205 (plist[ftype].param[j].AJ == at1) ) );
1214 fprintf(debug, " !bonded");
1225 fprintf(debug, " keeping");
1227 /* now copy the bond to the new array */
1228 ps->param[kept_i] = ps->param[i];
1231 else if (IS_CHEMBOND(cftype))
1233 srenew(plist[F_CONNBONDS].param, plist[F_CONNBONDS].nr+1);
1234 plist[F_CONNBONDS].param[plist[F_CONNBONDS].nr] = ps->param[i];
1235 plist[F_CONNBONDS].nr++;
1244 fprintf(debug, "\n");
1250 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1251 nremoved, interaction_function[cftype].longname, kept_i);
1255 fprintf(stderr, "Converted %4d %15ss with virtual sites to connections, %5d left\n",
1256 nconverted, interaction_function[cftype].longname, kept_i);
1260 fprintf(stderr, "Warning: removed %d %ss with vsite with %s construction\n"
1261 " This vsite construction does not guarantee constant "
1263 " If the constructions were generated by pdb2gmx ignore "
1265 nOut, interaction_function[cftype].longname,
1266 interaction_function[F_VSITE3OUT].longname );
1271 static void clean_vsite_angles(t_params *plist, t_pindex pindex[],
1272 int cftype, int vsite_type[],
1273 at2vsitecon_t *at2vc)
1275 int i, j, parnr, k, l, m, n, nvsite, kept_i, vsnral, vsitetype;
1276 atom_id atom, constr, at1, at2;
1277 atom_id vsiteatoms[MAXATOMLIST];
1278 gmx_bool bKeep, bUsed, bPresent, bAll3FAD, bFirstTwo;
1281 ps = &(plist[cftype]);
1284 for (i = 0; (i < ps->nr); i++) /* for all angles in the plist */
1288 /* check if all virtual sites are constructed from the same atoms */
1290 for (k = 0; (k < 3) && !bKeep; k++) /* for all atoms in the angle */
1292 atom = ps->param[i].a[k];
1293 if (vsite_type[atom] != NOTSET)
1296 bAll3FAD = bAll3FAD && (pindex[atom].ftype == F_VSITE3FAD);
1299 /* store construction atoms of first vsite */
1300 vsnral = NRAL(pindex[atom].ftype)-1;
1301 for (m = 0; (m < vsnral); m++)
1304 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1308 /* check if this vsite is constructed from the same atoms */
1309 if (vsnral == NRAL(pindex[atom].ftype)-1)
1311 for (m = 0; (m < vsnral) && !bKeep; m++)
1315 plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
1316 for (n = 0; (n < vsnral) && !bPresent; n++)
1318 if (constr == vsiteatoms[n])
1336 /* keep all angles with no virtual sites in them or
1337 with virtual sites with more than 3 constr. atoms */
1338 if (nvsite == 0 && vsnral > 3)
1343 /* check if all non-vsite atoms are used in construction: */
1345 for (k = 0; (k < 3) && !bKeep; k++) /* for all atoms in the angle */
1347 atom = ps->param[i].a[k];
1348 if (vsite_type[atom] == NOTSET)
1351 for (m = 0; (m < vsnral) && !bUsed; m++)
1353 if (atom == vsiteatoms[m])
1356 bFirstTwo = bFirstTwo && m < 2;
1366 if (!( bAll3FAD && bFirstTwo ) )
1368 /* check if all constructing atoms are constrained together */
1369 for (m = 0; m < vsnral && !bKeep; m++) /* all constr. atoms */
1371 at1 = vsiteatoms[m];
1372 at2 = vsiteatoms[(m+1) % vsnral];
1374 for (j = 0; j < at2vc[at1].nr; j++)
1376 if (at2vc[at1].aj[j] == at2)
1390 /* now copy the angle to the new array */
1391 ps->param[kept_i] = ps->param[i];
1396 if (ps->nr != kept_i)
1398 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1399 ps->nr-kept_i, interaction_function[cftype].longname, kept_i);
1404 static void clean_vsite_dihs(t_params *plist, t_pindex pindex[],
1405 int cftype, int vsite_type[])
1407 int ftype, i, parnr, k, l, m, n, nvsite, kept_i, vsnral;
1408 atom_id atom, constr;
1409 atom_id vsiteatoms[4];
1410 gmx_bool bKeep, bUsed, bPresent;
1413 ps = &(plist[cftype]);
1417 for (i = 0; (i < ps->nr); i++) /* for all dihedrals in the plist */
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)
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);
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])
1472 /* keep all dihedrals with no virtual sites in them */
1478 /* check if all atoms in dihedral are either virtual sites, or used in
1479 construction of virtual sites. If so, keep it, if not throw away: */
1480 for (k = 0; (k < 4) && !bKeep; k++) /* for all atoms in the dihedral */
1482 atom = ps->param[i].a[k];
1483 if (vsite_type[atom] == NOTSET)
1486 for (m = 0; (m < vsnral) && !bUsed; m++)
1488 if (atom == vsiteatoms[m])
1498 fprintf(debug, "unused atom in dih: %u\n", atom+1);
1506 ps->param[kept_i] = ps->param[i];
1511 if (ps->nr != kept_i)
1513 fprintf(stderr, "Removed %4d %15ss with virtual sites, %5d left\n",
1514 ps->nr-kept_i, interaction_function[cftype].longname, kept_i);
1519 void clean_vsite_bondeds(t_params *plist, int natoms, gmx_bool bRmVSiteBds)
1521 int i, k, nvsite, ftype, vsite, parnr;
1524 at2vsitecon_t *at2vc;
1526 pindex = 0; /* avoid warnings */
1527 /* make vsite_type array */
1528 snew(vsite_type, natoms);
1529 for (i = 0; i < natoms; i++)
1531 vsite_type[i] = NOTSET;
1534 for (ftype = 0; ftype < F_NRE; ftype++)
1536 if (interaction_function[ftype].flags & IF_VSITE)
1538 nvsite += plist[ftype].nr;
1540 while (i < plist[ftype].nr)
1542 vsite = plist[ftype].param[i].AI;
1543 if (vsite_type[vsite] == NOTSET)
1545 vsite_type[vsite] = ftype;
1549 gmx_fatal(FARGS, "multiple vsite constructions for atom %d", vsite+1);
1551 if (ftype == F_VSITEN)
1553 while (i < plist[ftype].nr && plist[ftype].param[i].AI == vsite)
1566 /* the rest only if we have virtual sites: */
1569 fprintf(stderr, "Cleaning up constraints %swith virtual sites\n",
1570 bRmVSiteBds ? "and constant bonded interactions " : "");
1572 /* Make a reverse list to avoid ninteractions^2 operations */
1573 at2vc = make_at2vsitecon(natoms, plist);
1575 snew(pindex, natoms);
1576 for (ftype = 0; ftype < F_NRE; ftype++)
1578 if ((interaction_function[ftype].flags & IF_VSITE) &&
1581 for (parnr = 0; (parnr < plist[ftype].nr); parnr++)
1583 k = plist[ftype].param[parnr].AI;
1584 pindex[k].ftype = ftype;
1585 pindex[k].parnr = parnr;
1592 for (i = 0; i < natoms; i++)
1594 fprintf(debug, "atom %d vsite_type %s\n", i,
1595 vsite_type[i] == NOTSET ? "NOTSET" :
1596 interaction_function[vsite_type[i]].name);
1600 /* remove things with vsite atoms */
1601 for (ftype = 0; ftype < F_NRE; ftype++)
1603 if ( ( ( interaction_function[ftype].flags & IF_BOND ) && bRmVSiteBds ) ||
1604 ( interaction_function[ftype].flags & IF_CONSTRAINT ) )
1606 if (interaction_function[ftype].flags & (IF_BTYPE | IF_CONSTRAINT) )
1608 clean_vsite_bonds (plist, pindex, ftype, vsite_type);
1610 else if (interaction_function[ftype].flags & IF_ATYPE)
1612 clean_vsite_angles(plist, pindex, ftype, vsite_type, at2vc);
1614 else if ( (ftype == F_PDIHS) || (ftype == F_IDIHS) )
1616 clean_vsite_dihs (plist, pindex, ftype, vsite_type);
1620 /* check if we have constraints left with virtual sites in them */
1621 for (ftype = 0; ftype < F_NRE; ftype++)
1623 if (interaction_function[ftype].flags & IF_CONSTRAINT)
1625 check_vsite_constraints(plist, ftype, vsite_type);
1629 done_at2vsitecon(natoms, at2vc);