Merge branch release-2016 into release-2018
[gromacs.git] / src / gromacs / pulling / pull.cpp
blob9a4720204103bf7a21045f4169f2b066e58a2336
1 /*
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,2016,2017, 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.
37 #include "gmxpre.h"
39 #include "pull.h"
41 #include "config.h"
43 #include <assert.h>
44 #include <stdio.h>
45 #include <stdlib.h>
46 #include <string.h>
48 #include <cmath>
50 #include <algorithm>
52 #include "gromacs/commandline/filenm.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/domdec/ga2la.h"
55 #include "gromacs/fileio/gmxfio.h"
56 #include "gromacs/fileio/xvgr.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/math/functions.h"
59 #include "gromacs/math/units.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/math/vectypes.h"
62 #include "gromacs/mdlib/gmx_omp_nthreads.h"
63 #include "gromacs/mdlib/mdrun.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/forceoutput.h"
66 #include "gromacs/mdtypes/inputrec.h"
67 #include "gromacs/mdtypes/md_enums.h"
68 #include "gromacs/mdtypes/mdatom.h"
69 #include "gromacs/pbcutil/pbc.h"
70 #include "gromacs/pulling/pull_internal.h"
71 #include "gromacs/topology/mtop_lookup.h"
72 #include "gromacs/topology/topology.h"
73 #include "gromacs/utility/cstringutil.h"
74 #include "gromacs/utility/exceptions.h"
75 #include "gromacs/utility/fatalerror.h"
76 #include "gromacs/utility/futil.h"
77 #include "gromacs/utility/gmxassert.h"
78 #include "gromacs/utility/mutex.h"
79 #include "gromacs/utility/pleasecite.h"
80 #include "gromacs/utility/real.h"
81 #include "gromacs/utility/smalloc.h"
83 bool pull_coordinate_is_angletype(const t_pull_coord *pcrd)
85 return (pcrd->eGeom == epullgANGLE ||
86 pcrd->eGeom == epullgDIHEDRAL ||
87 pcrd->eGeom == epullgANGLEAXIS);
90 static bool pull_coordinate_is_directional(const t_pull_coord *pcrd)
92 return (pcrd->eGeom == epullgDIR ||
93 pcrd->eGeom == epullgDIRPBC ||
94 pcrd->eGeom == epullgDIRRELATIVE ||
95 pcrd->eGeom == epullgCYL);
98 const char *pull_coordinate_units(const t_pull_coord *pcrd)
100 return pull_coordinate_is_angletype(pcrd) ? "deg" : "nm";
103 double pull_conversion_factor_userinput2internal(const t_pull_coord *pcrd)
105 if (pull_coordinate_is_angletype(pcrd))
107 return DEG2RAD;
109 else
111 return 1.0;
115 double pull_conversion_factor_internal2userinput(const t_pull_coord *pcrd)
117 if (pull_coordinate_is_angletype(pcrd))
119 return RAD2DEG;
121 else
123 return 1.0;
127 static std::string append_before_extension(std::string pathname,
128 std::string to_append)
130 /* Appends to_append before last '.' in pathname */
131 size_t extPos = pathname.find_last_of('.');
132 if (extPos == std::string::npos)
134 return pathname + to_append;
136 else
138 return pathname.substr(0, extPos) + to_append +
139 pathname.substr(extPos, std::string::npos);
143 static void pull_print_group_x(FILE *out, ivec dim,
144 const pull_group_work_t *pgrp)
146 int m;
148 for (m = 0; m < DIM; m++)
150 if (dim[m])
152 fprintf(out, "\t%g", pgrp->x[m]);
157 static void pull_print_coord_dr_components(FILE *out, const ivec dim, const dvec dr)
159 for (int m = 0; m < DIM; m++)
161 if (dim[m])
163 fprintf(out, "\t%g", dr[m]);
168 static void pull_print_coord_dr(FILE *out, const pull_coord_work_t *pcrd,
169 gmx_bool bPrintRefValue,
170 gmx_bool bPrintComponents)
172 double unit_factor = pull_conversion_factor_internal2userinput(&pcrd->params);
174 fprintf(out, "\t%g", pcrd->value*unit_factor);
176 if (bPrintRefValue)
178 fprintf(out, "\t%g", pcrd->value_ref*unit_factor);
181 if (bPrintComponents)
183 pull_print_coord_dr_components(out, pcrd->params.dim, pcrd->dr01);
184 if (pcrd->params.ngroup >= 4)
186 pull_print_coord_dr_components(out, pcrd->params.dim, pcrd->dr23);
188 if (pcrd->params.ngroup >= 6)
190 pull_print_coord_dr_components(out, pcrd->params.dim, pcrd->dr45);
195 static void pull_print_x(FILE *out, struct pull_t *pull, double t)
197 int c;
199 fprintf(out, "%.4f", t);
201 for (c = 0; c < pull->ncoord; c++)
203 pull_coord_work_t *pcrd;
205 pcrd = &pull->coord[c];
207 pull_print_coord_dr(out, pcrd,
208 pull->params.bPrintRefValue && pcrd->params.eType != epullEXTERNAL,
209 pull->params.bPrintComp);
211 if (pull->params.bPrintCOM)
213 if (pcrd->params.eGeom == epullgCYL)
215 pull_print_group_x(out, pcrd->params.dim, &pull->dyna[c]);
217 else
219 pull_print_group_x(out, pcrd->params.dim,
220 &pull->group[pcrd->params.group[0]]);
222 for (int g = 1; g < pcrd->params.ngroup; g++)
224 pull_print_group_x(out, pcrd->params.dim, &pull->group[pcrd->params.group[g]]);
228 fprintf(out, "\n");
231 static void pull_print_f(FILE *out, struct pull_t *pull, double t)
233 int c;
235 fprintf(out, "%.4f", t);
237 for (c = 0; c < pull->ncoord; c++)
239 fprintf(out, "\t%g", pull->coord[c].f_scal);
241 fprintf(out, "\n");
244 void pull_print_output(struct pull_t *pull, gmx_int64_t step, double time)
246 GMX_ASSERT(pull->numExternalPotentialsStillToBeAppliedThisStep == 0, "pull_print_output called before all external pull potentials have been applied");
248 if ((pull->params.nstxout != 0) && (step % pull->params.nstxout == 0))
250 pull_print_x(pull->out_x, pull, time);
253 if ((pull->params.nstfout != 0) && (step % pull->params.nstfout == 0))
255 pull_print_f(pull->out_f, pull, time);
259 static void set_legend_for_coord_components(const pull_coord_work_t *pcrd, int coord_index, char **setname, int *nsets_ptr)
261 /* Loop over the distance vectors and print their components. Each vector is made up of two consecutive groups. */
262 for (int g = 0; g < pcrd->params.ngroup; g += 2)
264 /* Loop over the components */
265 for (int m = 0; m < DIM; m++)
267 if (pcrd->params.dim[m])
269 char legend[STRLEN];
271 if (g == 0 && pcrd->params.ngroup <= 2)
273 /* For the simplest case we print a simplified legend without group indices, just the cooordinate index
274 and which dimensional component it is. */
275 sprintf(legend, "%d d%c", coord_index + 1, 'X' + m);
277 else
279 /* Otherwise, print also the group indices (relative to the coordinate, not the global ones). */
280 sprintf(legend, "%d g %d-%d d%c", coord_index + 1, g + 1, g + 2, 'X' + m);
283 setname[*nsets_ptr] = gmx_strdup(legend);
284 (*nsets_ptr)++;
290 static FILE *open_pull_out(const char *fn, struct pull_t *pull,
291 const gmx_output_env_t *oenv,
292 gmx_bool bCoord,
293 const ContinuationOptions &continuationOptions)
295 FILE *fp;
296 int nsets, c, m;
297 char **setname, buf[50];
299 if (continuationOptions.appendFiles)
301 fp = gmx_fio_fopen(fn, "a+");
303 else
305 fp = gmx_fio_fopen(fn, "w+");
306 if (bCoord)
308 sprintf(buf, "Position (nm%s)", pull->bAngle ? ", deg" : "");
309 xvgr_header(fp, "Pull COM", "Time (ps)", buf,
310 exvggtXNY, oenv);
312 else
314 sprintf(buf, "Force (kJ/mol/nm%s)", pull->bAngle ? ", kJ/mol/rad" : "");
315 xvgr_header(fp, "Pull force", "Time (ps)", buf,
316 exvggtXNY, oenv);
319 /* With default mdp options only the actual coordinate value is printed (1),
320 * but optionally the reference value (+ 1),
321 * the group COMs for all the groups (+ ngroups_max*DIM)
322 * and the components of the distance vectors can be printed (+ (ngroups_max/2)*DIM).
324 snew(setname, pull->ncoord*(1 + 1 + c_pullCoordNgroupMax*DIM + c_pullCoordNgroupMax/2*DIM));
326 nsets = 0;
327 for (c = 0; c < pull->ncoord; c++)
329 if (bCoord)
331 /* The order of this legend should match the order of printing
332 * the data in print_pull_x above.
335 /* The pull coord distance */
336 sprintf(buf, "%d", c+1);
337 setname[nsets] = gmx_strdup(buf);
338 nsets++;
339 if (pull->params.bPrintRefValue &&
340 pull->coord[c].params.eType != epullEXTERNAL)
342 sprintf(buf, "%d ref", c+1);
343 setname[nsets] = gmx_strdup(buf);
344 nsets++;
346 if (pull->params.bPrintComp)
348 set_legend_for_coord_components(&pull->coord[c], c, setname, &nsets);
351 if (pull->params.bPrintCOM)
353 for (int g = 0; g < pull->coord[c].params.ngroup; g++)
355 /* Legend for reference group position */
356 for (m = 0; m < DIM; m++)
358 if (pull->coord[c].params.dim[m])
360 sprintf(buf, "%d g %d %c", c+1, g + 1, 'X'+m);
361 setname[nsets] = gmx_strdup(buf);
362 nsets++;
368 else
370 /* For the pull force we always only use one scalar */
371 sprintf(buf, "%d", c+1);
372 setname[nsets] = gmx_strdup(buf);
373 nsets++;
376 if (nsets > 1)
378 xvgr_legend(fp, nsets, (const char**)setname, oenv);
380 for (c = 0; c < nsets; c++)
382 sfree(setname[c]);
384 sfree(setname);
387 return fp;
390 /* Apply forces in a mass weighted fashion for part of the pull group */
391 static void apply_forces_grp_part(const pull_group_work_t *pgrp,
392 int ind_start, int ind_end,
393 const t_mdatoms *md,
394 const dvec f_pull, int sign, rvec *f)
396 double inv_wm = pgrp->mwscale;
398 for (int i = ind_start; i < ind_end; i++)
400 int ii = pgrp->ind_loc[i];
401 double wmass = md->massT[ii];
402 if (pgrp->weight_loc)
404 wmass *= pgrp->weight_loc[i];
407 for (int d = 0; d < DIM; d++)
409 f[ii][d] += sign * wmass * f_pull[d] * inv_wm;
414 /* Apply forces in a mass weighted fashion */
415 static void apply_forces_grp(const pull_group_work_t *pgrp,
416 const t_mdatoms *md,
417 const dvec f_pull, int sign, rvec *f,
418 int nthreads)
420 if (pgrp->params.nat == 1 && pgrp->nat_loc == 1)
422 /* Only one atom and our rank has this atom: we can skip
423 * the mass weighting, which means that this code also works
424 * for mass=0, e.g. with a virtual site.
426 for (int d = 0; d < DIM; d++)
428 f[pgrp->ind_loc[0]][d] += sign*f_pull[d];
431 else
433 if (pgrp->nat_loc <= c_pullMaxNumLocalAtomsSingleThreaded)
435 apply_forces_grp_part(pgrp, 0, pgrp->nat_loc, md, f_pull, sign, f);
437 else
439 #pragma omp parallel for num_threads(nthreads) schedule(static)
440 for (int th = 0; th < nthreads; th++)
442 int ind_start = (pgrp->nat_loc*(th + 0))/nthreads;
443 int ind_end = (pgrp->nat_loc*(th + 1))/nthreads;
444 apply_forces_grp_part(pgrp, ind_start, ind_end,
445 md, f_pull, sign, f);
451 /* Apply forces in a mass weighted fashion to a cylinder group */
452 static void apply_forces_cyl_grp(const pull_group_work_t *pgrp,
453 const double dv_corr,
454 const t_mdatoms *md,
455 const dvec f_pull, double f_scal,
456 int sign, rvec *f,
457 int gmx_unused nthreads)
459 double inv_wm = pgrp->mwscale;
461 /* The cylinder group is always a slab in the system, thus large.
462 * Therefore we always thread-parallelize this group.
464 #pragma omp parallel for num_threads(nthreads) schedule(static)
465 for (int i = 0; i < pgrp->nat_loc; i++)
467 int ii = pgrp->ind_loc[i];
468 double mass = md->massT[ii];
469 double weight = pgrp->weight_loc[i];
470 /* The stored axial distance from the cylinder center (dv) needs
471 * to be corrected for an offset (dv_corr), which was unknown when
472 * we calculated dv.
474 double dv_com = pgrp->dv[i] + dv_corr;
476 /* Here we not only add the pull force working along vec (f_pull),
477 * but also a radial component, due to the dependence of the weights
478 * on the radial distance.
480 for (int m = 0; m < DIM; m++)
482 f[ii][m] += sign*inv_wm*(mass*weight*f_pull[m] +
483 pgrp->mdw[i][m]*dv_com*f_scal);
488 /* Apply torque forces in a mass weighted fashion to the groups that define
489 * the pull vector direction for pull coordinate pcrd.
491 static void apply_forces_vec_torque(const struct pull_t *pull,
492 const pull_coord_work_t *pcrd,
493 const t_mdatoms *md,
494 rvec *f)
496 double inpr;
497 int m;
498 dvec f_perp;
500 /* The component inpr along the pull vector is accounted for in the usual
501 * way. Here we account for the component perpendicular to vec.
503 inpr = 0;
504 for (m = 0; m < DIM; m++)
506 inpr += pcrd->dr01[m]*pcrd->vec[m];
508 /* The torque force works along the component of the distance vector
509 * of between the two "usual" pull groups that is perpendicular to
510 * the pull vector. The magnitude of this force is the "usual" scale force
511 * multiplied with the ratio of the distance between the two "usual" pull
512 * groups and the distance between the two groups that define the vector.
514 for (m = 0; m < DIM; m++)
516 f_perp[m] = (pcrd->dr01[m] - inpr*pcrd->vec[m])/pcrd->vec_len*pcrd->f_scal;
519 /* Apply the force to the groups defining the vector using opposite signs */
520 apply_forces_grp(&pull->group[pcrd->params.group[2]], md,
521 f_perp, -1, f, pull->nthreads);
522 apply_forces_grp(&pull->group[pcrd->params.group[3]], md,
523 f_perp, 1, f, pull->nthreads);
526 /* Apply forces in a mass weighted fashion */
527 static void apply_forces_coord(struct pull_t * pull, int coord,
528 const t_mdatoms * md,
529 rvec *f)
531 /* Here it would be more efficient to use one large thread-parallel
532 * region instead of potential parallel regions within apply_forces_grp.
533 * But there could be overlap between pull groups and this would lead
534 * to data races.
537 const pull_coord_work_t *pcrd = &pull->coord[coord];
539 if (pcrd->params.eGeom == epullgCYL)
541 apply_forces_cyl_grp(&pull->dyna[coord], pcrd->cyl_dev, md,
542 pcrd->f01, pcrd->f_scal, -1, f,
543 pull->nthreads);
545 /* Sum the force along the vector and the radial force */
546 dvec f_tot;
547 for (int m = 0; m < DIM; m++)
549 f_tot[m] = pcrd->f01[m] + pcrd->f_scal*pcrd->ffrad[m];
551 apply_forces_grp(&pull->group[pcrd->params.group[1]], md,
552 f_tot, 1, f, pull->nthreads);
554 else
556 if (pcrd->params.eGeom == epullgDIRRELATIVE)
558 /* We need to apply the torque forces to the pull groups
559 * that define the pull vector.
561 apply_forces_vec_torque(pull, pcrd, md, f);
564 if (pull->group[pcrd->params.group[0]].params.nat > 0)
566 apply_forces_grp(&pull->group[pcrd->params.group[0]], md,
567 pcrd->f01, -1, f, pull->nthreads);
569 apply_forces_grp(&pull->group[pcrd->params.group[1]], md,
570 pcrd->f01, 1, f, pull->nthreads);
572 if (pcrd->params.ngroup >= 4)
574 apply_forces_grp(&pull->group[pcrd->params.group[2]], md,
575 pcrd->f23, -1, f, pull->nthreads);
576 apply_forces_grp(&pull->group[pcrd->params.group[3]], md,
577 pcrd->f23, 1, f, pull->nthreads);
579 if (pcrd->params.ngroup >= 6)
581 apply_forces_grp(&pull->group[pcrd->params.group[4]], md,
582 pcrd->f45, -1, f, pull->nthreads);
583 apply_forces_grp(&pull->group[pcrd->params.group[5]], md,
584 pcrd->f45, 1, f, pull->nthreads);
589 real max_pull_distance2(const pull_coord_work_t *pcrd,
590 const t_pbc *pbc)
592 /* Note that this maximum distance calculation is more complex than
593 * most other cases in GROMACS, since here we have to take care of
594 * distance calculations that don't involve all three dimensions.
595 * For example, we can use distances that are larger than the
596 * box X and Y dimensions for a box that is elongated along Z.
599 real max_d2 = GMX_REAL_MAX;
601 if (pull_coordinate_is_directional(&pcrd->params))
603 /* Directional pulling along along direction pcrd->vec.
604 * Calculating the exact maximum distance is complex and bug-prone.
605 * So we take a safe approach by not allowing distances that
606 * are larger than half the distance between unit cell faces
607 * along dimensions involved in pcrd->vec.
609 for (int m = 0; m < DIM; m++)
611 if (m < pbc->ndim_ePBC && pcrd->vec[m] != 0)
613 real imageDistance2 = gmx::square(pbc->box[m][m]);
614 for (int d = m + 1; d < DIM; d++)
616 imageDistance2 -= gmx::square(pbc->box[d][m]);
618 max_d2 = std::min(max_d2, imageDistance2);
622 else
624 /* Distance pulling along dimensions with pcrd->params.dim[d]==1.
625 * We use half the minimum box vector length of the dimensions involved.
626 * This is correct for all cases, except for corner cases with
627 * triclinic boxes where e.g. dim[XX]=1 and dim[YY]=0 with
628 * box[YY][XX]!=0 and box[YY][YY] < box[XX][XX]. But since even
629 * in such corner cases the user could get correct results,
630 * depending on the details of the setup, we avoid further
631 * code complications.
633 for (int m = 0; m < DIM; m++)
635 if (m < pbc->ndim_ePBC && pcrd->params.dim[m] != 0)
637 real imageDistance2 = gmx::square(pbc->box[m][m]);
638 for (int d = 0; d < m; d++)
640 if (pcrd->params.dim[d] != 0)
642 imageDistance2 += gmx::square(pbc->box[m][d]);
645 max_d2 = std::min(max_d2, imageDistance2);
650 return 0.25*max_d2;
653 /* This function returns the distance based on coordinates xg and xref.
654 * Note that the pull coordinate struct pcrd is not modified.
656 static void low_get_pull_coord_dr(const struct pull_t *pull,
657 const pull_coord_work_t *pcrd,
658 const t_pbc *pbc,
659 dvec xg, dvec xref, double max_dist2,
660 dvec dr)
662 const pull_group_work_t *pgrp0 = &pull->group[pcrd->params.group[0]];
664 /* Only the first group can be an absolute reference, in that case nat=0 */
665 if (pgrp0->params.nat == 0)
667 for (int m = 0; m < DIM; m++)
669 xref[m] = pcrd->params.origin[m];
673 dvec xrefr;
674 copy_dvec(xref, xrefr);
676 dvec dref = {0, 0, 0};
677 if (pcrd->params.eGeom == epullgDIRPBC)
679 for (int m = 0; m < DIM; m++)
681 dref[m] = pcrd->value_ref*pcrd->vec[m];
683 /* Add the reference position, so we use the correct periodic image */
684 dvec_inc(xrefr, dref);
687 pbc_dx_d(pbc, xg, xrefr, dr);
689 bool directional = pull_coordinate_is_directional(&pcrd->params);
690 double dr2 = 0;
691 for (int m = 0; m < DIM; m++)
693 dr[m] *= pcrd->params.dim[m];
694 if (pcrd->params.dim[m] && !(directional && pcrd->vec[m] == 0))
696 dr2 += dr[m]*dr[m];
699 /* Check if we are close to switching to another periodic image */
700 if (max_dist2 > 0 && dr2 > 0.98*0.98*max_dist2)
702 /* Note that technically there is no issue with switching periodic
703 * image, as pbc_dx_d returns the distance to the closest periodic
704 * image. However in all cases where periodic image switches occur,
705 * the pull results will be useless in practice.
707 gmx_fatal(FARGS, "Distance between pull groups %d and %d (%f nm) is larger than 0.49 times the box size (%f).\n%s",
708 pcrd->params.group[0], pcrd->params.group[1],
709 sqrt(dr2), sqrt(0.98*0.98*max_dist2),
710 pcrd->params.eGeom == epullgDIR ? "You might want to consider using \"pull-geometry = direction-periodic\" instead.\n" : "");
713 if (pcrd->params.eGeom == epullgDIRPBC)
715 dvec_inc(dr, dref);
719 /* This function returns the distance based on the contents of the pull struct.
720 * pull->coord[coord_ind].dr, and potentially vector, are updated.
722 static void get_pull_coord_dr(struct pull_t *pull,
723 int coord_ind,
724 const t_pbc *pbc)
726 double md2;
727 pull_coord_work_t *pcrd;
728 pull_group_work_t *pgrp0, *pgrp1;
730 pcrd = &pull->coord[coord_ind];
732 if (pcrd->params.eGeom == epullgDIRPBC)
734 md2 = -1;
736 else
738 md2 = static_cast<double>(max_pull_distance2(pcrd, pbc));
741 if (pcrd->params.eGeom == epullgDIRRELATIVE)
743 /* We need to determine the pull vector */
744 const pull_group_work_t *pgrp2, *pgrp3;
745 dvec vec;
746 int m;
748 pgrp2 = &pull->group[pcrd->params.group[2]];
749 pgrp3 = &pull->group[pcrd->params.group[3]];
751 pbc_dx_d(pbc, pgrp3->x, pgrp2->x, vec);
753 for (m = 0; m < DIM; m++)
755 vec[m] *= pcrd->params.dim[m];
757 pcrd->vec_len = dnorm(vec);
758 for (m = 0; m < DIM; m++)
760 pcrd->vec[m] = vec[m]/pcrd->vec_len;
762 if (debug)
764 fprintf(debug, "pull coord %d vector: %6.3f %6.3f %6.3f normalized: %6.3f %6.3f %6.3f\n",
765 coord_ind,
766 vec[XX], vec[YY], vec[ZZ],
767 pcrd->vec[XX], pcrd->vec[YY], pcrd->vec[ZZ]);
771 pgrp0 = &pull->group[pcrd->params.group[0]];
772 pgrp1 = &pull->group[pcrd->params.group[1]];
774 low_get_pull_coord_dr(pull, pcrd, pbc,
775 pgrp1->x,
776 pcrd->params.eGeom == epullgCYL ? pull->dyna[coord_ind].x : pgrp0->x,
777 md2,
778 pcrd->dr01);
780 if (pcrd->params.ngroup >= 4)
782 pull_group_work_t *pgrp2, *pgrp3;
783 pgrp2 = &pull->group[pcrd->params.group[2]];
784 pgrp3 = &pull->group[pcrd->params.group[3]];
786 low_get_pull_coord_dr(pull, pcrd, pbc,
787 pgrp3->x,
788 pgrp2->x,
789 md2,
790 pcrd->dr23);
792 if (pcrd->params.ngroup >= 6)
794 pull_group_work_t *pgrp4, *pgrp5;
795 pgrp4 = &pull->group[pcrd->params.group[4]];
796 pgrp5 = &pull->group[pcrd->params.group[5]];
798 low_get_pull_coord_dr(pull, pcrd, pbc,
799 pgrp5->x,
800 pgrp4->x,
801 md2,
802 pcrd->dr45);
806 /* Modify x so that it is periodic in [-pi, pi)
807 * It is assumed that x is in [-3pi, 3pi) so that x
808 * needs to be shifted by at most one period.
810 static void make_periodic_2pi(double *x)
812 if (*x >= M_PI)
814 *x -= M_2PI;
816 else if (*x < -M_PI)
818 *x += M_2PI;
822 /* This function should always be used to modify pcrd->value_ref */
823 static void low_set_pull_coord_reference_value(pull_coord_work_t *pcrd,
824 int coord_ind,
825 double value_ref)
827 GMX_ASSERT(pcrd->params.eType != epullEXTERNAL, "The pull coord reference value should not be used with type external-potential");
829 if (pcrd->params.eGeom == epullgDIST)
831 if (value_ref < 0)
833 gmx_fatal(FARGS, "Pull reference distance for coordinate %d (%f) needs to be non-negative", coord_ind + 1, value_ref);
836 else if (pcrd->params.eGeom == epullgANGLE || pcrd->params.eGeom == epullgANGLEAXIS)
838 if (value_ref < 0 || value_ref > M_PI)
840 gmx_fatal(FARGS, "Pull reference angle for coordinate %d (%f) needs to be in the allowed interval [0,180] deg", coord_ind + 1, value_ref*pull_conversion_factor_internal2userinput(&pcrd->params));
843 else if (pcrd->params.eGeom == epullgDIHEDRAL)
845 /* Allow pulling to be periodic for dihedral angles by remapping the reference value to the interval [-pi, pi). */
846 make_periodic_2pi(&value_ref);
849 pcrd->value_ref = value_ref;
852 static void update_pull_coord_reference_value(pull_coord_work_t *pcrd, int coord_ind, double t)
854 /* With zero rate the reference value is set initially and doesn't change */
855 if (pcrd->params.rate != 0)
857 double value_ref = (pcrd->params.init + pcrd->params.rate*t)*pull_conversion_factor_userinput2internal(&pcrd->params);
858 low_set_pull_coord_reference_value(pcrd, coord_ind, value_ref);
862 /* Returns the dihedral angle. Updates the plane normal vectors m, n. */
863 static double get_dihedral_angle_coord(pull_coord_work_t *pcrd)
865 double phi, sign;
866 dvec dr32; /* store instead of dr23? */
868 dsvmul(-1, pcrd->dr23, dr32);
869 dcprod(pcrd->dr01, dr32, pcrd->planevec_m); /* Normal of first plane */
870 dcprod(dr32, pcrd->dr45, pcrd->planevec_n); /* Normal of second plane */
871 phi = gmx_angle_between_dvecs(pcrd->planevec_m, pcrd->planevec_n);
873 /* Note 1: the sign below is opposite of that in the bondeds or Bekker 1994
874 * because there r_ij = ri - rj, while here dr01 = r_1 - r_0
875 * Note 2: the angle between the plane normal vectors equals pi only when
876 * both planes coincide. Thus, when phi = pi dr01 will lie in both planes and
877 * we get a positive sign below. Thus, the range of the dihedral angle is (-180, 180].
879 sign = (diprod(pcrd->dr01, pcrd->planevec_n) < 0.0) ? 1.0 : -1.0;
880 return sign*phi;
883 /* Calculates pull->coord[coord_ind].value.
884 * This function also updates pull->coord[coord_ind].dr.
886 static void get_pull_coord_distance(struct pull_t *pull,
887 int coord_ind,
888 const t_pbc *pbc)
890 pull_coord_work_t *pcrd;
891 int m;
893 pcrd = &pull->coord[coord_ind];
895 get_pull_coord_dr(pull, coord_ind, pbc);
897 switch (pcrd->params.eGeom)
899 case epullgDIST:
900 /* Pull along the vector between the com's */
901 pcrd->value = dnorm(pcrd->dr01);
902 break;
903 case epullgDIR:
904 case epullgDIRPBC:
905 case epullgDIRRELATIVE:
906 case epullgCYL:
907 /* Pull along vec */
908 pcrd->value = 0;
909 for (m = 0; m < DIM; m++)
911 pcrd->value += pcrd->vec[m]*pcrd->dr01[m];
913 break;
914 case epullgANGLE:
915 pcrd->value = gmx_angle_between_dvecs(pcrd->dr01, pcrd->dr23);
916 break;
917 case epullgDIHEDRAL:
918 pcrd->value = get_dihedral_angle_coord(pcrd);
919 break;
920 case epullgANGLEAXIS:
921 pcrd->value = gmx_angle_between_dvecs(pcrd->dr01, pcrd->vec);
922 break;
923 default:
924 gmx_incons("Unsupported pull type in get_pull_coord_distance");
928 /* Returns the deviation from the reference value.
929 * Updates pull->coord[coord_ind].dr, .value and .value_ref.
931 static double get_pull_coord_deviation(struct pull_t *pull,
932 int coord_ind,
933 const t_pbc *pbc,
934 double t)
936 pull_coord_work_t *pcrd;
937 double dev = 0;
939 pcrd = &pull->coord[coord_ind];
941 get_pull_coord_distance(pull, coord_ind, pbc);
943 update_pull_coord_reference_value(pcrd, coord_ind, t);
945 /* Determine the deviation */
946 dev = pcrd->value - pcrd->value_ref;
948 /* Check that values are allowed */
949 if (pcrd->params.eGeom == epullgDIST && pcrd->value == 0)
951 /* With no vector we can not determine the direction for the force,
952 * so we set the force to zero.
954 dev = 0;
956 else if (pcrd->params.eGeom == epullgDIHEDRAL)
958 /* The reference value is in [-pi, pi). The coordinate value is in (-pi, pi].
959 Thus, the unwrapped deviation is here in (-2pi, 2pi].
960 After making it periodic, the deviation will be in [-pi, pi). */
961 make_periodic_2pi(&dev);
964 return dev;
967 double get_pull_coord_value(struct pull_t *pull,
968 int coord_ind,
969 const t_pbc *pbc)
971 get_pull_coord_distance(pull, coord_ind, pbc);
973 return pull->coord[coord_ind].value;
976 static void clear_pull_forces_coord(pull_coord_work_t *pcrd)
978 clear_dvec(pcrd->f01);
979 clear_dvec(pcrd->f23);
980 clear_dvec(pcrd->f45);
981 pcrd->f_scal = 0;
984 void clear_pull_forces(struct pull_t *pull)
986 int i;
988 /* Zeroing the forces is only required for constraint pulling.
989 * It can happen that multiple constraint steps need to be applied
990 * and therefore the constraint forces need to be accumulated.
992 for (i = 0; i < pull->ncoord; i++)
994 clear_pull_forces_coord(&pull->coord[i]);
998 /* Apply constraint using SHAKE */
999 static void do_constraint(struct pull_t *pull, t_pbc *pbc,
1000 rvec *x, rvec *v,
1001 gmx_bool bMaster, tensor vir,
1002 double dt, double t)
1005 dvec *r_ij; /* x[i] com of i in prev. step. Obeys constr. -> r_ij[i] */
1006 dvec unc_ij; /* xp[i] com of i this step, before constr. -> unc_ij */
1007 dvec *rnew; /* current 'new' positions of the groups */
1008 double *dr_tot; /* the total update of the coords */
1009 dvec vec;
1010 double inpr;
1011 double lambda, rm, invdt = 0;
1012 gmx_bool bConverged_all, bConverged = FALSE;
1013 int niter = 0, g, c, ii, j, m, max_iter = 100;
1014 double a;
1015 dvec tmp, tmp3;
1017 snew(r_ij, pull->ncoord);
1018 snew(dr_tot, pull->ncoord);
1020 snew(rnew, pull->ngroup);
1022 /* copy the current unconstrained positions for use in iterations. We
1023 iterate until rinew[i] and rjnew[j] obey the constraints. Then
1024 rinew - pull.x_unc[i] is the correction dr to group i */
1025 for (g = 0; g < pull->ngroup; g++)
1027 copy_dvec(pull->group[g].xp, rnew[g]);
1030 /* Determine the constraint directions from the old positions */
1031 for (c = 0; c < pull->ncoord; c++)
1033 pull_coord_work_t *pcrd;
1035 pcrd = &pull->coord[c];
1037 if (pcrd->params.eType != epullCONSTRAINT)
1039 continue;
1042 /* Note that get_pull_coord_distance also sets pcrd->dr and pcrd->value.
1043 * We don't modify dr and value anymore, so these values are also used
1044 * for printing.
1046 get_pull_coord_distance(pull, c, pbc);
1048 if (debug)
1050 fprintf(debug, "Pull coord %d dr %f %f %f\n",
1051 c, pcrd->dr01[XX], pcrd->dr01[YY], pcrd->dr01[ZZ]);
1054 if (pcrd->params.eGeom == epullgDIR ||
1055 pcrd->params.eGeom == epullgDIRPBC)
1057 /* Select the component along vec */
1058 a = 0;
1059 for (m = 0; m < DIM; m++)
1061 a += pcrd->vec[m]*pcrd->dr01[m];
1063 for (m = 0; m < DIM; m++)
1065 r_ij[c][m] = a*pcrd->vec[m];
1068 else
1070 copy_dvec(pcrd->dr01, r_ij[c]);
1073 if (dnorm2(r_ij[c]) == 0)
1075 gmx_fatal(FARGS, "Distance for pull coordinate %d is zero with constraint pulling, which is not allowed.", c + 1);
1079 bConverged_all = FALSE;
1080 while (!bConverged_all && niter < max_iter)
1082 bConverged_all = TRUE;
1084 /* loop over all constraints */
1085 for (c = 0; c < pull->ncoord; c++)
1087 pull_coord_work_t *pcrd;
1088 pull_group_work_t *pgrp0, *pgrp1;
1089 dvec dr0, dr1;
1091 pcrd = &pull->coord[c];
1093 if (pcrd->params.eType != epullCONSTRAINT)
1095 continue;
1098 update_pull_coord_reference_value(pcrd, c, t);
1100 pgrp0 = &pull->group[pcrd->params.group[0]];
1101 pgrp1 = &pull->group[pcrd->params.group[1]];
1103 /* Get the current difference vector */
1104 low_get_pull_coord_dr(pull, pcrd, pbc,
1105 rnew[pcrd->params.group[1]],
1106 rnew[pcrd->params.group[0]],
1107 -1, unc_ij);
1109 if (debug)
1111 fprintf(debug, "Pull coord %d, iteration %d\n", c, niter);
1114 rm = 1.0/(pgrp0->invtm + pgrp1->invtm);
1116 switch (pcrd->params.eGeom)
1118 case epullgDIST:
1119 if (pcrd->value_ref <= 0)
1121 gmx_fatal(FARGS, "The pull constraint reference distance for group %d is <= 0 (%f)", c, pcrd->value_ref);
1125 double q, c_a, c_b, c_c;
1127 c_a = diprod(r_ij[c], r_ij[c]);
1128 c_b = diprod(unc_ij, r_ij[c])*2;
1129 c_c = diprod(unc_ij, unc_ij) - gmx::square(pcrd->value_ref);
1131 if (c_b < 0)
1133 q = -0.5*(c_b - std::sqrt(c_b*c_b - 4*c_a*c_c));
1134 lambda = -q/c_a;
1136 else
1138 q = -0.5*(c_b + std::sqrt(c_b*c_b - 4*c_a*c_c));
1139 lambda = -c_c/q;
1142 if (debug)
1144 fprintf(debug,
1145 "Pull ax^2+bx+c=0: a=%e b=%e c=%e lambda=%e\n",
1146 c_a, c_b, c_c, lambda);
1150 /* The position corrections dr due to the constraints */
1151 dsvmul(-lambda*rm*pgrp1->invtm, r_ij[c], dr1);
1152 dsvmul( lambda*rm*pgrp0->invtm, r_ij[c], dr0);
1153 dr_tot[c] += -lambda*dnorm(r_ij[c]);
1154 break;
1155 case epullgDIR:
1156 case epullgDIRPBC:
1157 case epullgCYL:
1158 /* A 1-dimensional constraint along a vector */
1159 a = 0;
1160 for (m = 0; m < DIM; m++)
1162 vec[m] = pcrd->vec[m];
1163 a += unc_ij[m]*vec[m];
1165 /* Select only the component along the vector */
1166 dsvmul(a, vec, unc_ij);
1167 lambda = a - pcrd->value_ref;
1168 if (debug)
1170 fprintf(debug, "Pull inpr %e lambda: %e\n", a, lambda);
1173 /* The position corrections dr due to the constraints */
1174 dsvmul(-lambda*rm*pgrp1->invtm, vec, dr1);
1175 dsvmul( lambda*rm*pgrp0->invtm, vec, dr0);
1176 dr_tot[c] += -lambda;
1177 break;
1178 default:
1179 gmx_incons("Invalid enumeration value for eGeom");
1180 /* Keep static analyzer happy */
1181 clear_dvec(dr0);
1182 clear_dvec(dr1);
1185 /* DEBUG */
1186 if (debug)
1188 int g0, g1;
1190 g0 = pcrd->params.group[0];
1191 g1 = pcrd->params.group[1];
1192 low_get_pull_coord_dr(pull, pcrd, pbc, rnew[g1], rnew[g0], -1, tmp);
1193 low_get_pull_coord_dr(pull, pcrd, pbc, dr1, dr0, -1, tmp3);
1194 fprintf(debug,
1195 "Pull cur %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
1196 rnew[g0][0], rnew[g0][1], rnew[g0][2],
1197 rnew[g1][0], rnew[g1][1], rnew[g1][2], dnorm(tmp));
1198 fprintf(debug,
1199 "Pull ref %8s %8s %8s %8s %8s %8s d: %8.5f\n",
1200 "", "", "", "", "", "", pcrd->value_ref);
1201 fprintf(debug,
1202 "Pull cor %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
1203 dr0[0], dr0[1], dr0[2],
1204 dr1[0], dr1[1], dr1[2],
1205 dnorm(tmp3));
1206 } /* END DEBUG */
1208 /* Update the COMs with dr */
1209 dvec_inc(rnew[pcrd->params.group[1]], dr1);
1210 dvec_inc(rnew[pcrd->params.group[0]], dr0);
1213 /* Check if all constraints are fullfilled now */
1214 for (c = 0; c < pull->ncoord; c++)
1216 pull_coord_work_t *pcrd;
1218 pcrd = &pull->coord[c];
1220 if (pcrd->params.eType != epullCONSTRAINT)
1222 continue;
1225 low_get_pull_coord_dr(pull, pcrd, pbc,
1226 rnew[pcrd->params.group[1]],
1227 rnew[pcrd->params.group[0]],
1228 -1, unc_ij);
1230 switch (pcrd->params.eGeom)
1232 case epullgDIST:
1233 bConverged =
1234 fabs(dnorm(unc_ij) - pcrd->value_ref) < pull->params.constr_tol;
1235 break;
1236 case epullgDIR:
1237 case epullgDIRPBC:
1238 case epullgCYL:
1239 for (m = 0; m < DIM; m++)
1241 vec[m] = pcrd->vec[m];
1243 inpr = diprod(unc_ij, vec);
1244 dsvmul(inpr, vec, unc_ij);
1245 bConverged =
1246 fabs(diprod(unc_ij, vec) - pcrd->value_ref) < pull->params.constr_tol;
1247 break;
1250 if (!bConverged)
1252 if (debug)
1254 fprintf(debug, "NOT CONVERGED YET: Group %d:"
1255 "d_ref = %f, current d = %f\n",
1256 g, pcrd->value_ref, dnorm(unc_ij));
1259 bConverged_all = FALSE;
1263 niter++;
1264 /* if after all constraints are dealt with and bConverged is still TRUE
1265 we're finished, if not we do another iteration */
1267 if (niter > max_iter)
1269 gmx_fatal(FARGS, "Too many iterations for constraint run: %d", niter);
1272 /* DONE ITERATING, NOW UPDATE COORDINATES AND CALC. CONSTRAINT FORCES */
1274 if (v)
1276 invdt = 1/dt;
1279 /* update atoms in the groups */
1280 for (g = 0; g < pull->ngroup; g++)
1282 const pull_group_work_t *pgrp;
1283 dvec dr;
1285 pgrp = &pull->group[g];
1287 /* get the final constraint displacement dr for group g */
1288 dvec_sub(rnew[g], pgrp->xp, dr);
1290 if (dnorm2(dr) == 0)
1292 /* No displacement, no update necessary */
1293 continue;
1296 /* update the atom positions */
1297 copy_dvec(dr, tmp);
1298 for (j = 0; j < pgrp->nat_loc; j++)
1300 ii = pgrp->ind_loc[j];
1301 if (pgrp->weight_loc)
1303 dsvmul(pgrp->wscale*pgrp->weight_loc[j], dr, tmp);
1305 for (m = 0; m < DIM; m++)
1307 x[ii][m] += tmp[m];
1309 if (v)
1311 for (m = 0; m < DIM; m++)
1313 v[ii][m] += invdt*tmp[m];
1319 /* calculate the constraint forces, used for output and virial only */
1320 for (c = 0; c < pull->ncoord; c++)
1322 pull_coord_work_t *pcrd;
1324 pcrd = &pull->coord[c];
1326 if (pcrd->params.eType != epullCONSTRAINT)
1328 continue;
1331 /* Accumulate the forces, in case we have multiple constraint steps */
1332 pcrd->f_scal += dr_tot[c]/((pull->group[pcrd->params.group[0]].invtm + pull->group[pcrd->params.group[1]].invtm)*dt*dt);
1334 if (vir != nullptr && pcrd->params.eGeom != epullgDIRPBC && bMaster)
1336 double f_invr;
1338 /* Add the pull contribution to the virial */
1339 /* We have already checked above that r_ij[c] != 0 */
1340 f_invr = pcrd->f_scal/dnorm(r_ij[c]);
1342 for (j = 0; j < DIM; j++)
1344 for (m = 0; m < DIM; m++)
1346 vir[j][m] -= 0.5*f_invr*r_ij[c][j]*r_ij[c][m];
1352 /* finished! I hope. Give back some memory */
1353 sfree(r_ij);
1354 sfree(dr_tot);
1355 sfree(rnew);
1358 static void add_virial_coord_dr(tensor vir, const dvec dr, const dvec f)
1360 for (int j = 0; j < DIM; j++)
1362 for (int m = 0; m < DIM; m++)
1364 vir[j][m] -= 0.5*f[j]*dr[m];
1369 /* Adds the pull contribution to the virial */
1370 static void add_virial_coord(tensor vir, const pull_coord_work_t *pcrd)
1372 if (vir != nullptr && pcrd->params.eGeom != epullgDIRPBC)
1374 /* Add the pull contribution for each distance vector to the virial. */
1375 add_virial_coord_dr(vir, pcrd->dr01, pcrd->f01);
1376 if (pcrd->params.ngroup >= 4)
1378 add_virial_coord_dr(vir, pcrd->dr23, pcrd->f23);
1380 if (pcrd->params.ngroup >= 6)
1382 add_virial_coord_dr(vir, pcrd->dr45, pcrd->f45);
1387 static void calc_pull_coord_scalar_force_and_potential(pull_coord_work_t *pcrd,
1388 double dev, real lambda,
1389 real *V, real *dVdl)
1391 real k, dkdl;
1393 k = (1.0 - lambda)*pcrd->params.k + lambda*pcrd->params.kB;
1394 dkdl = pcrd->params.kB - pcrd->params.k;
1396 switch (pcrd->params.eType)
1398 case epullUMBRELLA:
1399 case epullFLATBOTTOM:
1400 case epullFLATBOTTOMHIGH:
1401 /* The only difference between an umbrella and a flat-bottom
1402 * potential is that a flat-bottom is zero above or below
1403 the reference value.
1405 if ((pcrd->params.eType == epullFLATBOTTOM && dev < 0) ||
1406 (pcrd->params.eType == epullFLATBOTTOMHIGH && dev > 0))
1408 dev = 0;
1411 pcrd->f_scal = -k*dev;
1412 *V += 0.5* k*gmx::square(dev);
1413 *dVdl += 0.5*dkdl*gmx::square(dev);
1414 break;
1415 case epullCONST_F:
1416 pcrd->f_scal = -k;
1417 *V += k*pcrd->value;
1418 *dVdl += dkdl*pcrd->value;
1419 break;
1420 case epullEXTERNAL:
1421 gmx_incons("the scalar pull force should not be calculated internally for pull type external");
1422 break;
1423 default:
1424 gmx_incons("Unsupported pull type in do_pull_pot");
1428 static void calc_pull_coord_vector_force(pull_coord_work_t *pcrd)
1430 /* The geometry of the coordinate determines how the scalar force relates to the force on each group */
1432 if (pcrd->params.eGeom == epullgDIST)
1434 double invdr01 = pcrd->value > 0 ? 1./pcrd->value : 0.;
1435 for (int m = 0; m < DIM; m++)
1437 pcrd->f01[m] = pcrd->f_scal*pcrd->dr01[m]*invdr01;
1440 else if (pcrd->params.eGeom == epullgANGLE)
1443 double cos_theta, cos_theta2;
1445 cos_theta = cos(pcrd->value);
1446 cos_theta2 = gmx::square(cos_theta);
1448 /* The force at theta = 0, pi is undefined so we should not apply any force.
1449 * cos(theta) = 1 for theta = 0, pi so this is what we check for below.
1450 * Note that dr01 = 0 or dr23 = 0 evaluates to theta = 0 so we do not
1451 * have to check for this before dividing by their norm below.
1453 if (cos_theta2 < 1)
1455 /* The forces f01 and f23 can be decomposed into one vector parallel to dr01
1456 * and another vector parallel to dr23:
1457 * f01 = -dV/dtheta*(a*dr23/|dr23| - b*dr01*|dr01|)/|dr01|,
1458 * f23 = -dV/dtheta*(a*dr01/|dr01| - b*dr23*|dr23|)/|dr23|,
1460 double a = -gmx::invsqrt(1 - cos_theta2); /* comes from d/dx acos(x) */
1461 double b = a*cos_theta;
1462 double invdr01 = 1./dnorm(pcrd->dr01);
1463 double invdr23 = 1./dnorm(pcrd->dr23);
1464 dvec normalized_dr01, normalized_dr23;
1465 dsvmul(invdr01, pcrd->dr01, normalized_dr01);
1466 dsvmul(invdr23, pcrd->dr23, normalized_dr23);
1468 for (int m = 0; m < DIM; m++)
1470 /* Here, f_scal is -dV/dtheta */
1471 pcrd->f01[m] = pcrd->f_scal*invdr01*(a*normalized_dr23[m] - b*normalized_dr01[m]);
1472 pcrd->f23[m] = pcrd->f_scal*invdr23*(a*normalized_dr01[m] - b*normalized_dr23[m]);
1475 else
1477 /* No forces to apply for ill-defined cases*/
1478 clear_pull_forces_coord(pcrd);
1481 else if (pcrd->params.eGeom == epullgANGLEAXIS)
1483 double cos_theta, cos_theta2;
1485 /* The angle-axis force is exactly the same as the angle force (above) except that in
1486 this case the second vector (dr23) is replaced by the pull vector. */
1487 cos_theta = cos(pcrd->value);
1488 cos_theta2 = gmx::square(cos_theta);
1490 if (cos_theta2 < 1)
1492 double a, b;
1493 double invdr01;
1494 dvec normalized_dr01;
1496 invdr01 = 1./dnorm(pcrd->dr01);
1497 dsvmul(invdr01, pcrd->dr01, normalized_dr01);
1498 a = -gmx::invsqrt(1 - cos_theta2); /* comes from d/dx acos(x) */
1499 b = a*cos_theta;
1501 for (int m = 0; m < DIM; m++)
1503 pcrd->f01[m] = pcrd->f_scal*invdr01*(a*pcrd->vec[m] - b*normalized_dr01[m]);
1506 else
1508 clear_pull_forces_coord(pcrd);
1511 else if (pcrd->params.eGeom == epullgDIHEDRAL)
1513 double m2, n2, tol, sqrdist_32;
1514 dvec dr32;
1515 /* Note: there is a small difference here compared to the
1516 dihedral force calculations in the bondeds (ref: Bekker 1994).
1517 There rij = ri - rj, while here dr01 = r1 - r0.
1518 However, all distance vectors occur in form of cross or inner products
1519 so that two signs cancel and we end up with the same expressions.
1520 Also, we treat the more general case of 6 groups (0..5) instead of 4 (i, j, k, l).
1522 m2 = diprod(pcrd->planevec_m, pcrd->planevec_m);
1523 n2 = diprod(pcrd->planevec_n, pcrd->planevec_n);
1524 dsvmul(-1, pcrd->dr23, dr32);
1525 sqrdist_32 = diprod(dr32, dr32);
1526 tol = sqrdist_32*GMX_REAL_EPS; /* Avoid tiny angles */
1527 if ((m2 > tol) && (n2 > tol))
1529 double a_01, a_23_01, a_23_45, a_45;
1530 double inv_dist_32, inv_sqrdist_32, dist_32;
1531 dvec u, v;
1532 inv_dist_32 = gmx::invsqrt(sqrdist_32);
1533 inv_sqrdist_32 = inv_dist_32*inv_dist_32;
1534 dist_32 = sqrdist_32*inv_dist_32;
1536 /* Forces on groups 0, 1 */
1537 a_01 = pcrd->f_scal*dist_32/m2; /* f_scal is -dV/dphi */
1538 dsvmul(-a_01, pcrd->planevec_m, pcrd->f01); /* added sign to get force on group 1, not 0 */
1540 /* Forces on groups 4, 5 */
1541 a_45 = -pcrd->f_scal*dist_32/n2;
1542 dsvmul(a_45, pcrd->planevec_n, pcrd->f45); /* force on group 5 */
1544 /* Force on groups 2, 3 (defining the axis) */
1545 a_23_01 = -diprod(pcrd->dr01, dr32)*inv_sqrdist_32;
1546 a_23_45 = -diprod(pcrd->dr45, dr32)*inv_sqrdist_32;
1547 dsvmul(-a_23_01, pcrd->f01, u); /* added sign to get force from group 0, not 1 */
1548 dsvmul(a_23_45, pcrd->f45, v);
1549 dvec_sub(u, v, pcrd->f23); /* force on group 3 */
1551 else
1553 /* No force to apply for ill-defined cases */
1554 clear_pull_forces_coord(pcrd);
1557 else
1559 for (int m = 0; m < DIM; m++)
1561 pcrd->f01[m] = pcrd->f_scal*pcrd->vec[m];
1567 /* We use a global mutex for locking access to the pull data structure
1568 * during registration of external pull potential providers.
1569 * We could use a different, local mutex for each pull object, but the overhead
1570 * is extremely small here and registration is only done during initialization.
1572 static gmx::Mutex registrationMutex;
1574 using Lock = gmx::lock_guard<gmx::Mutex>;
1576 void register_external_pull_potential(struct pull_t *pull,
1577 int coord_index,
1578 const char *provider)
1580 GMX_RELEASE_ASSERT(pull != nullptr, "register_external_pull_potential called before init_pull");
1581 GMX_RELEASE_ASSERT(provider != nullptr, "register_external_pull_potential called with NULL as provider name");
1583 if (coord_index < 0 || coord_index > pull->ncoord - 1)
1585 gmx_fatal(FARGS, "Module '%s' attempted to register an external potential for pull coordinate %d which is out of the pull coordinate range %d - %d\n",
1586 provider, coord_index + 1, 1, pull->ncoord);
1589 pull_coord_work_t *pcrd = &pull->coord[coord_index];
1591 if (pcrd->params.eType != epullEXTERNAL)
1593 gmx_fatal(FARGS, "Module '%s' attempted to register an external potential for pull coordinate %d which of type '%s', whereas external potentials are only supported with type '%s'",
1594 provider, coord_index + 1, epull_names[pcrd->params.eType], epull_names[epullEXTERNAL]);
1597 GMX_RELEASE_ASSERT(pcrd->params.externalPotentialProvider != nullptr, "The external potential provider string for a pull coordinate is NULL");
1599 if (gmx_strcasecmp(provider, pcrd->params.externalPotentialProvider) != 0)
1601 gmx_fatal(FARGS, "Module '%s' attempted to register an external potential for pull coordinate %d which expects the external potential to be provided by a module named '%s'",
1602 provider, coord_index + 1, pcrd->params.externalPotentialProvider);
1605 /* Lock to avoid (extremely unlikely) simultaneous reading and writing of
1606 * pcrd->bExternalPotentialProviderHasBeenRegistered and
1607 * pull->numUnregisteredExternalPotentials.
1609 Lock registrationLock(registrationMutex);
1611 if (pcrd->bExternalPotentialProviderHasBeenRegistered)
1613 gmx_fatal(FARGS, "Module '%s' attempted to register an external potential for pull coordinate %d more than once",
1614 provider, coord_index + 1);
1617 pcrd->bExternalPotentialProviderHasBeenRegistered = true;
1618 pull->numUnregisteredExternalPotentials--;
1620 GMX_RELEASE_ASSERT(pull->numUnregisteredExternalPotentials >= 0, "Negative unregistered potentials, the pull code is inconsistent");
1624 static void check_external_potential_registration(const struct pull_t *pull)
1626 if (pull->numUnregisteredExternalPotentials > 0)
1628 int c;
1629 for (c = 0; c < pull->ncoord; c++)
1631 if (pull->coord[c].params.eType == epullEXTERNAL &&
1632 !pull->coord[c].bExternalPotentialProviderHasBeenRegistered)
1634 break;
1638 GMX_RELEASE_ASSERT(c < pull->ncoord, "Internal inconsistency in the pull potential provider counting");
1640 gmx_fatal(FARGS, "No external provider for external pull potentials have been provided for %d pull coordinates. The first coordinate without provider is number %d, which expects a module named '%s' to provide the external potential.",
1641 pull->numUnregisteredExternalPotentials,
1642 c + 1, pull->coord[c].params.externalPotentialProvider);
1647 /* Pull takes care of adding the forces of the external potential.
1648 * The external potential module has to make sure that the corresponding
1649 * potential energy is added either to the pull term or to a term
1650 * specific to the external module.
1652 void apply_external_pull_coord_force(struct pull_t *pull,
1653 int coord_index,
1654 double coord_force,
1655 const t_mdatoms *mdatoms,
1656 gmx::ForceWithVirial *forceWithVirial)
1658 pull_coord_work_t *pcrd;
1660 GMX_ASSERT(coord_index >= 0 && coord_index < pull->ncoord, "apply_external_pull_coord_force called with coord_index out of range");
1662 if (pull->comm.bParticipate)
1664 pcrd = &pull->coord[coord_index];
1666 GMX_RELEASE_ASSERT(pcrd->params.eType == epullEXTERNAL, "The pull force can only be set externally on pull coordinates of external type");
1668 GMX_ASSERT(pcrd->bExternalPotentialProviderHasBeenRegistered, "apply_external_pull_coord_force called for an unregistered pull coordinate");
1670 /* Set the force */
1671 pcrd->f_scal = coord_force;
1673 /* Calculate the forces on the pull groups */
1674 calc_pull_coord_vector_force(pcrd);
1676 /* Add the forces for this coordinate to the total virial and force */
1677 if (forceWithVirial->computeVirial_)
1679 matrix virial = { { 0 } };
1680 add_virial_coord(virial, pcrd);
1681 forceWithVirial->addVirialContribution(virial);
1684 apply_forces_coord(pull, coord_index, mdatoms,
1685 as_rvec_array(forceWithVirial->force_.data()));
1688 pull->numExternalPotentialsStillToBeAppliedThisStep--;
1691 /* Calculate the pull potential and scalar force for a pull coordinate */
1692 static void do_pull_pot_coord(struct pull_t *pull, int coord_ind, t_pbc *pbc,
1693 double t, real lambda,
1694 real *V, tensor vir, real *dVdl)
1696 pull_coord_work_t *pcrd;
1697 double dev;
1699 pcrd = &pull->coord[coord_ind];
1701 assert(pcrd->params.eType != epullCONSTRAINT);
1703 dev = get_pull_coord_deviation(pull, coord_ind, pbc, t);
1705 calc_pull_coord_scalar_force_and_potential(pcrd, dev, lambda, V, dVdl);
1707 calc_pull_coord_vector_force(pcrd);
1709 add_virial_coord(vir, pcrd);
1712 real pull_potential(struct pull_t *pull, t_mdatoms *md, t_pbc *pbc,
1713 t_commrec *cr, double t, real lambda,
1714 rvec *x, gmx::ForceWithVirial *force, real *dvdlambda)
1716 real V = 0;
1718 assert(pull != NULL);
1720 /* Ideally we should check external potential registration only during
1721 * the initialization phase, but that requires another function call
1722 * that should be done exactly in the right order. So we check here.
1724 check_external_potential_registration(pull);
1726 if (pull->comm.bParticipate)
1728 real dVdl = 0;
1730 pull_calc_coms(cr, pull, md, pbc, t, x, nullptr);
1732 rvec *f = as_rvec_array(force->force_.data());
1733 matrix virial = { { 0 } };
1734 const bool computeVirial = (force->computeVirial_ && MASTER(cr));
1735 for (int c = 0; c < pull->ncoord; c++)
1737 /* For external potential the force is assumed to be given by an external module by a call to
1738 apply_pull_coord_external_force */
1739 if (pull->coord[c].params.eType == epullCONSTRAINT || pull->coord[c].params.eType == epullEXTERNAL)
1741 continue;
1744 do_pull_pot_coord(pull, c, pbc, t, lambda,
1746 computeVirial ? virial : nullptr,
1747 &dVdl);
1749 /* Distribute the force over the atoms in the pulled groups */
1750 apply_forces_coord(pull, c, md, f);
1753 if (MASTER(cr))
1755 force->addVirialContribution(virial);
1756 *dvdlambda += dVdl;
1760 GMX_ASSERT(pull->numExternalPotentialsStillToBeAppliedThisStep == 0, "Too few or too many external pull potentials have been applied the previous step");
1761 /* All external pull potentials still need to be applied */
1762 pull->numExternalPotentialsStillToBeAppliedThisStep =
1763 pull->numCoordinatesWithExternalPotential;
1765 return (MASTER(cr) ? V : 0.0);
1768 void pull_constraint(struct pull_t *pull, t_mdatoms *md, t_pbc *pbc,
1769 t_commrec *cr, double dt, double t,
1770 rvec *x, rvec *xp, rvec *v, tensor vir)
1772 assert(pull != NULL);
1774 if (pull->comm.bParticipate)
1776 pull_calc_coms(cr, pull, md, pbc, t, x, xp);
1778 do_constraint(pull, pbc, xp, v, MASTER(cr), vir, dt, t);
1782 static void make_local_pull_group(gmx_ga2la_t *ga2la,
1783 pull_group_work_t *pg, int start, int end)
1785 int i, ii;
1787 pg->nat_loc = 0;
1788 for (i = 0; i < pg->params.nat; i++)
1790 ii = pg->params.ind[i];
1791 if (ga2la)
1793 if (!ga2la_get_home(ga2la, ii, &ii))
1795 ii = -1;
1798 if (ii >= start && ii < end)
1800 /* This is a home atom, add it to the local pull group */
1801 if (pg->nat_loc >= pg->nalloc_loc)
1803 pg->nalloc_loc = over_alloc_dd(pg->nat_loc+1);
1804 srenew(pg->ind_loc, pg->nalloc_loc);
1805 if (pg->epgrppbc == epgrppbcCOS || pg->params.weight != nullptr)
1807 srenew(pg->weight_loc, pg->nalloc_loc);
1810 pg->ind_loc[pg->nat_loc] = ii;
1811 if (pg->params.weight != nullptr)
1813 pg->weight_loc[pg->nat_loc] = pg->params.weight[i];
1815 pg->nat_loc++;
1820 void dd_make_local_pull_groups(t_commrec *cr, struct pull_t *pull, t_mdatoms *md)
1822 gmx_domdec_t *dd;
1823 pull_comm_t *comm;
1824 gmx_ga2la_t *ga2la;
1825 gmx_bool bMustParticipate;
1826 int g;
1828 dd = cr->dd;
1830 comm = &pull->comm;
1832 if (dd)
1834 ga2la = dd->ga2la;
1836 else
1838 ga2la = nullptr;
1841 /* We always make the master node participate, such that it can do i/o
1842 * and to simplify MC type extensions people might have.
1844 bMustParticipate = (comm->bParticipateAll || dd == nullptr || DDMASTER(dd));
1846 for (g = 0; g < pull->ngroup; g++)
1848 int a;
1850 make_local_pull_group(ga2la, &pull->group[g],
1851 0, md->homenr);
1853 /* We should participate if we have pull or pbc atoms */
1854 if (!bMustParticipate &&
1855 (pull->group[g].nat_loc > 0 ||
1856 (pull->group[g].params.pbcatom >= 0 &&
1857 ga2la_get_home(dd->ga2la, pull->group[g].params.pbcatom, &a))))
1859 bMustParticipate = TRUE;
1863 if (!comm->bParticipateAll)
1865 /* Keep currently not required ranks in the communicator
1866 * if they needed to participate up to 20 decompositions ago.
1867 * This avoids frequent rebuilds due to atoms jumping back and forth.
1869 const gmx_int64_t history_count = 20;
1870 gmx_bool bWillParticipate;
1871 int count[2];
1873 /* Increase the decomposition counter for the current call */
1874 comm->setup_count++;
1876 if (bMustParticipate)
1878 comm->must_count = comm->setup_count;
1881 bWillParticipate =
1882 bMustParticipate ||
1883 (comm->bParticipate &&
1884 comm->must_count >= comm->setup_count - history_count);
1886 if (debug && dd != nullptr)
1888 fprintf(debug, "Our DD rank (%3d) pull #atoms>0 or master: %d, will be part %d\n",
1889 dd->rank, bMustParticipate, bWillParticipate);
1892 if (bWillParticipate)
1894 /* Count the number of ranks that we want to have participating */
1895 count[0] = 1;
1896 /* Count the number of ranks that need to be added */
1897 count[1] = comm->bParticipate ? 0 : 1;
1899 else
1901 count[0] = 0;
1902 count[1] = 0;
1905 /* The cost of this global operation will be less that the cost
1906 * of the extra MPI_Comm_split calls that we can avoid.
1908 gmx_sumi(2, count, cr);
1910 /* If we are missing ranks or if we have 20% more ranks than needed
1911 * we make a new sub-communicator.
1913 if (count[1] > 0 || 6*count[0] < 5*comm->nparticipate)
1915 if (debug)
1917 fprintf(debug, "Creating new pull subcommunicator of size %d\n",
1918 count[0]);
1921 #if GMX_MPI
1922 if (comm->mpi_comm_com != MPI_COMM_NULL)
1924 MPI_Comm_free(&comm->mpi_comm_com);
1927 /* This might be an extremely expensive operation, so we try
1928 * to avoid this splitting as much as possible.
1930 assert(dd != NULL);
1931 MPI_Comm_split(dd->mpi_comm_all, bWillParticipate ? 0 : 1, dd->rank,
1932 &comm->mpi_comm_com);
1933 #endif
1935 comm->bParticipate = bWillParticipate;
1936 comm->nparticipate = count[0];
1940 /* Since the PBC of atoms might have changed, we need to update the PBC */
1941 pull->bSetPBCatoms = TRUE;
1944 static void init_pull_group_index(FILE *fplog, t_commrec *cr,
1945 int g, pull_group_work_t *pg,
1946 gmx_bool bConstraint, ivec pulldim_con,
1947 const gmx_mtop_t *mtop,
1948 const t_inputrec *ir, real lambda)
1950 if (EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD)
1952 /* There are no masses in the integrator.
1953 * But we still want to have the correct mass-weighted COMs.
1954 * So we store the real masses in the weights.
1955 * We do not set nweight, so these weights do not end up in the tpx file.
1957 if (pg->params.nweight == 0)
1959 snew(pg->params.weight, pg->params.nat);
1963 if (cr && PAR(cr))
1965 pg->nat_loc = 0;
1966 pg->nalloc_loc = 0;
1967 pg->ind_loc = nullptr;
1968 pg->weight_loc = nullptr;
1970 else
1972 pg->nat_loc = pg->params.nat;
1973 pg->ind_loc = pg->params.ind;
1974 if (pg->epgrppbc == epgrppbcCOS)
1976 snew(pg->weight_loc, pg->params.nat);
1978 else
1980 pg->weight_loc = pg->params.weight;
1984 const gmx_groups_t *groups = &mtop->groups;
1986 /* Count frozen dimensions and (weighted) mass */
1987 int nfrozen = 0;
1988 double tmass = 0;
1989 double wmass = 0;
1990 double wwmass = 0;
1991 int molb = 0;
1992 for (int i = 0; i < pg->params.nat; i++)
1994 int ii = pg->params.ind[i];
1995 if (bConstraint && ir->opts.nFreeze)
1997 for (int d = 0; d < DIM; d++)
1999 if (pulldim_con[d] == 1 &&
2000 ir->opts.nFreeze[ggrpnr(groups, egcFREEZE, ii)][d])
2002 nfrozen++;
2006 const t_atom &atom = mtopGetAtomParameters(mtop, ii, &molb);
2007 real m;
2008 if (ir->efep == efepNO)
2010 m = atom.m;
2012 else
2014 m = (1 - lambda)*atom.m + lambda*atom.mB;
2016 real w;
2017 if (pg->params.nweight > 0)
2019 w = pg->params.weight[i];
2021 else
2023 w = 1;
2025 if (EI_ENERGY_MINIMIZATION(ir->eI))
2027 /* Move the mass to the weight */
2028 w *= m;
2029 m = 1;
2030 pg->params.weight[i] = w;
2032 else if (ir->eI == eiBD)
2034 real mbd;
2035 if (ir->bd_fric)
2037 mbd = ir->bd_fric*ir->delta_t;
2039 else
2041 if (groups->grpnr[egcTC] == nullptr)
2043 mbd = ir->delta_t/ir->opts.tau_t[0];
2045 else
2047 mbd = ir->delta_t/ir->opts.tau_t[groups->grpnr[egcTC][ii]];
2050 w *= m/mbd;
2051 m = mbd;
2052 pg->params.weight[i] = w;
2054 tmass += m;
2055 wmass += m*w;
2056 wwmass += m*w*w;
2059 if (wmass == 0)
2061 /* We can have single atom groups with zero mass with potential pulling
2062 * without cosine weighting.
2064 if (pg->params.nat == 1 && !bConstraint && pg->epgrppbc != epgrppbcCOS)
2066 /* With one atom the mass doesn't matter */
2067 wwmass = 1;
2069 else
2071 gmx_fatal(FARGS, "The total%s mass of pull group %d is zero",
2072 pg->params.weight ? " weighted" : "", g);
2075 if (fplog)
2077 fprintf(fplog,
2078 "Pull group %d: %5d atoms, mass %9.3f",
2079 g, pg->params.nat, tmass);
2080 if (pg->params.weight ||
2081 EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD)
2083 fprintf(fplog, ", weighted mass %9.3f", wmass*wmass/wwmass);
2085 if (pg->epgrppbc == epgrppbcCOS)
2087 fprintf(fplog, ", cosine weighting will be used");
2089 fprintf(fplog, "\n");
2092 if (nfrozen == 0)
2094 /* A value != 0 signals not frozen, it is updated later */
2095 pg->invtm = -1.0;
2097 else
2099 int ndim = 0;
2100 for (int d = 0; d < DIM; d++)
2102 ndim += pulldim_con[d]*pg->params.nat;
2104 if (fplog && nfrozen > 0 && nfrozen < ndim)
2106 fprintf(fplog,
2107 "\nWARNING: In pull group %d some, but not all of the degrees of freedom\n"
2108 " that are subject to pulling are frozen.\n"
2109 " For constraint pulling the whole group will be frozen.\n\n",
2112 pg->invtm = 0.0;
2113 pg->wscale = 1.0;
2117 struct pull_t *
2118 init_pull(FILE *fplog, const pull_params_t *pull_params, const t_inputrec *ir,
2119 int nfile, const t_filenm fnm[],
2120 const gmx_mtop_t *mtop, t_commrec *cr,
2121 const gmx_output_env_t *oenv, real lambda,
2122 gmx_bool bOutFile,
2123 const ContinuationOptions &continuationOptions)
2125 struct pull_t *pull;
2126 pull_comm_t *comm;
2127 int g, c, m;
2129 snew(pull, 1);
2131 /* Copy the pull parameters */
2132 pull->params = *pull_params;
2133 /* Avoid pointer copies */
2134 pull->params.group = nullptr;
2135 pull->params.coord = nullptr;
2137 pull->ncoord = pull_params->ncoord;
2138 pull->ngroup = pull_params->ngroup;
2139 snew(pull->coord, pull->ncoord);
2140 snew(pull->group, pull->ngroup);
2142 pull->bPotential = FALSE;
2143 pull->bConstraint = FALSE;
2144 pull->bCylinder = FALSE;
2145 pull->bAngle = FALSE;
2147 for (g = 0; g < pull->ngroup; g++)
2149 pull_group_work_t *pgrp;
2150 int i;
2152 pgrp = &pull->group[g];
2154 /* Copy the pull group parameters */
2155 pgrp->params = pull_params->group[g];
2157 /* Avoid pointer copies by allocating and copying arrays */
2158 snew(pgrp->params.ind, pgrp->params.nat);
2159 for (i = 0; i < pgrp->params.nat; i++)
2161 pgrp->params.ind[i] = pull_params->group[g].ind[i];
2163 if (pgrp->params.nweight > 0)
2165 snew(pgrp->params.weight, pgrp->params.nweight);
2166 for (i = 0; i < pgrp->params.nweight; i++)
2168 pgrp->params.weight[i] = pull_params->group[g].weight[i];
2173 pull->numCoordinatesWithExternalPotential = 0;
2175 for (c = 0; c < pull->ncoord; c++)
2177 pull_coord_work_t *pcrd;
2178 int calc_com_start, calc_com_end, g;
2180 pcrd = &pull->coord[c];
2182 /* Copy all pull coordinate parameters */
2183 pcrd->params = pull_params->coord[c];
2185 switch (pcrd->params.eGeom)
2187 case epullgDIST:
2188 case epullgDIRRELATIVE: /* Direction vector is determined at each step */
2189 case epullgANGLE:
2190 case epullgDIHEDRAL:
2191 break;
2192 case epullgDIR:
2193 case epullgDIRPBC:
2194 case epullgCYL:
2195 case epullgANGLEAXIS:
2196 copy_rvec_to_dvec(pull_params->coord[c].vec, pcrd->vec);
2197 break;
2198 default:
2199 /* We allow reading of newer tpx files with new pull geometries,
2200 * but with the same tpx format, with old code. A new geometry
2201 * only adds a new enum value, which will be out of range for
2202 * old code. The first place we need to generate an error is
2203 * here, since the pull code can't handle this.
2204 * The enum can be used before arriving here only for printing
2205 * the string corresponding to the geometry, which will result
2206 * in printing "UNDEFINED".
2208 gmx_fatal(FARGS, "Pull geometry not supported for pull coordinate %d. The geometry enum %d in the input is larger than that supported by the code (up to %d). You are probably reading a tpr file generated with a newer version of Gromacs with an binary from an older version of Gromacs.",
2209 c + 1, pcrd->params.eGeom, epullgNR - 1);
2212 if (pcrd->params.eType == epullCONSTRAINT)
2214 /* Check restrictions of the constraint pull code */
2215 if (pcrd->params.eGeom == epullgCYL ||
2216 pcrd->params.eGeom == epullgDIRRELATIVE ||
2217 pcrd->params.eGeom == epullgANGLE ||
2218 pcrd->params.eGeom == epullgDIHEDRAL ||
2219 pcrd->params.eGeom == epullgANGLEAXIS)
2221 gmx_fatal(FARGS, "Pulling of type %s can not be combined with geometry %s. Consider using pull type %s.",
2222 epull_names[pcrd->params.eType],
2223 epullg_names[pcrd->params.eGeom],
2224 epull_names[epullUMBRELLA]);
2227 pull->bConstraint = TRUE;
2229 else
2231 pull->bPotential = TRUE;
2234 if (pcrd->params.eGeom == epullgCYL)
2236 pull->bCylinder = TRUE;
2238 else if (pcrd->params.eGeom == epullgANGLE || pcrd->params.eGeom == epullgDIHEDRAL || pcrd->params.eGeom == epullgANGLEAXIS)
2240 pull->bAngle = TRUE;
2243 /* We only need to calculate the plain COM of a group
2244 * when it is not only used as a cylinder group.
2246 calc_com_start = (pcrd->params.eGeom == epullgCYL ? 1 : 0);
2247 calc_com_end = pcrd->params.ngroup;
2249 for (g = calc_com_start; g <= calc_com_end; g++)
2251 pull->group[pcrd->params.group[g]].bCalcCOM = TRUE;
2254 /* With non-zero rate the reference value is set at every step */
2255 if (pcrd->params.rate == 0)
2257 /* Initialize the constant reference value */
2258 if (pcrd->params.eType != epullEXTERNAL)
2260 low_set_pull_coord_reference_value(pcrd, c, pcrd->params.init*pull_conversion_factor_userinput2internal(&pcrd->params));
2262 else
2264 /* With an external pull potential, the reference value loses
2265 * it's meaning and should not be used. Setting it to zero
2266 * makes any terms dependent on it disappear.
2267 * The only issue this causes is that with cylinder pulling
2268 * no atoms of the cylinder group within the cylinder radius
2269 * should be more than half a box length away from the COM of
2270 * the pull group along the axial direction.
2272 pcrd->value_ref = 0.0;
2276 if (pcrd->params.eType == epullEXTERNAL)
2278 GMX_RELEASE_ASSERT(pcrd->params.rate == 0, "With an external potential, a pull coordinate should have rate = 0");
2280 /* This external potential needs to be registered later */
2281 pull->numCoordinatesWithExternalPotential++;
2283 pcrd->bExternalPotentialProviderHasBeenRegistered = false;
2286 pull->numUnregisteredExternalPotentials =
2287 pull->numCoordinatesWithExternalPotential;
2288 pull->numExternalPotentialsStillToBeAppliedThisStep = 0;
2290 pull->ePBC = ir->ePBC;
2291 switch (pull->ePBC)
2293 case epbcNONE: pull->npbcdim = 0; break;
2294 case epbcXY: pull->npbcdim = 2; break;
2295 default: pull->npbcdim = 3; break;
2298 if (fplog)
2300 gmx_bool bAbs, bCos;
2302 bAbs = FALSE;
2303 for (c = 0; c < pull->ncoord; c++)
2305 if (pull->group[pull->coord[c].params.group[0]].params.nat == 0 ||
2306 pull->group[pull->coord[c].params.group[1]].params.nat == 0)
2308 bAbs = TRUE;
2312 fprintf(fplog, "\n");
2313 if (pull->bPotential)
2315 fprintf(fplog, "Will apply potential COM pulling\n");
2317 if (pull->bConstraint)
2319 fprintf(fplog, "Will apply constraint COM pulling\n");
2321 fprintf(fplog, "with %d pull coordinate%s and %d group%s\n",
2322 pull->ncoord, pull->ncoord == 1 ? "" : "s",
2323 pull->ngroup, pull->ngroup == 1 ? "" : "s");
2324 if (bAbs)
2326 fprintf(fplog, "with an absolute reference\n");
2328 bCos = FALSE;
2329 for (g = 0; g < pull->ngroup; g++)
2331 if (pull->group[g].params.nat > 1 &&
2332 pull->group[g].params.pbcatom < 0)
2334 /* We are using cosine weighting */
2335 fprintf(fplog, "Cosine weighting is used for group %d\n", g);
2336 bCos = TRUE;
2339 if (bCos)
2341 please_cite(fplog, "Engin2010");
2345 pull->bRefAt = FALSE;
2346 pull->cosdim = -1;
2347 for (g = 0; g < pull->ngroup; g++)
2349 pull_group_work_t *pgrp;
2351 pgrp = &pull->group[g];
2352 pgrp->epgrppbc = epgrppbcNONE;
2353 if (pgrp->params.nat > 0)
2355 /* There is an issue when a group is used in multiple coordinates
2356 * and constraints are applied in different dimensions with atoms
2357 * frozen in some, but not all dimensions.
2358 * Since there is only one mass array per group, we can't have
2359 * frozen/non-frozen atoms for different coords at the same time.
2360 * But since this is a very exotic case, we don't check for this.
2361 * A warning is printed in init_pull_group_index.
2364 gmx_bool bConstraint;
2365 ivec pulldim, pulldim_con;
2367 /* Loop over all pull coordinates to see along which dimensions
2368 * this group is pulled and if it is involved in constraints.
2370 bConstraint = FALSE;
2371 clear_ivec(pulldim);
2372 clear_ivec(pulldim_con);
2373 for (c = 0; c < pull->ncoord; c++)
2375 pull_coord_work_t *pcrd;
2376 int gi;
2377 gmx_bool bGroupUsed;
2379 pcrd = &pull->coord[c];
2381 bGroupUsed = FALSE;
2382 for (gi = 0; gi < pcrd->params.ngroup; gi++)
2384 if (pcrd->params.group[gi] == g)
2386 bGroupUsed = TRUE;
2390 if (bGroupUsed)
2392 for (m = 0; m < DIM; m++)
2394 if (pcrd->params.dim[m] == 1)
2396 pulldim[m] = 1;
2398 if (pcrd->params.eType == epullCONSTRAINT)
2400 bConstraint = TRUE;
2401 pulldim_con[m] = 1;
2408 /* Determine if we need to take PBC into account for calculating
2409 * the COM's of the pull groups.
2411 GMX_RELEASE_ASSERT(pull->npbcdim <= DIM, "npbcdim must be <= the number of dimensions");
2412 for (m = 0; m < pull->npbcdim; m++)
2414 GMX_RELEASE_ASSERT(m <= DIM, "npbcdim must be <= the number of dimensions");
2415 if (pulldim[m] == 1 && pgrp->params.nat > 1)
2417 if (pgrp->params.pbcatom >= 0)
2419 pgrp->epgrppbc = epgrppbcREFAT;
2420 pull->bRefAt = TRUE;
2422 else
2424 if (pgrp->params.weight != nullptr)
2426 gmx_fatal(FARGS, "Pull groups can not have relative weights and cosine weighting at same time");
2428 pgrp->epgrppbc = epgrppbcCOS;
2429 if (pull->cosdim >= 0 && pull->cosdim != m)
2431 gmx_fatal(FARGS, "Can only use cosine weighting with pulling in one dimension (use mdp option pull-coord?-dim)");
2433 pull->cosdim = m;
2437 /* Set the indices */
2438 init_pull_group_index(fplog, cr, g, pgrp,
2439 bConstraint, pulldim_con,
2440 mtop, ir, lambda);
2442 else
2444 /* Absolute reference, set the inverse mass to zero.
2445 * This is only relevant (and used) with constraint pulling.
2447 pgrp->invtm = 0;
2448 pgrp->wscale = 1;
2452 /* If we use cylinder coordinates, do some initialising for them */
2453 if (pull->bCylinder)
2455 snew(pull->dyna, pull->ncoord);
2457 for (c = 0; c < pull->ncoord; c++)
2459 const pull_coord_work_t *pcrd;
2461 pcrd = &pull->coord[c];
2463 if (pcrd->params.eGeom == epullgCYL)
2465 if (pull->group[pcrd->params.group[0]].params.nat == 0)
2467 gmx_fatal(FARGS, "A cylinder pull group is not supported when using absolute reference!\n");
2473 /* The gmx_omp_nthreads module might not be initialized here, so max(1,) */
2474 pull->nthreads = std::max(1, gmx_omp_nthreads_get(emntDefault));
2475 snew(pull->sum_com, pull->nthreads);
2477 comm = &pull->comm;
2479 #if GMX_MPI
2480 /* Use a sub-communicator when we have more than 32 ranks */
2481 comm->bParticipateAll = (cr == nullptr || !DOMAINDECOMP(cr) ||
2482 cr->dd->nnodes <= 32 ||
2483 getenv("GMX_PULL_PARTICIPATE_ALL") != nullptr);
2484 /* This sub-commicator is not used with comm->bParticipateAll,
2485 * so we can always initialize it to NULL.
2487 comm->mpi_comm_com = MPI_COMM_NULL;
2488 comm->nparticipate = 0;
2489 #else
2490 /* No MPI: 1 rank: all ranks pull */
2491 comm->bParticipateAll = TRUE;
2492 #endif
2493 comm->bParticipate = comm->bParticipateAll;
2494 comm->setup_count = 0;
2495 comm->must_count = 0;
2497 if (!comm->bParticipateAll && fplog != nullptr)
2499 fprintf(fplog, "Will use a sub-communicator for pull communication\n");
2502 comm->rbuf = nullptr;
2503 comm->dbuf = nullptr;
2504 comm->dbuf_cyl = nullptr;
2506 /* We still need to initialize the PBC reference coordinates */
2507 pull->bSetPBCatoms = TRUE;
2509 /* Only do I/O when we are doing dynamics and if we are the MASTER */
2510 pull->out_x = nullptr;
2511 pull->out_f = nullptr;
2512 if (bOutFile)
2514 /* Check for px and pf filename collision, if we are writing
2515 both files */
2516 std::string px_filename, pf_filename;
2517 std::string px_appended, pf_appended;
2520 px_filename = std::string(opt2fn("-px", nfile, fnm));
2521 pf_filename = std::string(opt2fn("-pf", nfile, fnm));
2523 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2526 if ((pull->params.nstxout != 0) &&
2527 (pull->params.nstfout != 0) &&
2528 (px_filename == pf_filename))
2530 if (!opt2bSet("-px", nfile, fnm) && !opt2bSet("-pf", nfile, fnm))
2532 /* We are writing both pull files but neither set directly. */
2535 px_appended = append_before_extension(px_filename, "_pullx");
2536 pf_appended = append_before_extension(pf_filename, "_pullf");
2538 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2539 pull->out_x = open_pull_out(px_appended.c_str(), pull, oenv,
2540 TRUE, continuationOptions);
2541 pull->out_f = open_pull_out(pf_appended.c_str(), pull, oenv,
2542 FALSE, continuationOptions);
2543 return pull;
2545 else
2547 /* If one of -px and -pf is set but the filenames are identical: */
2548 gmx_fatal(FARGS, "Identical pull_x and pull_f output filenames %s",
2549 px_filename.c_str());
2552 if (pull->params.nstxout != 0)
2554 pull->out_x = open_pull_out(opt2fn("-px", nfile, fnm), pull, oenv,
2555 TRUE, continuationOptions);
2557 if (pull->params.nstfout != 0)
2559 pull->out_f = open_pull_out(opt2fn("-pf", nfile, fnm), pull, oenv,
2560 FALSE, continuationOptions);
2564 return pull;
2567 static void destroy_pull_group(pull_group_work_t *pgrp)
2569 if (pgrp->ind_loc != pgrp->params.ind)
2571 sfree(pgrp->ind_loc);
2573 if (pgrp->weight_loc != pgrp->params.weight)
2575 sfree(pgrp->weight_loc);
2577 sfree(pgrp->mdw);
2578 sfree(pgrp->dv);
2580 sfree(pgrp->params.ind);
2581 sfree(pgrp->params.weight);
2584 static void destroy_pull(struct pull_t *pull)
2586 int i;
2588 for (i = 0; i < pull->ngroup; i++)
2590 destroy_pull_group(&pull->group[i]);
2592 for (i = 0; i < pull->ncoord; i++)
2594 if (pull->coord[i].params.eGeom == epullgCYL)
2596 destroy_pull_group(&(pull->dyna[i]));
2599 sfree(pull->group);
2600 sfree(pull->dyna);
2601 sfree(pull->coord);
2603 #if GMX_MPI
2604 if (pull->comm.mpi_comm_com != MPI_COMM_NULL)
2606 MPI_Comm_free(&pull->comm.mpi_comm_com);
2608 #endif
2609 sfree(pull->comm.rbuf);
2610 sfree(pull->comm.dbuf);
2611 sfree(pull->comm.dbuf_cyl);
2613 sfree(pull);
2616 void finish_pull(struct pull_t *pull)
2618 check_external_potential_registration(pull);
2620 if (pull->out_x)
2622 gmx_fio_fclose(pull->out_x);
2624 if (pull->out_f)
2626 gmx_fio_fclose(pull->out_f);
2629 destroy_pull(pull);
2632 gmx_bool pull_have_potential(const struct pull_t *pull)
2634 return pull->bPotential;
2637 gmx_bool pull_have_constraint(const struct pull_t *pull)
2639 return pull->bConstraint;