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.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/gmxlib/network.h"
47 #include "gromacs/math/functions.h"
48 #include "gromacs/math/utilities.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/mdtypes/commrec.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/mdtypes/mdatom.h"
54 #include "gromacs/mdtypes/state.h"
55 #include "gromacs/pbcutil/pbc.h"
56 #include "gromacs/pulling/pull.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/real.h"
61 #include "gromacs/utility/smalloc.h"
63 #include "pull_internal.h"
67 // Helper function to deduce MPI datatype from the type of data
68 gmx_unused static MPI_Datatype mpiDatatype(const float gmx_unused* data)
73 // Helper function to deduce MPI datatype from the type of data
74 gmx_unused static MPI_Datatype mpiDatatype(const double gmx_unused* data)
82 // Helper function; note that gmx_sum(d) should actually be templated
83 gmx_unused static void gmxAllReduce(int n, real* data, const t_commrec* cr)
89 // Helper function; note that gmx_sum(d) should actually be templated
90 gmx_unused static void gmxAllReduce(int n, double* data, const t_commrec* cr)
92 gmx_sumd(n, data, cr);
95 // Reduce data of n elements over all ranks currently participating in pull
97 static void pullAllReduce(const t_commrec* cr, pull_comm_t* comm, int n, T* data)
99 if (cr != nullptr && PAR(cr))
101 if (comm->bParticipateAll)
103 /* Sum the contributions over all DD ranks */
104 gmxAllReduce(n, data, cr);
108 /* Separate branch because gmx_sum uses cr->mpi_comm_mygroup */
110 MPI_Allreduce(MPI_IN_PLACE, data, n, mpiDatatype(data), MPI_SUM, comm->mpi_comm_com);
112 gmx_incons("comm->bParticipateAll=FALSE without GMX_MPI");
118 /* Copies the coordinates of the PBC atom of pgrp to x_pbc.
119 * When those coordinates are not available on this rank, clears x_pbc.
121 static void setPbcAtomCoords(const pull_group_work_t& pgrp, const rvec* x, rvec x_pbc)
123 if (pgrp.pbcAtomSet != nullptr)
125 if (pgrp.pbcAtomSet->numAtomsLocal() > 0)
127 /* We have the atom locally, copy its coordinates */
128 copy_rvec(x[pgrp.pbcAtomSet->localIndex()[0]], x_pbc);
132 /* Another rank has it, clear the coordinates for MPI_Allreduce */
138 copy_rvec(x[pgrp.params.pbcatom], x_pbc);
142 static void pull_set_pbcatoms(const t_commrec* cr, struct pull_t* pull, const rvec* x, gmx::ArrayRef<gmx::RVec> x_pbc)
145 for (size_t g = 0; g < pull->group.size(); g++)
147 const pull_group_work_t& group = pull->group[g];
148 if (group.needToCalcCom && (group.epgrppbc == epgrppbcREFAT || group.epgrppbc == epgrppbcPREVSTEPCOM))
150 setPbcAtomCoords(pull->group[g], x, x_pbc[g]);
155 clear_rvec(x_pbc[g]);
159 if (cr && PAR(cr) && numPbcAtoms > 0)
161 /* Sum over participating ranks to get x_pbc from the home ranks.
162 * This can be very expensive at high parallelization, so we only
163 * do this after each DD repartitioning.
165 pullAllReduce(cr, &pull->comm, pull->group.size() * DIM, static_cast<real*>(x_pbc[0]));
170 make_cyl_refgrps(const t_commrec* cr, pull_t* pull, const real* masses, t_pbc* pbc, double t, const rvec* x)
172 pull_comm_t* comm = &pull->comm;
174 GMX_ASSERT(comm->cylinderBuffer.size() == pull->coord.size() * c_cylinderBufferStride,
175 "cylinderBuffer should have the correct size");
177 double inv_cyl_r2 = 1.0 / gmx::square(pull->params.cylinder_r);
179 /* loop over all groups to make a reference group for each*/
180 for (size_t c = 0; c < pull->coord.size(); c++)
182 pull_coord_work_t* pcrd;
183 double sum_a, wmass, wwmass;
184 dvec radf_fac0, radf_fac1;
186 pcrd = &pull->coord[c];
191 clear_dvec(radf_fac0);
192 clear_dvec(radf_fac1);
194 if (pcrd->params.eGeom == epullgCYL)
196 /* pref will be the same group for all pull coordinates */
197 const pull_group_work_t& pref = pull->group[pcrd->params.group[0]];
198 const pull_group_work_t& pgrp = pull->group[pcrd->params.group[1]];
199 pull_group_work_t& pdyna = pull->dyna[c];
201 copy_dvec_to_rvec(pcrd->spatialData.vec, direction);
203 /* Since we have not calculated the COM of the cylinder group yet,
204 * we calculate distances with respect to location of the pull
205 * group minus the reference position along the vector.
206 * here we already have the COM of the pull group. This resolves
207 * any PBC issues and we don't need to use a PBC-atom here.
209 if (pcrd->params.rate != 0)
211 /* With rate=0, value_ref is set initially */
212 pcrd->value_ref = pcrd->params.init + pcrd->params.rate * t;
215 for (int m = 0; m < DIM; m++)
217 reference[m] = pgrp.x[m] - pcrd->spatialData.vec[m] * pcrd->value_ref;
220 auto localAtomIndices = pref.atomSet.localIndex();
222 /* This actually only needs to be done at init or DD time,
223 * but resizing with the same size does not cause much overhead.
225 pdyna.localWeights.resize(localAtomIndices.size());
226 pdyna.mdw.resize(localAtomIndices.size());
227 pdyna.dv.resize(localAtomIndices.size());
229 /* loop over all atoms in the main ref group */
230 for (gmx::index indexInSet = 0; indexInSet < localAtomIndices.ssize(); indexInSet++)
232 int atomIndex = localAtomIndices[indexInSet];
234 pbc_dx_aiuc(pbc, x[atomIndex], reference, dx);
235 double axialLocation = iprod(direction, dx);
238 for (int m = 0; m < DIM; m++)
240 /* Determine the radial components */
241 radialLocation[m] = dx[m] - axialLocation * direction[m];
242 dr2 += gmx::square(radialLocation[m]);
244 double dr2_rel = dr2 * inv_cyl_r2;
248 /* add atom to sum of COM and to weight array */
250 double mass = masses[atomIndex];
251 /* The radial weight function is 1-2x^2+x^4,
252 * where x=r/cylinder_r. Since this function depends
253 * on the radial component, we also get radial forces
256 double weight = 1 + (-2 + dr2_rel) * dr2_rel;
257 double dweight_r = (-4 + 4 * dr2_rel) * inv_cyl_r2;
258 pdyna.localWeights[indexInSet] = weight;
259 sum_a += mass * weight * axialLocation;
260 wmass += mass * weight;
261 wwmass += mass * weight * weight;
263 dsvmul(mass * dweight_r, radialLocation, mdw);
264 copy_dvec(mdw, pdyna.mdw[indexInSet]);
265 /* Currently we only have the axial component of the
266 * offset from the cylinder COM up to an unkown offset.
267 * We add this offset after the reduction needed
268 * for determining the COM of the cylinder group.
270 pdyna.dv[indexInSet] = axialLocation;
271 for (int m = 0; m < DIM; m++)
273 radf_fac0[m] += mdw[m];
274 radf_fac1[m] += mdw[m] * axialLocation;
279 pdyna.localWeights[indexInSet] = 0;
284 auto buffer = gmx::arrayRefFromArray(
285 comm->cylinderBuffer.data() + c * c_cylinderBufferStride, c_cylinderBufferStride);
291 buffer[3] = radf_fac0[XX];
292 buffer[4] = radf_fac0[YY];
293 buffer[5] = radf_fac0[ZZ];
295 buffer[6] = radf_fac1[XX];
296 buffer[7] = radf_fac1[YY];
297 buffer[8] = radf_fac1[ZZ];
300 if (cr != nullptr && PAR(cr))
302 /* Sum the contributions over the ranks */
303 pullAllReduce(cr, comm, pull->coord.size() * c_cylinderBufferStride, comm->cylinderBuffer.data());
306 for (size_t c = 0; c < pull->coord.size(); c++)
308 pull_coord_work_t* pcrd;
310 pcrd = &pull->coord[c];
312 if (pcrd->params.eGeom == epullgCYL)
314 pull_group_work_t* pdyna = &pull->dyna[c];
315 pull_group_work_t* pgrp = &pull->group[pcrd->params.group[1]];
316 PullCoordSpatialData& spatialData = pcrd->spatialData;
318 auto buffer = gmx::constArrayRefFromArray(
319 comm->cylinderBuffer.data() + c * c_cylinderBufferStride, c_cylinderBufferStride);
320 double wmass = buffer[0];
321 double wwmass = buffer[1];
322 pdyna->mwscale = 1.0 / wmass;
323 /* Cylinder pulling can't be used with constraints, but we set
324 * wscale and invtm anyhow, in case someone would like to use them.
326 pdyna->wscale = wmass / wwmass;
327 pdyna->invtm = wwmass / (wmass * wmass);
329 /* We store the deviation of the COM from the reference location
330 * used above, since we need it when we apply the radial forces
331 * to the atoms in the cylinder group.
333 spatialData.cyl_dev = 0;
334 for (int m = 0; m < DIM; m++)
336 double reference = pgrp->x[m] - spatialData.vec[m] * pcrd->value_ref;
337 double dist = -spatialData.vec[m] * buffer[2] * pdyna->mwscale;
338 pdyna->x[m] = reference - dist;
339 spatialData.cyl_dev += dist;
341 /* Now we know the exact COM of the cylinder reference group,
342 * we can determine the radial force factor (ffrad) that when
343 * multiplied with the axial pull force will give the radial
344 * force on the pulled (non-cylinder) group.
346 for (int m = 0; m < DIM; m++)
348 spatialData.ffrad[m] = (buffer[6 + m] + buffer[3 + m] * spatialData.cyl_dev) / wmass;
354 "Pull cylinder group %zu:%8.3f%8.3f%8.3f m:%8.3f\n",
361 "ffrad %8.3f %8.3f %8.3f\n",
362 spatialData.ffrad[XX],
363 spatialData.ffrad[YY],
364 spatialData.ffrad[ZZ]);
370 static double atan2_0_2pi(double y, double x)
382 static void sum_com_part(const pull_group_work_t* pgrp,
394 dvec sum_wmx = { 0, 0, 0 };
395 dvec sum_wmxp = { 0, 0, 0 };
397 auto localAtomIndices = pgrp->atomSet.localIndex();
398 for (int i = ind_start; i < ind_end; i++)
400 int ii = localAtomIndices[i];
402 if (pgrp->localWeights.empty())
411 w = pgrp->localWeights[i];
416 if (pgrp->epgrppbc == epgrppbcNONE)
418 /* Plain COM: sum the coordinates */
419 for (int d = 0; d < DIM; d++)
421 sum_wmx[d] += wm * x[ii][d];
425 for (int d = 0; d < DIM; d++)
427 sum_wmxp[d] += wm * xp[ii][d];
435 /* Sum the difference with the reference atom */
436 pbc_dx(pbc, x[ii], x_pbc, dx);
437 for (int d = 0; d < DIM; d++)
439 sum_wmx[d] += wm * dx[d];
443 /* For xp add the difference between xp and x to dx,
444 * such that we use the same periodic image,
445 * also when xp has a large displacement.
447 for (int d = 0; d < DIM; d++)
449 sum_wmxp[d] += wm * (dx[d] + xp[ii][d] - x[ii][d]);
455 sum_com->sum_wm = sum_wm;
456 sum_com->sum_wwm = sum_wwm;
457 copy_dvec(sum_wmx, sum_com->sum_wmx);
460 copy_dvec(sum_wmxp, sum_com->sum_wmxp);
464 static void sum_com_part_cosweight(const pull_group_work_t* pgrp,
474 /* Cosine weighting geometry */
483 auto localAtomIndices = pgrp->atomSet.localIndex();
485 for (int i = ind_start; i < ind_end; i++)
487 int ii = localAtomIndices[i];
489 /* Determine cos and sin sums */
490 real cw = std::cos(x[ii][cosdim] * twopi_box);
491 real sw = std::sin(x[ii][cosdim] * twopi_box);
492 sum_cm += static_cast<double>(cw * m);
493 sum_sm += static_cast<double>(sw * m);
494 sum_ccm += static_cast<double>(cw * cw * m);
495 sum_csm += static_cast<double>(cw * sw * m);
496 sum_ssm += static_cast<double>(sw * sw * m);
500 real cw = std::cos(xp[ii][cosdim] * twopi_box);
501 real sw = std::sin(xp[ii][cosdim] * twopi_box);
502 sum_cmp += static_cast<double>(cw * m);
503 sum_smp += static_cast<double>(sw * m);
507 sum_com->sum_cm = sum_cm;
508 sum_com->sum_sm = sum_sm;
509 sum_com->sum_ccm = sum_ccm;
510 sum_com->sum_csm = sum_csm;
511 sum_com->sum_ssm = sum_ssm;
512 sum_com->sum_cmp = sum_cmp;
513 sum_com->sum_smp = sum_smp;
516 /* calculates center of mass of selection index from all coordinates x */
517 void pull_calc_coms(const t_commrec* cr, pull_t* pull, const real* masses, t_pbc* pbc, double t, const rvec x[], rvec* xp)
524 GMX_ASSERT(comm->pbcAtomBuffer.size() == pull->group.size(),
525 "pbcAtomBuffer should have size number of groups");
526 GMX_ASSERT(comm->comBuffer.size() == pull->group.size() * c_comBufferStride,
527 "comBuffer should have size #group*c_comBufferStride");
529 if (pull->bRefAt && pull->bSetPBCatoms)
531 pull_set_pbcatoms(cr, pull, x, comm->pbcAtomBuffer);
533 if (cr != nullptr && DOMAINDECOMP(cr))
535 /* We can keep these PBC reference coordinates fixed for nstlist
536 * steps, since atoms won't jump over PBC.
537 * This avoids a global reduction at the next nstlist-1 steps.
538 * Note that the exact values of the pbc reference coordinates
539 * are irrelevant, as long all atoms in the group are within
540 * half a box distance of the reference coordinate.
542 pull->bSetPBCatoms = FALSE;
546 if (pull->cosdim >= 0)
550 assert(pull->npbcdim <= DIM);
552 for (m = pull->cosdim + 1; m < pull->npbcdim; m++)
554 if (pbc->box[m][pull->cosdim] != 0)
556 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
559 twopi_box = 2.0 * M_PI / pbc->box[pull->cosdim][pull->cosdim];
562 for (size_t g = 0; g < pull->group.size(); g++)
564 pull_group_work_t* pgrp = &pull->group[g];
566 /* Cosine-weighted COMs behave different from all other weighted COMs
567 * in the sense that the weights depend on instantaneous coordinates,
568 * not on pre-set weights. Thus we resize the local weight buffer here.
570 if (pgrp->epgrppbc == epgrppbcCOS)
572 pgrp->localWeights.resize(pgrp->atomSet.localIndex().size());
575 auto comBuffer = gmx::arrayRefFromArray(comm->comBuffer.data() + g * c_comBufferStride,
578 if (pgrp->needToCalcCom)
580 if (pgrp->epgrppbc != epgrppbcCOS)
582 rvec x_pbc = { 0, 0, 0 };
584 switch (pgrp->epgrppbc)
587 /* Set the pbc atom */
588 copy_rvec(comm->pbcAtomBuffer[g], x_pbc);
590 case epgrppbcPREVSTEPCOM:
591 /* Set the pbc reference to the COM of the group of the last step */
592 copy_dvec_to_rvec(pgrp->x_prev_step, comm->pbcAtomBuffer[g]);
593 copy_dvec_to_rvec(pgrp->x_prev_step, x_pbc);
596 /* The final sums should end up in comSums[0] */
597 ComSums& comSumsTotal = pull->comSums[0];
599 /* If we have a single-atom group the mass is irrelevant, so
600 * we can remove the mass factor to avoid division by zero.
601 * Note that with constraint pulling the mass does matter, but
602 * in that case a check group mass != 0 has been done before.
604 if (pgrp->params.ind.size() == 1 && pgrp->atomSet.numAtomsLocal() == 1
605 && masses[pgrp->atomSet.localIndex()[0]] == 0)
607 GMX_ASSERT(xp == nullptr,
608 "We should not have groups with zero mass with constraints, i.e. "
611 /* Copy the single atom coordinate */
612 for (int d = 0; d < DIM; d++)
614 comSumsTotal.sum_wmx[d] = x[pgrp->atomSet.localIndex()[0]][d];
616 /* Set all mass factors to 1 to get the correct COM */
617 comSumsTotal.sum_wm = 1;
618 comSumsTotal.sum_wwm = 1;
620 else if (pgrp->atomSet.numAtomsLocal() <= c_pullMaxNumLocalAtomsSingleThreaded)
622 sum_com_part(pgrp, 0, pgrp->atomSet.numAtomsLocal(), x, xp, masses, pbc, x_pbc, &comSumsTotal);
626 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
627 for (int t = 0; t < pull->nthreads; t++)
629 int ind_start = (pgrp->atomSet.numAtomsLocal() * (t + 0)) / pull->nthreads;
630 int ind_end = (pgrp->atomSet.numAtomsLocal() * (t + 1)) / pull->nthreads;
632 pgrp, ind_start, ind_end, x, xp, masses, pbc, x_pbc, &pull->comSums[t]);
635 /* Reduce the thread contributions to sum_com[0] */
636 for (int t = 1; t < pull->nthreads; t++)
638 comSumsTotal.sum_wm += pull->comSums[t].sum_wm;
639 comSumsTotal.sum_wwm += pull->comSums[t].sum_wwm;
640 dvec_inc(comSumsTotal.sum_wmx, pull->comSums[t].sum_wmx);
641 dvec_inc(comSumsTotal.sum_wmxp, pull->comSums[t].sum_wmxp);
645 if (pgrp->localWeights.empty())
647 comSumsTotal.sum_wwm = comSumsTotal.sum_wm;
650 /* Copy local sums to a buffer for global summing */
651 copy_dvec(comSumsTotal.sum_wmx, comBuffer[0]);
653 copy_dvec(comSumsTotal.sum_wmxp, comBuffer[1]);
655 comBuffer[2][0] = comSumsTotal.sum_wm;
656 comBuffer[2][1] = comSumsTotal.sum_wwm;
661 /* Cosine weighting geometry.
662 * This uses a slab of the system, thus we always have many
663 * atoms in the pull groups. Therefore, always use threads.
665 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
666 for (int t = 0; t < pull->nthreads; t++)
668 int ind_start = (pgrp->atomSet.numAtomsLocal() * (t + 0)) / pull->nthreads;
669 int ind_end = (pgrp->atomSet.numAtomsLocal() * (t + 1)) / pull->nthreads;
670 sum_com_part_cosweight(
671 pgrp, ind_start, ind_end, pull->cosdim, twopi_box, x, xp, masses, &pull->comSums[t]);
674 /* Reduce the thread contributions to comSums[0] */
675 ComSums& comSumsTotal = pull->comSums[0];
676 for (int t = 1; t < pull->nthreads; t++)
678 comSumsTotal.sum_cm += pull->comSums[t].sum_cm;
679 comSumsTotal.sum_sm += pull->comSums[t].sum_sm;
680 comSumsTotal.sum_ccm += pull->comSums[t].sum_ccm;
681 comSumsTotal.sum_csm += pull->comSums[t].sum_csm;
682 comSumsTotal.sum_ssm += pull->comSums[t].sum_ssm;
683 comSumsTotal.sum_cmp += pull->comSums[t].sum_cmp;
684 comSumsTotal.sum_smp += pull->comSums[t].sum_smp;
687 /* Copy local sums to a buffer for global summing */
688 comBuffer[0][0] = comSumsTotal.sum_cm;
689 comBuffer[0][1] = comSumsTotal.sum_sm;
691 comBuffer[1][0] = comSumsTotal.sum_ccm;
692 comBuffer[1][1] = comSumsTotal.sum_csm;
693 comBuffer[1][2] = comSumsTotal.sum_ssm;
694 comBuffer[2][0] = comSumsTotal.sum_cmp;
695 comBuffer[2][1] = comSumsTotal.sum_smp;
701 clear_dvec(comBuffer[0]);
702 clear_dvec(comBuffer[1]);
703 clear_dvec(comBuffer[2]);
709 pull->group.size() * c_comBufferStride * DIM,
710 static_cast<double*>(comm->comBuffer[0]));
712 for (size_t g = 0; g < pull->group.size(); g++)
714 pull_group_work_t* pgrp;
716 pgrp = &pull->group[g];
717 if (pgrp->needToCalcCom)
719 GMX_ASSERT(!pgrp->params.ind.empty(),
720 "Normal pull groups should have atoms, only group 0, which should have "
721 "bCalcCom=FALSE has nat=0");
723 const auto comBuffer = gmx::constArrayRefFromArray(
724 comm->comBuffer.data() + g * c_comBufferStride, c_comBufferStride);
726 if (pgrp->epgrppbc != epgrppbcCOS)
728 double wmass, wwmass;
731 /* Determine the inverse mass */
732 wmass = comBuffer[2][0];
733 wwmass = comBuffer[2][1];
734 pgrp->mwscale = 1.0 / wmass;
735 /* invtm==0 signals a frozen group, so then we should keep it zero */
736 if (pgrp->invtm != 0)
738 pgrp->wscale = wmass / wwmass;
739 pgrp->invtm = wwmass / (wmass * wmass);
741 /* Divide by the total mass */
742 for (m = 0; m < DIM; m++)
744 pgrp->x[m] = comBuffer[0][m] * pgrp->mwscale;
747 pgrp->xp[m] = comBuffer[1][m] * pgrp->mwscale;
749 if (pgrp->epgrppbc == epgrppbcREFAT || pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
751 pgrp->x[m] += comm->pbcAtomBuffer[g][m];
754 pgrp->xp[m] += comm->pbcAtomBuffer[g][m];
761 /* Cosine weighting geometry */
762 double csw, snw, wmass, wwmass;
764 /* Determine the optimal location of the cosine weight */
765 csw = comBuffer[0][0];
766 snw = comBuffer[0][1];
767 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw) / twopi_box;
768 /* Set the weights for the local atoms */
769 wmass = sqrt(csw * csw + snw * snw);
770 wwmass = (comBuffer[1][0] * csw * csw + comBuffer[1][1] * csw * snw
771 + comBuffer[1][2] * snw * snw)
774 pgrp->mwscale = 1.0 / wmass;
775 pgrp->wscale = wmass / wwmass;
776 pgrp->invtm = wwmass / (wmass * wmass);
777 /* Set the weights for the local atoms */
780 for (size_t i = 0; i < pgrp->atomSet.numAtomsLocal(); i++)
782 int ii = pgrp->atomSet.localIndex()[i];
783 pgrp->localWeights[i] = csw * std::cos(twopi_box * x[ii][pull->cosdim])
784 + snw * std::sin(twopi_box * x[ii][pull->cosdim]);
788 csw = comBuffer[2][0];
789 snw = comBuffer[2][1];
790 pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw) / twopi_box;
795 fprintf(debug, "Pull group %zu wmass %f invtm %f\n", g, 1.0 / pgrp->mwscale, pgrp->invtm);
802 /* Calculate the COMs for the cyclinder reference groups */
803 make_cyl_refgrps(cr, pull, masses, pbc, t, x);
807 using BoolVec = gmx::BasicVector<bool>;
809 /* Returns whether the pull group obeys the PBC restrictions */
810 static bool pullGroupObeysPbcRestrictions(const pull_group_work_t& group,
811 const BoolVec& dimUsed,
814 const gmx::RVec& x_pbc,
815 const real pbcMargin)
817 /* Determine which dimensions are relevant for PBC */
818 BoolVec dimUsesPbc = { false, false, false };
819 bool pbcIsRectangular = true;
820 for (int d = 0; d < pbc.ndim_ePBC; d++)
824 dimUsesPbc[d] = true;
825 /* All non-zero dimensions of vector v are involved in PBC */
826 for (int d2 = d + 1; d2 < pbc.ndim_ePBC; d2++)
829 if (pbc.box[d2][d] != 0)
831 dimUsesPbc[d2] = true;
832 pbcIsRectangular = false;
838 rvec marginPerDim = {};
839 real marginDistance2 = 0;
840 if (pbcIsRectangular)
842 /* Use margins for dimensions independently */
843 for (int d = 0; d < pbc.ndim_ePBC; d++)
845 marginPerDim[d] = pbcMargin * pbc.hbox_diag[d];
850 /* Check the total distance along the relevant dimensions */
851 for (int d = 0; d < pbc.ndim_ePBC; d++)
855 marginDistance2 += pbcMargin * gmx::square(0.5) * norm2(pbc.box[d]);
860 auto localAtomIndices = group.atomSet.localIndex();
861 for (gmx::index indexInSet = 0; indexInSet < localAtomIndices.ssize(); indexInSet++)
864 pbc_dx(&pbc, x[localAtomIndices[indexInSet]], x_pbc, dx);
866 bool atomIsTooFar = false;
867 if (pbcIsRectangular)
869 for (int d = 0; d < pbc.ndim_ePBC; d++)
871 if (dimUsesPbc[d] && (dx[d] < -marginPerDim[d] || dx[d] > marginPerDim[d]))
879 real pbcDistance2 = 0;
880 for (int d = 0; d < pbc.ndim_ePBC; d++)
884 pbcDistance2 += gmx::square(dx[d]);
887 atomIsTooFar = (pbcDistance2 > marginDistance2);
898 int pullCheckPbcWithinGroups(const pull_t& pull, const rvec* x, const t_pbc& pbc, real pbcMargin)
900 if (pbc.pbcType == PbcType::No)
905 /* Determine what dimensions are used for each group by pull coordinates */
906 std::vector<BoolVec> dimUsed(pull.group.size(), { false, false, false });
907 for (size_t c = 0; c < pull.coord.size(); c++)
909 const t_pull_coord& coordParams = pull.coord[c].params;
910 for (int groupIndex = 0; groupIndex < coordParams.ngroup; groupIndex++)
912 for (int d = 0; d < DIM; d++)
914 if (coordParams.dim[d] && !(coordParams.eGeom == epullgCYL && groupIndex == 0))
916 dimUsed[coordParams.group[groupIndex]][d] = true;
922 /* Check PBC for every group that uses a PBC reference atom treatment */
923 for (size_t g = 0; g < pull.group.size(); g++)
925 const pull_group_work_t& group = pull.group[g];
926 if ((group.epgrppbc == epgrppbcREFAT || group.epgrppbc == epgrppbcPREVSTEPCOM)
927 && !pullGroupObeysPbcRestrictions(group, dimUsed[g], x, pbc, pull.comm.pbcAtomBuffer[g], pbcMargin))
936 bool pullCheckPbcWithinGroup(const pull_t& pull,
937 gmx::ArrayRef<const gmx::RVec> x,
942 if (pbc.pbcType == PbcType::No)
946 GMX_ASSERT(groupNr < gmx::ssize(pull.group), "groupNr is out of range");
948 /* Check PBC if the group uses a PBC reference atom treatment. */
949 const pull_group_work_t& group = pull.group[groupNr];
950 if (group.epgrppbc != epgrppbcREFAT && group.epgrppbc != epgrppbcPREVSTEPCOM)
955 /* Determine what dimensions are used for each group by pull coordinates */
956 BoolVec dimUsed = { false, false, false };
957 for (size_t c = 0; c < pull.coord.size(); c++)
959 const t_pull_coord& coordParams = pull.coord[c].params;
960 for (int groupIndex = 0; groupIndex < coordParams.ngroup; groupIndex++)
962 if (coordParams.group[groupIndex] == groupNr)
964 for (int d = 0; d < DIM; d++)
966 if (coordParams.dim[d] && !(coordParams.eGeom == epullgCYL && groupIndex == 0))
975 return (pullGroupObeysPbcRestrictions(
976 group, dimUsed, as_rvec_array(x.data()), pbc, pull.comm.pbcAtomBuffer[groupNr], pbcMargin));
979 void setPrevStepPullComFromState(struct pull_t* pull, const t_state* state)
981 for (size_t g = 0; g < pull->group.size(); g++)
983 for (int j = 0; j < DIM; j++)
985 pull->group[g].x_prev_step[j] = state->pull_com_prev_step[g * DIM + j];
990 void updatePrevStepPullCom(struct pull_t* pull, t_state* state)
992 for (size_t g = 0; g < pull->group.size(); g++)
994 if (pull->group[g].needToCalcCom)
996 for (int j = 0; j < DIM; j++)
998 pull->group[g].x_prev_step[j] = pull->group[g].x[j];
999 state->pull_com_prev_step[g * DIM + j] = pull->group[g].x[j];
1005 void allocStatePrevStepPullCom(t_state* state, const pull_t* pull)
1009 state->pull_com_prev_step.clear();
1012 size_t ngroup = pull->group.size();
1013 if (state->pull_com_prev_step.size() / DIM != ngroup)
1015 state->pull_com_prev_step.resize(ngroup * DIM, NAN);
1019 void initPullComFromPrevStep(const t_commrec* cr, pull_t* pull, const real* masses, t_pbc* pbc, const rvec x[])
1021 pull_comm_t* comm = &pull->comm;
1022 size_t ngroup = pull->group.size();
1024 if (!comm->bParticipate)
1029 GMX_ASSERT(comm->pbcAtomBuffer.size() == pull->group.size(),
1030 "pbcAtomBuffer should have size number of groups");
1031 GMX_ASSERT(comm->comBuffer.size() == pull->group.size() * c_comBufferStride,
1032 "comBuffer should have size #group*c_comBufferStride");
1034 pull_set_pbcatoms(cr, pull, x, comm->pbcAtomBuffer);
1036 for (size_t g = 0; g < ngroup; g++)
1038 pull_group_work_t* pgrp;
1040 pgrp = &pull->group[g];
1042 if (pgrp->needToCalcCom && pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
1044 GMX_ASSERT(pgrp->params.ind.size() > 1,
1045 "Groups with no atoms, or only one atom, should not "
1046 "use the COM from the previous step as reference.");
1048 rvec x_pbc = { 0, 0, 0 };
1049 copy_rvec(comm->pbcAtomBuffer[g], x_pbc);
1053 fprintf(debug, "Initialising prev step COM of pull group %zu. x_pbc =", g);
1054 for (int m = 0; m < DIM; m++)
1056 fprintf(debug, " %f", x_pbc[m]);
1058 fprintf(debug, "\n");
1061 /* The following is to a large extent similar to pull_calc_coms() */
1063 /* The final sums should end up in sum_com[0] */
1064 ComSums& comSumsTotal = pull->comSums[0];
1066 if (pgrp->atomSet.numAtomsLocal() <= c_pullMaxNumLocalAtomsSingleThreaded)
1068 sum_com_part(pgrp, 0, pgrp->atomSet.numAtomsLocal(), x, nullptr, masses, pbc, x_pbc, &comSumsTotal);
1072 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
1073 for (int t = 0; t < pull->nthreads; t++)
1075 int ind_start = (pgrp->atomSet.numAtomsLocal() * (t + 0)) / pull->nthreads;
1076 int ind_end = (pgrp->atomSet.numAtomsLocal() * (t + 1)) / pull->nthreads;
1078 pgrp, ind_start, ind_end, x, nullptr, masses, pbc, x_pbc, &pull->comSums[t]);
1081 /* Reduce the thread contributions to sum_com[0] */
1082 for (int t = 1; t < pull->nthreads; t++)
1084 comSumsTotal.sum_wm += pull->comSums[t].sum_wm;
1085 comSumsTotal.sum_wwm += pull->comSums[t].sum_wwm;
1086 dvec_inc(comSumsTotal.sum_wmx, pull->comSums[t].sum_wmx);
1087 dvec_inc(comSumsTotal.sum_wmxp, pull->comSums[t].sum_wmxp);
1091 if (pgrp->localWeights.empty())
1093 comSumsTotal.sum_wwm = comSumsTotal.sum_wm;
1096 /* Copy local sums to a buffer for global summing */
1097 auto localSums = gmx::arrayRefFromArray(comm->comBuffer.data() + g * c_comBufferStride,
1100 localSums[0] = comSumsTotal.sum_wmx;
1101 localSums[1] = comSumsTotal.sum_wmxp;
1102 localSums[2][0] = comSumsTotal.sum_wm;
1103 localSums[2][1] = comSumsTotal.sum_wwm;
1104 localSums[2][2] = 0;
1108 pullAllReduce(cr, comm, ngroup * c_comBufferStride * DIM, static_cast<double*>(comm->comBuffer[0]));
1110 for (size_t g = 0; g < ngroup; g++)
1112 pull_group_work_t* pgrp;
1114 pgrp = &pull->group[g];
1115 if (pgrp->needToCalcCom)
1117 if (pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
1119 auto localSums = gmx::arrayRefFromArray(
1120 comm->comBuffer.data() + g * c_comBufferStride, c_comBufferStride);
1121 double wmass, wwmass;
1123 /* Determine the inverse mass */
1124 wmass = localSums[2][0];
1125 wwmass = localSums[2][1];
1126 pgrp->mwscale = 1.0 / wmass;
1127 /* invtm==0 signals a frozen group, so then we should keep it zero */
1128 if (pgrp->invtm != 0)
1130 pgrp->wscale = wmass / wwmass;
1131 pgrp->invtm = wwmass / (wmass * wmass);
1133 /* Divide by the total mass */
1134 for (int m = 0; m < DIM; m++)
1136 pgrp->x[m] = localSums[0][m] * pgrp->mwscale;
1137 pgrp->x[m] += comm->pbcAtomBuffer[g][m];
1141 fprintf(debug, "Pull group %zu wmass %f invtm %f\n", g, 1.0 / pgrp->mwscale, pgrp->invtm);
1142 fprintf(debug, "Initialising prev step COM of pull group %zu to", g);
1143 for (int m = 0; m < DIM; m++)
1145 fprintf(debug, " %f", pgrp->x[m]);
1147 fprintf(debug, "\n");
1149 copy_dvec(pgrp->x, pgrp->x_prev_step);