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/units.h"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/mdtypes/commrec.h"
53 #include "gromacs/mdtypes/inputrec.h"
54 #include "gromacs/mdtypes/md_enums.h"
55 #include "gromacs/mdtypes/state.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/pulling/pull.h"
58 #include "gromacs/utility/arrayref.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/futil.h"
61 #include "gromacs/utility/gmxassert.h"
62 #include "gromacs/utility/real.h"
63 #include "gromacs/utility/smalloc.h"
65 #include "pull_internal.h"
69 // Helper function to deduce MPI datatype from the type of data
70 gmx_unused static MPI_Datatype mpiDatatype(const float gmx_unused* data)
75 // Helper function to deduce MPI datatype from the type of data
76 gmx_unused static MPI_Datatype mpiDatatype(const double gmx_unused* data)
84 // Helper function; note that gmx_sum(d) should actually be templated
85 gmx_unused static void gmxAllReduce(int n, real* data, const t_commrec* cr)
91 // Helper function; note that gmx_sum(d) should actually be templated
92 gmx_unused static void gmxAllReduce(int n, double* data, const t_commrec* cr)
94 gmx_sumd(n, data, cr);
97 // Reduce data of n elements over all ranks currently participating in pull
99 static void pullAllReduce(const t_commrec* cr, pull_comm_t* comm, int n, T* data)
101 if (cr != nullptr && PAR(cr))
103 if (comm->bParticipateAll)
105 /* Sum the contributions over all DD ranks */
106 gmxAllReduce(n, data, cr);
110 /* Separate branch because gmx_sum uses cr->mpi_comm_mygroup */
112 MPI_Allreduce(MPI_IN_PLACE, data, n, mpiDatatype(data), MPI_SUM, comm->mpi_comm_com);
114 gmx_incons("comm->bParticipateAll=FALSE without GMX_MPI");
120 /* Copies the coordinates of the PBC atom of pgrp to x_pbc.
121 * When those coordinates are not available on this rank, clears x_pbc.
123 static void setPbcAtomCoords(const pull_group_work_t& pgrp, gmx::ArrayRef<const gmx::RVec> x, rvec x_pbc)
125 if (pgrp.pbcAtomSet != nullptr)
127 if (pgrp.pbcAtomSet->numAtomsLocal() > 0)
129 /* We have the atom locally, copy its coordinates */
130 copy_rvec(x[pgrp.pbcAtomSet->localIndex()[0]], x_pbc);
134 /* Another rank has it, clear the coordinates for MPI_Allreduce */
140 copy_rvec(x[pgrp.params.pbcatom], x_pbc);
144 static void pull_set_pbcatoms(const t_commrec* cr,
146 gmx::ArrayRef<const gmx::RVec> x,
147 gmx::ArrayRef<gmx::RVec> x_pbc)
150 for (size_t g = 0; g < pull->group.size(); g++)
152 const pull_group_work_t& group = pull->group[g];
153 if (group.needToCalcCom && (group.epgrppbc == epgrppbcREFAT || group.epgrppbc == epgrppbcPREVSTEPCOM))
155 setPbcAtomCoords(pull->group[g], x, x_pbc[g]);
160 clear_rvec(x_pbc[g]);
164 if (cr && PAR(cr) && numPbcAtoms > 0)
166 /* Sum over participating ranks to get x_pbc from the home ranks.
167 * This can be very expensive at high parallelization, so we only
168 * do this after each DD repartitioning.
170 pullAllReduce(cr, &pull->comm, pull->group.size() * DIM, static_cast<real*>(x_pbc[0]));
174 static void make_cyl_refgrps(const t_commrec* cr,
176 gmx::ArrayRef<const real> masses,
179 gmx::ArrayRef<const gmx::RVec> x)
181 pull_comm_t* comm = &pull->comm;
183 GMX_ASSERT(comm->cylinderBuffer.size() == pull->coord.size() * c_cylinderBufferStride,
184 "cylinderBuffer should have the correct size");
186 double inv_cyl_r2 = 1.0 / gmx::square(pull->params.cylinder_r);
188 /* loop over all groups to make a reference group for each*/
189 for (size_t c = 0; c < pull->coord.size(); c++)
191 pull_coord_work_t* pcrd;
192 double sum_a, wmass, wwmass;
193 dvec radf_fac0, radf_fac1;
195 pcrd = &pull->coord[c];
200 clear_dvec(radf_fac0);
201 clear_dvec(radf_fac1);
203 if (pcrd->params.eGeom == PullGroupGeometry::Cylinder)
205 /* pref will be the same group for all pull coordinates */
206 const pull_group_work_t& pref = pull->group[pcrd->params.group[0]];
207 const pull_group_work_t& pgrp = pull->group[pcrd->params.group[1]];
208 pull_group_work_t& pdyna = pull->dyna[c];
210 copy_dvec_to_rvec(pcrd->spatialData.vec, direction);
212 /* Since we have not calculated the COM of the cylinder group yet,
213 * we calculate distances with respect to location of the pull
214 * group minus the reference position along the vector.
215 * here we already have the COM of the pull group. This resolves
216 * any PBC issues and we don't need to use a PBC-atom here.
218 if (pcrd->params.rate != 0)
220 /* With rate=0, value_ref is set initially */
221 pcrd->value_ref = pcrd->params.init + pcrd->params.rate * t;
224 for (int m = 0; m < DIM; m++)
226 reference[m] = pgrp.x[m] - pcrd->spatialData.vec[m] * pcrd->value_ref;
229 auto localAtomIndices = pref.atomSet.localIndex();
231 /* This actually only needs to be done at init or DD time,
232 * but resizing with the same size does not cause much overhead.
234 pdyna.localWeights.resize(localAtomIndices.size());
235 pdyna.mdw.resize(localAtomIndices.size());
236 pdyna.dv.resize(localAtomIndices.size());
238 /* loop over all atoms in the main ref group */
239 for (gmx::index indexInSet = 0; indexInSet < localAtomIndices.ssize(); indexInSet++)
241 int atomIndex = localAtomIndices[indexInSet];
243 pbc_dx_aiuc(pbc, x[atomIndex], reference, dx);
244 double axialLocation = iprod(direction, dx);
247 for (int m = 0; m < DIM; m++)
249 /* Determine the radial components */
250 radialLocation[m] = dx[m] - axialLocation * direction[m];
251 dr2 += gmx::square(radialLocation[m]);
253 double dr2_rel = dr2 * inv_cyl_r2;
257 /* add atom to sum of COM and to weight array */
259 double mass = masses[atomIndex];
260 /* The radial weight function is 1-2x^2+x^4,
261 * where x=r/cylinder_r. Since this function depends
262 * on the radial component, we also get radial forces
265 double weight = 1 + (-2 + dr2_rel) * dr2_rel;
266 double dweight_r = (-4 + 4 * dr2_rel) * inv_cyl_r2;
267 pdyna.localWeights[indexInSet] = weight;
268 sum_a += mass * weight * axialLocation;
269 wmass += mass * weight;
270 wwmass += mass * weight * weight;
272 dsvmul(mass * dweight_r, radialLocation, mdw);
273 copy_dvec(mdw, pdyna.mdw[indexInSet]);
274 /* Currently we only have the axial component of the
275 * offset from the cylinder COM up to an unkown offset.
276 * We add this offset after the reduction needed
277 * for determining the COM of the cylinder group.
279 pdyna.dv[indexInSet] = axialLocation;
280 for (int m = 0; m < DIM; m++)
282 radf_fac0[m] += mdw[m];
283 radf_fac1[m] += mdw[m] * axialLocation;
288 pdyna.localWeights[indexInSet] = 0;
293 auto buffer = gmx::arrayRefFromArray(
294 comm->cylinderBuffer.data() + c * c_cylinderBufferStride, c_cylinderBufferStride);
300 buffer[3] = radf_fac0[XX];
301 buffer[4] = radf_fac0[YY];
302 buffer[5] = radf_fac0[ZZ];
304 buffer[6] = radf_fac1[XX];
305 buffer[7] = radf_fac1[YY];
306 buffer[8] = radf_fac1[ZZ];
309 if (cr != nullptr && PAR(cr))
311 /* Sum the contributions over the ranks */
312 pullAllReduce(cr, comm, pull->coord.size() * c_cylinderBufferStride, comm->cylinderBuffer.data());
315 for (size_t c = 0; c < pull->coord.size(); c++)
317 pull_coord_work_t* pcrd;
319 pcrd = &pull->coord[c];
321 if (pcrd->params.eGeom == PullGroupGeometry::Cylinder)
323 pull_group_work_t* pdyna = &pull->dyna[c];
324 pull_group_work_t* pgrp = &pull->group[pcrd->params.group[1]];
325 PullCoordSpatialData& spatialData = pcrd->spatialData;
327 auto buffer = gmx::constArrayRefFromArray(
328 comm->cylinderBuffer.data() + c * c_cylinderBufferStride, c_cylinderBufferStride);
329 double wmass = buffer[0];
330 double wwmass = buffer[1];
331 pdyna->mwscale = 1.0 / wmass;
332 /* Cylinder pulling can't be used with constraints, but we set
333 * wscale and invtm anyhow, in case someone would like to use them.
335 pdyna->wscale = wmass / wwmass;
336 pdyna->invtm = wwmass / (wmass * wmass);
338 /* We store the deviation of the COM from the reference location
339 * used above, since we need it when we apply the radial forces
340 * to the atoms in the cylinder group.
342 spatialData.cyl_dev = 0;
343 for (int m = 0; m < DIM; m++)
345 double reference = pgrp->x[m] - spatialData.vec[m] * pcrd->value_ref;
346 double dist = -spatialData.vec[m] * buffer[2] * pdyna->mwscale;
347 pdyna->x[m] = reference - dist;
348 spatialData.cyl_dev += dist;
350 /* Now we know the exact COM of the cylinder reference group,
351 * we can determine the radial force factor (ffrad) that when
352 * multiplied with the axial pull force will give the radial
353 * force on the pulled (non-cylinder) group.
355 for (int m = 0; m < DIM; m++)
357 spatialData.ffrad[m] = (buffer[6 + m] + buffer[3 + m] * spatialData.cyl_dev) / wmass;
363 "Pull cylinder group %zu:%8.3f%8.3f%8.3f m:%8.3f\n",
370 "ffrad %8.3f %8.3f %8.3f\n",
371 spatialData.ffrad[XX],
372 spatialData.ffrad[YY],
373 spatialData.ffrad[ZZ]);
379 static double atan2_0_2pi(double y, double x)
391 static void sum_com_part(const pull_group_work_t* pgrp,
394 gmx::ArrayRef<const gmx::RVec> x,
395 gmx::ArrayRef<const gmx::RVec> xp,
396 gmx::ArrayRef<const real> mass,
403 dvec sum_wmx = { 0, 0, 0 };
404 dvec sum_wmxp = { 0, 0, 0 };
406 auto localAtomIndices = pgrp->atomSet.localIndex();
407 for (int i = ind_start; i < ind_end; i++)
409 int ii = localAtomIndices[i];
411 if (pgrp->localWeights.empty())
420 w = pgrp->localWeights[i];
425 if (pgrp->epgrppbc == epgrppbcNONE)
427 /* Plain COM: sum the coordinates */
428 for (int d = 0; d < DIM; d++)
430 sum_wmx[d] += wm * x[ii][d];
434 for (int d = 0; d < DIM; d++)
436 sum_wmxp[d] += wm * xp[ii][d];
444 /* Sum the difference with the reference atom */
445 pbc_dx(pbc, x[ii], x_pbc, dx);
446 for (int d = 0; d < DIM; d++)
448 sum_wmx[d] += wm * dx[d];
452 /* For xp add the difference between xp and x to dx,
453 * such that we use the same periodic image,
454 * also when xp has a large displacement.
456 for (int d = 0; d < DIM; d++)
458 sum_wmxp[d] += wm * (dx[d] + xp[ii][d] - x[ii][d]);
464 sum_com->sum_wm = sum_wm;
465 sum_com->sum_wwm = sum_wwm;
466 copy_dvec(sum_wmx, sum_com->sum_wmx);
469 copy_dvec(sum_wmxp, sum_com->sum_wmxp);
473 static void sum_com_part_cosweight(const pull_group_work_t* pgrp,
478 gmx::ArrayRef<const gmx::RVec> x,
479 gmx::ArrayRef<const gmx::RVec> xp,
480 gmx::ArrayRef<const real> mass,
483 /* Cosine weighting geometry */
492 auto localAtomIndices = pgrp->atomSet.localIndex();
494 for (int i = ind_start; i < ind_end; i++)
496 int ii = localAtomIndices[i];
498 /* Determine cos and sin sums */
499 real cw = std::cos(x[ii][cosdim] * twopi_box);
500 real sw = std::sin(x[ii][cosdim] * twopi_box);
501 sum_cm += static_cast<double>(cw * m);
502 sum_sm += static_cast<double>(sw * m);
503 sum_ccm += static_cast<double>(cw * cw * m);
504 sum_csm += static_cast<double>(cw * sw * m);
505 sum_ssm += static_cast<double>(sw * sw * m);
509 real cw = std::cos(xp[ii][cosdim] * twopi_box);
510 real sw = std::sin(xp[ii][cosdim] * twopi_box);
511 sum_cmp += static_cast<double>(cw * m);
512 sum_smp += static_cast<double>(sw * m);
516 sum_com->sum_cm = sum_cm;
517 sum_com->sum_sm = sum_sm;
518 sum_com->sum_ccm = sum_ccm;
519 sum_com->sum_csm = sum_csm;
520 sum_com->sum_ssm = sum_ssm;
521 sum_com->sum_cmp = sum_cmp;
522 sum_com->sum_smp = sum_smp;
525 /* calculates center of mass of selection index from all coordinates x */
526 void pull_calc_coms(const t_commrec* cr,
528 gmx::ArrayRef<const real> masses,
531 gmx::ArrayRef<const gmx::RVec> x,
532 gmx::ArrayRef<gmx::RVec> xp)
539 GMX_ASSERT(comm->pbcAtomBuffer.size() == pull->group.size(),
540 "pbcAtomBuffer should have size number of groups");
541 GMX_ASSERT(comm->comBuffer.size() == pull->group.size() * c_comBufferStride,
542 "comBuffer should have size #group*c_comBufferStride");
544 if (pull->bRefAt && pull->bSetPBCatoms)
546 pull_set_pbcatoms(cr, pull, x, comm->pbcAtomBuffer);
548 if (cr != nullptr && DOMAINDECOMP(cr))
550 /* We can keep these PBC reference coordinates fixed for nstlist
551 * steps, since atoms won't jump over PBC.
552 * This avoids a global reduction at the next nstlist-1 steps.
553 * Note that the exact values of the pbc reference coordinates
554 * are irrelevant, as long all atoms in the group are within
555 * half a box distance of the reference coordinate.
557 pull->bSetPBCatoms = FALSE;
561 if (pull->cosdim >= 0)
565 assert(pull->npbcdim <= DIM);
567 for (m = pull->cosdim + 1; m < pull->npbcdim; m++)
569 if (pbc->box[m][pull->cosdim] != 0)
571 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
574 twopi_box = 2.0 * M_PI / pbc->box[pull->cosdim][pull->cosdim];
577 for (size_t g = 0; g < pull->group.size(); g++)
579 pull_group_work_t* pgrp = &pull->group[g];
581 /* Cosine-weighted COMs behave different from all other weighted COMs
582 * in the sense that the weights depend on instantaneous coordinates,
583 * not on pre-set weights. Thus we resize the local weight buffer here.
585 if (pgrp->epgrppbc == epgrppbcCOS)
587 pgrp->localWeights.resize(pgrp->atomSet.localIndex().size());
590 auto comBuffer = gmx::arrayRefFromArray(comm->comBuffer.data() + g * c_comBufferStride,
593 if (pgrp->needToCalcCom)
595 if (pgrp->epgrppbc != epgrppbcCOS)
597 rvec x_pbc = { 0, 0, 0 };
599 switch (pgrp->epgrppbc)
602 /* Set the pbc atom */
603 copy_rvec(comm->pbcAtomBuffer[g], x_pbc);
605 case epgrppbcPREVSTEPCOM:
606 /* Set the pbc reference to the COM of the group of the last step */
607 copy_dvec_to_rvec(pgrp->x_prev_step, comm->pbcAtomBuffer[g]);
608 copy_dvec_to_rvec(pgrp->x_prev_step, x_pbc);
611 /* The final sums should end up in comSums[0] */
612 ComSums& comSumsTotal = pull->comSums[0];
614 /* If we have a single-atom group the mass is irrelevant, so
615 * we can remove the mass factor to avoid division by zero.
616 * Note that with constraint pulling the mass does matter, but
617 * in that case a check group mass != 0 has been done before.
619 if (pgrp->params.ind.size() == 1 && pgrp->atomSet.numAtomsLocal() == 1
620 && masses[pgrp->atomSet.localIndex()[0]] == 0)
622 GMX_ASSERT(xp.empty(),
623 "We should not have groups with zero mass with constraints, i.e. "
626 /* Copy the single atom coordinate */
627 for (int d = 0; d < DIM; d++)
629 comSumsTotal.sum_wmx[d] = x[pgrp->atomSet.localIndex()[0]][d];
631 /* Set all mass factors to 1 to get the correct COM */
632 comSumsTotal.sum_wm = 1;
633 comSumsTotal.sum_wwm = 1;
635 else if (pgrp->atomSet.numAtomsLocal() <= c_pullMaxNumLocalAtomsSingleThreaded)
637 sum_com_part(pgrp, 0, pgrp->atomSet.numAtomsLocal(), x, xp, masses, pbc, x_pbc, &comSumsTotal);
641 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
642 for (int t = 0; t < pull->nthreads; t++)
644 int ind_start = (pgrp->atomSet.numAtomsLocal() * (t + 0)) / pull->nthreads;
645 int ind_end = (pgrp->atomSet.numAtomsLocal() * (t + 1)) / pull->nthreads;
647 pgrp, ind_start, ind_end, x, xp, masses, pbc, x_pbc, &pull->comSums[t]);
650 /* Reduce the thread contributions to sum_com[0] */
651 for (int t = 1; t < pull->nthreads; t++)
653 comSumsTotal.sum_wm += pull->comSums[t].sum_wm;
654 comSumsTotal.sum_wwm += pull->comSums[t].sum_wwm;
655 dvec_inc(comSumsTotal.sum_wmx, pull->comSums[t].sum_wmx);
656 dvec_inc(comSumsTotal.sum_wmxp, pull->comSums[t].sum_wmxp);
660 if (pgrp->localWeights.empty())
662 comSumsTotal.sum_wwm = comSumsTotal.sum_wm;
665 /* Copy local sums to a buffer for global summing */
666 copy_dvec(comSumsTotal.sum_wmx, comBuffer[0]);
668 copy_dvec(comSumsTotal.sum_wmxp, comBuffer[1]);
670 comBuffer[2][0] = comSumsTotal.sum_wm;
671 comBuffer[2][1] = comSumsTotal.sum_wwm;
676 /* Cosine weighting geometry.
677 * This uses a slab of the system, thus we always have many
678 * atoms in the pull groups. Therefore, always use threads.
680 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
681 for (int t = 0; t < pull->nthreads; t++)
683 int ind_start = (pgrp->atomSet.numAtomsLocal() * (t + 0)) / pull->nthreads;
684 int ind_end = (pgrp->atomSet.numAtomsLocal() * (t + 1)) / pull->nthreads;
685 sum_com_part_cosweight(
686 pgrp, ind_start, ind_end, pull->cosdim, twopi_box, x, xp, masses, &pull->comSums[t]);
689 /* Reduce the thread contributions to comSums[0] */
690 ComSums& comSumsTotal = pull->comSums[0];
691 for (int t = 1; t < pull->nthreads; t++)
693 comSumsTotal.sum_cm += pull->comSums[t].sum_cm;
694 comSumsTotal.sum_sm += pull->comSums[t].sum_sm;
695 comSumsTotal.sum_ccm += pull->comSums[t].sum_ccm;
696 comSumsTotal.sum_csm += pull->comSums[t].sum_csm;
697 comSumsTotal.sum_ssm += pull->comSums[t].sum_ssm;
698 comSumsTotal.sum_cmp += pull->comSums[t].sum_cmp;
699 comSumsTotal.sum_smp += pull->comSums[t].sum_smp;
702 /* Copy local sums to a buffer for global summing */
703 comBuffer[0][0] = comSumsTotal.sum_cm;
704 comBuffer[0][1] = comSumsTotal.sum_sm;
706 comBuffer[1][0] = comSumsTotal.sum_ccm;
707 comBuffer[1][1] = comSumsTotal.sum_csm;
708 comBuffer[1][2] = comSumsTotal.sum_ssm;
709 comBuffer[2][0] = comSumsTotal.sum_cmp;
710 comBuffer[2][1] = comSumsTotal.sum_smp;
716 clear_dvec(comBuffer[0]);
717 clear_dvec(comBuffer[1]);
718 clear_dvec(comBuffer[2]);
724 pull->group.size() * c_comBufferStride * DIM,
725 static_cast<double*>(comm->comBuffer[0]));
727 for (size_t g = 0; g < pull->group.size(); g++)
729 pull_group_work_t* pgrp;
731 pgrp = &pull->group[g];
732 if (pgrp->needToCalcCom)
734 GMX_ASSERT(!pgrp->params.ind.empty(),
735 "Normal pull groups should have atoms, only group 0, which should have "
736 "bCalcCom=FALSE has nat=0");
738 const auto comBuffer = gmx::constArrayRefFromArray(
739 comm->comBuffer.data() + g * c_comBufferStride, c_comBufferStride);
741 if (pgrp->epgrppbc != epgrppbcCOS)
743 double wmass, wwmass;
746 /* Determine the inverse mass */
747 wmass = comBuffer[2][0];
748 wwmass = comBuffer[2][1];
749 pgrp->mwscale = 1.0 / wmass;
750 /* invtm==0 signals a frozen group, so then we should keep it zero */
751 if (pgrp->invtm != 0)
753 pgrp->wscale = wmass / wwmass;
754 pgrp->invtm = wwmass / (wmass * wmass);
756 /* Divide by the total mass */
757 for (m = 0; m < DIM; m++)
759 pgrp->x[m] = comBuffer[0][m] * pgrp->mwscale;
762 pgrp->xp[m] = comBuffer[1][m] * pgrp->mwscale;
764 if (pgrp->epgrppbc == epgrppbcREFAT || pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
766 pgrp->x[m] += comm->pbcAtomBuffer[g][m];
769 pgrp->xp[m] += comm->pbcAtomBuffer[g][m];
776 /* Cosine weighting geometry */
777 double csw, snw, wmass, wwmass;
779 /* Determine the optimal location of the cosine weight */
780 csw = comBuffer[0][0];
781 snw = comBuffer[0][1];
782 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw) / twopi_box;
783 /* Set the weights for the local atoms */
784 wmass = sqrt(csw * csw + snw * snw);
785 wwmass = (comBuffer[1][0] * csw * csw + comBuffer[1][1] * csw * snw
786 + comBuffer[1][2] * snw * snw)
789 pgrp->mwscale = 1.0 / wmass;
790 pgrp->wscale = wmass / wwmass;
791 pgrp->invtm = wwmass / (wmass * wmass);
792 /* Set the weights for the local atoms */
795 for (size_t i = 0; i < pgrp->atomSet.numAtomsLocal(); i++)
797 int ii = pgrp->atomSet.localIndex()[i];
798 pgrp->localWeights[i] = csw * std::cos(twopi_box * x[ii][pull->cosdim])
799 + snw * std::sin(twopi_box * x[ii][pull->cosdim]);
803 csw = comBuffer[2][0];
804 snw = comBuffer[2][1];
805 pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw) / twopi_box;
810 fprintf(debug, "Pull group %zu wmass %f invtm %f\n", g, 1.0 / pgrp->mwscale, pgrp->invtm);
817 /* Calculate the COMs for the cyclinder reference groups */
818 make_cyl_refgrps(cr, pull, masses, pbc, t, x);
822 using BoolVec = gmx::BasicVector<bool>;
824 /* Returns whether the pull group obeys the PBC restrictions */
825 static bool pullGroupObeysPbcRestrictions(const pull_group_work_t& group,
826 const BoolVec& dimUsed,
827 gmx::ArrayRef<const gmx::RVec> x,
829 const gmx::RVec& x_pbc,
830 const real pbcMargin)
832 /* Determine which dimensions are relevant for PBC */
833 BoolVec dimUsesPbc = { false, false, false };
834 bool pbcIsRectangular = true;
835 for (int d = 0; d < pbc.ndim_ePBC; d++)
839 dimUsesPbc[d] = true;
840 /* All non-zero dimensions of vector v are involved in PBC */
841 for (int d2 = d + 1; d2 < pbc.ndim_ePBC; d2++)
844 if (pbc.box[d2][d] != 0)
846 dimUsesPbc[d2] = true;
847 pbcIsRectangular = false;
853 rvec marginPerDim = {};
854 real marginDistance2 = 0;
855 if (pbcIsRectangular)
857 /* Use margins for dimensions independently */
858 for (int d = 0; d < pbc.ndim_ePBC; d++)
860 marginPerDim[d] = pbcMargin * pbc.hbox_diag[d];
865 /* Check the total distance along the relevant dimensions */
866 for (int d = 0; d < pbc.ndim_ePBC; d++)
870 marginDistance2 += pbcMargin * gmx::square(0.5) * norm2(pbc.box[d]);
875 auto localAtomIndices = group.atomSet.localIndex();
876 for (gmx::index indexInSet = 0; indexInSet < localAtomIndices.ssize(); indexInSet++)
879 pbc_dx(&pbc, x[localAtomIndices[indexInSet]], x_pbc, dx);
881 bool atomIsTooFar = false;
882 if (pbcIsRectangular)
884 for (int d = 0; d < pbc.ndim_ePBC; d++)
886 if (dimUsesPbc[d] && (dx[d] < -marginPerDim[d] || dx[d] > marginPerDim[d]))
894 real pbcDistance2 = 0;
895 for (int d = 0; d < pbc.ndim_ePBC; d++)
899 pbcDistance2 += gmx::square(dx[d]);
902 atomIsTooFar = (pbcDistance2 > marginDistance2);
913 int pullCheckPbcWithinGroups(const pull_t& pull, gmx::ArrayRef<const gmx::RVec> x, const t_pbc& pbc, real pbcMargin)
915 if (pbc.pbcType == PbcType::No)
920 /* Determine what dimensions are used for each group by pull coordinates */
921 std::vector<BoolVec> dimUsed(pull.group.size(), { false, false, false });
922 for (size_t c = 0; c < pull.coord.size(); c++)
924 const t_pull_coord& coordParams = pull.coord[c].params;
925 for (int groupIndex = 0; groupIndex < coordParams.ngroup; groupIndex++)
927 for (int d = 0; d < DIM; d++)
929 if (coordParams.dim[d]
930 && !(coordParams.eGeom == PullGroupGeometry::Cylinder && groupIndex == 0))
932 dimUsed[coordParams.group[groupIndex]][d] = true;
938 /* Check PBC for every group that uses a PBC reference atom treatment */
939 for (size_t g = 0; g < pull.group.size(); g++)
941 const pull_group_work_t& group = pull.group[g];
942 if ((group.epgrppbc == epgrppbcREFAT || group.epgrppbc == epgrppbcPREVSTEPCOM)
943 && !pullGroupObeysPbcRestrictions(group, dimUsed[g], x, pbc, pull.comm.pbcAtomBuffer[g], pbcMargin))
952 bool pullCheckPbcWithinGroup(const pull_t& pull,
953 gmx::ArrayRef<const gmx::RVec> x,
958 if (pbc.pbcType == PbcType::No)
962 GMX_ASSERT(groupNr < gmx::ssize(pull.group), "groupNr is out of range");
964 /* Check PBC if the group uses a PBC reference atom treatment. */
965 const pull_group_work_t& group = pull.group[groupNr];
966 if (group.epgrppbc != epgrppbcREFAT && group.epgrppbc != epgrppbcPREVSTEPCOM)
971 /* Determine what dimensions are used for each group by pull coordinates */
972 BoolVec dimUsed = { false, false, false };
973 for (size_t c = 0; c < pull.coord.size(); c++)
975 const t_pull_coord& coordParams = pull.coord[c].params;
976 for (int groupIndex = 0; groupIndex < coordParams.ngroup; groupIndex++)
978 if (coordParams.group[groupIndex] == groupNr)
980 for (int d = 0; d < DIM; d++)
982 if (coordParams.dim[d]
983 && !(coordParams.eGeom == PullGroupGeometry::Cylinder && groupIndex == 0))
992 return (pullGroupObeysPbcRestrictions(
993 group, dimUsed, x, pbc, pull.comm.pbcAtomBuffer[groupNr], pbcMargin));
996 void setPrevStepPullComFromState(struct pull_t* pull, const t_state* state)
998 for (size_t g = 0; g < pull->group.size(); g++)
1000 for (int j = 0; j < DIM; j++)
1002 pull->group[g].x_prev_step[j] = state->pull_com_prev_step[g * DIM + j];
1007 void updatePrevStepPullCom(struct pull_t* pull, t_state* state)
1009 for (size_t g = 0; g < pull->group.size(); g++)
1011 if (pull->group[g].needToCalcCom)
1013 for (int j = 0; j < DIM; j++)
1015 pull->group[g].x_prev_step[j] = pull->group[g].x[j];
1016 state->pull_com_prev_step[g * DIM + j] = pull->group[g].x[j];
1022 void allocStatePrevStepPullCom(t_state* state, const pull_t* pull)
1026 state->pull_com_prev_step.clear();
1029 size_t ngroup = pull->group.size();
1030 if (state->pull_com_prev_step.size() / DIM != ngroup)
1032 state->pull_com_prev_step.resize(ngroup * DIM, NAN);
1036 void initPullComFromPrevStep(const t_commrec* cr,
1038 gmx::ArrayRef<const real> masses,
1040 gmx::ArrayRef<const gmx::RVec> x)
1042 pull_comm_t* comm = &pull->comm;
1043 size_t ngroup = pull->group.size();
1045 if (!comm->bParticipate)
1050 GMX_ASSERT(comm->pbcAtomBuffer.size() == pull->group.size(),
1051 "pbcAtomBuffer should have size number of groups");
1052 GMX_ASSERT(comm->comBuffer.size() == pull->group.size() * c_comBufferStride,
1053 "comBuffer should have size #group*c_comBufferStride");
1055 pull_set_pbcatoms(cr, pull, x, comm->pbcAtomBuffer);
1057 for (size_t g = 0; g < ngroup; g++)
1059 pull_group_work_t* pgrp;
1061 pgrp = &pull->group[g];
1063 if (pgrp->needToCalcCom && pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
1065 GMX_ASSERT(pgrp->params.ind.size() > 1,
1066 "Groups with no atoms, or only one atom, should not "
1067 "use the COM from the previous step as reference.");
1069 gmx::RVec x_pbc = { 0, 0, 0 };
1070 copy_rvec(comm->pbcAtomBuffer[g], x_pbc);
1074 fprintf(debug, "Initialising prev step COM of pull group %zu. x_pbc =", g);
1075 for (int m = 0; m < DIM; m++)
1077 fprintf(debug, " %f", x_pbc[m]);
1079 fprintf(debug, "\n");
1082 /* The following is to a large extent similar to pull_calc_coms() */
1084 /* The final sums should end up in sum_com[0] */
1085 ComSums& comSumsTotal = pull->comSums[0];
1087 if (pgrp->atomSet.numAtomsLocal() <= c_pullMaxNumLocalAtomsSingleThreaded)
1089 sum_com_part(pgrp, 0, pgrp->atomSet.numAtomsLocal(), x, {}, masses, pbc, x_pbc, &comSumsTotal);
1093 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
1094 for (int t = 0; t < pull->nthreads; t++)
1096 int ind_start = (pgrp->atomSet.numAtomsLocal() * (t + 0)) / pull->nthreads;
1097 int ind_end = (pgrp->atomSet.numAtomsLocal() * (t + 1)) / pull->nthreads;
1098 sum_com_part(pgrp, ind_start, ind_end, x, {}, masses, pbc, x_pbc, &pull->comSums[t]);
1101 /* Reduce the thread contributions to sum_com[0] */
1102 for (int t = 1; t < pull->nthreads; t++)
1104 comSumsTotal.sum_wm += pull->comSums[t].sum_wm;
1105 comSumsTotal.sum_wwm += pull->comSums[t].sum_wwm;
1106 dvec_inc(comSumsTotal.sum_wmx, pull->comSums[t].sum_wmx);
1107 dvec_inc(comSumsTotal.sum_wmxp, pull->comSums[t].sum_wmxp);
1111 if (pgrp->localWeights.empty())
1113 comSumsTotal.sum_wwm = comSumsTotal.sum_wm;
1116 /* Copy local sums to a buffer for global summing */
1117 auto localSums = gmx::arrayRefFromArray(comm->comBuffer.data() + g * c_comBufferStride,
1120 localSums[0] = comSumsTotal.sum_wmx;
1121 localSums[1] = comSumsTotal.sum_wmxp;
1122 localSums[2][0] = comSumsTotal.sum_wm;
1123 localSums[2][1] = comSumsTotal.sum_wwm;
1124 localSums[2][2] = 0;
1128 pullAllReduce(cr, comm, ngroup * c_comBufferStride * DIM, static_cast<double*>(comm->comBuffer[0]));
1130 for (size_t g = 0; g < ngroup; g++)
1132 pull_group_work_t* pgrp;
1134 pgrp = &pull->group[g];
1135 if (pgrp->needToCalcCom)
1137 if (pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
1139 auto localSums = gmx::arrayRefFromArray(
1140 comm->comBuffer.data() + g * c_comBufferStride, c_comBufferStride);
1141 double wmass, wwmass;
1143 /* Determine the inverse mass */
1144 wmass = localSums[2][0];
1145 wwmass = localSums[2][1];
1146 pgrp->mwscale = 1.0 / wmass;
1147 /* invtm==0 signals a frozen group, so then we should keep it zero */
1148 if (pgrp->invtm != 0)
1150 pgrp->wscale = wmass / wwmass;
1151 pgrp->invtm = wwmass / (wmass * wmass);
1153 /* Divide by the total mass */
1154 for (int m = 0; m < DIM; m++)
1156 pgrp->x[m] = localSums[0][m] * pgrp->mwscale;
1157 pgrp->x[m] += comm->pbcAtomBuffer[g][m];
1161 fprintf(debug, "Pull group %zu wmass %f invtm %f\n", g, 1.0 / pgrp->mwscale, pgrp->invtm);
1162 fprintf(debug, "Initialising prev step COM of pull group %zu to", g);
1163 for (int m = 0; m < DIM; m++)
1165 fprintf(debug, " %f", pgrp->x[m]);
1167 fprintf(debug, "\n");
1169 copy_dvec(pgrp->x, pgrp->x_prev_step);