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,2015, 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, const t_pull_coord *pcrd,
77 gmx_bool bPrintRefValue, double t,
78 gmx_bool bPrintComponents)
80 double distance, distance2;
83 if (pcrd->eGeom == epullgDIST)
85 /* Geometry=distance: print the length of the distance vector */
87 for (m = 0; m < DIM; m++)
91 distance2 += pcrd->dr[m]*pcrd->dr[m];
94 distance = sqrt(distance2);
98 /* Geometry=direction-like: print distance along the pull vector */
100 for (m = 0; m < DIM; m++)
104 distance += pcrd->dr[m]*pcrd->vec[m];
108 fprintf(out, "\t%g", distance);
112 fprintf(out, "\t%g", pcrd->init + pcrd->rate*t);
115 if (bPrintComponents)
117 for (m = 0; m < DIM; m++)
121 fprintf(out, "\t%g", pcrd->dr[m]);
127 static void pull_print_x(FILE *out, t_pull *pull, double t)
132 fprintf(out, "%.4f", t);
134 for (c = 0; c < pull->ncoord; c++)
136 pcrd = &pull->coord[c];
138 if (pull->bPrintCOM1)
140 if (pcrd->eGeom == epullgCYL)
142 pull_print_group_x(out, pcrd->dim, &pull->dyna[c]);
146 pull_print_group_x(out, pcrd->dim, &pull->group[pcrd->group[0]]);
149 if (pull->bPrintCOM2)
151 pull_print_group_x(out, pcrd->dim, &pull->group[pcrd->group[1]]);
153 pull_print_coord_dr(out, pcrd,
154 pull->bPrintRefValue, t, pull->bPrintComp);
159 static void pull_print_f(FILE *out, t_pull *pull, double t)
163 fprintf(out, "%.4f", t);
165 for (c = 0; c < pull->ncoord; c++)
167 fprintf(out, "\t%g", pull->coord[c].f_scal);
172 void pull_print_output(t_pull *pull, gmx_int64_t step, double time)
174 if ((pull->nstxout != 0) && (step % pull->nstxout == 0))
176 pull_print_x(pull->out_x, pull, time);
179 if ((pull->nstfout != 0) && (step % pull->nstfout == 0))
181 pull_print_f(pull->out_f, pull, time);
185 static FILE *open_pull_out(const char *fn, t_pull *pull, const output_env_t oenv,
186 gmx_bool bCoord, unsigned long Flags)
190 char **setname, buf[10];
192 if (Flags & MD_APPENDFILES)
194 fp = gmx_fio_fopen(fn, "a+");
198 fp = gmx_fio_fopen(fn, "w+");
201 xvgr_header(fp, "Pull COM", "Time (ps)", "Position (nm)",
206 xvgr_header(fp, "Pull force", "Time (ps)", "Force (kJ/mol/nm)",
210 /* With default mdp options only the actual distance is printed,
211 * but optionally 2 COMs, the reference distance and distance
212 * components can also be printed.
214 snew(setname, pull->ncoord*(DIM + DIM + 1 + 1 + DIM));
216 for (c = 0; c < pull->ncoord; c++)
220 /* The order of this legend should match the order of printing
221 * the data in print_pull_x above.
224 if (pull->bPrintCOM1)
226 /* Legend for reference group position */
227 for (m = 0; m < DIM; m++)
229 if (pull->coord[c].dim[m])
231 sprintf(buf, "%d %s%d%c", c+1, "c", 1, 'X'+m);
232 setname[nsets] = gmx_strdup(buf);
237 if (pull->bPrintCOM2)
239 /* Legend for reference group position */
240 for (m = 0; m < DIM; m++)
242 if (pull->coord[c].dim[m])
244 sprintf(buf, "%d %s%d%c", c+1, "c", 2, 'X'+m);
245 setname[nsets] = gmx_strdup(buf);
250 /* The pull coord distance */
251 sprintf(buf, "%d", c+1);
252 setname[nsets] = gmx_strdup(buf);
254 if (pull->bPrintRefValue)
256 sprintf(buf, "%c%d", 'r', c+1);
257 setname[nsets] = gmx_strdup(buf);
260 if (pull->bPrintComp)
262 /* The pull coord distance components */
263 for (m = 0; m < DIM; m++)
265 if (pull->coord[c].dim[m])
267 sprintf(buf, "%d %s%c", c+1, "d", 'X'+m);
268 setname[nsets] = gmx_strdup(buf);
276 /* For the pull force we always only use one scalar */
277 sprintf(buf, "%d", c+1);
278 setname[nsets] = gmx_strdup(buf);
284 xvgr_legend(fp, nsets, (const char**)setname, oenv);
286 for (c = 0; c < nsets; c++)
296 /* Apply forces in a mass weighted fashion */
297 static void apply_forces_grp(const t_pull_group *pgrp, const t_mdatoms *md,
298 const dvec f_pull, int sign, rvec *f)
301 double wmass, inv_wm;
303 inv_wm = pgrp->mwscale;
305 if (pgrp->nat == 1 && pgrp->nat_loc == 1)
307 /* Only one atom and our rank has this atom: we can skip
308 * the mass weighting, which means that this code also works
309 * for mass=0, e.g. with a virtual site.
311 for (m = 0; m < DIM; m++)
313 f[pgrp->ind_loc[0]][m] += sign*f_pull[m];
318 for (i = 0; i < pgrp->nat_loc; i++)
320 ii = pgrp->ind_loc[i];
321 wmass = md->massT[ii];
322 if (pgrp->weight_loc)
324 wmass *= pgrp->weight_loc[i];
327 for (m = 0; m < DIM; m++)
329 f[ii][m] += sign * wmass * f_pull[m] * inv_wm;
335 /* Apply forces in a mass weighted fashion to a cylinder group */
336 static void apply_forces_cyl_grp(const t_pull_group *pgrp, double dv_corr,
338 const dvec f_pull, double f_scal,
342 double mass, weight, inv_wm, dv_com;
344 inv_wm = pgrp->mwscale;
346 for (i = 0; i < pgrp->nat_loc; i++)
348 ii = pgrp->ind_loc[i];
349 mass = md->massT[ii];
350 weight = pgrp->weight_loc[i];
351 /* The stored axial distance from the cylinder center (dv) needs
352 * to be corrected for an offset (dv_corr), which was unknown when
355 dv_com = pgrp->dv[i] + dv_corr;
357 /* Here we not only add the pull force working along vec (f_pull),
358 * but also a radial component, due to the dependence of the weights
359 * on the radial distance.
361 for (m = 0; m < DIM; m++)
363 f[ii][m] += sign*inv_wm*(mass*weight*f_pull[m] +
364 pgrp->mdw[i][m]*dv_com*f_scal);
369 /* Apply forces in a mass weighted fashion */
370 static void apply_forces(t_pull * pull, t_mdatoms * md, rvec *f)
373 const t_pull_coord *pcrd;
375 for (c = 0; c < pull->ncoord; c++)
377 pcrd = &pull->coord[c];
379 if (pcrd->eGeom == epullgCYL)
384 apply_forces_cyl_grp(&pull->dyna[c], pcrd->cyl_dev, md,
385 pcrd->f, pcrd->f_scal, -1, f);
387 /* Sum the force along the vector and the radial force */
388 for (m = 0; m < DIM; m++)
390 f_tot[m] = pcrd->f[m] + pcrd->f_scal*pcrd->ffrad[m];
392 apply_forces_grp(&pull->group[pcrd->group[1]], md, f_tot, 1, f);
396 if (pull->group[pcrd->group[0]].nat > 0)
398 apply_forces_grp(&pull->group[pcrd->group[0]], md, pcrd->f, -1, f);
400 apply_forces_grp(&pull->group[pcrd->group[1]], md, pcrd->f, 1, f);
405 static double max_pull_distance2(const t_pull_coord *pcrd, const t_pbc *pbc)
410 max_d2 = GMX_DOUBLE_MAX;
412 for (m = 0; m < pbc->ndim_ePBC; m++)
414 if (pcrd->dim[m] != 0)
416 max_d2 = min(max_d2, norm2(pbc->box[m]));
423 static void low_get_pull_coord_dr(const t_pull *pull,
424 const t_pull_coord *pcrd,
425 const t_pbc *pbc, double t,
426 dvec xg, dvec xref, double max_dist2,
429 const t_pull_group *pgrp0, *pgrp1;
431 dvec xrefr, dref = {0, 0, 0};
434 pgrp0 = &pull->group[pcrd->group[0]];
435 pgrp1 = &pull->group[pcrd->group[1]];
437 /* Only the first group can be an absolute reference, in that case nat=0 */
440 for (m = 0; m < DIM; m++)
442 xref[m] = pcrd->origin[m];
446 copy_dvec(xref, xrefr);
448 if (pcrd->eGeom == epullgDIRPBC)
450 for (m = 0; m < DIM; m++)
452 dref[m] = (pcrd->init + pcrd->rate*t)*pcrd->vec[m];
454 /* Add the reference position, so we use the correct periodic image */
455 dvec_inc(xrefr, dref);
458 pbc_dx_d(pbc, xg, xrefr, dr);
460 for (m = 0; m < DIM; m++)
462 dr[m] *= pcrd->dim[m];
465 if (max_dist2 >= 0 && dr2 > 0.98*0.98*max_dist2)
467 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",
468 pcrd->group[0], pcrd->group[1], sqrt(dr2), sqrt(max_dist2));
471 if (pcrd->eGeom == epullgDIRPBC)
477 static void get_pull_coord_dr(const t_pull *pull,
479 const t_pbc *pbc, double t,
483 const t_pull_coord *pcrd;
485 pcrd = &pull->coord[coord_ind];
487 if (pcrd->eGeom == epullgDIRPBC)
493 md2 = max_pull_distance2(pcrd, pbc);
496 low_get_pull_coord_dr(pull, pcrd, pbc, t,
497 pull->group[pcrd->group[1]].x,
498 pcrd->eGeom == epullgCYL ? pull->dyna[coord_ind].x : pull->group[pcrd->group[0]].x,
503 void get_pull_coord_distance(const t_pull *pull,
505 const t_pbc *pbc, double t,
506 dvec dr, double *dev)
508 static gmx_bool bWarned = FALSE; /* TODO: this should be fixed for thread-safety,
509 but is fairly benign */
510 const t_pull_coord *pcrd;
512 double ref, drs, inpr;
514 pcrd = &pull->coord[coord_ind];
516 get_pull_coord_dr(pull, coord_ind, pbc, t, dr);
518 ref = pcrd->init + pcrd->rate*t;
523 /* Pull along the vector between the com's */
524 if (ref < 0 && !bWarned)
526 fprintf(stderr, "\nPull reference distance for coordinate %d is negative (%f)\n", coord_ind+1, ref);
532 /* With no vector we can not determine the direction for the force,
533 * so we set the force to zero.
539 /* Determine the deviation */
548 for (m = 0; m < DIM; m++)
550 inpr += pcrd->vec[m]*dr[m];
557 void clear_pull_forces(t_pull *pull)
561 /* Zeroing the forces is only required for constraint pulling.
562 * It can happen that multiple constraint steps need to be applied
563 * and therefore the constraint forces need to be accumulated.
565 for (i = 0; i < pull->ncoord; i++)
567 clear_dvec(pull->coord[i].f);
568 pull->coord[i].f_scal = 0;
572 /* Apply constraint using SHAKE */
573 static void do_constraint(t_pull *pull, t_pbc *pbc,
575 gmx_bool bMaster, tensor vir,
579 dvec *r_ij; /* x[i] com of i in prev. step. Obeys constr. -> r_ij[i] */
580 dvec unc_ij; /* xp[i] com of i this step, before constr. -> unc_ij */
581 dvec *rnew; /* current 'new' positions of the groups */
582 double *dr_tot; /* the total update of the coords */
586 double lambda, rm, invdt = 0;
587 gmx_bool bConverged_all, bConverged = FALSE;
588 int niter = 0, g, c, ii, j, m, max_iter = 100;
590 dvec f; /* the pull force */
592 t_pull_group *pgrp0, *pgrp1;
595 snew(r_ij, pull->ncoord);
596 snew(dr_tot, pull->ncoord);
598 snew(rnew, pull->ngroup);
600 /* copy the current unconstrained positions for use in iterations. We
601 iterate until rinew[i] and rjnew[j] obey the constraints. Then
602 rinew - pull.x_unc[i] is the correction dr to group i */
603 for (g = 0; g < pull->ngroup; g++)
605 copy_dvec(pull->group[g].xp, rnew[g]);
608 /* Determine the constraint directions from the old positions */
609 for (c = 0; c < pull->ncoord; c++)
611 pcrd = &pull->coord[c];
613 if (pcrd->eType != epullCONSTRAINT)
618 get_pull_coord_dr(pull, c, pbc, t, r_ij[c]);
619 /* Store the difference vector at time t for printing */
620 copy_dvec(r_ij[c], pcrd->dr);
623 fprintf(debug, "Pull coord %d dr %f %f %f\n",
624 c, r_ij[c][XX], r_ij[c][YY], r_ij[c][ZZ]);
627 if (pcrd->eGeom == epullgDIR || pcrd->eGeom == epullgDIRPBC)
629 /* Select the component along vec */
631 for (m = 0; m < DIM; m++)
633 a += pcrd->vec[m]*r_ij[c][m];
635 for (m = 0; m < DIM; m++)
637 r_ij[c][m] = a*pcrd->vec[m];
641 if (dnorm2(r_ij[c]) == 0)
643 gmx_fatal(FARGS, "Distance for pull coordinate %d is zero with constraint pulling, which is not allowed.", c + 1);
647 bConverged_all = FALSE;
648 while (!bConverged_all && niter < max_iter)
650 bConverged_all = TRUE;
652 /* loop over all constraints */
653 for (c = 0; c < pull->ncoord; c++)
657 pcrd = &pull->coord[c];
659 if (pcrd->eType != epullCONSTRAINT)
664 pgrp0 = &pull->group[pcrd->group[0]];
665 pgrp1 = &pull->group[pcrd->group[1]];
667 /* Get the current difference vector */
668 low_get_pull_coord_dr(pull, pcrd, pbc, t,
669 rnew[pcrd->group[1]],
670 rnew[pcrd->group[0]],
673 ref = pcrd->init + pcrd->rate*t;
677 fprintf(debug, "Pull coord %d, iteration %d\n", c, niter);
680 rm = 1.0/(pgrp0->invtm + pgrp1->invtm);
687 gmx_fatal(FARGS, "The pull constraint reference distance for group %d is <= 0 (%f)", c, ref);
691 double q, c_a, c_b, c_c;
693 c_a = diprod(r_ij[c], r_ij[c]);
694 c_b = diprod(unc_ij, r_ij[c])*2;
695 c_c = diprod(unc_ij, unc_ij) - dsqr(ref);
699 q = -0.5*(c_b - sqrt(c_b*c_b - 4*c_a*c_c));
704 q = -0.5*(c_b + sqrt(c_b*c_b - 4*c_a*c_c));
711 "Pull ax^2+bx+c=0: a=%e b=%e c=%e lambda=%e\n",
712 c_a, c_b, c_c, lambda);
716 /* The position corrections dr due to the constraints */
717 dsvmul(-lambda*rm*pgrp1->invtm, r_ij[c], dr1);
718 dsvmul( lambda*rm*pgrp0->invtm, r_ij[c], dr0);
719 dr_tot[c] += -lambda*dnorm(r_ij[c]);
724 /* A 1-dimensional constraint along a vector */
726 for (m = 0; m < DIM; m++)
728 vec[m] = pcrd->vec[m];
729 a += unc_ij[m]*vec[m];
731 /* Select only the component along the vector */
732 dsvmul(a, vec, unc_ij);
736 fprintf(debug, "Pull inpr %e lambda: %e\n", a, lambda);
739 /* The position corrections dr due to the constraints */
740 dsvmul(-lambda*rm*pgrp1->invtm, vec, dr1);
741 dsvmul( lambda*rm*pgrp0->invtm, vec, dr0);
742 dr_tot[c] += -lambda;
753 low_get_pull_coord_dr(pull, pcrd, pbc, t, rnew[g1], rnew[g0], -1, tmp);
754 low_get_pull_coord_dr(pull, pcrd, pbc, t, dr1, dr0, -1, tmp3);
756 "Pull cur %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
757 rnew[g0][0], rnew[g0][1], rnew[g0][2],
758 rnew[g1][0], rnew[g1][1], rnew[g1][2], dnorm(tmp));
760 "Pull ref %8s %8s %8s %8s %8s %8s d: %8.5f\n",
761 "", "", "", "", "", "", ref);
763 "Pull cor %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
764 dr0[0], dr0[1], dr0[2],
765 dr1[0], dr1[1], dr1[2],
769 /* Update the COMs with dr */
770 dvec_inc(rnew[pcrd->group[1]], dr1);
771 dvec_inc(rnew[pcrd->group[0]], dr0);
774 /* Check if all constraints are fullfilled now */
775 for (c = 0; c < pull->ncoord; c++)
777 pcrd = &pull->coord[c];
779 if (pcrd->eType != epullCONSTRAINT)
784 ref = pcrd->init + pcrd->rate*t;
786 low_get_pull_coord_dr(pull, pcrd, pbc, t,
787 rnew[pcrd->group[1]],
788 rnew[pcrd->group[0]],
794 bConverged = fabs(dnorm(unc_ij) - ref) < pull->constr_tol;
799 for (m = 0; m < DIM; m++)
801 vec[m] = pcrd->vec[m];
803 inpr = diprod(unc_ij, vec);
804 dsvmul(inpr, vec, unc_ij);
806 fabs(diprod(unc_ij, vec) - ref) < pull->constr_tol;
814 fprintf(debug, "NOT CONVERGED YET: Group %d:"
815 "d_ref = %f, current d = %f\n",
816 g, ref, dnorm(unc_ij));
819 bConverged_all = FALSE;
824 /* if after all constraints are dealt with and bConverged is still TRUE
825 we're finished, if not we do another iteration */
827 if (niter > max_iter)
829 gmx_fatal(FARGS, "Too many iterations for constraint run: %d", niter);
832 /* DONE ITERATING, NOW UPDATE COORDINATES AND CALC. CONSTRAINT FORCES */
839 /* update atoms in the groups */
840 for (g = 0; g < pull->ngroup; g++)
842 const t_pull_group *pgrp;
845 pgrp = &pull->group[g];
847 /* get the final constraint displacement dr for group g */
848 dvec_sub(rnew[g], pgrp->xp, dr);
852 /* No displacement, no update necessary */
856 /* update the atom positions */
858 for (j = 0; j < pgrp->nat_loc; j++)
860 ii = pgrp->ind_loc[j];
861 if (pgrp->weight_loc)
863 dsvmul(pgrp->wscale*pgrp->weight_loc[j], dr, tmp);
865 for (m = 0; m < DIM; m++)
871 for (m = 0; m < DIM; m++)
873 v[ii][m] += invdt*tmp[m];
879 /* calculate the constraint forces, used for output and virial only */
880 for (c = 0; c < pull->ncoord; c++)
882 pcrd = &pull->coord[c];
884 if (pcrd->eType != epullCONSTRAINT)
889 pcrd->f_scal = dr_tot[c]/((pull->group[pcrd->group[0]].invtm + pull->group[pcrd->group[1]].invtm)*dt*dt);
891 if (vir != NULL && pcrd->eGeom != epullgDIRPBC && bMaster)
895 /* Add the pull contribution to the virial */
896 /* We have already checked above that r_ij[c] != 0 */
897 f_invr = pcrd->f_scal/dnorm(r_ij[c]);
899 for (j = 0; j < DIM; j++)
901 for (m = 0; m < DIM; m++)
903 vir[j][m] -= 0.5*f_invr*r_ij[c][j]*r_ij[c][m];
909 /* finished! I hope. Give back some memory */
915 /* Pulling with a harmonic umbrella potential or constant force */
916 static void do_pull_pot(t_pull *pull, t_pbc *pbc, double t, real lambda,
917 real *V, tensor vir, real *dVdl)
920 double dev, ndr, invdr = 0;
924 /* loop over the pull coordinates */
927 for (c = 0; c < pull->ncoord; c++)
929 pcrd = &pull->coord[c];
931 if (pcrd->eType == epullCONSTRAINT)
936 get_pull_coord_distance(pull, c, pbc, t, pcrd->dr, &dev);
938 k = (1.0 - lambda)*pcrd->k + lambda*pcrd->kB;
939 dkdl = pcrd->kB - pcrd->k;
941 if (pcrd->eGeom == epullgDIST)
943 ndr = dnorm(pcrd->dr);
950 /* With an harmonic umbrella, the force is 0 at r=0,
951 * so we can set invdr to any value.
952 * With a constant force, the force at r=0 is not defined,
953 * so we zero it (this is anyhow a very rare event).
961 for (m = 0; m < DIM; m++)
963 ndr += pcrd->vec[m]*pcrd->dr[m];
970 case epullFLATBOTTOM:
971 /* The only difference between an umbrella and a flat-bottom
972 * potential is that a flat-bottom is zero below zero.
974 if (pcrd->eType == epullFLATBOTTOM && dev < 0)
979 pcrd->f_scal = -k*dev;
980 *V += 0.5* k*dsqr(dev);
981 *dVdl += 0.5*dkdl*dsqr(dev);
989 gmx_incons("Unsupported pull type in do_pull_pot");
992 if (pcrd->eGeom == epullgDIST)
994 for (m = 0; m < DIM; m++)
996 pcrd->f[m] = pcrd->f_scal*pcrd->dr[m]*invdr;
1001 for (m = 0; m < DIM; m++)
1003 pcrd->f[m] = pcrd->f_scal*pcrd->vec[m];
1007 if (vir != NULL && pcrd->eGeom != epullgDIRPBC)
1009 /* Add the pull contribution to the virial */
1010 for (j = 0; j < DIM; j++)
1012 for (m = 0; m < DIM; m++)
1014 vir[j][m] -= 0.5*pcrd->f[j]*pcrd->dr[m];
1021 real pull_potential(t_pull *pull, t_mdatoms *md, t_pbc *pbc,
1022 t_commrec *cr, double t, real lambda,
1023 rvec *x, rvec *f, tensor vir, real *dvdlambda)
1027 pull_calc_coms(cr, pull, md, pbc, t, x, NULL);
1029 do_pull_pot(pull, pbc, t, lambda,
1030 &V, MASTER(cr) ? vir : NULL, &dVdl);
1032 /* Distribute forces over pulled groups */
1033 apply_forces(pull, md, f);
1040 return (MASTER(cr) ? V : 0.0);
1043 void pull_constraint(t_pull *pull, t_mdatoms *md, t_pbc *pbc,
1044 t_commrec *cr, double dt, double t,
1045 rvec *x, rvec *xp, rvec *v, tensor vir)
1047 pull_calc_coms(cr, pull, md, pbc, t, x, xp);
1049 do_constraint(pull, pbc, xp, v, MASTER(cr), vir, dt, t);
1052 static void make_local_pull_group(gmx_ga2la_t ga2la,
1053 t_pull_group *pg, int start, int end)
1058 for (i = 0; i < pg->nat; i++)
1063 if (!ga2la_get_home(ga2la, ii, &ii))
1068 if (ii >= start && ii < end)
1070 /* This is a home atom, add it to the local pull group */
1071 if (pg->nat_loc >= pg->nalloc_loc)
1073 pg->nalloc_loc = over_alloc_dd(pg->nat_loc+1);
1074 srenew(pg->ind_loc, pg->nalloc_loc);
1075 if (pg->epgrppbc == epgrppbcCOS || pg->weight)
1077 srenew(pg->weight_loc, pg->nalloc_loc);
1080 pg->ind_loc[pg->nat_loc] = ii;
1083 pg->weight_loc[pg->nat_loc] = pg->weight[i];
1090 void dd_make_local_pull_groups(gmx_domdec_t *dd, t_pull *pull, t_mdatoms *md)
1104 for (g = 0; g < pull->ngroup; g++)
1106 make_local_pull_group(ga2la, &pull->group[g],
1110 /* Since the PBC of atoms might have changed, we need to update the PBC */
1111 pull->bSetPBCatoms = TRUE;
1114 static void init_pull_group_index(FILE *fplog, t_commrec *cr,
1115 int g, t_pull_group *pg,
1116 gmx_bool bConstraint, ivec pulldim_con,
1117 gmx_mtop_t *mtop, t_inputrec *ir, real lambda)
1119 int i, ii, d, nfrozen, ndim;
1121 double tmass, wmass, wwmass;
1122 gmx_groups_t *groups;
1123 gmx_mtop_atomlookup_t alook;
1126 if (EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD)
1128 /* There are no masses in the integrator.
1129 * But we still want to have the correct mass-weighted COMs.
1130 * So we store the real masses in the weights.
1131 * We do not set nweight, so these weights do not end up in the tpx file.
1133 if (pg->nweight == 0)
1135 snew(pg->weight, pg->nat);
1144 pg->weight_loc = NULL;
1148 pg->nat_loc = pg->nat;
1149 pg->ind_loc = pg->ind;
1150 if (pg->epgrppbc == epgrppbcCOS)
1152 snew(pg->weight_loc, pg->nat);
1156 pg->weight_loc = pg->weight;
1160 groups = &mtop->groups;
1162 alook = gmx_mtop_atomlookup_init(mtop);
1168 for (i = 0; i < pg->nat; i++)
1171 gmx_mtop_atomnr_to_atom(alook, ii, &atom);
1172 if (bConstraint && ir->opts.nFreeze)
1174 for (d = 0; d < DIM; d++)
1176 if (pulldim_con[d] == 1 &&
1177 ir->opts.nFreeze[ggrpnr(groups, egcFREEZE, ii)][d])
1183 if (ir->efep == efepNO)
1189 m = (1 - lambda)*atom->m + lambda*atom->mB;
1191 if (pg->nweight > 0)
1199 if (EI_ENERGY_MINIMIZATION(ir->eI))
1201 /* Move the mass to the weight */
1206 else if (ir->eI == eiBD)
1210 mbd = ir->bd_fric*ir->delta_t;
1214 if (groups->grpnr[egcTC] == NULL)
1216 mbd = ir->delta_t/ir->opts.tau_t[0];
1220 mbd = ir->delta_t/ir->opts.tau_t[groups->grpnr[egcTC][ii]];
1232 gmx_mtop_atomlookup_destroy(alook);
1236 /* We can have single atom groups with zero mass with potential pulling
1237 * without cosine weighting.
1239 if (pg->nat == 1 && !bConstraint && pg->epgrppbc != epgrppbcCOS)
1241 /* With one atom the mass doesn't matter */
1246 gmx_fatal(FARGS, "The total%s mass of pull group %d is zero",
1247 pg->weight ? " weighted" : "", g);
1253 "Pull group %d: %5d atoms, mass %9.3f", g, pg->nat, tmass);
1254 if (pg->weight || EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD)
1256 fprintf(fplog, ", weighted mass %9.3f", wmass*wmass/wwmass);
1258 if (pg->epgrppbc == epgrppbcCOS)
1260 fprintf(fplog, ", cosine weighting will be used");
1262 fprintf(fplog, "\n");
1267 /* A value != 0 signals not frozen, it is updated later */
1273 for (d = 0; d < DIM; d++)
1275 ndim += pulldim_con[d]*pg->nat;
1277 if (fplog && nfrozen > 0 && nfrozen < ndim)
1280 "\nWARNING: In pull group %d some, but not all of the degrees of freedom\n"
1281 " that are subject to pulling are frozen.\n"
1282 " For constraint pulling the whole group will be frozen.\n\n",
1290 void init_pull(FILE *fplog, t_inputrec *ir, int nfile, const t_filenm fnm[],
1291 gmx_mtop_t *mtop, t_commrec *cr, const output_env_t oenv, real lambda,
1292 gmx_bool bOutFile, unsigned long Flags)
1300 pull->bPotential = FALSE;
1301 pull->bConstraint = FALSE;
1302 pull->bCylinder = FALSE;
1303 for (c = 0; c < pull->ncoord; c++)
1307 pcrd = &pull->coord[c];
1309 if (pcrd->eType == epullCONSTRAINT)
1311 pull->bConstraint = TRUE;
1315 pull->bPotential = TRUE;
1318 if (pcrd->eGeom == epullgCYL)
1320 pull->bCylinder = TRUE;
1322 if (pcrd->eType == epullCONSTRAINT)
1324 gmx_fatal(FARGS, "Pulling of type %s can not be combined with geometry %s. Consider using pull type %s.",
1325 epull_names[pcrd->eType],
1326 epullg_names[pcrd->eGeom],
1327 epull_names[epullUMBRELLA]);
1332 /* We only need to calculate the plain COM of a group
1333 * when it is not only used as a cylinder group.
1335 if (pull->group[pcrd->group[0]].nat > 0)
1337 pull->group[pcrd->group[0]].bCalcCOM = TRUE;
1340 if (pull->group[pcrd->group[1]].nat > 0)
1342 pull->group[pcrd->group[1]].bCalcCOM = TRUE;
1346 pull->ePBC = ir->ePBC;
1349 case epbcNONE: pull->npbcdim = 0; break;
1350 case epbcXY: pull->npbcdim = 2; break;
1351 default: pull->npbcdim = 3; break;
1356 gmx_bool bAbs, bCos;
1359 for (c = 0; c < pull->ncoord; c++)
1361 if (pull->group[pull->coord[c].group[0]].nat == 0 ||
1362 pull->group[pull->coord[c].group[1]].nat == 0)
1368 fprintf(fplog, "\n");
1369 if (pull->bPotential)
1371 fprintf(fplog, "Will apply potential COM pulling\n");
1373 if (pull->bConstraint)
1375 fprintf(fplog, "Will apply constraint COM pulling\n");
1377 fprintf(fplog, "with %d pull coordinate%s and %d group%s\n",
1378 pull->ncoord, pull->ncoord == 1 ? "" : "s",
1379 pull->ngroup, pull->ngroup == 1 ? "" : "s");
1382 fprintf(fplog, "with an absolute reference\n");
1385 for (g = 0; g < pull->ngroup; g++)
1387 if (pull->group[g].nat > 1 &&
1388 pull->group[g].pbcatom < 0)
1390 /* We are using cosine weighting */
1391 fprintf(fplog, "Cosine weighting is used for group %d\n", g);
1397 please_cite(fplog, "Engin2010");
1403 pull->dbuf_cyl = NULL;
1404 pull->bRefAt = FALSE;
1406 for (g = 0; g < pull->ngroup; g++)
1408 pgrp = &pull->group[g];
1409 pgrp->epgrppbc = epgrppbcNONE;
1412 gmx_bool bConstraint;
1413 ivec pulldim, pulldim_con;
1415 /* There is an issue when a group is used in multiple coordinates
1416 * and constraints are applied in different dimensions with atoms
1417 * frozen in some, but not all dimensions.
1418 * Since there is only one mass array per group, we can't have
1419 * frozen/non-frozen atoms for different coords at the same time.
1420 * But since this is a very exotic case, we don't check for this.
1421 * A warning is printed in init_pull_group_index.
1423 bConstraint = FALSE;
1424 clear_ivec(pulldim);
1425 clear_ivec(pulldim_con);
1426 for (c = 0; c < pull->ncoord; c++)
1428 if (pull->coord[c].group[0] == g ||
1429 pull->coord[c].group[1] == g)
1431 for (m = 0; m < DIM; m++)
1433 if (pull->coord[c].dim[m] == 1)
1437 if (pull->coord[c].eType == epullCONSTRAINT)
1447 /* Determine if we need to take PBC into account for calculating
1448 * the COM's of the pull groups.
1450 for (m = 0; m < pull->npbcdim; m++)
1452 if (pulldim[m] == 1 && pgrp->nat > 1)
1454 if (pgrp->pbcatom >= 0)
1456 pgrp->epgrppbc = epgrppbcREFAT;
1457 pull->bRefAt = TRUE;
1463 gmx_fatal(FARGS, "Pull groups can not have relative weights and cosine weighting at same time");
1465 pgrp->epgrppbc = epgrppbcCOS;
1466 if (pull->cosdim >= 0 && pull->cosdim != m)
1468 gmx_fatal(FARGS, "Can only use cosine weighting with pulling in one dimension (use mdp option pull-coord?-dim)");
1474 /* Set the indices */
1475 init_pull_group_index(fplog, cr, g, pgrp,
1476 bConstraint, pulldim_con,
1481 /* Absolute reference, set the inverse mass to zero.
1482 * This is only relevant (and used) with constraint pulling.
1489 /* If we use cylinder coordinates, do some initialising for them */
1490 if (pull->bCylinder)
1492 snew(pull->dyna, pull->ncoord);
1494 for (c = 0; c < pull->ncoord; c++)
1496 const t_pull_coord *pcrd;
1498 pcrd = &pull->coord[c];
1500 if (pcrd->eGeom == epullgCYL)
1502 if (pull->group[pcrd->group[0]].nat == 0)
1504 gmx_fatal(FARGS, "A cylinder pull group is not supported when using absolute reference!\n");
1510 /* We still need to initialize the PBC reference coordinates */
1511 pull->bSetPBCatoms = TRUE;
1513 /* Only do I/O when we are doing dynamics and if we are the MASTER */
1518 if (pull->nstxout > 0)
1520 pull->out_x = open_pull_out(opt2fn("-px", nfile, fnm), pull, oenv,
1523 if (pull->nstfout > 0)
1525 pull->out_f = open_pull_out(opt2fn("-pf", nfile, fnm), pull, oenv,
1531 void finish_pull(t_pull *pull)
1535 gmx_fio_fclose(pull->out_x);
1539 gmx_fio_fclose(pull->out_f);