2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/fileio/futil.h"
47 #include "gromacs/fileio/gmxfio.h"
50 #include "types/commrec.h"
52 #include "gromacs/fileio/filenm.h"
54 #include "gromacs/utility/smalloc.h"
59 #include "mtop_util.h"
61 #include "gmx_ga2la.h"
66 static void pull_print_group_x(FILE *out, ivec dim, const t_pull_group *pgrp)
70 for (m = 0; m < DIM; m++)
74 fprintf(out, "\t%g", pgrp->x[m]);
79 static void pull_print_coord_dr(FILE *out, ivec dim, const t_pull_coord *pcrd)
83 for (m = 0; m < DIM; m++)
87 fprintf(out, "\t%g", pcrd->dr[m]);
92 static void pull_print_x(FILE *out, t_pull *pull, double t)
95 const t_pull_coord *pcrd;
97 fprintf(out, "%.4f", t);
99 for (c = 0; c < pull->ncoord; c++)
101 pcrd = &pull->coord[c];
107 pull_print_group_x(out, pull->dim, &pull->dyna[c]);
111 pull_print_group_x(out, pull->dim, &pull->group[pcrd->group[0]]);
114 pull_print_coord_dr(out, pull->dim, pcrd);
119 static void pull_print_f(FILE *out, t_pull *pull, double t)
123 fprintf(out, "%.4f", t);
125 for (c = 0; c < pull->ncoord; c++)
127 fprintf(out, "\t%g", pull->coord[c].f_scal);
132 void pull_print_output(t_pull *pull, gmx_int64_t step, double time)
134 if ((pull->nstxout != 0) && (step % pull->nstxout == 0))
136 pull_print_x(pull->out_x, pull, time);
139 if ((pull->nstfout != 0) && (step % pull->nstfout == 0))
141 pull_print_f(pull->out_f, pull, time);
145 static FILE *open_pull_out(const char *fn, t_pull *pull, const output_env_t oenv,
146 gmx_bool bCoord, unsigned long Flags)
150 char **setname, buf[10];
152 if (Flags & MD_APPENDFILES)
154 fp = gmx_fio_fopen(fn, "a+");
158 fp = gmx_fio_fopen(fn, "w+");
161 xvgr_header(fp, "Pull COM", "Time (ps)", "Position (nm)",
166 xvgr_header(fp, "Pull force", "Time (ps)", "Force (kJ/mol/nm)",
170 snew(setname, 2*pull->ncoord*DIM);
172 for (c = 0; c < pull->ncoord; c++)
178 for (m = 0; m < DIM; m++)
182 sprintf(buf, "%d %s%c", c+1, "c", 'X'+m);
183 setname[nsets] = strdup(buf);
188 for (m = 0; m < DIM; m++)
192 sprintf(buf, "%d %s%c", c+1, "d", 'X'+m);
193 setname[nsets] = strdup(buf);
200 sprintf(buf, "%d", c+1);
201 setname[nsets] = strdup(buf);
207 xvgr_legend(fp, nsets, (const char**)setname, oenv);
209 for (c = 0; c < nsets; c++)
219 /* Apply forces in a mass weighted fashion */
220 static void apply_forces_grp(const t_pull_group *pgrp, const t_mdatoms *md,
221 const dvec f_pull, int sign, rvec *f)
224 double wmass, inv_wm;
226 inv_wm = pgrp->wscale*pgrp->invtm;
228 for (i = 0; i < pgrp->nat_loc; i++)
230 ii = pgrp->ind_loc[i];
231 wmass = md->massT[ii];
232 if (pgrp->weight_loc)
234 wmass *= pgrp->weight_loc[i];
237 for (m = 0; m < DIM; m++)
239 f[ii][m] += sign * wmass * f_pull[m] * inv_wm;
244 /* Apply forces in a mass weighted fashion */
245 static void apply_forces(t_pull * pull, t_mdatoms * md, rvec *f)
248 const t_pull_coord *pcrd;
250 for (c = 0; c < pull->ncoord; c++)
252 pcrd = &pull->coord[c];
256 apply_forces_grp(&pull->dyna[c], md, pcrd->f, -1, f);
260 if (pull->group[pcrd->group[0]].nat > 0)
262 apply_forces_grp(&pull->group[pcrd->group[0]], md, pcrd->f, -1, f);
265 apply_forces_grp(&pull->group[pcrd->group[1]], md, pcrd->f, 1, f);
269 static double max_pull_distance2(const t_pull *pull, const t_pbc *pbc)
274 max_d2 = GMX_DOUBLE_MAX;
276 if (pull->eGeom != epullgDIRPBC)
278 for (m = 0; m < pbc->ndim_ePBC; m++)
280 if (pull->dim[m] != 0)
282 max_d2 = min(max_d2, norm2(pbc->box[m]));
290 static void low_get_pull_coord_dr(const t_pull *pull,
291 const t_pull_coord *pcrd,
292 const t_pbc *pbc, double t,
293 dvec xg, dvec xref, double max_dist2,
296 const t_pull_group *pgrp0, *pgrp1;
298 dvec xrefr, dref = {0, 0, 0};
301 pgrp0 = &pull->group[pcrd->group[0]];
302 pgrp1 = &pull->group[pcrd->group[1]];
304 /* Only the first group can be an absolute reference, in that case nat=0 */
307 for (m = 0; m < DIM; m++)
309 xref[m] = pcrd->origin[m];
313 copy_dvec(xref, xrefr);
315 if (pull->eGeom == epullgDIRPBC)
317 for (m = 0; m < DIM; m++)
319 dref[m] = (pcrd->init + pcrd->rate*t)*pcrd->vec[m];
321 /* Add the reference position, so we use the correct periodic image */
322 dvec_inc(xrefr, dref);
325 pbc_dx_d(pbc, xg, xrefr, dr);
327 for (m = 0; m < DIM; m++)
329 dr[m] *= pull->dim[m];
332 if (max_dist2 >= 0 && dr2 > 0.98*0.98*max_dist2)
334 gmx_fatal(FARGS, "Distance between pull groups %d and %d (%f nm) is larger than 0.49 times the box size (%f)",
335 pcrd->group[0], pcrd->group[1], sqrt(dr2), sqrt(max_dist2));
338 if (pull->eGeom == epullgDIRPBC)
344 static void get_pull_coord_dr(const t_pull *pull,
346 const t_pbc *pbc, double t,
350 const t_pull_coord *pcrd;
352 if (pull->eGeom == epullgDIRPBC)
358 md2 = max_pull_distance2(pull, pbc);
361 pcrd = &pull->coord[coord_ind];
363 low_get_pull_coord_dr(pull, pcrd, pbc, t,
364 pull->group[pcrd->group[1]].x,
365 PULL_CYL(pull) ? pull->dyna[coord_ind].x : pull->group[pcrd->group[0]].x,
370 void get_pull_coord_distance(const t_pull *pull,
372 const t_pbc *pbc, double t,
373 dvec dr, double *dev)
375 static gmx_bool bWarned = FALSE; /* TODO: this should be fixed for thread-safety,
376 but is fairly benign */
377 const t_pull_coord *pcrd;
379 double ref, drs, inpr;
381 pcrd = &pull->coord[coord_ind];
383 get_pull_coord_dr(pull, coord_ind, pbc, t, dr);
385 ref = pcrd->init + pcrd->rate*t;
390 /* Pull along the vector between the com's */
391 if (ref < 0 && !bWarned)
393 fprintf(stderr, "\nPull reference distance for coordinate %d is negative (%f)\n", coord_ind+1, ref);
399 /* With no vector we can not determine the direction for the force,
400 * so we set the force to zero.
406 /* Determine the deviation */
415 for (m = 0; m < DIM; m++)
417 inpr += pcrd->vec[m]*dr[m];
424 void clear_pull_forces(t_pull *pull)
428 /* Zeroing the forces is only required for constraint pulling.
429 * It can happen that multiple constraint steps need to be applied
430 * and therefore the constraint forces need to be accumulated.
432 for (i = 0; i < pull->ncoord; i++)
434 clear_dvec(pull->coord[i].f);
435 pull->coord[i].f_scal = 0;
439 /* Apply constraint using SHAKE */
440 static void do_constraint(t_pull *pull, t_pbc *pbc,
442 gmx_bool bMaster, tensor vir,
446 dvec *r_ij; /* x[i] com of i in prev. step. Obeys constr. -> r_ij[i] */
447 dvec unc_ij; /* xp[i] com of i this step, before constr. -> unc_ij */
448 dvec *rnew; /* current 'new' positions of the groups */
449 double *dr_tot; /* the total update of the coords */
453 double lambda, rm, mass, invdt = 0;
454 gmx_bool bConverged_all, bConverged = FALSE;
455 int niter = 0, g, c, ii, j, m, max_iter = 100;
457 dvec f; /* the pull force */
459 t_pull_group *pdyna, *pgrp0, *pgrp1;
462 snew(r_ij, pull->ncoord);
463 snew(dr_tot, pull->ncoord);
465 snew(rnew, pull->ngroup);
467 /* copy the current unconstrained positions for use in iterations. We
468 iterate until rinew[i] and rjnew[j] obey the constraints. Then
469 rinew - pull.x_unc[i] is the correction dr to group i */
470 for (g = 0; g < pull->ngroup; g++)
472 copy_dvec(pull->group[g].xp, rnew[g]);
476 /* There is only one pull coordinate and reference group */
477 copy_dvec(pull->dyna[0].xp, rnew[pull->coord[0].group[0]]);
480 /* Determine the constraint directions from the old positions */
481 for (c = 0; c < pull->ncoord; c++)
483 get_pull_coord_dr(pull, c, pbc, t, r_ij[c]);
484 /* Store the difference vector at time t for printing */
485 copy_dvec(r_ij[c], pull->coord[c].dr);
488 fprintf(debug, "Pull coord %d dr %f %f %f\n",
489 c, r_ij[c][XX], r_ij[c][YY], r_ij[c][ZZ]);
492 if (pull->eGeom == epullgDIR || pull->eGeom == epullgDIRPBC)
494 /* Select the component along vec */
496 for (m = 0; m < DIM; m++)
498 a += pull->coord[c].vec[m]*r_ij[c][m];
500 for (m = 0; m < DIM; m++)
502 r_ij[c][m] = a*pull->coord[c].vec[m];
507 bConverged_all = FALSE;
508 while (!bConverged_all && niter < max_iter)
510 bConverged_all = TRUE;
512 /* loop over all constraints */
513 for (c = 0; c < pull->ncoord; c++)
517 pcrd = &pull->coord[c];
518 pgrp0 = &pull->group[pcrd->group[0]];
519 pgrp1 = &pull->group[pcrd->group[1]];
521 /* Get the current difference vector */
522 low_get_pull_coord_dr(pull, pcrd, pbc, t,
523 rnew[pcrd->group[1]],
524 rnew[pcrd->group[0]],
527 ref = pcrd->init + pcrd->rate*t;
531 fprintf(debug, "Pull coord %d, iteration %d\n", c, niter);
534 rm = 1.0/(pgrp0->invtm + pgrp1->invtm);
541 gmx_fatal(FARGS, "The pull constraint reference distance for group %d is <= 0 (%f)", c, ref);
545 double q, c_a, c_b, c_c;
547 c_a = diprod(r_ij[c], r_ij[c]);
548 c_b = diprod(unc_ij, r_ij[c])*2;
549 c_c = diprod(unc_ij, unc_ij) - dsqr(ref);
553 q = -0.5*(c_b - sqrt(c_b*c_b - 4*c_a*c_c));
558 q = -0.5*(c_b + sqrt(c_b*c_b - 4*c_a*c_c));
565 "Pull ax^2+bx+c=0: a=%e b=%e c=%e lambda=%e\n",
566 c_a, c_b, c_c, lambda);
570 /* The position corrections dr due to the constraints */
571 dsvmul(-lambda*rm*pgrp1->invtm, r_ij[c], dr1);
572 dsvmul( lambda*rm*pgrp0->invtm, r_ij[c], dr0);
573 dr_tot[c] += -lambda*dnorm(r_ij[c]);
578 /* A 1-dimensional constraint along a vector */
580 for (m = 0; m < DIM; m++)
582 vec[m] = pcrd->vec[m];
583 a += unc_ij[m]*vec[m];
585 /* Select only the component along the vector */
586 dsvmul(a, vec, unc_ij);
590 fprintf(debug, "Pull inpr %e lambda: %e\n", a, lambda);
593 /* The position corrections dr due to the constraints */
594 dsvmul(-lambda*rm*pgrp1->invtm, vec, dr1);
595 dsvmul( lambda*rm*pgrp0->invtm, vec, dr0);
596 dr_tot[c] += -lambda;
607 low_get_pull_coord_dr(pull, pcrd, pbc, t, rnew[g1], rnew[g0], -1, tmp);
608 low_get_pull_coord_dr(pull, pcrd, pbc, t, dr1, dr0, -1, tmp3);
610 "Pull cur %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
611 rnew[g0][0], rnew[g0][1], rnew[g0][2],
612 rnew[g1][0], rnew[g1][1], rnew[g1][2], dnorm(tmp));
614 "Pull ref %8s %8s %8s %8s %8s %8s d: %8.5f\n",
615 "", "", "", "", "", "", ref);
617 "Pull cor %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
618 dr0[0], dr0[1], dr0[2],
619 dr1[0], dr1[1], dr1[2],
623 /* Update the COMs with dr */
624 dvec_inc(rnew[pcrd->group[1]], dr1);
625 dvec_inc(rnew[pcrd->group[0]], dr0);
628 /* Check if all constraints are fullfilled now */
629 for (c = 0; c < pull->ncoord; c++)
631 pcrd = &pull->coord[c];
633 low_get_pull_coord_dr(pull, pcrd, pbc, t,
634 rnew[pcrd->group[1]],
635 rnew[pcrd->group[0]],
641 bConverged = fabs(dnorm(unc_ij) - ref) < pull->constr_tol;
646 for (m = 0; m < DIM; m++)
648 vec[m] = pcrd->vec[m];
650 inpr = diprod(unc_ij, vec);
651 dsvmul(inpr, vec, unc_ij);
653 fabs(diprod(unc_ij, vec) - ref) < pull->constr_tol;
661 fprintf(debug, "NOT CONVERGED YET: Group %d:"
662 "d_ref = %f, current d = %f\n",
663 g, ref, dnorm(unc_ij));
666 bConverged_all = FALSE;
671 /* if after all constraints are dealt with and bConverged is still TRUE
672 we're finished, if not we do another iteration */
674 if (niter > max_iter)
676 gmx_fatal(FARGS, "Too many iterations for constraint run: %d", niter);
679 /* DONE ITERATING, NOW UPDATE COORDINATES AND CALC. CONSTRAINT FORCES */
686 /* update atoms in the groups */
687 for (g = 0; g < pull->ngroup; g++)
689 const t_pull_group *pgrp;
692 if (PULL_CYL(pull) && g == pull->coord[0].group[0])
694 pgrp = &pull->dyna[0];
698 pgrp = &pull->group[g];
701 /* get the final constraint displacement dr for group g */
702 dvec_sub(rnew[g], pgrp->xp, dr);
703 /* select components of dr */
704 for (m = 0; m < DIM; m++)
706 dr[m] *= pull->dim[m];
709 /* update the atom positions */
711 for (j = 0; j < pgrp->nat_loc; j++)
713 ii = pgrp->ind_loc[j];
714 if (pgrp->weight_loc)
716 dsvmul(pgrp->wscale*pgrp->weight_loc[j], dr, tmp);
718 for (m = 0; m < DIM; m++)
724 for (m = 0; m < DIM; m++)
726 v[ii][m] += invdt*tmp[m];
732 /* calculate the constraint forces, used for output and virial only */
733 for (c = 0; c < pull->ncoord; c++)
735 pcrd = &pull->coord[c];
736 pcrd->f_scal = dr_tot[c]/((pull->group[pcrd->group[0]].invtm + pull->group[pcrd->group[1]].invtm)*dt*dt);
742 /* Add the pull contribution to the virial */
743 f_invr = pcrd->f_scal/dnorm(r_ij[c]);
745 for (j = 0; j < DIM; j++)
747 for (m = 0; m < DIM; m++)
749 vir[j][m] -= 0.5*f_invr*r_ij[c][j]*r_ij[c][m];
755 /* finished! I hope. Give back some memory */
761 /* Pulling with a harmonic umbrella potential or constant force */
762 static void do_pull_pot(int ePull,
763 t_pull *pull, t_pbc *pbc, double t, real lambda,
764 real *V, tensor vir, real *dVdl)
767 double dev, ndr, invdr;
771 /* loop over the pull coordinates */
774 for (c = 0; c < pull->ncoord; c++)
776 pcrd = &pull->coord[c];
778 get_pull_coord_distance(pull, c, pbc, t, pcrd->dr, &dev);
780 k = (1.0 - lambda)*pcrd->k + lambda*pcrd->kB;
781 dkdl = pcrd->kB - pcrd->k;
786 ndr = dnorm(pcrd->dr);
788 if (ePull == epullUMBRELLA)
790 pcrd->f_scal = -k*dev;
791 *V += 0.5* k*dsqr(dev);
792 *dVdl += 0.5*dkdl*dsqr(dev);
800 for (m = 0; m < DIM; m++)
802 pcrd->f[m] = pcrd->f_scal*pcrd->dr[m]*invdr;
808 if (ePull == epullUMBRELLA)
810 pcrd->f_scal = -k*dev;
811 *V += 0.5* k*dsqr(dev);
812 *dVdl += 0.5*dkdl*dsqr(dev);
817 for (m = 0; m < DIM; m++)
819 ndr += pcrd->vec[m]*pcrd->dr[m];
825 for (m = 0; m < DIM; m++)
827 pcrd->f[m] = pcrd->f_scal*pcrd->vec[m];
834 /* Add the pull contribution to the virial */
835 for (j = 0; j < DIM; j++)
837 for (m = 0; m < DIM; m++)
839 vir[j][m] -= 0.5*pcrd->f[j]*pcrd->dr[m];
846 real pull_potential(int ePull, t_pull *pull, t_mdatoms *md, t_pbc *pbc,
847 t_commrec *cr, double t, real lambda,
848 rvec *x, rvec *f, tensor vir, real *dvdlambda)
852 pull_calc_coms(cr, pull, md, pbc, t, x, NULL);
854 do_pull_pot(ePull, pull, pbc, t, lambda,
855 &V, pull->bVirial && MASTER(cr) ? vir : NULL, &dVdl);
857 /* Distribute forces over pulled groups */
858 apply_forces(pull, md, f);
865 return (MASTER(cr) ? V : 0.0);
868 void pull_constraint(t_pull *pull, t_mdatoms *md, t_pbc *pbc,
869 t_commrec *cr, double dt, double t,
870 rvec *x, rvec *xp, rvec *v, tensor vir)
872 pull_calc_coms(cr, pull, md, pbc, t, x, xp);
874 do_constraint(pull, pbc, xp, v, pull->bVirial && MASTER(cr), vir, dt, t);
877 static void make_local_pull_group(gmx_ga2la_t ga2la,
878 t_pull_group *pg, int start, int end)
883 for (i = 0; i < pg->nat; i++)
888 if (!ga2la_get_home(ga2la, ii, &ii))
893 if (ii >= start && ii < end)
895 /* This is a home atom, add it to the local pull group */
896 if (pg->nat_loc >= pg->nalloc_loc)
898 pg->nalloc_loc = over_alloc_dd(pg->nat_loc+1);
899 srenew(pg->ind_loc, pg->nalloc_loc);
900 if (pg->epgrppbc == epgrppbcCOS || pg->weight)
902 srenew(pg->weight_loc, pg->nalloc_loc);
905 pg->ind_loc[pg->nat_loc] = ii;
908 pg->weight_loc[pg->nat_loc] = pg->weight[i];
915 void dd_make_local_pull_groups(gmx_domdec_t *dd, t_pull *pull, t_mdatoms *md)
929 for (g = 0; g < pull->ngroup; g++)
931 make_local_pull_group(ga2la, &pull->group[g],
936 static void init_pull_group_index(FILE *fplog, t_commrec *cr,
937 int g, t_pull_group *pg, ivec pulldims,
938 gmx_mtop_t *mtop, t_inputrec *ir, real lambda)
940 int i, ii, d, nfrozen, ndim;
942 double tmass, wmass, wwmass;
943 gmx_groups_t *groups;
944 gmx_mtop_atomlookup_t alook;
947 if (EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD)
949 /* There are no masses in the integrator.
950 * But we still want to have the correct mass-weighted COMs.
951 * So we store the real masses in the weights.
952 * We do not set nweight, so these weights do not end up in the tpx file.
954 if (pg->nweight == 0)
956 snew(pg->weight, pg->nat);
965 pg->weight_loc = NULL;
969 pg->nat_loc = pg->nat;
970 pg->ind_loc = pg->ind;
971 if (pg->epgrppbc == epgrppbcCOS)
973 snew(pg->weight_loc, pg->nat);
977 pg->weight_loc = pg->weight;
981 groups = &mtop->groups;
983 alook = gmx_mtop_atomlookup_init(mtop);
989 for (i = 0; i < pg->nat; i++)
992 gmx_mtop_atomnr_to_atom(alook, ii, &atom);
993 if (ir->opts.nFreeze)
995 for (d = 0; d < DIM; d++)
997 if (pulldims[d] && ir->opts.nFreeze[ggrpnr(groups, egcFREEZE, ii)][d])
1003 if (ir->efep == efepNO)
1009 m = (1 - lambda)*atom->m + lambda*atom->mB;
1011 if (pg->nweight > 0)
1019 if (EI_ENERGY_MINIMIZATION(ir->eI))
1021 /* Move the mass to the weight */
1026 else if (ir->eI == eiBD)
1030 mbd = ir->bd_fric*ir->delta_t;
1034 if (groups->grpnr[egcTC] == NULL)
1036 mbd = ir->delta_t/ir->opts.tau_t[0];
1040 mbd = ir->delta_t/ir->opts.tau_t[groups->grpnr[egcTC][ii]];
1052 gmx_mtop_atomlookup_destroy(alook);
1056 gmx_fatal(FARGS, "The total%s mass of pull group %d is zero",
1057 pg->weight ? " weighted" : "", g);
1062 "Pull group %d: %5d atoms, mass %9.3f", g, pg->nat, tmass);
1063 if (pg->weight || EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD)
1065 fprintf(fplog, ", weighted mass %9.3f", wmass*wmass/wwmass);
1067 if (pg->epgrppbc == epgrppbcCOS)
1069 fprintf(fplog, ", cosine weighting will be used");
1071 fprintf(fplog, "\n");
1076 /* A value > 0 signals not frozen, it is updated later */
1082 for (d = 0; d < DIM; d++)
1084 ndim += pulldims[d]*pg->nat;
1086 if (fplog && nfrozen > 0 && nfrozen < ndim)
1089 "\nWARNING: In pull group %d some, but not all of the degrees of freedom\n"
1090 " that are subject to pulling are frozen.\n"
1091 " For pulling the whole group will be frozen.\n\n",
1099 void init_pull(FILE *fplog, t_inputrec *ir, int nfile, const t_filenm fnm[],
1100 gmx_mtop_t *mtop, t_commrec *cr, const output_env_t oenv, real lambda,
1101 gmx_bool bOutFile, unsigned long Flags)
1105 int c, g, start = 0, end = 0, m;
1109 pull->ePBC = ir->ePBC;
1112 case epbcNONE: pull->npbcdim = 0; break;
1113 case epbcXY: pull->npbcdim = 2; break;
1114 default: pull->npbcdim = 3; break;
1119 gmx_bool bAbs, bCos;
1122 for (c = 0; c < pull->ncoord; c++)
1124 if (pull->group[pull->coord[c].group[0]].nat == 0 ||
1125 pull->group[pull->coord[c].group[1]].nat == 0)
1131 fprintf(fplog, "\nWill apply %s COM pulling in geometry '%s'\n",
1132 EPULLTYPE(ir->ePull), EPULLGEOM(pull->eGeom));
1133 fprintf(fplog, "with %d pull coordinate%s and %d group%s\n",
1134 pull->ncoord, pull->ncoord == 1 ? "" : "s",
1135 pull->ngroup, pull->ngroup == 1 ? "" : "s");
1138 fprintf(fplog, "with an absolute reference\n");
1141 for (g = 0; g < pull->ngroup; g++)
1143 if (pull->group[g].nat > 1 &&
1144 pull->group[g].pbcatom < 0)
1146 /* We are using cosine weighting */
1147 fprintf(fplog, "Cosine weighting is used for group %d\n", g);
1153 please_cite(fplog, "Engin2010");
1157 /* We always add the virial contribution,
1158 * except for geometry = direction_periodic where this is impossible.
1160 pull->bVirial = (pull->eGeom != epullgDIRPBC);
1161 if (getenv("GMX_NO_PULLVIR") != NULL)
1165 fprintf(fplog, "Found env. var., will not add the virial contribution of the COM pull forces\n");
1167 pull->bVirial = FALSE;
1172 pull->dbuf_cyl = NULL;
1173 pull->bRefAt = FALSE;
1175 for (g = 0; g < pull->ngroup; g++)
1177 pgrp = &pull->group[g];
1178 pgrp->epgrppbc = epgrppbcNONE;
1181 /* Determine if we need to take PBC into account for calculating
1182 * the COM's of the pull groups.
1184 for (m = 0; m < pull->npbcdim; m++)
1186 if (pull->dim[m] && pgrp->nat > 1)
1188 if (pgrp->pbcatom >= 0)
1190 pgrp->epgrppbc = epgrppbcREFAT;
1191 pull->bRefAt = TRUE;
1197 gmx_fatal(FARGS, "Pull groups can not have relative weights and cosine weighting at same time");
1199 pgrp->epgrppbc = epgrppbcCOS;
1200 if (pull->cosdim >= 0 && pull->cosdim != m)
1202 gmx_fatal(FARGS, "Can only use cosine weighting with pulling in one dimension (use mdp option pull_dim)");
1208 /* Set the indices */
1209 init_pull_group_index(fplog, cr, g, pgrp, pull->dim, mtop, ir, lambda);
1210 if (PULL_CYL(pull) && pgrp->invtm == 0)
1212 gmx_fatal(FARGS, "Can not have frozen atoms in a cylinder pull group");
1217 /* Absolute reference, set the inverse mass to zero */
1223 /* if we use dynamic reference groups, do some initialising for them */
1226 if (ir->ePull == epullCONSTRAINT && pull->ncoord > 1)
1228 /* We can't easily update the single reference group with multiple
1229 * constraints. This would require recalculating COMs.
1231 gmx_fatal(FARGS, "Constraint COM pulling supports only one coordinate with geometry=cylinder, you can use umbrella pulling with multiple coordinates");
1234 for (c = 0; c < pull->ncoord; c++)
1236 if (pull->group[pull->coord[c].group[0]].nat == 0)
1238 gmx_fatal(FARGS, "Dynamic reference groups are not supported when using absolute reference!\n");
1242 snew(pull->dyna, pull->ncoord);
1245 /* Only do I/O when we are doing dynamics and if we are the MASTER */
1250 if (pull->nstxout > 0)
1252 pull->out_x = open_pull_out(opt2fn("-px", nfile, fnm), pull, oenv, TRUE, Flags);
1254 if (pull->nstfout > 0)
1256 pull->out_f = open_pull_out(opt2fn("-pf", nfile, fnm), pull, oenv,
1262 void finish_pull(t_pull *pull)
1266 gmx_fio_fclose(pull->out_x);
1270 gmx_fio_fclose(pull->out_f);