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,2018,2019, 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.
51 #include "gromacs/commandline/filenm.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/domdec/localatomset.h"
54 #include "gromacs/domdec/localatomsetmanager.h"
55 #include "gromacs/fileio/gmxfio.h"
56 #include "gromacs/gmxlib/network.h"
57 #include "gromacs/math/functions.h"
58 #include "gromacs/math/units.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/math/vectypes.h"
61 #include "gromacs/mdlib/gmx_omp_nthreads.h"
62 #include "gromacs/mdtypes/commrec.h"
63 #include "gromacs/mdtypes/forceoutput.h"
64 #include "gromacs/mdtypes/inputrec.h"
65 #include "gromacs/mdtypes/md_enums.h"
66 #include "gromacs/mdtypes/mdatom.h"
67 #include "gromacs/mdtypes/state.h"
68 #include "gromacs/pbcutil/pbc.h"
69 #include "gromacs/topology/mtop_lookup.h"
70 #include "gromacs/topology/topology.h"
71 #include "gromacs/utility/cstringutil.h"
72 #include "gromacs/utility/exceptions.h"
73 #include "gromacs/utility/fatalerror.h"
74 #include "gromacs/utility/futil.h"
75 #include "gromacs/utility/gmxassert.h"
76 #include "gromacs/utility/mutex.h"
77 #include "gromacs/utility/pleasecite.h"
78 #include "gromacs/utility/real.h"
79 #include "gromacs/utility/smalloc.h"
80 #include "gromacs/utility/strconvert.h"
82 #include "pull_internal.h"
84 static int groupPbcFromParams(const t_pull_group
¶ms
, bool setPbcRefToPrevStepCOM
)
91 else if (params
.pbcatom
>= 0)
93 if (setPbcRefToPrevStepCOM
)
95 return epgrppbcPREVSTEPCOM
;
108 /* NOTE: The params initialization currently copies pointers.
109 * So the lifetime of the source, currently always inputrec,
110 * should not end before that of this object.
111 * This will be fixed when the pointers are replacted by std::vector.
113 pull_group_work_t::pull_group_work_t(const t_pull_group
¶ms
,
114 gmx::LocalAtomSet atomSet
,
115 bool bSetPbcRefToPrevStepCOM
) :
117 epgrppbc(groupPbcFromParams(params
, bSetPbcRefToPrevStepCOM
)),
118 needToCalcCom(false),
128 bool pull_coordinate_is_angletype(const t_pull_coord
*pcrd
)
130 return (pcrd
->eGeom
== epullgANGLE
||
131 pcrd
->eGeom
== epullgDIHEDRAL
||
132 pcrd
->eGeom
== epullgANGLEAXIS
);
135 static bool pull_coordinate_is_directional(const t_pull_coord
*pcrd
)
137 return (pcrd
->eGeom
== epullgDIR
||
138 pcrd
->eGeom
== epullgDIRPBC
||
139 pcrd
->eGeom
== epullgDIRRELATIVE
||
140 pcrd
->eGeom
== epullgCYL
);
143 const char *pull_coordinate_units(const t_pull_coord
*pcrd
)
145 return pull_coordinate_is_angletype(pcrd
) ? "deg" : "nm";
148 double pull_conversion_factor_userinput2internal(const t_pull_coord
*pcrd
)
150 if (pull_coordinate_is_angletype(pcrd
))
160 double pull_conversion_factor_internal2userinput(const t_pull_coord
*pcrd
)
162 if (pull_coordinate_is_angletype(pcrd
))
172 /* Apply forces in a mass weighted fashion for part of the pull group */
173 static void apply_forces_grp_part(const pull_group_work_t
*pgrp
,
174 int ind_start
, int ind_end
,
176 const dvec f_pull
, int sign
, rvec
*f
)
178 double inv_wm
= pgrp
->mwscale
;
180 auto localAtomIndices
= pgrp
->atomSet
.localIndex();
181 for (int i
= ind_start
; i
< ind_end
; i
++)
183 int ii
= localAtomIndices
[i
];
184 double wmass
= md
->massT
[ii
];
185 if (!pgrp
->localWeights
.empty())
187 wmass
*= pgrp
->localWeights
[i
];
190 for (int d
= 0; d
< DIM
; d
++)
192 f
[ii
][d
] += sign
* wmass
* f_pull
[d
] * inv_wm
;
197 /* Apply forces in a mass weighted fashion */
198 static void apply_forces_grp(const pull_group_work_t
*pgrp
,
200 const dvec f_pull
, int sign
, rvec
*f
,
203 auto localAtomIndices
= pgrp
->atomSet
.localIndex();
205 if (pgrp
->params
.nat
== 1 && pgrp
->atomSet
.numAtomsLocal() == 1)
207 /* Only one atom and our rank has this atom: we can skip
208 * the mass weighting, which means that this code also works
209 * for mass=0, e.g. with a virtual site.
211 for (int d
= 0; d
< DIM
; d
++)
213 f
[localAtomIndices
[0]][d
] += sign
*f_pull
[d
];
218 if (localAtomIndices
.size() <= c_pullMaxNumLocalAtomsSingleThreaded
)
220 apply_forces_grp_part(pgrp
, 0, localAtomIndices
.size(), md
, f_pull
, sign
, f
);
224 #pragma omp parallel for num_threads(nthreads) schedule(static)
225 for (int th
= 0; th
< nthreads
; th
++)
227 int ind_start
= (localAtomIndices
.size()*(th
+ 0))/nthreads
;
228 int ind_end
= (localAtomIndices
.size()*(th
+ 1))/nthreads
;
229 apply_forces_grp_part(pgrp
, ind_start
, ind_end
,
230 md
, f_pull
, sign
, f
);
236 /* Apply forces in a mass weighted fashion to a cylinder group */
237 static void apply_forces_cyl_grp(const pull_group_work_t
*pgrp
,
238 const double dv_corr
,
240 const dvec f_pull
, double f_scal
,
242 int gmx_unused nthreads
)
244 double inv_wm
= pgrp
->mwscale
;
246 auto localAtomIndices
= pgrp
->atomSet
.localIndex();
248 /* The cylinder group is always a slab in the system, thus large.
249 * Therefore we always thread-parallelize this group.
251 int numAtomsLocal
= localAtomIndices
.size();
252 #pragma omp parallel for num_threads(nthreads) schedule(static)
253 for (int i
= 0; i
< numAtomsLocal
; i
++)
255 double weight
= pgrp
->localWeights
[i
];
260 int ii
= localAtomIndices
[i
];
261 double mass
= md
->massT
[ii
];
262 /* The stored axial distance from the cylinder center (dv) needs
263 * to be corrected for an offset (dv_corr), which was unknown when
266 double dv_com
= pgrp
->dv
[i
] + dv_corr
;
268 /* Here we not only add the pull force working along vec (f_pull),
269 * but also a radial component, due to the dependence of the weights
270 * on the radial distance.
272 for (int m
= 0; m
< DIM
; m
++)
274 f
[ii
][m
] += sign
*inv_wm
*(mass
*weight
*f_pull
[m
] +
275 pgrp
->mdw
[i
][m
]*dv_com
*f_scal
);
280 /* Apply torque forces in a mass weighted fashion to the groups that define
281 * the pull vector direction for pull coordinate pcrd.
283 static void apply_forces_vec_torque(const struct pull_t
*pull
,
284 const pull_coord_work_t
*pcrd
,
288 const PullCoordSpatialData
&spatialData
= pcrd
->spatialData
;
290 /* The component inpr along the pull vector is accounted for in the usual
291 * way. Here we account for the component perpendicular to vec.
294 for (int m
= 0; m
< DIM
; m
++)
296 inpr
+= spatialData
.dr01
[m
]*spatialData
.vec
[m
];
298 /* The torque force works along the component of the distance vector
299 * of between the two "usual" pull groups that is perpendicular to
300 * the pull vector. The magnitude of this force is the "usual" scale force
301 * multiplied with the ratio of the distance between the two "usual" pull
302 * groups and the distance between the two groups that define the vector.
305 for (int m
= 0; m
< DIM
; m
++)
307 f_perp
[m
] = (spatialData
.dr01
[m
] - inpr
*spatialData
.vec
[m
])/spatialData
.vec_len
*pcrd
->scalarForce
;
310 /* Apply the force to the groups defining the vector using opposite signs */
311 apply_forces_grp(&pull
->group
[pcrd
->params
.group
[2]], md
,
312 f_perp
, -1, f
, pull
->nthreads
);
313 apply_forces_grp(&pull
->group
[pcrd
->params
.group
[3]], md
,
314 f_perp
, 1, f
, pull
->nthreads
);
317 /* Apply forces in a mass weighted fashion */
318 static void apply_forces_coord(struct pull_t
* pull
, int coord
,
319 const PullCoordVectorForces
&forces
,
320 const t_mdatoms
* md
,
323 /* Here it would be more efficient to use one large thread-parallel
324 * region instead of potential parallel regions within apply_forces_grp.
325 * But there could be overlap between pull groups and this would lead
329 const pull_coord_work_t
&pcrd
= pull
->coord
[coord
];
331 if (pcrd
.params
.eGeom
== epullgCYL
)
333 apply_forces_cyl_grp(&pull
->dyna
[coord
], pcrd
.spatialData
.cyl_dev
, md
,
334 forces
.force01
, pcrd
.scalarForce
, -1, f
,
337 /* Sum the force along the vector and the radial force */
339 for (int m
= 0; m
< DIM
; m
++)
341 f_tot
[m
] = forces
.force01
[m
] + pcrd
.scalarForce
*pcrd
.spatialData
.ffrad
[m
];
343 apply_forces_grp(&pull
->group
[pcrd
.params
.group
[1]], md
,
344 f_tot
, 1, f
, pull
->nthreads
);
348 if (pcrd
.params
.eGeom
== epullgDIRRELATIVE
)
350 /* We need to apply the torque forces to the pull groups
351 * that define the pull vector.
353 apply_forces_vec_torque(pull
, &pcrd
, md
, f
);
356 if (pull
->group
[pcrd
.params
.group
[0]].params
.nat
> 0)
358 apply_forces_grp(&pull
->group
[pcrd
.params
.group
[0]], md
,
359 forces
.force01
, -1, f
, pull
->nthreads
);
361 apply_forces_grp(&pull
->group
[pcrd
.params
.group
[1]], md
,
362 forces
.force01
, 1, f
, pull
->nthreads
);
364 if (pcrd
.params
.ngroup
>= 4)
366 apply_forces_grp(&pull
->group
[pcrd
.params
.group
[2]], md
,
367 forces
.force23
, -1, f
, pull
->nthreads
);
368 apply_forces_grp(&pull
->group
[pcrd
.params
.group
[3]], md
,
369 forces
.force23
, 1, f
, pull
->nthreads
);
371 if (pcrd
.params
.ngroup
>= 6)
373 apply_forces_grp(&pull
->group
[pcrd
.params
.group
[4]], md
,
374 forces
.force45
, -1, f
, pull
->nthreads
);
375 apply_forces_grp(&pull
->group
[pcrd
.params
.group
[5]], md
,
376 forces
.force45
, 1, f
, pull
->nthreads
);
381 real
max_pull_distance2(const pull_coord_work_t
*pcrd
,
384 /* Note that this maximum distance calculation is more complex than
385 * most other cases in GROMACS, since here we have to take care of
386 * distance calculations that don't involve all three dimensions.
387 * For example, we can use distances that are larger than the
388 * box X and Y dimensions for a box that is elongated along Z.
391 real max_d2
= GMX_REAL_MAX
;
393 if (pull_coordinate_is_directional(&pcrd
->params
))
395 /* Directional pulling along along direction pcrd->vec.
396 * Calculating the exact maximum distance is complex and bug-prone.
397 * So we take a safe approach by not allowing distances that
398 * are larger than half the distance between unit cell faces
399 * along dimensions involved in pcrd->vec.
401 for (int m
= 0; m
< DIM
; m
++)
403 if (m
< pbc
->ndim_ePBC
&& pcrd
->spatialData
.vec
[m
] != 0)
405 real imageDistance2
= gmx::square(pbc
->box
[m
][m
]);
406 for (int d
= m
+ 1; d
< DIM
; d
++)
408 imageDistance2
-= gmx::square(pbc
->box
[d
][m
]);
410 max_d2
= std::min(max_d2
, imageDistance2
);
416 /* Distance pulling along dimensions with pcrd->params.dim[d]==1.
417 * We use half the minimum box vector length of the dimensions involved.
418 * This is correct for all cases, except for corner cases with
419 * triclinic boxes where e.g. dim[XX]=1 and dim[YY]=0 with
420 * box[YY][XX]!=0 and box[YY][YY] < box[XX][XX]. But since even
421 * in such corner cases the user could get correct results,
422 * depending on the details of the setup, we avoid further
423 * code complications.
425 for (int m
= 0; m
< DIM
; m
++)
427 if (m
< pbc
->ndim_ePBC
&& pcrd
->params
.dim
[m
] != 0)
429 real imageDistance2
= gmx::square(pbc
->box
[m
][m
]);
430 for (int d
= 0; d
< m
; d
++)
432 if (pcrd
->params
.dim
[d
] != 0)
434 imageDistance2
+= gmx::square(pbc
->box
[m
][d
]);
437 max_d2
= std::min(max_d2
, imageDistance2
);
445 /* This function returns the distance based on coordinates xg and xref.
446 * Note that the pull coordinate struct pcrd is not modified.
448 static void low_get_pull_coord_dr(const struct pull_t
*pull
,
449 const pull_coord_work_t
*pcrd
,
451 dvec xg
, dvec xref
, double max_dist2
,
454 const pull_group_work_t
*pgrp0
= &pull
->group
[pcrd
->params
.group
[0]];
456 /* Only the first group can be an absolute reference, in that case nat=0 */
457 if (pgrp0
->params
.nat
== 0)
459 for (int m
= 0; m
< DIM
; m
++)
461 xref
[m
] = pcrd
->params
.origin
[m
];
466 copy_dvec(xref
, xrefr
);
468 dvec dref
= {0, 0, 0};
469 if (pcrd
->params
.eGeom
== epullgDIRPBC
)
471 for (int m
= 0; m
< DIM
; m
++)
473 dref
[m
] = pcrd
->value_ref
*pcrd
->spatialData
.vec
[m
];
475 /* Add the reference position, so we use the correct periodic image */
476 dvec_inc(xrefr
, dref
);
479 pbc_dx_d(pbc
, xg
, xrefr
, dr
);
481 bool directional
= pull_coordinate_is_directional(&pcrd
->params
);
483 for (int m
= 0; m
< DIM
; m
++)
485 dr
[m
] *= pcrd
->params
.dim
[m
];
486 if (pcrd
->params
.dim
[m
] && !(directional
&& pcrd
->spatialData
.vec
[m
] == 0))
491 /* Check if we are close to switching to another periodic image */
492 if (max_dist2
> 0 && dr2
> 0.98*0.98*max_dist2
)
494 /* Note that technically there is no issue with switching periodic
495 * image, as pbc_dx_d returns the distance to the closest periodic
496 * image. However in all cases where periodic image switches occur,
497 * the pull results will be useless in practice.
499 gmx_fatal(FARGS
, "Distance between pull groups %d and %d (%f nm) is larger than 0.49 times the box size (%f).\n%s",
500 pcrd
->params
.group
[0], pcrd
->params
.group
[1],
501 sqrt(dr2
), sqrt(0.98*0.98*max_dist2
),
502 pcrd
->params
.eGeom
== epullgDIR
? "You might want to consider using \"pull-geometry = direction-periodic\" instead.\n" : "");
505 if (pcrd
->params
.eGeom
== epullgDIRPBC
)
511 /* This function returns the distance based on the contents of the pull struct.
512 * pull->coord[coord_ind].dr, and potentially vector, are updated.
514 static void get_pull_coord_dr(struct pull_t
*pull
,
518 pull_coord_work_t
*pcrd
= &pull
->coord
[coord_ind
];
519 PullCoordSpatialData
&spatialData
= pcrd
->spatialData
;
522 if (pcrd
->params
.eGeom
== epullgDIRPBC
)
528 md2
= static_cast<double>(max_pull_distance2(pcrd
, pbc
));
531 if (pcrd
->params
.eGeom
== epullgDIRRELATIVE
)
533 /* We need to determine the pull vector */
537 const pull_group_work_t
&pgrp2
= pull
->group
[pcrd
->params
.group
[2]];
538 const pull_group_work_t
&pgrp3
= pull
->group
[pcrd
->params
.group
[3]];
540 pbc_dx_d(pbc
, pgrp3
.x
, pgrp2
.x
, vec
);
542 for (m
= 0; m
< DIM
; m
++)
544 vec
[m
] *= pcrd
->params
.dim
[m
];
546 spatialData
.vec_len
= dnorm(vec
);
547 for (m
= 0; m
< DIM
; m
++)
549 spatialData
.vec
[m
] = vec
[m
]/spatialData
.vec_len
;
553 fprintf(debug
, "pull coord %d vector: %6.3f %6.3f %6.3f normalized: %6.3f %6.3f %6.3f\n",
555 vec
[XX
], vec
[YY
], vec
[ZZ
],
556 spatialData
.vec
[XX
], spatialData
.vec
[YY
], spatialData
.vec
[ZZ
]);
560 pull_group_work_t
*pgrp0
= &pull
->group
[pcrd
->params
.group
[0]];
561 pull_group_work_t
*pgrp1
= &pull
->group
[pcrd
->params
.group
[1]];
563 low_get_pull_coord_dr(pull
, pcrd
, pbc
,
565 pcrd
->params
.eGeom
== epullgCYL
? pull
->dyna
[coord_ind
].x
: pgrp0
->x
,
569 if (pcrd
->params
.ngroup
>= 4)
571 pull_group_work_t
*pgrp2
, *pgrp3
;
572 pgrp2
= &pull
->group
[pcrd
->params
.group
[2]];
573 pgrp3
= &pull
->group
[pcrd
->params
.group
[3]];
575 low_get_pull_coord_dr(pull
, pcrd
, pbc
,
581 if (pcrd
->params
.ngroup
>= 6)
583 pull_group_work_t
*pgrp4
, *pgrp5
;
584 pgrp4
= &pull
->group
[pcrd
->params
.group
[4]];
585 pgrp5
= &pull
->group
[pcrd
->params
.group
[5]];
587 low_get_pull_coord_dr(pull
, pcrd
, pbc
,
595 /* Modify x so that it is periodic in [-pi, pi)
596 * It is assumed that x is in [-3pi, 3pi) so that x
597 * needs to be shifted by at most one period.
599 static void make_periodic_2pi(double *x
)
611 /* This function should always be used to modify pcrd->value_ref */
612 static void low_set_pull_coord_reference_value(pull_coord_work_t
*pcrd
,
616 GMX_ASSERT(pcrd
->params
.eType
!= epullEXTERNAL
, "The pull coord reference value should not be used with type external-potential");
618 if (pcrd
->params
.eGeom
== epullgDIST
)
622 gmx_fatal(FARGS
, "Pull reference distance for coordinate %d (%f) needs to be non-negative", coord_ind
+ 1, value_ref
);
625 else if (pcrd
->params
.eGeom
== epullgANGLE
|| pcrd
->params
.eGeom
== epullgANGLEAXIS
)
627 if (value_ref
< 0 || value_ref
> M_PI
)
629 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
));
632 else if (pcrd
->params
.eGeom
== epullgDIHEDRAL
)
634 /* Allow pulling to be periodic for dihedral angles by remapping the reference value to the interval [-pi, pi). */
635 make_periodic_2pi(&value_ref
);
638 pcrd
->value_ref
= value_ref
;
641 static void update_pull_coord_reference_value(pull_coord_work_t
*pcrd
, int coord_ind
, double t
)
643 /* With zero rate the reference value is set initially and doesn't change */
644 if (pcrd
->params
.rate
!= 0)
646 double value_ref
= (pcrd
->params
.init
+ pcrd
->params
.rate
*t
)*pull_conversion_factor_userinput2internal(&pcrd
->params
);
647 low_set_pull_coord_reference_value(pcrd
, coord_ind
, value_ref
);
651 /* Returns the dihedral angle. Updates the plane normal vectors m, n. */
652 static double get_dihedral_angle_coord(PullCoordSpatialData
*spatialData
)
655 dvec dr32
; /* store instead of dr23? */
657 dsvmul(-1, spatialData
->dr23
, dr32
);
658 dcprod(spatialData
->dr01
, dr32
, spatialData
->planevec_m
); /* Normal of first plane */
659 dcprod(dr32
, spatialData
->dr45
, spatialData
->planevec_n
); /* Normal of second plane */
660 phi
= gmx_angle_between_dvecs(spatialData
->planevec_m
, spatialData
->planevec_n
);
662 /* Note 1: the sign below is opposite of that in the bondeds or Bekker 1994
663 * because there r_ij = ri - rj, while here dr01 = r_1 - r_0
664 * Note 2: the angle between the plane normal vectors equals pi only when
665 * both planes coincide. Thus, when phi = pi dr01 will lie in both planes and
666 * we get a positive sign below. Thus, the range of the dihedral angle is (-180, 180].
668 sign
= (diprod(spatialData
->dr01
, spatialData
->planevec_n
) < 0.0) ? 1.0 : -1.0;
672 /* Calculates pull->coord[coord_ind].value.
673 * This function also updates pull->coord[coord_ind].dr.
675 static void get_pull_coord_distance(struct pull_t
*pull
,
679 pull_coord_work_t
*pcrd
;
682 pcrd
= &pull
->coord
[coord_ind
];
684 get_pull_coord_dr(pull
, coord_ind
, pbc
);
686 PullCoordSpatialData
&spatialData
= pcrd
->spatialData
;
688 switch (pcrd
->params
.eGeom
)
691 /* Pull along the vector between the com's */
692 spatialData
.value
= dnorm(spatialData
.dr01
);
696 case epullgDIRRELATIVE
:
699 spatialData
.value
= 0;
700 for (m
= 0; m
< DIM
; m
++)
702 spatialData
.value
+= spatialData
.vec
[m
]*spatialData
.dr01
[m
];
706 spatialData
.value
= gmx_angle_between_dvecs(spatialData
.dr01
, spatialData
.dr23
);
709 spatialData
.value
= get_dihedral_angle_coord(&spatialData
);
711 case epullgANGLEAXIS
:
712 spatialData
.value
= gmx_angle_between_dvecs(spatialData
.dr01
, spatialData
.vec
);
715 gmx_incons("Unsupported pull type in get_pull_coord_distance");
719 /* Returns the deviation from the reference value.
720 * Updates pull->coord[coord_ind].dr, .value and .value_ref.
722 static double get_pull_coord_deviation(struct pull_t
*pull
,
727 pull_coord_work_t
*pcrd
;
730 pcrd
= &pull
->coord
[coord_ind
];
732 /* Update the reference value before computing the distance,
733 * since it is used in the distance computation with periodic pulling.
735 update_pull_coord_reference_value(pcrd
, coord_ind
, t
);
737 get_pull_coord_distance(pull
, coord_ind
, pbc
);
739 get_pull_coord_distance(pull
, coord_ind
, pbc
);
741 /* Determine the deviation */
742 dev
= pcrd
->spatialData
.value
- pcrd
->value_ref
;
744 /* Check that values are allowed */
745 if (pcrd
->params
.eGeom
== epullgDIST
&& pcrd
->spatialData
.value
== 0)
747 /* With no vector we can not determine the direction for the force,
748 * so we set the force to zero.
752 else if (pcrd
->params
.eGeom
== epullgDIHEDRAL
)
754 /* The reference value is in [-pi, pi). The coordinate value is in (-pi, pi].
755 Thus, the unwrapped deviation is here in (-2pi, 2pi].
756 After making it periodic, the deviation will be in [-pi, pi). */
757 make_periodic_2pi(&dev
);
763 double get_pull_coord_value(struct pull_t
*pull
,
767 get_pull_coord_distance(pull
, coord_ind
, pbc
);
769 return pull
->coord
[coord_ind
].spatialData
.value
;
772 void clear_pull_forces(pull_t
*pull
)
774 /* Zeroing the forces is only required for constraint pulling.
775 * It can happen that multiple constraint steps need to be applied
776 * and therefore the constraint forces need to be accumulated.
778 for (pull_coord_work_t
&coord
: pull
->coord
)
780 coord
.scalarForce
= 0;
784 /* Apply constraint using SHAKE */
785 static void do_constraint(struct pull_t
*pull
, t_pbc
*pbc
,
787 gmx_bool bMaster
, tensor vir
,
791 dvec
*r_ij
; /* x[i] com of i in prev. step. Obeys constr. -> r_ij[i] */
792 dvec unc_ij
; /* xp[i] com of i this step, before constr. -> unc_ij */
793 dvec
*rnew
; /* current 'new' positions of the groups */
794 double *dr_tot
; /* the total update of the coords */
797 double lambda
, rm
, invdt
= 0;
798 gmx_bool bConverged_all
, bConverged
= FALSE
;
799 int niter
= 0, ii
, j
, m
, max_iter
= 100;
803 snew(r_ij
, pull
->coord
.size());
804 snew(dr_tot
, pull
->coord
.size());
806 snew(rnew
, pull
->group
.size());
808 /* copy the current unconstrained positions for use in iterations. We
809 iterate until rinew[i] and rjnew[j] obey the constraints. Then
810 rinew - pull.x_unc[i] is the correction dr to group i */
811 for (size_t g
= 0; g
< pull
->group
.size(); g
++)
813 copy_dvec(pull
->group
[g
].xp
, rnew
[g
]);
816 /* Determine the constraint directions from the old positions */
817 for (size_t c
= 0; c
< pull
->coord
.size(); c
++)
819 pull_coord_work_t
*pcrd
;
821 pcrd
= &pull
->coord
[c
];
823 if (pcrd
->params
.eType
!= epullCONSTRAINT
)
828 /* Note that get_pull_coord_distance also sets pcrd->dr and pcrd->value.
829 * We don't modify dr and value anymore, so these values are also used
832 get_pull_coord_distance(pull
, c
, pbc
);
834 const PullCoordSpatialData
&spatialData
= pcrd
->spatialData
;
837 fprintf(debug
, "Pull coord %zu dr %f %f %f\n",
838 c
, spatialData
.dr01
[XX
], spatialData
.dr01
[YY
], spatialData
.dr01
[ZZ
]);
841 if (pcrd
->params
.eGeom
== epullgDIR
||
842 pcrd
->params
.eGeom
== epullgDIRPBC
)
844 /* Select the component along vec */
846 for (m
= 0; m
< DIM
; m
++)
848 a
+= spatialData
.vec
[m
]*spatialData
.dr01
[m
];
850 for (m
= 0; m
< DIM
; m
++)
852 r_ij
[c
][m
] = a
*spatialData
.vec
[m
];
857 copy_dvec(spatialData
.dr01
, r_ij
[c
]);
860 if (dnorm2(r_ij
[c
]) == 0)
862 gmx_fatal(FARGS
, "Distance for pull coordinate %zu is zero with constraint pulling, which is not allowed.", c
+ 1);
866 bConverged_all
= FALSE
;
867 while (!bConverged_all
&& niter
< max_iter
)
869 bConverged_all
= TRUE
;
871 /* loop over all constraints */
872 for (size_t c
= 0; c
< pull
->coord
.size(); c
++)
874 pull_coord_work_t
*pcrd
;
875 pull_group_work_t
*pgrp0
, *pgrp1
;
878 pcrd
= &pull
->coord
[c
];
880 if (pcrd
->params
.eType
!= epullCONSTRAINT
)
885 update_pull_coord_reference_value(pcrd
, c
, t
);
887 pgrp0
= &pull
->group
[pcrd
->params
.group
[0]];
888 pgrp1
= &pull
->group
[pcrd
->params
.group
[1]];
890 /* Get the current difference vector */
891 low_get_pull_coord_dr(pull
, pcrd
, pbc
,
892 rnew
[pcrd
->params
.group
[1]],
893 rnew
[pcrd
->params
.group
[0]],
898 fprintf(debug
, "Pull coord %zu, iteration %d\n", c
, niter
);
901 rm
= 1.0/(pgrp0
->invtm
+ pgrp1
->invtm
);
903 switch (pcrd
->params
.eGeom
)
906 if (pcrd
->value_ref
<= 0)
908 gmx_fatal(FARGS
, "The pull constraint reference distance for group %zu is <= 0 (%f)", c
, pcrd
->value_ref
);
912 double q
, c_a
, c_b
, c_c
;
914 c_a
= diprod(r_ij
[c
], r_ij
[c
]);
915 c_b
= diprod(unc_ij
, r_ij
[c
])*2;
916 c_c
= diprod(unc_ij
, unc_ij
) - gmx::square(pcrd
->value_ref
);
920 q
= -0.5*(c_b
- std::sqrt(c_b
*c_b
- 4*c_a
*c_c
));
925 q
= -0.5*(c_b
+ std::sqrt(c_b
*c_b
- 4*c_a
*c_c
));
932 "Pull ax^2+bx+c=0: a=%e b=%e c=%e lambda=%e\n",
933 c_a
, c_b
, c_c
, lambda
);
937 /* The position corrections dr due to the constraints */
938 dsvmul(-lambda
*rm
*pgrp1
->invtm
, r_ij
[c
], dr1
);
939 dsvmul( lambda
*rm
*pgrp0
->invtm
, r_ij
[c
], dr0
);
940 dr_tot
[c
] += -lambda
*dnorm(r_ij
[c
]);
945 /* A 1-dimensional constraint along a vector */
947 for (m
= 0; m
< DIM
; m
++)
949 vec
[m
] = pcrd
->spatialData
.vec
[m
];
950 a
+= unc_ij
[m
]*vec
[m
];
952 /* Select only the component along the vector */
953 dsvmul(a
, vec
, unc_ij
);
954 lambda
= a
- pcrd
->value_ref
;
957 fprintf(debug
, "Pull inpr %e lambda: %e\n", a
, lambda
);
960 /* The position corrections dr due to the constraints */
961 dsvmul(-lambda
*rm
*pgrp1
->invtm
, vec
, dr1
);
962 dsvmul( lambda
*rm
*pgrp0
->invtm
, vec
, dr0
);
963 dr_tot
[c
] += -lambda
;
966 gmx_incons("Invalid enumeration value for eGeom");
974 g0
= pcrd
->params
.group
[0];
975 g1
= pcrd
->params
.group
[1];
976 low_get_pull_coord_dr(pull
, pcrd
, pbc
, rnew
[g1
], rnew
[g0
], -1, tmp
);
977 low_get_pull_coord_dr(pull
, pcrd
, pbc
, dr1
, dr0
, -1, tmp3
);
979 "Pull cur %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
980 rnew
[g0
][0], rnew
[g0
][1], rnew
[g0
][2],
981 rnew
[g1
][0], rnew
[g1
][1], rnew
[g1
][2], dnorm(tmp
));
983 "Pull ref %8s %8s %8s %8s %8s %8s d: %8.5f\n",
984 "", "", "", "", "", "", pcrd
->value_ref
);
986 "Pull cor %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
987 dr0
[0], dr0
[1], dr0
[2],
988 dr1
[0], dr1
[1], dr1
[2],
992 /* Update the COMs with dr */
993 dvec_inc(rnew
[pcrd
->params
.group
[1]], dr1
);
994 dvec_inc(rnew
[pcrd
->params
.group
[0]], dr0
);
997 /* Check if all constraints are fullfilled now */
998 for (pull_coord_work_t
&coord
: pull
->coord
)
1000 if (coord
.params
.eType
!= epullCONSTRAINT
)
1005 low_get_pull_coord_dr(pull
, &coord
, pbc
,
1006 rnew
[coord
.params
.group
[1]],
1007 rnew
[coord
.params
.group
[0]],
1010 switch (coord
.params
.eGeom
)
1014 fabs(dnorm(unc_ij
) - coord
.value_ref
) < pull
->params
.constr_tol
;
1019 for (m
= 0; m
< DIM
; m
++)
1021 vec
[m
] = coord
.spatialData
.vec
[m
];
1023 inpr
= diprod(unc_ij
, vec
);
1024 dsvmul(inpr
, vec
, unc_ij
);
1026 fabs(diprod(unc_ij
, vec
) - coord
.value_ref
) < pull
->params
.constr_tol
;
1034 fprintf(debug
, "Pull constraint not converged: "
1036 "d_ref = %f, current d = %f\n",
1037 coord
.params
.group
[0], coord
.params
.group
[1],
1038 coord
.value_ref
, dnorm(unc_ij
));
1041 bConverged_all
= FALSE
;
1046 /* if after all constraints are dealt with and bConverged is still TRUE
1047 we're finished, if not we do another iteration */
1049 if (niter
> max_iter
)
1051 gmx_fatal(FARGS
, "Too many iterations for constraint run: %d", niter
);
1054 /* DONE ITERATING, NOW UPDATE COORDINATES AND CALC. CONSTRAINT FORCES */
1061 /* update atoms in the groups */
1062 for (size_t g
= 0; g
< pull
->group
.size(); g
++)
1064 const pull_group_work_t
*pgrp
;
1067 pgrp
= &pull
->group
[g
];
1069 /* get the final constraint displacement dr for group g */
1070 dvec_sub(rnew
[g
], pgrp
->xp
, dr
);
1072 if (dnorm2(dr
) == 0)
1074 /* No displacement, no update necessary */
1078 /* update the atom positions */
1079 auto localAtomIndices
= pgrp
->atomSet
.localIndex();
1081 for (gmx::index j
= 0; j
< localAtomIndices
.ssize(); j
++)
1083 ii
= localAtomIndices
[j
];
1084 if (!pgrp
->localWeights
.empty())
1086 dsvmul(pgrp
->wscale
*pgrp
->localWeights
[j
], dr
, tmp
);
1088 for (m
= 0; m
< DIM
; m
++)
1094 for (m
= 0; m
< DIM
; m
++)
1096 v
[ii
][m
] += invdt
*tmp
[m
];
1102 /* calculate the constraint forces, used for output and virial only */
1103 for (size_t c
= 0; c
< pull
->coord
.size(); c
++)
1105 pull_coord_work_t
*pcrd
;
1107 pcrd
= &pull
->coord
[c
];
1109 if (pcrd
->params
.eType
!= epullCONSTRAINT
)
1114 /* Accumulate the forces, in case we have multiple constraint steps */
1115 double force
= dr_tot
[c
]/((pull
->group
[pcrd
->params
.group
[0]].invtm
+ pull
->group
[pcrd
->params
.group
[1]].invtm
)*dt
*dt
);
1116 pcrd
->scalarForce
+= force
;
1118 if (vir
!= nullptr && pcrd
->params
.eGeom
!= epullgDIRPBC
&& bMaster
)
1122 /* Add the pull contribution to the virial */
1123 /* We have already checked above that r_ij[c] != 0 */
1124 f_invr
= pcrd
->scalarForce
/dnorm(r_ij
[c
]);
1126 for (j
= 0; j
< DIM
; j
++)
1128 for (m
= 0; m
< DIM
; m
++)
1130 vir
[j
][m
] -= 0.5*f_invr
*r_ij
[c
][j
]*r_ij
[c
][m
];
1136 /* finished! I hope. Give back some memory */
1142 static void add_virial_coord_dr(tensor vir
, const dvec dr
, const dvec f
)
1144 for (int j
= 0; j
< DIM
; j
++)
1146 for (int m
= 0; m
< DIM
; m
++)
1148 vir
[j
][m
] -= 0.5*f
[j
]*dr
[m
];
1153 /* Adds the pull contribution to the virial */
1154 static void add_virial_coord(tensor vir
,
1155 const pull_coord_work_t
&pcrd
,
1156 const PullCoordVectorForces
&forces
)
1158 if (vir
!= nullptr && pcrd
.params
.eGeom
!= epullgDIRPBC
)
1160 /* Add the pull contribution for each distance vector to the virial. */
1161 add_virial_coord_dr(vir
, pcrd
.spatialData
.dr01
, forces
.force01
);
1162 if (pcrd
.params
.ngroup
>= 4)
1164 add_virial_coord_dr(vir
, pcrd
.spatialData
.dr23
, forces
.force23
);
1166 if (pcrd
.params
.ngroup
>= 6)
1168 add_virial_coord_dr(vir
, pcrd
.spatialData
.dr45
, forces
.force45
);
1173 static void calc_pull_coord_scalar_force_and_potential(pull_coord_work_t
*pcrd
,
1174 double dev
, real lambda
,
1175 real
*V
, real
*dVdl
)
1179 k
= (1.0 - lambda
)*pcrd
->params
.k
+ lambda
*pcrd
->params
.kB
;
1180 dkdl
= pcrd
->params
.kB
- pcrd
->params
.k
;
1182 switch (pcrd
->params
.eType
)
1185 case epullFLATBOTTOM
:
1186 case epullFLATBOTTOMHIGH
:
1187 /* The only difference between an umbrella and a flat-bottom
1188 * potential is that a flat-bottom is zero above or below
1189 the reference value.
1191 if ((pcrd
->params
.eType
== epullFLATBOTTOM
&& dev
< 0) ||
1192 (pcrd
->params
.eType
== epullFLATBOTTOMHIGH
&& dev
> 0))
1197 pcrd
->scalarForce
= -k
*dev
;
1198 *V
+= 0.5* k
*gmx::square(dev
);
1199 *dVdl
+= 0.5*dkdl
*gmx::square(dev
);
1202 pcrd
->scalarForce
= -k
;
1203 *V
+= k
*pcrd
->spatialData
.value
;
1204 *dVdl
+= dkdl
*pcrd
->spatialData
.value
;
1207 gmx_incons("the scalar pull force should not be calculated internally for pull type external");
1209 gmx_incons("Unsupported pull type in do_pull_pot");
1213 static PullCoordVectorForces
1214 calculateVectorForces(const pull_coord_work_t
&pcrd
)
1216 const t_pull_coord
¶ms
= pcrd
.params
;
1217 const PullCoordSpatialData
&spatialData
= pcrd
.spatialData
;
1219 /* The geometry of the coordinate determines how the scalar force relates to the force on each group */
1220 PullCoordVectorForces forces
;
1222 if (params
.eGeom
== epullgDIST
)
1224 double invdr01
= spatialData
.value
> 0 ? 1./spatialData
.value
: 0.;
1225 for (int m
= 0; m
< DIM
; m
++)
1227 forces
.force01
[m
] = pcrd
.scalarForce
*spatialData
.dr01
[m
]*invdr01
;
1230 else if (params
.eGeom
== epullgANGLE
)
1233 double cos_theta
, cos_theta2
;
1235 cos_theta
= cos(spatialData
.value
);
1236 cos_theta2
= gmx::square(cos_theta
);
1238 /* The force at theta = 0, pi is undefined so we should not apply any force.
1239 * cos(theta) = 1 for theta = 0, pi so this is what we check for below.
1240 * Note that dr01 = 0 or dr23 = 0 evaluates to theta = 0 so we do not
1241 * have to check for this before dividing by their norm below.
1245 /* The forces f01 and f23 can be decomposed into one vector parallel to dr01
1246 * and another vector parallel to dr23:
1247 * f01 = -dV/dtheta*(a*dr23/|dr23| - b*dr01*|dr01|)/|dr01|,
1248 * f23 = -dV/dtheta*(a*dr01/|dr01| - b*dr23*|dr23|)/|dr23|,
1250 double a
= -gmx::invsqrt(1 - cos_theta2
); /* comes from d/dx acos(x) */
1251 double b
= a
*cos_theta
;
1252 double invdr01
= 1./dnorm(spatialData
.dr01
);
1253 double invdr23
= 1./dnorm(spatialData
.dr23
);
1254 dvec normalized_dr01
, normalized_dr23
;
1255 dsvmul(invdr01
, spatialData
.dr01
, normalized_dr01
);
1256 dsvmul(invdr23
, spatialData
.dr23
, normalized_dr23
);
1258 for (int m
= 0; m
< DIM
; m
++)
1260 /* Here, f_scal is -dV/dtheta */
1261 forces
.force01
[m
] = pcrd
.scalarForce
*invdr01
*(a
*normalized_dr23
[m
] - b
*normalized_dr01
[m
]);
1262 forces
.force23
[m
] = pcrd
.scalarForce
*invdr23
*(a
*normalized_dr01
[m
] - b
*normalized_dr23
[m
]);
1267 /* No forces to apply for ill-defined cases*/
1268 clear_dvec(forces
.force01
);
1269 clear_dvec(forces
.force23
);
1272 else if (params
.eGeom
== epullgANGLEAXIS
)
1274 double cos_theta
, cos_theta2
;
1276 /* The angle-axis force is exactly the same as the angle force (above) except that in
1277 this case the second vector (dr23) is replaced by the pull vector. */
1278 cos_theta
= cos(spatialData
.value
);
1279 cos_theta2
= gmx::square(cos_theta
);
1285 dvec normalized_dr01
;
1287 invdr01
= 1./dnorm(spatialData
.dr01
);
1288 dsvmul(invdr01
, spatialData
.dr01
, normalized_dr01
);
1289 a
= -gmx::invsqrt(1 - cos_theta2
); /* comes from d/dx acos(x) */
1292 for (int m
= 0; m
< DIM
; m
++)
1294 forces
.force01
[m
] = pcrd
.scalarForce
*invdr01
*(a
*spatialData
.vec
[m
] - b
*normalized_dr01
[m
]);
1299 clear_dvec(forces
.force01
);
1302 else if (params
.eGeom
== epullgDIHEDRAL
)
1304 double m2
, n2
, tol
, sqrdist_32
;
1306 /* Note: there is a small difference here compared to the
1307 dihedral force calculations in the bondeds (ref: Bekker 1994).
1308 There rij = ri - rj, while here dr01 = r1 - r0.
1309 However, all distance vectors occur in form of cross or inner products
1310 so that two signs cancel and we end up with the same expressions.
1311 Also, we treat the more general case of 6 groups (0..5) instead of 4 (i, j, k, l).
1313 m2
= diprod(spatialData
.planevec_m
, spatialData
.planevec_m
);
1314 n2
= diprod(spatialData
.planevec_n
, spatialData
.planevec_n
);
1315 dsvmul(-1, spatialData
.dr23
, dr32
);
1316 sqrdist_32
= diprod(dr32
, dr32
);
1317 tol
= sqrdist_32
*GMX_REAL_EPS
; /* Avoid tiny angles */
1318 if ((m2
> tol
) && (n2
> tol
))
1320 double a_01
, a_23_01
, a_23_45
, a_45
;
1321 double inv_dist_32
, inv_sqrdist_32
, dist_32
;
1323 inv_dist_32
= gmx::invsqrt(sqrdist_32
);
1324 inv_sqrdist_32
= inv_dist_32
*inv_dist_32
;
1325 dist_32
= sqrdist_32
*inv_dist_32
;
1327 /* Forces on groups 0, 1 */
1328 a_01
= pcrd
.scalarForce
*dist_32
/m2
; /* scalarForce is -dV/dphi */
1329 dsvmul(-a_01
, spatialData
.planevec_m
, forces
.force01
); /* added sign to get force on group 1, not 0 */
1331 /* Forces on groups 4, 5 */
1332 a_45
= -pcrd
.scalarForce
*dist_32
/n2
;
1333 dsvmul(a_45
, spatialData
.planevec_n
, forces
.force45
); /* force on group 5 */
1335 /* Force on groups 2, 3 (defining the axis) */
1336 a_23_01
= -diprod(spatialData
.dr01
, dr32
)*inv_sqrdist_32
;
1337 a_23_45
= -diprod(spatialData
.dr45
, dr32
)*inv_sqrdist_32
;
1338 dsvmul(-a_23_01
, forces
.force01
, u
); /* added sign to get force from group 0, not 1 */
1339 dsvmul(a_23_45
, forces
.force45
, v
);
1340 dvec_sub(u
, v
, forces
.force23
); /* force on group 3 */
1344 /* No force to apply for ill-defined cases */
1345 clear_dvec(forces
.force01
);
1346 clear_dvec(forces
.force23
);
1347 clear_dvec(forces
.force45
);
1352 for (int m
= 0; m
< DIM
; m
++)
1354 forces
.force01
[m
] = pcrd
.scalarForce
*spatialData
.vec
[m
];
1362 /* We use a global mutex for locking access to the pull data structure
1363 * during registration of external pull potential providers.
1364 * We could use a different, local mutex for each pull object, but the overhead
1365 * is extremely small here and registration is only done during initialization.
1367 static gmx::Mutex registrationMutex
;
1369 using Lock
= gmx::lock_guard
<gmx::Mutex
>;
1371 void register_external_pull_potential(struct pull_t
*pull
,
1373 const char *provider
)
1375 GMX_RELEASE_ASSERT(pull
!= nullptr, "register_external_pull_potential called before init_pull");
1376 GMX_RELEASE_ASSERT(provider
!= nullptr, "register_external_pull_potential called with NULL as provider name");
1378 if (coord_index
< 0 || coord_index
>= gmx::ssize(pull
->coord
))
1380 gmx_fatal(FARGS
, "Module '%s' attempted to register an external potential for pull coordinate %d which is out of the pull coordinate range %d - %zu\n",
1381 provider
, coord_index
+ 1, 1, pull
->coord
.size());
1384 pull_coord_work_t
*pcrd
= &pull
->coord
[coord_index
];
1386 if (pcrd
->params
.eType
!= epullEXTERNAL
)
1388 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'",
1389 provider
, coord_index
+ 1, epull_names
[pcrd
->params
.eType
], epull_names
[epullEXTERNAL
]);
1392 GMX_RELEASE_ASSERT(pcrd
->params
.externalPotentialProvider
!= nullptr, "The external potential provider string for a pull coordinate is NULL");
1394 if (gmx_strcasecmp(provider
, pcrd
->params
.externalPotentialProvider
) != 0)
1396 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'",
1397 provider
, coord_index
+ 1, pcrd
->params
.externalPotentialProvider
);
1400 /* Lock to avoid (extremely unlikely) simultaneous reading and writing of
1401 * pcrd->bExternalPotentialProviderHasBeenRegistered and
1402 * pull->numUnregisteredExternalPotentials.
1404 Lock
registrationLock(registrationMutex
);
1406 if (pcrd
->bExternalPotentialProviderHasBeenRegistered
)
1408 gmx_fatal(FARGS
, "Module '%s' attempted to register an external potential for pull coordinate %d more than once",
1409 provider
, coord_index
+ 1);
1412 pcrd
->bExternalPotentialProviderHasBeenRegistered
= true;
1413 pull
->numUnregisteredExternalPotentials
--;
1415 GMX_RELEASE_ASSERT(pull
->numUnregisteredExternalPotentials
>= 0, "Negative unregistered potentials, the pull code is inconsistent");
1419 static void check_external_potential_registration(const struct pull_t
*pull
)
1421 if (pull
->numUnregisteredExternalPotentials
> 0)
1424 for (c
= 0; c
< pull
->coord
.size(); c
++)
1426 if (pull
->coord
[c
].params
.eType
== epullEXTERNAL
&&
1427 !pull
->coord
[c
].bExternalPotentialProviderHasBeenRegistered
)
1433 GMX_RELEASE_ASSERT(c
< pull
->coord
.size(), "Internal inconsistency in the pull potential provider counting");
1435 gmx_fatal(FARGS
, "No external provider for external pull potentials have been provided for %d pull coordinates. The first coordinate without provider is number %zu, which expects a module named '%s' to provide the external potential.",
1436 pull
->numUnregisteredExternalPotentials
,
1437 c
+ 1, pull
->coord
[c
].params
.externalPotentialProvider
);
1442 /* Pull takes care of adding the forces of the external potential.
1443 * The external potential module has to make sure that the corresponding
1444 * potential energy is added either to the pull term or to a term
1445 * specific to the external module.
1447 void apply_external_pull_coord_force(struct pull_t
*pull
,
1450 const t_mdatoms
*mdatoms
,
1451 gmx::ForceWithVirial
*forceWithVirial
)
1453 pull_coord_work_t
*pcrd
;
1455 GMX_ASSERT(coord_index
>= 0 && coord_index
< gmx::ssize(pull
->coord
), "apply_external_pull_coord_force called with coord_index out of range");
1457 if (pull
->comm
.bParticipate
)
1459 pcrd
= &pull
->coord
[coord_index
];
1461 GMX_RELEASE_ASSERT(pcrd
->params
.eType
== epullEXTERNAL
, "The pull force can only be set externally on pull coordinates of external type");
1463 GMX_ASSERT(pcrd
->bExternalPotentialProviderHasBeenRegistered
, "apply_external_pull_coord_force called for an unregistered pull coordinate");
1466 pcrd
->scalarForce
= coord_force
;
1468 /* Calculate the forces on the pull groups */
1469 PullCoordVectorForces pullCoordForces
= calculateVectorForces(*pcrd
);
1471 /* Add the forces for this coordinate to the total virial and force */
1472 if (forceWithVirial
->computeVirial_
&& pull
->comm
.isMasterRank
)
1474 matrix virial
= { { 0 } };
1475 add_virial_coord(virial
, *pcrd
, pullCoordForces
);
1476 forceWithVirial
->addVirialContribution(virial
);
1479 apply_forces_coord(pull
, coord_index
, pullCoordForces
, mdatoms
,
1480 as_rvec_array(forceWithVirial
->force_
.data()));
1483 pull
->numExternalPotentialsStillToBeAppliedThisStep
--;
1486 /* Calculate the pull potential and scalar force for a pull coordinate.
1487 * Returns the vector forces for the pull coordinate.
1489 static PullCoordVectorForces
1490 do_pull_pot_coord(struct pull_t
*pull
, int coord_ind
, t_pbc
*pbc
,
1491 double t
, real lambda
,
1492 real
*V
, tensor vir
, real
*dVdl
)
1494 pull_coord_work_t
&pcrd
= pull
->coord
[coord_ind
];
1496 assert(pcrd
.params
.eType
!= epullCONSTRAINT
);
1498 double dev
= get_pull_coord_deviation(pull
, coord_ind
, pbc
, t
);
1500 calc_pull_coord_scalar_force_and_potential(&pcrd
, dev
, lambda
, V
, dVdl
);
1502 PullCoordVectorForces pullCoordForces
= calculateVectorForces(pcrd
);
1504 add_virial_coord(vir
, pcrd
, pullCoordForces
);
1506 return pullCoordForces
;
1509 real
pull_potential(struct pull_t
*pull
, const t_mdatoms
*md
, t_pbc
*pbc
,
1510 const t_commrec
*cr
, double t
, real lambda
,
1511 const rvec
*x
, gmx::ForceWithVirial
*force
, real
*dvdlambda
)
1515 assert(pull
!= nullptr);
1517 /* Ideally we should check external potential registration only during
1518 * the initialization phase, but that requires another function call
1519 * that should be done exactly in the right order. So we check here.
1521 check_external_potential_registration(pull
);
1523 if (pull
->comm
.bParticipate
)
1527 pull_calc_coms(cr
, pull
, md
, pbc
, t
, x
, nullptr);
1529 rvec
*f
= as_rvec_array(force
->force_
.data());
1530 matrix virial
= { { 0 } };
1531 const bool computeVirial
= (force
->computeVirial_
&& MASTER(cr
));
1532 for (size_t c
= 0; c
< pull
->coord
.size(); c
++)
1534 pull_coord_work_t
*pcrd
;
1535 pcrd
= &pull
->coord
[c
];
1537 /* For external potential the force is assumed to be given by an external module by a call to
1538 apply_pull_coord_external_force */
1539 if (pcrd
->params
.eType
== epullCONSTRAINT
|| pcrd
->params
.eType
== epullEXTERNAL
)
1544 PullCoordVectorForces pullCoordForces
=
1545 do_pull_pot_coord(pull
, c
, pbc
, t
, lambda
,
1547 computeVirial
? virial
: nullptr,
1550 /* Distribute the force over the atoms in the pulled groups */
1551 apply_forces_coord(pull
, c
, pullCoordForces
, md
, f
);
1556 force
->addVirialContribution(virial
);
1561 GMX_ASSERT(pull
->numExternalPotentialsStillToBeAppliedThisStep
== 0, "Too few or too many external pull potentials have been applied the previous step");
1562 /* All external pull potentials still need to be applied */
1563 pull
->numExternalPotentialsStillToBeAppliedThisStep
=
1564 pull
->numCoordinatesWithExternalPotential
;
1566 return (MASTER(cr
) ? V
: 0.0);
1569 void pull_constraint(struct pull_t
*pull
, const t_mdatoms
*md
, t_pbc
*pbc
,
1570 const t_commrec
*cr
, double dt
, double t
,
1571 rvec
*x
, rvec
*xp
, rvec
*v
, tensor vir
)
1573 assert(pull
!= nullptr);
1575 if (pull
->comm
.bParticipate
)
1577 pull_calc_coms(cr
, pull
, md
, pbc
, t
, x
, xp
);
1579 do_constraint(pull
, pbc
, xp
, v
, MASTER(cr
), vir
, dt
, t
);
1583 void dd_make_local_pull_groups(const t_commrec
*cr
, struct pull_t
*pull
)
1587 gmx_bool bMustParticipate
;
1593 /* We always make the master node participate, such that it can do i/o,
1594 * add the virial and to simplify MC type extensions people might have.
1596 bMustParticipate
= (comm
->bParticipateAll
|| comm
->isMasterRank
);
1598 for (pull_group_work_t
&group
: pull
->group
)
1600 if (!group
.globalWeights
.empty())
1602 group
.localWeights
.resize(group
.atomSet
.numAtomsLocal());
1603 for (size_t i
= 0; i
< group
.atomSet
.numAtomsLocal(); ++i
)
1605 group
.localWeights
[i
] = group
.globalWeights
[group
.atomSet
.collectiveIndex()[i
]];
1609 GMX_ASSERT(bMustParticipate
|| dd
!= nullptr, "Either all ranks (including this rank) participate, or we use DD and need to have access to dd here");
1611 /* We should participate if we have pull or pbc atoms */
1612 if (!bMustParticipate
&&
1613 (group
.atomSet
.numAtomsLocal() > 0 || (group
.epgrppbc
== epgrppbcREFAT
&&
1614 group
.pbcAtomSet
->numAtomsLocal() > 0)))
1616 bMustParticipate
= TRUE
;
1620 if (!comm
->bParticipateAll
)
1622 /* Keep currently not required ranks in the communicator
1623 * if they needed to participate up to 20 decompositions ago.
1624 * This avoids frequent rebuilds due to atoms jumping back and forth.
1626 const int64_t history_count
= 20;
1627 gmx_bool bWillParticipate
;
1630 /* Increase the decomposition counter for the current call */
1631 comm
->setup_count
++;
1633 if (bMustParticipate
)
1635 comm
->must_count
= comm
->setup_count
;
1640 (comm
->bParticipate
&&
1641 comm
->must_count
>= comm
->setup_count
- history_count
);
1643 if (debug
&& dd
!= nullptr)
1645 fprintf(debug
, "Our DD rank (%3d) pull #atoms>0 or master: %s, will be part %s\n",
1646 dd
->rank
, gmx::boolToString(bMustParticipate
), gmx::boolToString(bWillParticipate
));
1649 if (bWillParticipate
)
1651 /* Count the number of ranks that we want to have participating */
1653 /* Count the number of ranks that need to be added */
1654 count
[1] = comm
->bParticipate
? 0 : 1;
1662 /* The cost of this global operation will be less that the cost
1663 * of the extra MPI_Comm_split calls that we can avoid.
1665 gmx_sumi(2, count
, cr
);
1667 /* If we are missing ranks or if we have 20% more ranks than needed
1668 * we make a new sub-communicator.
1670 if (count
[1] > 0 || 6*count
[0] < 5*comm
->nparticipate
)
1674 fprintf(debug
, "Creating new pull subcommunicator of size %d\n",
1679 if (comm
->mpi_comm_com
!= MPI_COMM_NULL
)
1681 MPI_Comm_free(&comm
->mpi_comm_com
);
1684 /* This might be an extremely expensive operation, so we try
1685 * to avoid this splitting as much as possible.
1687 assert(dd
!= nullptr);
1688 MPI_Comm_split(dd
->mpi_comm_all
, bWillParticipate
? 0 : 1, dd
->rank
,
1689 &comm
->mpi_comm_com
);
1692 comm
->bParticipate
= bWillParticipate
;
1693 comm
->nparticipate
= count
[0];
1695 /* When we use the previous COM for PBC, we need to broadcast
1696 * the previous COM to ranks that have joined the communicator.
1698 for (pull_group_work_t
&group
: pull
->group
)
1700 if (group
.epgrppbc
== epgrppbcPREVSTEPCOM
)
1702 GMX_ASSERT(comm
->bParticipate
|| !MASTER(cr
),
1703 "The master rank has to participate, as it should pass an up to date prev. COM "
1704 "to bcast here as well as to e.g. checkpointing");
1706 gmx_bcast(sizeof(dvec
), group
.x_prev_step
, cr
);
1712 /* Since the PBC of atoms might have changed, we need to update the PBC */
1713 pull
->bSetPBCatoms
= TRUE
;
1716 static void init_pull_group_index(FILE *fplog
, const t_commrec
*cr
,
1717 int g
, pull_group_work_t
*pg
,
1718 gmx_bool bConstraint
, const ivec pulldim_con
,
1719 const gmx_mtop_t
*mtop
,
1720 const t_inputrec
*ir
, real lambda
)
1722 /* With EM and BD there are no masses in the integrator.
1723 * But we still want to have the correct mass-weighted COMs.
1724 * So we store the real masses in the weights.
1726 const bool setWeights
= (pg
->params
.nweight
> 0 ||
1727 EI_ENERGY_MINIMIZATION(ir
->eI
) ||
1730 /* In parallel, store we need to extract localWeights from weights at DD time */
1731 std::vector
<real
> &weights
= ((cr
&& PAR(cr
)) ? pg
->globalWeights
: pg
->localWeights
);
1733 const SimulationGroups
&groups
= mtop
->groups
;
1735 /* Count frozen dimensions and (weighted) mass */
1741 for (int i
= 0; i
< pg
->params
.nat
; i
++)
1743 int ii
= pg
->params
.ind
[i
];
1744 if (bConstraint
&& ir
->opts
.nFreeze
)
1746 for (int d
= 0; d
< DIM
; d
++)
1748 if (pulldim_con
[d
] == 1 &&
1749 ir
->opts
.nFreeze
[getGroupType(groups
, SimulationAtomGroupType::Freeze
, ii
)][d
])
1755 const t_atom
&atom
= mtopGetAtomParameters(mtop
, ii
, &molb
);
1757 if (ir
->efep
== efepNO
)
1763 m
= (1 - lambda
)*atom
.m
+ lambda
*atom
.mB
;
1766 if (pg
->params
.nweight
> 0)
1768 w
= pg
->params
.weight
[i
];
1774 if (EI_ENERGY_MINIMIZATION(ir
->eI
))
1776 /* Move the mass to the weight */
1780 else if (ir
->eI
== eiBD
)
1783 if (ir
->bd_fric
!= 0.0F
)
1785 mbd
= ir
->bd_fric
*ir
->delta_t
;
1789 if (groups
.groupNumbers
[SimulationAtomGroupType::TemperatureCoupling
].empty())
1791 mbd
= ir
->delta_t
/ir
->opts
.tau_t
[0];
1795 mbd
= ir
->delta_t
/ir
->opts
.tau_t
[groups
.groupNumbers
[SimulationAtomGroupType::TemperatureCoupling
][ii
]];
1803 weights
.push_back(w
);
1812 /* We can have single atom groups with zero mass with potential pulling
1813 * without cosine weighting.
1815 if (pg
->params
.nat
== 1 && !bConstraint
&& pg
->epgrppbc
!= epgrppbcCOS
)
1817 /* With one atom the mass doesn't matter */
1822 gmx_fatal(FARGS
, "The total%s mass of pull group %d is zero",
1823 pg
->params
.weight
? " weighted" : "", g
);
1829 "Pull group %d: %5d atoms, mass %9.3f",
1830 g
, pg
->params
.nat
, tmass
);
1831 if (pg
->params
.weight
||
1832 EI_ENERGY_MINIMIZATION(ir
->eI
) || ir
->eI
== eiBD
)
1834 fprintf(fplog
, ", weighted mass %9.3f", wmass
*wmass
/wwmass
);
1836 if (pg
->epgrppbc
== epgrppbcCOS
)
1838 fprintf(fplog
, ", cosine weighting will be used");
1840 fprintf(fplog
, "\n");
1845 /* A value != 0 signals not frozen, it is updated later */
1851 for (int d
= 0; d
< DIM
; d
++)
1853 ndim
+= pulldim_con
[d
]*pg
->params
.nat
;
1855 if (fplog
&& nfrozen
> 0 && nfrozen
< ndim
)
1858 "\nWARNING: In pull group %d some, but not all of the degrees of freedom\n"
1859 " that are subject to pulling are frozen.\n"
1860 " For constraint pulling the whole group will be frozen.\n\n",
1869 init_pull(FILE *fplog
, const pull_params_t
*pull_params
, const t_inputrec
*ir
,
1870 const gmx_mtop_t
*mtop
, const t_commrec
*cr
, gmx::LocalAtomSetManager
*atomSets
,
1873 struct pull_t
*pull
;
1878 /* Copy the pull parameters */
1879 pull
->params
= *pull_params
;
1880 /* Avoid pointer copies */
1881 pull
->params
.group
= nullptr;
1882 pull
->params
.coord
= nullptr;
1884 for (int i
= 0; i
< pull_params
->ngroup
; ++i
)
1886 pull
->group
.emplace_back(pull_params
->group
[i
], atomSets
->add({pull_params
->group
[i
].ind
, pull_params
->group
[i
].ind
+pull_params
->group
[i
].nat
}),
1887 pull_params
->bSetPbcRefToPrevStepCOM
);
1890 if (cr
!= nullptr && DOMAINDECOMP(cr
))
1892 /* Set up the global to local atom mapping for PBC atoms */
1893 for (pull_group_work_t
&group
: pull
->group
)
1895 if (group
.epgrppbc
== epgrppbcREFAT
|| group
.epgrppbc
== epgrppbcPREVSTEPCOM
)
1897 /* pbcAtomSet consists of a single atom */
1898 group
.pbcAtomSet
= std::make_unique
<gmx::LocalAtomSet
>(atomSets
->add({&group
.params
.pbcatom
, &group
.params
.pbcatom
+ 1}));
1903 pull
->bPotential
= FALSE
;
1904 pull
->bConstraint
= FALSE
;
1905 pull
->bCylinder
= FALSE
;
1906 pull
->bAngle
= FALSE
;
1907 pull
->bXOutAverage
= pull_params
->bXOutAverage
;
1908 pull
->bFOutAverage
= pull_params
->bFOutAverage
;
1910 GMX_RELEASE_ASSERT(pull
->group
[0].params
.nat
== 0, "pull group 0 is an absolute reference group and should not contain atoms");
1912 pull
->numCoordinatesWithExternalPotential
= 0;
1914 for (int c
= 0; c
< pull
->params
.ncoord
; c
++)
1916 /* Construct a pull coordinate, copying all coordinate parameters */
1917 pull
->coord
.emplace_back(pull_params
->coord
[c
]);
1919 pull_coord_work_t
*pcrd
= &pull
->coord
.back();
1921 switch (pcrd
->params
.eGeom
)
1924 case epullgDIRRELATIVE
: /* Direction vector is determined at each step */
1926 case epullgDIHEDRAL
:
1931 case epullgANGLEAXIS
:
1932 copy_rvec_to_dvec(pull_params
->coord
[c
].vec
, pcrd
->spatialData
.vec
);
1935 /* We allow reading of newer tpx files with new pull geometries,
1936 * but with the same tpx format, with old code. A new geometry
1937 * only adds a new enum value, which will be out of range for
1938 * old code. The first place we need to generate an error is
1939 * here, since the pull code can't handle this.
1940 * The enum can be used before arriving here only for printing
1941 * the string corresponding to the geometry, which will result
1942 * in printing "UNDEFINED".
1944 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.",
1945 c
+ 1, pcrd
->params
.eGeom
, epullgNR
- 1);
1948 if (pcrd
->params
.eType
== epullCONSTRAINT
)
1950 /* Check restrictions of the constraint pull code */
1951 if (pcrd
->params
.eGeom
== epullgCYL
||
1952 pcrd
->params
.eGeom
== epullgDIRRELATIVE
||
1953 pcrd
->params
.eGeom
== epullgANGLE
||
1954 pcrd
->params
.eGeom
== epullgDIHEDRAL
||
1955 pcrd
->params
.eGeom
== epullgANGLEAXIS
)
1957 gmx_fatal(FARGS
, "Pulling of type %s can not be combined with geometry %s. Consider using pull type %s.",
1958 epull_names
[pcrd
->params
.eType
],
1959 epullg_names
[pcrd
->params
.eGeom
],
1960 epull_names
[epullUMBRELLA
]);
1963 pull
->bConstraint
= TRUE
;
1967 pull
->bPotential
= TRUE
;
1970 if (pcrd
->params
.eGeom
== epullgCYL
)
1972 pull
->bCylinder
= TRUE
;
1974 else if (pcrd
->params
.eGeom
== epullgANGLE
|| pcrd
->params
.eGeom
== epullgDIHEDRAL
|| pcrd
->params
.eGeom
== epullgANGLEAXIS
)
1976 pull
->bAngle
= TRUE
;
1979 /* We only need to calculate the plain COM of a group
1980 * when it is not only used as a cylinder group.
1981 * Also the absolute reference group 0 needs no COM computation.
1983 for (int i
= 0; i
< pcrd
->params
.ngroup
; i
++)
1985 int groupIndex
= pcrd
->params
.group
[i
];
1986 if (groupIndex
> 0 &&
1987 !(pcrd
->params
.eGeom
== epullgCYL
&& i
== 0))
1989 pull
->group
[groupIndex
].needToCalcCom
= true;
1993 /* With non-zero rate the reference value is set at every step */
1994 if (pcrd
->params
.rate
== 0)
1996 /* Initialize the constant reference value */
1997 if (pcrd
->params
.eType
!= epullEXTERNAL
)
1999 low_set_pull_coord_reference_value(pcrd
, c
, pcrd
->params
.init
*pull_conversion_factor_userinput2internal(&pcrd
->params
));
2003 /* With an external pull potential, the reference value loses
2004 * it's meaning and should not be used. Setting it to zero
2005 * makes any terms dependent on it disappear.
2006 * The only issue this causes is that with cylinder pulling
2007 * no atoms of the cylinder group within the cylinder radius
2008 * should be more than half a box length away from the COM of
2009 * the pull group along the axial direction.
2011 pcrd
->value_ref
= 0.0;
2015 if (pcrd
->params
.eType
== epullEXTERNAL
)
2017 GMX_RELEASE_ASSERT(pcrd
->params
.rate
== 0, "With an external potential, a pull coordinate should have rate = 0");
2019 /* This external potential needs to be registered later */
2020 pull
->numCoordinatesWithExternalPotential
++;
2022 pcrd
->bExternalPotentialProviderHasBeenRegistered
= false;
2025 pull
->numUnregisteredExternalPotentials
=
2026 pull
->numCoordinatesWithExternalPotential
;
2027 pull
->numExternalPotentialsStillToBeAppliedThisStep
= 0;
2029 pull
->ePBC
= ir
->ePBC
;
2032 case epbcNONE
: pull
->npbcdim
= 0; break;
2033 case epbcXY
: pull
->npbcdim
= 2; break;
2034 default: pull
->npbcdim
= 3; break;
2039 gmx_bool bAbs
, bCos
;
2042 for (const pull_coord_work_t
&coord
: pull
->coord
)
2044 if (pull
->group
[coord
.params
.group
[0]].params
.nat
== 0 ||
2045 pull
->group
[coord
.params
.group
[1]].params
.nat
== 0)
2051 fprintf(fplog
, "\n");
2052 if (pull
->bPotential
)
2054 fprintf(fplog
, "Will apply potential COM pulling\n");
2056 if (pull
->bConstraint
)
2058 fprintf(fplog
, "Will apply constraint COM pulling\n");
2060 // Don't include the reference group 0 in output, so we report ngroup-1
2061 int numRealGroups
= pull
->group
.size() - 1;
2062 GMX_RELEASE_ASSERT(numRealGroups
> 0, "The reference absolute position pull group should always be present");
2063 fprintf(fplog
, "with %zu pull coordinate%s and %d group%s\n",
2064 pull
->coord
.size(), pull
->coord
.size() == 1 ? "" : "s",
2065 numRealGroups
, numRealGroups
== 1 ? "" : "s");
2068 fprintf(fplog
, "with an absolute reference\n");
2071 // Don't include the reference group 0 in loop
2072 for (size_t g
= 1; g
< pull
->group
.size(); g
++)
2074 if (pull
->group
[g
].params
.nat
> 1 &&
2075 pull
->group
[g
].params
.pbcatom
< 0)
2077 /* We are using cosine weighting */
2078 fprintf(fplog
, "Cosine weighting is used for group %zu\n", g
);
2084 please_cite(fplog
, "Engin2010");
2088 pull
->bRefAt
= FALSE
;
2090 for (size_t g
= 0; g
< pull
->group
.size(); g
++)
2092 pull_group_work_t
*pgrp
;
2094 pgrp
= &pull
->group
[g
];
2095 if (pgrp
->params
.nat
> 0)
2097 /* There is an issue when a group is used in multiple coordinates
2098 * and constraints are applied in different dimensions with atoms
2099 * frozen in some, but not all dimensions.
2100 * Since there is only one mass array per group, we can't have
2101 * frozen/non-frozen atoms for different coords at the same time.
2102 * But since this is a very exotic case, we don't check for this.
2103 * A warning is printed in init_pull_group_index.
2106 gmx_bool bConstraint
;
2107 ivec pulldim
, pulldim_con
;
2109 /* Loop over all pull coordinates to see along which dimensions
2110 * this group is pulled and if it is involved in constraints.
2112 bConstraint
= FALSE
;
2113 clear_ivec(pulldim
);
2114 clear_ivec(pulldim_con
);
2115 for (const pull_coord_work_t
&coord
: pull
->coord
)
2117 gmx_bool bGroupUsed
= FALSE
;
2118 for (int gi
= 0; gi
< coord
.params
.ngroup
; gi
++)
2120 if (coord
.params
.group
[gi
] == static_cast<int>(g
))
2128 for (int m
= 0; m
< DIM
; m
++)
2130 if (coord
.params
.dim
[m
] == 1)
2134 if (coord
.params
.eType
== epullCONSTRAINT
)
2144 /* Determine if we need to take PBC into account for calculating
2145 * the COM's of the pull groups.
2147 switch (pgrp
->epgrppbc
)
2150 pull
->bRefAt
= TRUE
;
2153 if (pgrp
->params
.weight
!= nullptr)
2155 gmx_fatal(FARGS
, "Pull groups can not have relative weights and cosine weighting at same time");
2157 for (int m
= 0; m
< DIM
; m
++)
2159 if (m
< pull
->npbcdim
&& pulldim
[m
] == 1)
2161 if (pull
->cosdim
>= 0 && pull
->cosdim
!= m
)
2163 gmx_fatal(FARGS
, "Can only use cosine weighting with pulling in one dimension (use mdp option pull-coord?-dim)");
2173 /* Set the indices */
2174 init_pull_group_index(fplog
, cr
, g
, pgrp
,
2175 bConstraint
, pulldim_con
,
2180 /* Absolute reference, set the inverse mass to zero.
2181 * This is only relevant (and used) with constraint pulling.
2188 /* If we use cylinder coordinates, do some initialising for them */
2189 if (pull
->bCylinder
)
2191 for (const pull_coord_work_t
&coord
: pull
->coord
)
2193 if (coord
.params
.eGeom
== epullgCYL
)
2195 if (pull
->group
[coord
.params
.group
[0]].params
.nat
== 0)
2197 gmx_fatal(FARGS
, "A cylinder pull group is not supported when using absolute reference!\n");
2200 const auto &referenceGroup
= pull
->group
[coord
.params
.group
[0]];
2201 pull
->dyna
.emplace_back(referenceGroup
.params
, referenceGroup
.atomSet
, pull
->params
.bSetPbcRefToPrevStepCOM
);
2205 /* The gmx_omp_nthreads module might not be initialized here, so max(1,) */
2206 pull
->nthreads
= std::max(1, gmx_omp_nthreads_get(emntDefault
));
2207 pull
->comSums
.resize(pull
->nthreads
);
2212 /* Use a sub-communicator when we have more than 32 ranks, but not
2213 * when we have an external pull potential, since then the external
2214 * potential provider expects each rank to have the coordinate.
2216 comm
->bParticipateAll
= (cr
== nullptr || !DOMAINDECOMP(cr
) ||
2217 cr
->dd
->nnodes
<= 32 ||
2218 pull
->numCoordinatesWithExternalPotential
> 0 ||
2219 getenv("GMX_PULL_PARTICIPATE_ALL") != nullptr);
2220 /* This sub-commicator is not used with comm->bParticipateAll,
2221 * so we can always initialize it to NULL.
2223 comm
->mpi_comm_com
= MPI_COMM_NULL
;
2224 comm
->nparticipate
= 0;
2225 comm
->isMasterRank
= (cr
== nullptr || MASTER(cr
));
2227 /* No MPI: 1 rank: all ranks pull */
2228 comm
->bParticipateAll
= TRUE
;
2229 comm
->isMasterRank
= true;
2231 comm
->bParticipate
= comm
->bParticipateAll
;
2232 comm
->setup_count
= 0;
2233 comm
->must_count
= 0;
2235 if (!comm
->bParticipateAll
&& fplog
!= nullptr)
2237 fprintf(fplog
, "Will use a sub-communicator for pull communication\n");
2240 comm
->pbcAtomBuffer
.resize(pull
->group
.size());
2241 comm
->comBuffer
.resize(pull
->group
.size()*c_comBufferStride
);
2242 if (pull
->bCylinder
)
2244 comm
->cylinderBuffer
.resize(pull
->coord
.size()*c_cylinderBufferStride
);
2247 /* We still need to initialize the PBC reference coordinates */
2248 pull
->bSetPBCatoms
= TRUE
;
2250 pull
->out_x
= nullptr;
2251 pull
->out_f
= nullptr;
2256 static void destroy_pull(struct pull_t
*pull
)
2259 if (pull
->comm
.mpi_comm_com
!= MPI_COMM_NULL
)
2261 MPI_Comm_free(&pull
->comm
.mpi_comm_com
);
2268 void preparePrevStepPullCom(const t_inputrec
*ir
, pull_t
*pull_work
,
2269 const t_mdatoms
*md
, t_state
*state
, const t_state
*state_global
,
2270 const t_commrec
*cr
, bool startingFromCheckpoint
)
2272 if (!ir
->pull
|| !ir
->pull
->bSetPbcRefToPrevStepCOM
)
2276 allocStatePrevStepPullCom(state
, pull_work
);
2277 if (startingFromCheckpoint
)
2281 state
->pull_com_prev_step
= state_global
->pull_com_prev_step
;
2285 /* Only the master rank has the checkpointed COM from the previous step */
2286 gmx_bcast(sizeof(double) * state
->pull_com_prev_step
.size(), &state
->pull_com_prev_step
[0], cr
);
2288 setPrevStepPullComFromState(pull_work
, state
);
2293 set_pbc(&pbc
, ir
->ePBC
, state
->box
);
2294 initPullComFromPrevStep(cr
, pull_work
, md
, &pbc
, state
->x
.rvec_array());
2295 updatePrevStepPullCom(pull_work
, state
);
2299 void finish_pull(struct pull_t
*pull
)
2301 check_external_potential_registration(pull
);
2305 gmx_fio_fclose(pull
->out_x
);
2309 gmx_fio_fclose(pull
->out_f
);
2315 gmx_bool
pull_have_potential(const struct pull_t
*pull
)
2317 return pull
->bPotential
;
2320 gmx_bool
pull_have_constraint(const struct pull_t
*pull
)
2322 return pull
->bConstraint
;