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.
46 #include "gromacs/fileio/filenm.h"
47 #include "gromacs/fileio/gmxfio.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/legacyheaders/copyrite.h"
50 #include "gromacs/legacyheaders/gmx_ga2la.h"
51 #include "gromacs/legacyheaders/macros.h"
52 #include "gromacs/legacyheaders/mdrun.h"
53 #include "gromacs/legacyheaders/names.h"
54 #include "gromacs/legacyheaders/network.h"
55 #include "gromacs/legacyheaders/typedefs.h"
56 #include "gromacs/legacyheaders/types/commrec.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/topology/mtop_util.h"
60 #include "gromacs/utility/futil.h"
61 #include "gromacs/utility/smalloc.h"
63 static void pull_print_group_x(FILE *out, ivec dim, const t_pull_group *pgrp)
67 for (m = 0; m < DIM; m++)
71 fprintf(out, "\t%g", pgrp->x[m]);
76 static void pull_print_coord_dr(FILE *out, ivec dim, const t_pull_coord *pcrd)
80 for (m = 0; m < DIM; m++)
84 fprintf(out, "\t%g", pcrd->dr[m]);
89 static void pull_print_x(FILE *out, t_pull *pull, double t)
92 const t_pull_coord *pcrd;
94 fprintf(out, "%.4f", t);
96 for (c = 0; c < pull->ncoord; c++)
98 pcrd = &pull->coord[c];
104 pull_print_group_x(out, pull->dim, &pull->dyna[c]);
108 pull_print_group_x(out, pull->dim, &pull->group[pcrd->group[0]]);
111 pull_print_coord_dr(out, pull->dim, pcrd);
116 static void pull_print_f(FILE *out, t_pull *pull, double t)
120 fprintf(out, "%.4f", t);
122 for (c = 0; c < pull->ncoord; c++)
124 fprintf(out, "\t%g", pull->coord[c].f_scal);
129 void pull_print_output(t_pull *pull, gmx_int64_t step, double time)
131 if ((pull->nstxout != 0) && (step % pull->nstxout == 0))
133 pull_print_x(pull->out_x, pull, time);
136 if ((pull->nstfout != 0) && (step % pull->nstfout == 0))
138 pull_print_f(pull->out_f, pull, time);
142 static FILE *open_pull_out(const char *fn, t_pull *pull, const output_env_t oenv,
143 gmx_bool bCoord, unsigned long Flags)
147 char **setname, buf[10];
149 if (Flags & MD_APPENDFILES)
151 fp = gmx_fio_fopen(fn, "a+");
155 fp = gmx_fio_fopen(fn, "w+");
158 xvgr_header(fp, "Pull COM", "Time (ps)", "Position (nm)",
163 xvgr_header(fp, "Pull force", "Time (ps)", "Force (kJ/mol/nm)",
167 snew(setname, 2*pull->ncoord*DIM);
169 for (c = 0; c < pull->ncoord; c++)
175 for (m = 0; m < DIM; m++)
179 sprintf(buf, "%d %s%c", c+1, "c", 'X'+m);
180 setname[nsets] = gmx_strdup(buf);
185 for (m = 0; m < DIM; m++)
189 sprintf(buf, "%d %s%c", c+1, "d", 'X'+m);
190 setname[nsets] = gmx_strdup(buf);
197 sprintf(buf, "%d", c+1);
198 setname[nsets] = gmx_strdup(buf);
204 xvgr_legend(fp, nsets, (const char**)setname, oenv);
206 for (c = 0; c < nsets; c++)
216 /* Apply forces in a mass weighted fashion */
217 static void apply_forces_grp(const t_pull_group *pgrp, const t_mdatoms *md,
218 const dvec f_pull, int sign, rvec *f)
221 double wmass, inv_wm;
223 inv_wm = pgrp->wscale*pgrp->invtm;
225 for (i = 0; i < pgrp->nat_loc; i++)
227 ii = pgrp->ind_loc[i];
228 wmass = md->massT[ii];
229 if (pgrp->weight_loc)
231 wmass *= pgrp->weight_loc[i];
234 for (m = 0; m < DIM; m++)
236 f[ii][m] += sign * wmass * f_pull[m] * inv_wm;
241 /* Apply forces in a mass weighted fashion */
242 static void apply_forces(t_pull * pull, t_mdatoms * md, rvec *f)
245 const t_pull_coord *pcrd;
247 for (c = 0; c < pull->ncoord; c++)
249 pcrd = &pull->coord[c];
253 apply_forces_grp(&pull->dyna[c], md, pcrd->f, -1, f);
257 if (pull->group[pcrd->group[0]].nat > 0)
259 apply_forces_grp(&pull->group[pcrd->group[0]], md, pcrd->f, -1, f);
262 apply_forces_grp(&pull->group[pcrd->group[1]], md, pcrd->f, 1, f);
266 static double max_pull_distance2(const t_pull *pull, const t_pbc *pbc)
271 max_d2 = GMX_DOUBLE_MAX;
273 if (pull->eGeom != epullgDIRPBC)
275 for (m = 0; m < pbc->ndim_ePBC; m++)
277 if (pull->dim[m] != 0)
279 max_d2 = min(max_d2, norm2(pbc->box[m]));
287 static void low_get_pull_coord_dr(const t_pull *pull,
288 const t_pull_coord *pcrd,
289 const t_pbc *pbc, double t,
290 dvec xg, dvec xref, double max_dist2,
293 const t_pull_group *pgrp0, *pgrp1;
295 dvec xrefr, dref = {0, 0, 0};
298 pgrp0 = &pull->group[pcrd->group[0]];
299 pgrp1 = &pull->group[pcrd->group[1]];
301 /* Only the first group can be an absolute reference, in that case nat=0 */
304 for (m = 0; m < DIM; m++)
306 xref[m] = pcrd->origin[m];
310 copy_dvec(xref, xrefr);
312 if (pull->eGeom == epullgDIRPBC)
314 for (m = 0; m < DIM; m++)
316 dref[m] = (pcrd->init + pcrd->rate*t)*pcrd->vec[m];
318 /* Add the reference position, so we use the correct periodic image */
319 dvec_inc(xrefr, dref);
322 pbc_dx_d(pbc, xg, xrefr, dr);
324 for (m = 0; m < DIM; m++)
326 dr[m] *= pull->dim[m];
329 if (max_dist2 >= 0 && dr2 > 0.98*0.98*max_dist2)
331 gmx_fatal(FARGS, "Distance between pull groups %d and %d (%f nm) is larger than 0.49 times the box size (%f).\nYou might want to consider using \"pull-geometry = direction-periodic\" instead.\n",
332 pcrd->group[0], pcrd->group[1], sqrt(dr2), sqrt(max_dist2));
335 if (pull->eGeom == epullgDIRPBC)
341 static void get_pull_coord_dr(const t_pull *pull,
343 const t_pbc *pbc, double t,
347 const t_pull_coord *pcrd;
349 if (pull->eGeom == epullgDIRPBC)
355 md2 = max_pull_distance2(pull, pbc);
358 pcrd = &pull->coord[coord_ind];
360 low_get_pull_coord_dr(pull, pcrd, pbc, t,
361 pull->group[pcrd->group[1]].x,
362 PULL_CYL(pull) ? pull->dyna[coord_ind].x : pull->group[pcrd->group[0]].x,
367 void get_pull_coord_distance(const t_pull *pull,
369 const t_pbc *pbc, double t,
370 dvec dr, double *dev)
372 static gmx_bool bWarned = FALSE; /* TODO: this should be fixed for thread-safety,
373 but is fairly benign */
374 const t_pull_coord *pcrd;
376 double ref, drs, inpr;
378 pcrd = &pull->coord[coord_ind];
380 get_pull_coord_dr(pull, coord_ind, pbc, t, dr);
382 ref = pcrd->init + pcrd->rate*t;
387 /* Pull along the vector between the com's */
388 if (ref < 0 && !bWarned)
390 fprintf(stderr, "\nPull reference distance for coordinate %d is negative (%f)\n", coord_ind+1, ref);
396 /* With no vector we can not determine the direction for the force,
397 * so we set the force to zero.
403 /* Determine the deviation */
412 for (m = 0; m < DIM; m++)
414 inpr += pcrd->vec[m]*dr[m];
421 void clear_pull_forces(t_pull *pull)
425 /* Zeroing the forces is only required for constraint pulling.
426 * It can happen that multiple constraint steps need to be applied
427 * and therefore the constraint forces need to be accumulated.
429 for (i = 0; i < pull->ncoord; i++)
431 clear_dvec(pull->coord[i].f);
432 pull->coord[i].f_scal = 0;
436 /* Apply constraint using SHAKE */
437 static void do_constraint(t_pull *pull, t_pbc *pbc,
439 gmx_bool bMaster, tensor vir,
443 dvec *r_ij; /* x[i] com of i in prev. step. Obeys constr. -> r_ij[i] */
444 dvec unc_ij; /* xp[i] com of i this step, before constr. -> unc_ij */
445 dvec *rnew; /* current 'new' positions of the groups */
446 double *dr_tot; /* the total update of the coords */
450 double lambda, rm, mass, invdt = 0;
451 gmx_bool bConverged_all, bConverged = FALSE;
452 int niter = 0, g, c, ii, j, m, max_iter = 100;
454 dvec f; /* the pull force */
456 t_pull_group *pdyna, *pgrp0, *pgrp1;
459 snew(r_ij, pull->ncoord);
460 snew(dr_tot, pull->ncoord);
462 snew(rnew, pull->ngroup);
464 /* copy the current unconstrained positions for use in iterations. We
465 iterate until rinew[i] and rjnew[j] obey the constraints. Then
466 rinew - pull.x_unc[i] is the correction dr to group i */
467 for (g = 0; g < pull->ngroup; g++)
469 copy_dvec(pull->group[g].xp, rnew[g]);
473 /* There is only one pull coordinate and reference group */
474 copy_dvec(pull->dyna[0].xp, rnew[pull->coord[0].group[0]]);
477 /* Determine the constraint directions from the old positions */
478 for (c = 0; c < pull->ncoord; c++)
480 get_pull_coord_dr(pull, c, pbc, t, r_ij[c]);
481 /* Store the difference vector at time t for printing */
482 copy_dvec(r_ij[c], pull->coord[c].dr);
485 fprintf(debug, "Pull coord %d dr %f %f %f\n",
486 c, r_ij[c][XX], r_ij[c][YY], r_ij[c][ZZ]);
489 if (pull->eGeom == epullgDIR || pull->eGeom == epullgDIRPBC)
491 /* Select the component along vec */
493 for (m = 0; m < DIM; m++)
495 a += pull->coord[c].vec[m]*r_ij[c][m];
497 for (m = 0; m < DIM; m++)
499 r_ij[c][m] = a*pull->coord[c].vec[m];
503 if (dnorm2(r_ij[c]) == 0)
505 gmx_fatal(FARGS, "Distance for pull coordinate %d is zero with constraint pulling, which is not allowed.", c + 1);
509 bConverged_all = FALSE;
510 while (!bConverged_all && niter < max_iter)
512 bConverged_all = TRUE;
514 /* loop over all constraints */
515 for (c = 0; c < pull->ncoord; c++)
519 pcrd = &pull->coord[c];
520 pgrp0 = &pull->group[pcrd->group[0]];
521 pgrp1 = &pull->group[pcrd->group[1]];
523 /* Get the current difference vector */
524 low_get_pull_coord_dr(pull, pcrd, pbc, t,
525 rnew[pcrd->group[1]],
526 rnew[pcrd->group[0]],
529 ref = pcrd->init + pcrd->rate*t;
533 fprintf(debug, "Pull coord %d, iteration %d\n", c, niter);
536 rm = 1.0/(pgrp0->invtm + pgrp1->invtm);
543 gmx_fatal(FARGS, "The pull constraint reference distance for group %d is <= 0 (%f)", c, ref);
547 double q, c_a, c_b, c_c;
549 c_a = diprod(r_ij[c], r_ij[c]);
550 c_b = diprod(unc_ij, r_ij[c])*2;
551 c_c = diprod(unc_ij, unc_ij) - dsqr(ref);
555 q = -0.5*(c_b - sqrt(c_b*c_b - 4*c_a*c_c));
560 q = -0.5*(c_b + sqrt(c_b*c_b - 4*c_a*c_c));
567 "Pull ax^2+bx+c=0: a=%e b=%e c=%e lambda=%e\n",
568 c_a, c_b, c_c, lambda);
572 /* The position corrections dr due to the constraints */
573 dsvmul(-lambda*rm*pgrp1->invtm, r_ij[c], dr1);
574 dsvmul( lambda*rm*pgrp0->invtm, r_ij[c], dr0);
575 dr_tot[c] += -lambda*dnorm(r_ij[c]);
580 /* A 1-dimensional constraint along a vector */
582 for (m = 0; m < DIM; m++)
584 vec[m] = pcrd->vec[m];
585 a += unc_ij[m]*vec[m];
587 /* Select only the component along the vector */
588 dsvmul(a, vec, unc_ij);
592 fprintf(debug, "Pull inpr %e lambda: %e\n", a, lambda);
595 /* The position corrections dr due to the constraints */
596 dsvmul(-lambda*rm*pgrp1->invtm, vec, dr1);
597 dsvmul( lambda*rm*pgrp0->invtm, vec, dr0);
598 dr_tot[c] += -lambda;
609 low_get_pull_coord_dr(pull, pcrd, pbc, t, rnew[g1], rnew[g0], -1, tmp);
610 low_get_pull_coord_dr(pull, pcrd, pbc, t, dr1, dr0, -1, tmp3);
612 "Pull cur %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
613 rnew[g0][0], rnew[g0][1], rnew[g0][2],
614 rnew[g1][0], rnew[g1][1], rnew[g1][2], dnorm(tmp));
616 "Pull ref %8s %8s %8s %8s %8s %8s d: %8.5f\n",
617 "", "", "", "", "", "", ref);
619 "Pull cor %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
620 dr0[0], dr0[1], dr0[2],
621 dr1[0], dr1[1], dr1[2],
625 /* Update the COMs with dr */
626 dvec_inc(rnew[pcrd->group[1]], dr1);
627 dvec_inc(rnew[pcrd->group[0]], dr0);
630 /* Check if all constraints are fullfilled now */
631 for (c = 0; c < pull->ncoord; c++)
633 pcrd = &pull->coord[c];
634 ref = pcrd->init + pcrd->rate*t;
636 low_get_pull_coord_dr(pull, pcrd, pbc, t,
637 rnew[pcrd->group[1]],
638 rnew[pcrd->group[0]],
644 bConverged = fabs(dnorm(unc_ij) - ref) < pull->constr_tol;
649 for (m = 0; m < DIM; m++)
651 vec[m] = pcrd->vec[m];
653 inpr = diprod(unc_ij, vec);
654 dsvmul(inpr, vec, unc_ij);
656 fabs(diprod(unc_ij, vec) - ref) < pull->constr_tol;
664 fprintf(debug, "NOT CONVERGED YET: Group %d:"
665 "d_ref = %f, current d = %f\n",
666 g, ref, dnorm(unc_ij));
669 bConverged_all = FALSE;
674 /* if after all constraints are dealt with and bConverged is still TRUE
675 we're finished, if not we do another iteration */
677 if (niter > max_iter)
679 gmx_fatal(FARGS, "Too many iterations for constraint run: %d", niter);
682 /* DONE ITERATING, NOW UPDATE COORDINATES AND CALC. CONSTRAINT FORCES */
689 /* update atoms in the groups */
690 for (g = 0; g < pull->ngroup; g++)
692 const t_pull_group *pgrp;
695 if (PULL_CYL(pull) && g == pull->coord[0].group[0])
697 pgrp = &pull->dyna[0];
701 pgrp = &pull->group[g];
704 /* get the final constraint displacement dr for group g */
705 dvec_sub(rnew[g], pgrp->xp, dr);
706 /* select components of dr */
707 for (m = 0; m < DIM; m++)
709 dr[m] *= pull->dim[m];
712 /* update the atom positions */
714 for (j = 0; j < pgrp->nat_loc; j++)
716 ii = pgrp->ind_loc[j];
717 if (pgrp->weight_loc)
719 dsvmul(pgrp->wscale*pgrp->weight_loc[j], dr, tmp);
721 for (m = 0; m < DIM; m++)
727 for (m = 0; m < DIM; m++)
729 v[ii][m] += invdt*tmp[m];
735 /* calculate the constraint forces, used for output and virial only */
736 for (c = 0; c < pull->ncoord; c++)
738 pcrd = &pull->coord[c];
739 pcrd->f_scal = dr_tot[c]/((pull->group[pcrd->group[0]].invtm + pull->group[pcrd->group[1]].invtm)*dt*dt);
745 /* Add the pull contribution to the virial */
746 /* We have already checked above that r_ij[c] != 0 */
747 f_invr = pcrd->f_scal/dnorm(r_ij[c]);
749 for (j = 0; j < DIM; j++)
751 for (m = 0; m < DIM; m++)
753 vir[j][m] -= 0.5*f_invr*r_ij[c][j]*r_ij[c][m];
759 /* finished! I hope. Give back some memory */
765 /* Pulling with a harmonic umbrella potential or constant force */
766 static void do_pull_pot(int ePull,
767 t_pull *pull, t_pbc *pbc, double t, real lambda,
768 real *V, tensor vir, real *dVdl)
771 double dev, ndr, invdr;
775 /* loop over the pull coordinates */
778 for (c = 0; c < pull->ncoord; c++)
780 pcrd = &pull->coord[c];
782 get_pull_coord_distance(pull, c, pbc, t, pcrd->dr, &dev);
784 k = (1.0 - lambda)*pcrd->k + lambda*pcrd->kB;
785 dkdl = pcrd->kB - pcrd->k;
790 ndr = dnorm(pcrd->dr);
797 /* With an harmonic umbrella, the force is 0 at r=0,
798 * so we can set invdr to any value.
799 * With a constant force, the force at r=0 is not defined,
800 * so we zero it (this is anyhow a very rare event).
804 if (ePull == epullUMBRELLA)
806 pcrd->f_scal = -k*dev;
807 *V += 0.5* k*dsqr(dev);
808 *dVdl += 0.5*dkdl*dsqr(dev);
816 for (m = 0; m < DIM; m++)
818 pcrd->f[m] = pcrd->f_scal*pcrd->dr[m]*invdr;
824 if (ePull == epullUMBRELLA)
826 pcrd->f_scal = -k*dev;
827 *V += 0.5* k*dsqr(dev);
828 *dVdl += 0.5*dkdl*dsqr(dev);
833 for (m = 0; m < DIM; m++)
835 ndr += pcrd->vec[m]*pcrd->dr[m];
841 for (m = 0; m < DIM; m++)
843 pcrd->f[m] = pcrd->f_scal*pcrd->vec[m];
850 /* Add the pull contribution to the virial */
851 for (j = 0; j < DIM; j++)
853 for (m = 0; m < DIM; m++)
855 vir[j][m] -= 0.5*pcrd->f[j]*pcrd->dr[m];
862 real pull_potential(int ePull, t_pull *pull, t_mdatoms *md, t_pbc *pbc,
863 t_commrec *cr, double t, real lambda,
864 rvec *x, rvec *f, tensor vir, real *dvdlambda)
868 pull_calc_coms(cr, pull, md, pbc, t, x, NULL);
870 do_pull_pot(ePull, pull, pbc, t, lambda,
871 &V, pull->bVirial && MASTER(cr) ? vir : NULL, &dVdl);
873 /* Distribute forces over pulled groups */
874 apply_forces(pull, md, f);
881 return (MASTER(cr) ? V : 0.0);
884 void pull_constraint(t_pull *pull, t_mdatoms *md, t_pbc *pbc,
885 t_commrec *cr, double dt, double t,
886 rvec *x, rvec *xp, rvec *v, tensor vir)
888 pull_calc_coms(cr, pull, md, pbc, t, x, xp);
890 do_constraint(pull, pbc, xp, v, pull->bVirial && MASTER(cr), vir, dt, t);
893 static void make_local_pull_group(gmx_ga2la_t ga2la,
894 t_pull_group *pg, int start, int end)
899 for (i = 0; i < pg->nat; i++)
904 if (!ga2la_get_home(ga2la, ii, &ii))
909 if (ii >= start && ii < end)
911 /* This is a home atom, add it to the local pull group */
912 if (pg->nat_loc >= pg->nalloc_loc)
914 pg->nalloc_loc = over_alloc_dd(pg->nat_loc+1);
915 srenew(pg->ind_loc, pg->nalloc_loc);
916 if (pg->epgrppbc == epgrppbcCOS || pg->weight)
918 srenew(pg->weight_loc, pg->nalloc_loc);
921 pg->ind_loc[pg->nat_loc] = ii;
924 pg->weight_loc[pg->nat_loc] = pg->weight[i];
931 void dd_make_local_pull_groups(gmx_domdec_t *dd, t_pull *pull, t_mdatoms *md)
945 for (g = 0; g < pull->ngroup; g++)
947 make_local_pull_group(ga2la, &pull->group[g],
951 /* Since the PBC of atoms might have changed, we need to update the PBC */
952 pull->bSetPBCatoms = TRUE;
955 static void init_pull_group_index(FILE *fplog, t_commrec *cr,
956 int g, t_pull_group *pg, ivec pulldims,
957 gmx_mtop_t *mtop, t_inputrec *ir, real lambda)
959 int i, ii, d, nfrozen, ndim;
961 double tmass, wmass, wwmass;
962 gmx_groups_t *groups;
963 gmx_mtop_atomlookup_t alook;
966 if (EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD)
968 /* There are no masses in the integrator.
969 * But we still want to have the correct mass-weighted COMs.
970 * So we store the real masses in the weights.
971 * We do not set nweight, so these weights do not end up in the tpx file.
973 if (pg->nweight == 0)
975 snew(pg->weight, pg->nat);
984 pg->weight_loc = NULL;
988 pg->nat_loc = pg->nat;
989 pg->ind_loc = pg->ind;
990 if (pg->epgrppbc == epgrppbcCOS)
992 snew(pg->weight_loc, pg->nat);
996 pg->weight_loc = pg->weight;
1000 groups = &mtop->groups;
1002 alook = gmx_mtop_atomlookup_init(mtop);
1008 for (i = 0; i < pg->nat; i++)
1011 gmx_mtop_atomnr_to_atom(alook, ii, &atom);
1012 if (ir->opts.nFreeze)
1014 for (d = 0; d < DIM; d++)
1016 if (pulldims[d] && ir->opts.nFreeze[ggrpnr(groups, egcFREEZE, ii)][d])
1022 if (ir->efep == efepNO)
1028 m = (1 - lambda)*atom->m + lambda*atom->mB;
1030 if (pg->nweight > 0)
1038 if (EI_ENERGY_MINIMIZATION(ir->eI))
1040 /* Move the mass to the weight */
1045 else if (ir->eI == eiBD)
1049 mbd = ir->bd_fric*ir->delta_t;
1053 if (groups->grpnr[egcTC] == NULL)
1055 mbd = ir->delta_t/ir->opts.tau_t[0];
1059 mbd = ir->delta_t/ir->opts.tau_t[groups->grpnr[egcTC][ii]];
1071 gmx_mtop_atomlookup_destroy(alook);
1075 gmx_fatal(FARGS, "The total%s mass of pull group %d is zero",
1076 pg->weight ? " weighted" : "", g);
1081 "Pull group %d: %5d atoms, mass %9.3f", g, pg->nat, tmass);
1082 if (pg->weight || EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD)
1084 fprintf(fplog, ", weighted mass %9.3f", wmass*wmass/wwmass);
1086 if (pg->epgrppbc == epgrppbcCOS)
1088 fprintf(fplog, ", cosine weighting will be used");
1090 fprintf(fplog, "\n");
1095 /* A value > 0 signals not frozen, it is updated later */
1101 for (d = 0; d < DIM; d++)
1103 ndim += pulldims[d]*pg->nat;
1105 if (fplog && nfrozen > 0 && nfrozen < ndim)
1108 "\nWARNING: In pull group %d some, but not all of the degrees of freedom\n"
1109 " that are subject to pulling are frozen.\n"
1110 " For pulling the whole group will be frozen.\n\n",
1118 void init_pull(FILE *fplog, t_inputrec *ir, int nfile, const t_filenm fnm[],
1119 gmx_mtop_t *mtop, t_commrec *cr, const output_env_t oenv, real lambda,
1120 gmx_bool bOutFile, unsigned long Flags)
1124 int c, g, start = 0, end = 0, m;
1128 pull->ePBC = ir->ePBC;
1131 case epbcNONE: pull->npbcdim = 0; break;
1132 case epbcXY: pull->npbcdim = 2; break;
1133 default: pull->npbcdim = 3; break;
1138 gmx_bool bAbs, bCos;
1141 for (c = 0; c < pull->ncoord; c++)
1143 if (pull->group[pull->coord[c].group[0]].nat == 0 ||
1144 pull->group[pull->coord[c].group[1]].nat == 0)
1150 fprintf(fplog, "\nWill apply %s COM pulling in geometry '%s'\n",
1151 EPULLTYPE(ir->ePull), EPULLGEOM(pull->eGeom));
1152 fprintf(fplog, "with %d pull coordinate%s and %d group%s\n",
1153 pull->ncoord, pull->ncoord == 1 ? "" : "s",
1154 pull->ngroup, pull->ngroup == 1 ? "" : "s");
1157 fprintf(fplog, "with an absolute reference\n");
1160 for (g = 0; g < pull->ngroup; g++)
1162 if (pull->group[g].nat > 1 &&
1163 pull->group[g].pbcatom < 0)
1165 /* We are using cosine weighting */
1166 fprintf(fplog, "Cosine weighting is used for group %d\n", g);
1172 please_cite(fplog, "Engin2010");
1176 /* We always add the virial contribution,
1177 * except for geometry = direction_periodic where this is impossible.
1179 pull->bVirial = (pull->eGeom != epullgDIRPBC);
1180 if (getenv("GMX_NO_PULLVIR") != NULL)
1184 fprintf(fplog, "Found env. var., will not add the virial contribution of the COM pull forces\n");
1186 pull->bVirial = FALSE;
1191 pull->dbuf_cyl = NULL;
1192 pull->bRefAt = FALSE;
1194 for (g = 0; g < pull->ngroup; g++)
1196 pgrp = &pull->group[g];
1197 pgrp->epgrppbc = epgrppbcNONE;
1200 /* Determine if we need to take PBC into account for calculating
1201 * the COM's of the pull groups.
1203 for (m = 0; m < pull->npbcdim; m++)
1205 if (pull->dim[m] && pgrp->nat > 1)
1207 if (pgrp->pbcatom >= 0)
1209 pgrp->epgrppbc = epgrppbcREFAT;
1210 pull->bRefAt = TRUE;
1216 gmx_fatal(FARGS, "Pull groups can not have relative weights and cosine weighting at same time");
1218 pgrp->epgrppbc = epgrppbcCOS;
1219 if (pull->cosdim >= 0 && pull->cosdim != m)
1221 gmx_fatal(FARGS, "Can only use cosine weighting with pulling in one dimension (use mdp option pull_dim)");
1227 /* Set the indices */
1228 init_pull_group_index(fplog, cr, g, pgrp, pull->dim, mtop, ir, lambda);
1229 if (PULL_CYL(pull) && pgrp->invtm == 0)
1231 gmx_fatal(FARGS, "Can not have frozen atoms in a cylinder pull group");
1236 /* Absolute reference, set the inverse mass to zero */
1242 /* if we use dynamic reference groups, do some initialising for them */
1245 if (ir->ePull == epullCONSTRAINT && pull->ncoord > 1)
1247 /* We can't easily update the single reference group with multiple
1248 * constraints. This would require recalculating COMs.
1250 gmx_fatal(FARGS, "Constraint COM pulling supports only one coordinate with geometry=cylinder, you can use umbrella pulling with multiple coordinates");
1253 for (c = 0; c < pull->ncoord; c++)
1255 if (pull->group[pull->coord[c].group[0]].nat == 0)
1257 gmx_fatal(FARGS, "Dynamic reference groups are not supported when using absolute reference!\n");
1261 snew(pull->dyna, pull->ncoord);
1264 /* We still need to initialize the PBC reference coordinates */
1265 pull->bSetPBCatoms = TRUE;
1267 /* Only do I/O when we are doing dynamics and if we are the MASTER */
1272 if (pull->nstxout > 0)
1274 pull->out_x = open_pull_out(opt2fn("-px", nfile, fnm), pull, oenv, TRUE, Flags);
1276 if (pull->nstfout > 0)
1278 pull->out_f = open_pull_out(opt2fn("-pf", nfile, fnm), pull, oenv,
1284 void finish_pull(t_pull *pull)
1288 gmx_fio_fclose(pull->out_x);
1292 gmx_fio_fclose(pull->out_f);