fd8dbef6a6ea341d66cba3776c2678f1b39a9fea
[alexxy/gromacs.git] / src / gromacs / pulling / pullutil.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "config.h"
41
42 #include <cassert>
43 #include <cstdlib>
44
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"
62
63 #include "pull_internal.h"
64
65 #if GMX_MPI
66
67 // Helper function to deduce MPI datatype from the type of data
68 gmx_unused static MPI_Datatype mpiDatatype(const float gmx_unused* data)
69 {
70     return MPI_FLOAT;
71 }
72
73 // Helper function to deduce MPI datatype from the type of data
74 gmx_unused static MPI_Datatype mpiDatatype(const double gmx_unused* data)
75 {
76     return MPI_DOUBLE;
77 }
78
79 #endif // GMX_MPI
80
81 #if !GMX_DOUBLE
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)
84 {
85     gmx_sum(n, data, cr);
86 }
87 #endif
88
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)
91 {
92     gmx_sumd(n, data, cr);
93 }
94
95 // Reduce data of n elements over all ranks currently participating in pull
96 template<typename T>
97 static void pullAllReduce(const t_commrec* cr, pull_comm_t* comm, int n, T* data)
98 {
99     if (cr != nullptr && PAR(cr))
100     {
101         if (comm->bParticipateAll)
102         {
103             /* Sum the contributions over all DD ranks */
104             gmxAllReduce(n, data, cr);
105         }
106         else
107         {
108             /* Separate branch because gmx_sum uses cr->mpi_comm_mygroup */
109 #if GMX_MPI
110             MPI_Allreduce(MPI_IN_PLACE, data, n, mpiDatatype(data), MPI_SUM, comm->mpi_comm_com);
111 #else
112             gmx_incons("comm->bParticipateAll=FALSE without GMX_MPI");
113 #endif
114         }
115     }
116 }
117
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.
120  */
121 static void setPbcAtomCoords(const pull_group_work_t& pgrp, const rvec* x, rvec x_pbc)
122 {
123     if (pgrp.pbcAtomSet != nullptr)
124     {
125         if (pgrp.pbcAtomSet->numAtomsLocal() > 0)
126         {
127             /* We have the atom locally, copy its coordinates */
128             copy_rvec(x[pgrp.pbcAtomSet->localIndex()[0]], x_pbc);
129         }
130         else
131         {
132             /* Another rank has it, clear the coordinates for MPI_Allreduce */
133             clear_rvec(x_pbc);
134         }
135     }
136     else
137     {
138         copy_rvec(x[pgrp.params.pbcatom], x_pbc);
139     }
140 }
141
142 static void pull_set_pbcatoms(const t_commrec* cr, struct pull_t* pull, const rvec* x, gmx::ArrayRef<gmx::RVec> x_pbc)
143 {
144     int numPbcAtoms = 0;
145     for (size_t g = 0; g < pull->group.size(); g++)
146     {
147         const pull_group_work_t& group = pull->group[g];
148         if (group.needToCalcCom && (group.epgrppbc == epgrppbcREFAT || group.epgrppbc == epgrppbcPREVSTEPCOM))
149         {
150             setPbcAtomCoords(pull->group[g], x, x_pbc[g]);
151             numPbcAtoms++;
152         }
153         else
154         {
155             clear_rvec(x_pbc[g]);
156         }
157     }
158
159     if (cr && PAR(cr) && numPbcAtoms > 0)
160     {
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.
164          */
165         pullAllReduce(cr, &pull->comm, pull->group.size() * DIM, static_cast<real*>(x_pbc[0]));
166     }
167 }
168
169 static void
170 make_cyl_refgrps(const t_commrec* cr, pull_t* pull, const real* masses, t_pbc* pbc, double t, const rvec* x)
171 {
172     pull_comm_t* comm = &pull->comm;
173
174     GMX_ASSERT(comm->cylinderBuffer.size() == pull->coord.size() * c_cylinderBufferStride,
175                "cylinderBuffer should have the correct size");
176
177     double inv_cyl_r2 = 1.0 / gmx::square(pull->params.cylinder_r);
178
179     /* loop over all groups to make a reference group for each*/
180     for (size_t c = 0; c < pull->coord.size(); c++)
181     {
182         pull_coord_work_t* pcrd;
183         double             sum_a, wmass, wwmass;
184         dvec               radf_fac0, radf_fac1;
185
186         pcrd = &pull->coord[c];
187
188         sum_a  = 0;
189         wmass  = 0;
190         wwmass = 0;
191         clear_dvec(radf_fac0);
192         clear_dvec(radf_fac1);
193
194         if (pcrd->params.eGeom == epullgCYL)
195         {
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];
200             rvec                     direction;
201             copy_dvec_to_rvec(pcrd->spatialData.vec, direction);
202
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.
208              */
209             if (pcrd->params.rate != 0)
210             {
211                 /* With rate=0, value_ref is set initially */
212                 pcrd->value_ref = pcrd->params.init + pcrd->params.rate * t;
213             }
214             rvec reference;
215             for (int m = 0; m < DIM; m++)
216             {
217                 reference[m] = pgrp.x[m] - pcrd->spatialData.vec[m] * pcrd->value_ref;
218             }
219
220             auto localAtomIndices = pref.atomSet.localIndex();
221
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.
224              */
225             pdyna.localWeights.resize(localAtomIndices.size());
226             pdyna.mdw.resize(localAtomIndices.size());
227             pdyna.dv.resize(localAtomIndices.size());
228
229             /* loop over all atoms in the main ref group */
230             for (gmx::index indexInSet = 0; indexInSet < localAtomIndices.ssize(); indexInSet++)
231             {
232                 int  atomIndex = localAtomIndices[indexInSet];
233                 rvec dx;
234                 pbc_dx_aiuc(pbc, x[atomIndex], reference, dx);
235                 double axialLocation = iprod(direction, dx);
236                 dvec   radialLocation;
237                 double dr2 = 0;
238                 for (int m = 0; m < DIM; m++)
239                 {
240                     /* Determine the radial components */
241                     radialLocation[m] = dx[m] - axialLocation * direction[m];
242                     dr2 += gmx::square(radialLocation[m]);
243                 }
244                 double dr2_rel = dr2 * inv_cyl_r2;
245
246                 if (dr2_rel < 1)
247                 {
248                     /* add atom to sum of COM and to weight array */
249
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
254                      * on both groups.
255                      */
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;
262                     dvec mdw;
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.
269                      */
270                     pdyna.dv[indexInSet] = axialLocation;
271                     for (int m = 0; m < DIM; m++)
272                     {
273                         radf_fac0[m] += mdw[m];
274                         radf_fac1[m] += mdw[m] * axialLocation;
275                     }
276                 }
277                 else
278                 {
279                     pdyna.localWeights[indexInSet] = 0;
280                 }
281             }
282         }
283
284         auto buffer = gmx::arrayRefFromArray(
285                 comm->cylinderBuffer.data() + c * c_cylinderBufferStride, c_cylinderBufferStride);
286
287         buffer[0] = wmass;
288         buffer[1] = wwmass;
289         buffer[2] = sum_a;
290
291         buffer[3] = radf_fac0[XX];
292         buffer[4] = radf_fac0[YY];
293         buffer[5] = radf_fac0[ZZ];
294
295         buffer[6] = radf_fac1[XX];
296         buffer[7] = radf_fac1[YY];
297         buffer[8] = radf_fac1[ZZ];
298     }
299
300     if (cr != nullptr && PAR(cr))
301     {
302         /* Sum the contributions over the ranks */
303         pullAllReduce(cr, comm, pull->coord.size() * c_cylinderBufferStride, comm->cylinderBuffer.data());
304     }
305
306     for (size_t c = 0; c < pull->coord.size(); c++)
307     {
308         pull_coord_work_t* pcrd;
309
310         pcrd = &pull->coord[c];
311
312         if (pcrd->params.eGeom == epullgCYL)
313         {
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;
317
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.
325              */
326             pdyna->wscale = wmass / wwmass;
327             pdyna->invtm  = wwmass / (wmass * wmass);
328
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.
332              */
333             spatialData.cyl_dev = 0;
334             for (int m = 0; m < DIM; m++)
335             {
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;
340             }
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.
345              */
346             for (int m = 0; m < DIM; m++)
347             {
348                 spatialData.ffrad[m] = (buffer[6 + m] + buffer[3 + m] * spatialData.cyl_dev) / wmass;
349             }
350
351             if (debug)
352             {
353                 fprintf(debug,
354                         "Pull cylinder group %zu:%8.3f%8.3f%8.3f m:%8.3f\n",
355                         c,
356                         pdyna->x[0],
357                         pdyna->x[1],
358                         pdyna->x[2],
359                         1.0 / pdyna->invtm);
360                 fprintf(debug,
361                         "ffrad %8.3f %8.3f %8.3f\n",
362                         spatialData.ffrad[XX],
363                         spatialData.ffrad[YY],
364                         spatialData.ffrad[ZZ]);
365             }
366         }
367     }
368 }
369
370 static double atan2_0_2pi(double y, double x)
371 {
372     double a;
373
374     a = atan2(y, x);
375     if (a < 0)
376     {
377         a += 2.0 * M_PI;
378     }
379     return a;
380 }
381
382 static void sum_com_part(const pull_group_work_t* pgrp,
383                          int                      ind_start,
384                          int                      ind_end,
385                          const rvec*              x,
386                          const rvec*              xp,
387                          const real*              mass,
388                          const t_pbc*             pbc,
389                          const rvec               x_pbc,
390                          ComSums*                 sum_com)
391 {
392     double sum_wm   = 0;
393     double sum_wwm  = 0;
394     dvec   sum_wmx  = { 0, 0, 0 };
395     dvec   sum_wmxp = { 0, 0, 0 };
396
397     auto localAtomIndices = pgrp->atomSet.localIndex();
398     for (int i = ind_start; i < ind_end; i++)
399     {
400         int  ii = localAtomIndices[i];
401         real wm;
402         if (pgrp->localWeights.empty())
403         {
404             wm = mass[ii];
405             sum_wm += wm;
406         }
407         else
408         {
409             real w;
410
411             w  = pgrp->localWeights[i];
412             wm = w * mass[ii];
413             sum_wm += wm;
414             sum_wwm += wm * w;
415         }
416         if (pgrp->epgrppbc == epgrppbcNONE)
417         {
418             /* Plain COM: sum the coordinates */
419             for (int d = 0; d < DIM; d++)
420             {
421                 sum_wmx[d] += wm * x[ii][d];
422             }
423             if (xp)
424             {
425                 for (int d = 0; d < DIM; d++)
426                 {
427                     sum_wmxp[d] += wm * xp[ii][d];
428                 }
429             }
430         }
431         else
432         {
433             rvec dx;
434
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++)
438             {
439                 sum_wmx[d] += wm * dx[d];
440             }
441             if (xp)
442             {
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.
446                  */
447                 for (int d = 0; d < DIM; d++)
448                 {
449                     sum_wmxp[d] += wm * (dx[d] + xp[ii][d] - x[ii][d]);
450                 }
451             }
452         }
453     }
454
455     sum_com->sum_wm  = sum_wm;
456     sum_com->sum_wwm = sum_wwm;
457     copy_dvec(sum_wmx, sum_com->sum_wmx);
458     if (xp)
459     {
460         copy_dvec(sum_wmxp, sum_com->sum_wmxp);
461     }
462 }
463
464 static void sum_com_part_cosweight(const pull_group_work_t* pgrp,
465                                    int                      ind_start,
466                                    int                      ind_end,
467                                    int                      cosdim,
468                                    real                     twopi_box,
469                                    const rvec*              x,
470                                    const rvec*              xp,
471                                    const real*              mass,
472                                    ComSums*                 sum_com)
473 {
474     /* Cosine weighting geometry */
475     double sum_cm  = 0;
476     double sum_sm  = 0;
477     double sum_ccm = 0;
478     double sum_csm = 0;
479     double sum_ssm = 0;
480     double sum_cmp = 0;
481     double sum_smp = 0;
482
483     auto localAtomIndices = pgrp->atomSet.localIndex();
484
485     for (int i = ind_start; i < ind_end; i++)
486     {
487         int  ii = localAtomIndices[i];
488         real m  = mass[ii];
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);
497
498         if (xp != nullptr)
499         {
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);
504         }
505     }
506
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;
514 }
515
516 /* calculates center of mass of selection index from all coordinates x */
517 // Compiler segfault with 2019_update_5 and 2020_initial
518 #if defined(__INTEL_COMPILER)                                          \
519         && ((__INTEL_COMPILER == 1900 && __INTEL_COMPILER_UPDATE >= 5) \
520             || (__INTEL_COMPILER >= 1910 && __INTEL_COMPILER < 2021))
521 #    pragma intel optimization_level 2
522 #endif
523 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 {
525     real         twopi_box = 0;
526     pull_comm_t* comm;
527
528     comm = &pull->comm;
529
530     GMX_ASSERT(comm->pbcAtomBuffer.size() == pull->group.size(),
531                "pbcAtomBuffer should have size number of groups");
532     GMX_ASSERT(comm->comBuffer.size() == pull->group.size() * c_comBufferStride,
533                "comBuffer should have size #group*c_comBufferStride");
534
535     if (pull->bRefAt && pull->bSetPBCatoms)
536     {
537         pull_set_pbcatoms(cr, pull, x, comm->pbcAtomBuffer);
538
539         if (cr != nullptr && DOMAINDECOMP(cr))
540         {
541             /* We can keep these PBC reference coordinates fixed for nstlist
542              * steps, since atoms won't jump over PBC.
543              * This avoids a global reduction at the next nstlist-1 steps.
544              * Note that the exact values of the pbc reference coordinates
545              * are irrelevant, as long all atoms in the group are within
546              * half a box distance of the reference coordinate.
547              */
548             pull->bSetPBCatoms = FALSE;
549         }
550     }
551
552     if (pull->cosdim >= 0)
553     {
554         int m;
555
556         assert(pull->npbcdim <= DIM);
557
558         for (m = pull->cosdim + 1; m < pull->npbcdim; m++)
559         {
560             if (pbc->box[m][pull->cosdim] != 0)
561             {
562                 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
563             }
564         }
565         twopi_box = 2.0 * M_PI / pbc->box[pull->cosdim][pull->cosdim];
566     }
567
568     for (size_t g = 0; g < pull->group.size(); g++)
569     {
570         pull_group_work_t* pgrp = &pull->group[g];
571
572         /* Cosine-weighted COMs behave different from all other weighted COMs
573          * in the sense that the weights depend on instantaneous coordinates,
574          * not on pre-set weights. Thus we resize the local weight buffer here.
575          */
576         if (pgrp->epgrppbc == epgrppbcCOS)
577         {
578             pgrp->localWeights.resize(pgrp->atomSet.localIndex().size());
579         }
580
581         auto comBuffer = gmx::arrayRefFromArray(comm->comBuffer.data() + g * c_comBufferStride,
582                                                 c_comBufferStride);
583
584         if (pgrp->needToCalcCom)
585         {
586             if (pgrp->epgrppbc != epgrppbcCOS)
587             {
588                 rvec x_pbc = { 0, 0, 0 };
589
590                 switch (pgrp->epgrppbc)
591                 {
592                     case epgrppbcREFAT:
593                         /* Set the pbc atom */
594                         copy_rvec(comm->pbcAtomBuffer[g], x_pbc);
595                         break;
596                     case epgrppbcPREVSTEPCOM:
597                         /* Set the pbc reference to the COM of the group of the last step */
598                         copy_dvec_to_rvec(pgrp->x_prev_step, comm->pbcAtomBuffer[g]);
599                         copy_dvec_to_rvec(pgrp->x_prev_step, x_pbc);
600                 }
601
602                 /* The final sums should end up in comSums[0] */
603                 ComSums& comSumsTotal = pull->comSums[0];
604
605                 /* If we have a single-atom group the mass is irrelevant, so
606                  * we can remove the mass factor to avoid division by zero.
607                  * Note that with constraint pulling the mass does matter, but
608                  * in that case a check group mass != 0 has been done before.
609                  */
610                 if (pgrp->params.ind.size() == 1 && pgrp->atomSet.numAtomsLocal() == 1
611                     && masses[pgrp->atomSet.localIndex()[0]] == 0)
612                 {
613                     GMX_ASSERT(xp == nullptr,
614                                "We should not have groups with zero mass with constraints, i.e. "
615                                "xp!=NULL");
616
617                     /* Copy the single atom coordinate */
618                     for (int d = 0; d < DIM; d++)
619                     {
620                         comSumsTotal.sum_wmx[d] = x[pgrp->atomSet.localIndex()[0]][d];
621                     }
622                     /* Set all mass factors to 1 to get the correct COM */
623                     comSumsTotal.sum_wm  = 1;
624                     comSumsTotal.sum_wwm = 1;
625                 }
626                 else if (pgrp->atomSet.numAtomsLocal() <= c_pullMaxNumLocalAtomsSingleThreaded)
627                 {
628                     sum_com_part(pgrp, 0, pgrp->atomSet.numAtomsLocal(), x, xp, masses, pbc, x_pbc, &comSumsTotal);
629                 }
630                 else
631                 {
632 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
633                     for (int t = 0; t < pull->nthreads; t++)
634                     {
635                         int ind_start = (pgrp->atomSet.numAtomsLocal() * (t + 0)) / pull->nthreads;
636                         int ind_end   = (pgrp->atomSet.numAtomsLocal() * (t + 1)) / pull->nthreads;
637                         sum_com_part(
638                                 pgrp, ind_start, ind_end, x, xp, masses, pbc, x_pbc, &pull->comSums[t]);
639                     }
640
641                     /* Reduce the thread contributions to sum_com[0] */
642                     for (int t = 1; t < pull->nthreads; t++)
643                     {
644                         comSumsTotal.sum_wm += pull->comSums[t].sum_wm;
645                         comSumsTotal.sum_wwm += pull->comSums[t].sum_wwm;
646                         dvec_inc(comSumsTotal.sum_wmx, pull->comSums[t].sum_wmx);
647                         dvec_inc(comSumsTotal.sum_wmxp, pull->comSums[t].sum_wmxp);
648                     }
649                 }
650
651                 if (pgrp->localWeights.empty())
652                 {
653                     comSumsTotal.sum_wwm = comSumsTotal.sum_wm;
654                 }
655
656                 /* Copy local sums to a buffer for global summing */
657                 copy_dvec(comSumsTotal.sum_wmx, comBuffer[0]);
658
659                 copy_dvec(comSumsTotal.sum_wmxp, comBuffer[1]);
660
661                 comBuffer[2][0] = comSumsTotal.sum_wm;
662                 comBuffer[2][1] = comSumsTotal.sum_wwm;
663                 comBuffer[2][2] = 0;
664             }
665             else
666             {
667                 /* Cosine weighting geometry.
668                  * This uses a slab of the system, thus we always have many
669                  * atoms in the pull groups. Therefore, always use threads.
670                  */
671 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
672                 for (int t = 0; t < pull->nthreads; t++)
673                 {
674                     int ind_start = (pgrp->atomSet.numAtomsLocal() * (t + 0)) / pull->nthreads;
675                     int ind_end   = (pgrp->atomSet.numAtomsLocal() * (t + 1)) / pull->nthreads;
676                     sum_com_part_cosweight(
677                             pgrp, ind_start, ind_end, pull->cosdim, twopi_box, x, xp, masses, &pull->comSums[t]);
678                 }
679
680                 /* Reduce the thread contributions to comSums[0] */
681                 ComSums& comSumsTotal = pull->comSums[0];
682                 for (int t = 1; t < pull->nthreads; t++)
683                 {
684                     comSumsTotal.sum_cm += pull->comSums[t].sum_cm;
685                     comSumsTotal.sum_sm += pull->comSums[t].sum_sm;
686                     comSumsTotal.sum_ccm += pull->comSums[t].sum_ccm;
687                     comSumsTotal.sum_csm += pull->comSums[t].sum_csm;
688                     comSumsTotal.sum_ssm += pull->comSums[t].sum_ssm;
689                     comSumsTotal.sum_cmp += pull->comSums[t].sum_cmp;
690                     comSumsTotal.sum_smp += pull->comSums[t].sum_smp;
691                 }
692
693                 /* Copy local sums to a buffer for global summing */
694                 comBuffer[0][0] = comSumsTotal.sum_cm;
695                 comBuffer[0][1] = comSumsTotal.sum_sm;
696                 comBuffer[0][2] = 0;
697                 comBuffer[1][0] = comSumsTotal.sum_ccm;
698                 comBuffer[1][1] = comSumsTotal.sum_csm;
699                 comBuffer[1][2] = comSumsTotal.sum_ssm;
700                 comBuffer[2][0] = comSumsTotal.sum_cmp;
701                 comBuffer[2][1] = comSumsTotal.sum_smp;
702                 comBuffer[2][2] = 0;
703             }
704         }
705         else
706         {
707             clear_dvec(comBuffer[0]);
708             clear_dvec(comBuffer[1]);
709             clear_dvec(comBuffer[2]);
710         }
711     }
712
713     pullAllReduce(cr,
714                   comm,
715                   pull->group.size() * c_comBufferStride * DIM,
716                   static_cast<double*>(comm->comBuffer[0]));
717
718     for (size_t g = 0; g < pull->group.size(); g++)
719     {
720         pull_group_work_t* pgrp;
721
722         pgrp = &pull->group[g];
723         if (pgrp->needToCalcCom)
724         {
725             GMX_ASSERT(!pgrp->params.ind.empty(),
726                        "Normal pull groups should have atoms, only group 0, which should have "
727                        "bCalcCom=FALSE has nat=0");
728
729             const auto comBuffer = gmx::constArrayRefFromArray(
730                     comm->comBuffer.data() + g * c_comBufferStride, c_comBufferStride);
731
732             if (pgrp->epgrppbc != epgrppbcCOS)
733             {
734                 double wmass, wwmass;
735                 int    m;
736
737                 /* Determine the inverse mass */
738                 wmass         = comBuffer[2][0];
739                 wwmass        = comBuffer[2][1];
740                 pgrp->mwscale = 1.0 / wmass;
741                 /* invtm==0 signals a frozen group, so then we should keep it zero */
742                 if (pgrp->invtm != 0)
743                 {
744                     pgrp->wscale = wmass / wwmass;
745                     pgrp->invtm  = wwmass / (wmass * wmass);
746                 }
747                 /* Divide by the total mass */
748                 for (m = 0; m < DIM; m++)
749                 {
750                     pgrp->x[m] = comBuffer[0][m] * pgrp->mwscale;
751                     if (xp)
752                     {
753                         pgrp->xp[m] = comBuffer[1][m] * pgrp->mwscale;
754                     }
755                     if (pgrp->epgrppbc == epgrppbcREFAT || pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
756                     {
757                         pgrp->x[m] += comm->pbcAtomBuffer[g][m];
758                         if (xp)
759                         {
760                             pgrp->xp[m] += comm->pbcAtomBuffer[g][m];
761                         }
762                     }
763                 }
764             }
765             else
766             {
767                 /* Cosine weighting geometry */
768                 double csw, snw, wmass, wwmass;
769
770                 /* Determine the optimal location of the cosine weight */
771                 csw                   = comBuffer[0][0];
772                 snw                   = comBuffer[0][1];
773                 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw) / twopi_box;
774                 /* Set the weights for the local atoms */
775                 wmass  = sqrt(csw * csw + snw * snw);
776                 wwmass = (comBuffer[1][0] * csw * csw + comBuffer[1][1] * csw * snw
777                           + comBuffer[1][2] * snw * snw)
778                          / (wmass * wmass);
779
780                 pgrp->mwscale = 1.0 / wmass;
781                 pgrp->wscale  = wmass / wwmass;
782                 pgrp->invtm   = wwmass / (wmass * wmass);
783                 /* Set the weights for the local atoms */
784                 csw *= pgrp->invtm;
785                 snw *= pgrp->invtm;
786                 for (size_t i = 0; i < pgrp->atomSet.numAtomsLocal(); i++)
787                 {
788                     int ii                = pgrp->atomSet.localIndex()[i];
789                     pgrp->localWeights[i] = csw * std::cos(twopi_box * x[ii][pull->cosdim])
790                                             + snw * std::sin(twopi_box * x[ii][pull->cosdim]);
791                 }
792                 if (xp)
793                 {
794                     csw                    = comBuffer[2][0];
795                     snw                    = comBuffer[2][1];
796                     pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw) / twopi_box;
797                 }
798             }
799             if (debug)
800             {
801                 fprintf(debug, "Pull group %zu wmass %f invtm %f\n", g, 1.0 / pgrp->mwscale, pgrp->invtm);
802             }
803         }
804     }
805
806     if (pull->bCylinder)
807     {
808         /* Calculate the COMs for the cyclinder reference groups */
809         make_cyl_refgrps(cr, pull, masses, pbc, t, x);
810     }
811 }
812
813 using BoolVec = gmx::BasicVector<bool>;
814
815 /* Returns whether the pull group obeys the PBC restrictions */
816 static bool pullGroupObeysPbcRestrictions(const pull_group_work_t& group,
817                                           const BoolVec&           dimUsed,
818                                           const rvec*              x,
819                                           const t_pbc&             pbc,
820                                           const gmx::RVec&         x_pbc,
821                                           const real               pbcMargin)
822 {
823     /* Determine which dimensions are relevant for PBC */
824     BoolVec dimUsesPbc       = { false, false, false };
825     bool    pbcIsRectangular = true;
826     for (int d = 0; d < pbc.ndim_ePBC; d++)
827     {
828         if (dimUsed[d])
829         {
830             dimUsesPbc[d] = true;
831             /* All non-zero dimensions of vector v are involved in PBC */
832             for (int d2 = d + 1; d2 < pbc.ndim_ePBC; d2++)
833             {
834                 assert(d2 < DIM);
835                 if (pbc.box[d2][d] != 0)
836                 {
837                     dimUsesPbc[d2]   = true;
838                     pbcIsRectangular = false;
839                 }
840             }
841         }
842     }
843
844     rvec marginPerDim    = {};
845     real marginDistance2 = 0;
846     if (pbcIsRectangular)
847     {
848         /* Use margins for dimensions independently */
849         for (int d = 0; d < pbc.ndim_ePBC; d++)
850         {
851             marginPerDim[d] = pbcMargin * pbc.hbox_diag[d];
852         }
853     }
854     else
855     {
856         /* Check the total distance along the relevant dimensions */
857         for (int d = 0; d < pbc.ndim_ePBC; d++)
858         {
859             if (dimUsesPbc[d])
860             {
861                 marginDistance2 += pbcMargin * gmx::square(0.5) * norm2(pbc.box[d]);
862             }
863         }
864     }
865
866     auto localAtomIndices = group.atomSet.localIndex();
867     for (gmx::index indexInSet = 0; indexInSet < localAtomIndices.ssize(); indexInSet++)
868     {
869         rvec dx;
870         pbc_dx(&pbc, x[localAtomIndices[indexInSet]], x_pbc, dx);
871
872         bool atomIsTooFar = false;
873         if (pbcIsRectangular)
874         {
875             for (int d = 0; d < pbc.ndim_ePBC; d++)
876             {
877                 if (dimUsesPbc[d] && (dx[d] < -marginPerDim[d] || dx[d] > marginPerDim[d]))
878                 {
879                     atomIsTooFar = true;
880                 }
881             }
882         }
883         else
884         {
885             real pbcDistance2 = 0;
886             for (int d = 0; d < pbc.ndim_ePBC; d++)
887             {
888                 if (dimUsesPbc[d])
889                 {
890                     pbcDistance2 += gmx::square(dx[d]);
891                 }
892             }
893             atomIsTooFar = (pbcDistance2 > marginDistance2);
894         }
895         if (atomIsTooFar)
896         {
897             return false;
898         }
899     }
900
901     return true;
902 }
903
904 int pullCheckPbcWithinGroups(const pull_t& pull, const rvec* x, const t_pbc& pbc, real pbcMargin)
905 {
906     if (pbc.pbcType == PbcType::No)
907     {
908         return -1;
909     }
910
911     /* Determine what dimensions are used for each group by pull coordinates */
912     std::vector<BoolVec> dimUsed(pull.group.size(), { false, false, false });
913     for (size_t c = 0; c < pull.coord.size(); c++)
914     {
915         const t_pull_coord& coordParams = pull.coord[c].params;
916         for (int groupIndex = 0; groupIndex < coordParams.ngroup; groupIndex++)
917         {
918             for (int d = 0; d < DIM; d++)
919             {
920                 if (coordParams.dim[d] && !(coordParams.eGeom == epullgCYL && groupIndex == 0))
921                 {
922                     dimUsed[coordParams.group[groupIndex]][d] = true;
923                 }
924             }
925         }
926     }
927
928     /* Check PBC for every group that uses a PBC reference atom treatment */
929     for (size_t g = 0; g < pull.group.size(); g++)
930     {
931         const pull_group_work_t& group = pull.group[g];
932         if ((group.epgrppbc == epgrppbcREFAT || group.epgrppbc == epgrppbcPREVSTEPCOM)
933             && !pullGroupObeysPbcRestrictions(group, dimUsed[g], x, pbc, pull.comm.pbcAtomBuffer[g], pbcMargin))
934         {
935             return g;
936         }
937     }
938
939     return -1;
940 }
941
942 bool pullCheckPbcWithinGroup(const pull_t&                  pull,
943                              gmx::ArrayRef<const gmx::RVec> x,
944                              const t_pbc&                   pbc,
945                              int                            groupNr,
946                              real                           pbcMargin)
947 {
948     if (pbc.pbcType == PbcType::No)
949     {
950         return true;
951     }
952     GMX_ASSERT(groupNr < gmx::ssize(pull.group), "groupNr is out of range");
953
954     /* Check PBC if the group uses a PBC reference atom treatment. */
955     const pull_group_work_t& group = pull.group[groupNr];
956     if (group.epgrppbc != epgrppbcREFAT && group.epgrppbc != epgrppbcPREVSTEPCOM)
957     {
958         return true;
959     }
960
961     /* Determine what dimensions are used for each group by pull coordinates */
962     BoolVec dimUsed = { false, false, false };
963     for (size_t c = 0; c < pull.coord.size(); c++)
964     {
965         const t_pull_coord& coordParams = pull.coord[c].params;
966         for (int groupIndex = 0; groupIndex < coordParams.ngroup; groupIndex++)
967         {
968             if (coordParams.group[groupIndex] == groupNr)
969             {
970                 for (int d = 0; d < DIM; d++)
971                 {
972                     if (coordParams.dim[d] && !(coordParams.eGeom == epullgCYL && groupIndex == 0))
973                     {
974                         dimUsed[d] = true;
975                     }
976                 }
977             }
978         }
979     }
980
981     return (pullGroupObeysPbcRestrictions(
982             group, dimUsed, as_rvec_array(x.data()), pbc, pull.comm.pbcAtomBuffer[groupNr], pbcMargin));
983 }
984
985 void setPrevStepPullComFromState(struct pull_t* pull, const t_state* state)
986 {
987     for (size_t g = 0; g < pull->group.size(); g++)
988     {
989         for (int j = 0; j < DIM; j++)
990         {
991             pull->group[g].x_prev_step[j] = state->pull_com_prev_step[g * DIM + j];
992         }
993     }
994 }
995
996 void updatePrevStepPullCom(struct pull_t* pull, t_state* state)
997 {
998     for (size_t g = 0; g < pull->group.size(); g++)
999     {
1000         if (pull->group[g].needToCalcCom)
1001         {
1002             for (int j = 0; j < DIM; j++)
1003             {
1004                 pull->group[g].x_prev_step[j]          = pull->group[g].x[j];
1005                 state->pull_com_prev_step[g * DIM + j] = pull->group[g].x[j];
1006             }
1007         }
1008     }
1009 }
1010
1011 void allocStatePrevStepPullCom(t_state* state, const pull_t* pull)
1012 {
1013     if (!pull)
1014     {
1015         state->pull_com_prev_step.clear();
1016         return;
1017     }
1018     size_t ngroup = pull->group.size();
1019     if (state->pull_com_prev_step.size() / DIM != ngroup)
1020     {
1021         state->pull_com_prev_step.resize(ngroup * DIM, NAN);
1022     }
1023 }
1024
1025 void initPullComFromPrevStep(const t_commrec* cr, pull_t* pull, const real* masses, t_pbc* pbc, const rvec x[])
1026 {
1027     pull_comm_t* comm   = &pull->comm;
1028     size_t       ngroup = pull->group.size();
1029
1030     if (!comm->bParticipate)
1031     {
1032         return;
1033     }
1034
1035     GMX_ASSERT(comm->pbcAtomBuffer.size() == pull->group.size(),
1036                "pbcAtomBuffer should have size number of groups");
1037     GMX_ASSERT(comm->comBuffer.size() == pull->group.size() * c_comBufferStride,
1038                "comBuffer should have size #group*c_comBufferStride");
1039
1040     pull_set_pbcatoms(cr, pull, x, comm->pbcAtomBuffer);
1041
1042     for (size_t g = 0; g < ngroup; g++)
1043     {
1044         pull_group_work_t* pgrp;
1045
1046         pgrp = &pull->group[g];
1047
1048         if (pgrp->needToCalcCom && pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
1049         {
1050             GMX_ASSERT(pgrp->params.ind.size() > 1,
1051                        "Groups with no atoms, or only one atom, should not "
1052                        "use the COM from the previous step as reference.");
1053
1054             rvec x_pbc = { 0, 0, 0 };
1055             copy_rvec(comm->pbcAtomBuffer[g], x_pbc);
1056
1057             if (debug)
1058             {
1059                 fprintf(debug, "Initialising prev step COM of pull group %zu. x_pbc =", g);
1060                 for (int m = 0; m < DIM; m++)
1061                 {
1062                     fprintf(debug, " %f", x_pbc[m]);
1063                 }
1064                 fprintf(debug, "\n");
1065             }
1066
1067             /* The following is to a large extent similar to pull_calc_coms() */
1068
1069             /* The final sums should end up in sum_com[0] */
1070             ComSums& comSumsTotal = pull->comSums[0];
1071
1072             if (pgrp->atomSet.numAtomsLocal() <= c_pullMaxNumLocalAtomsSingleThreaded)
1073             {
1074                 sum_com_part(pgrp, 0, pgrp->atomSet.numAtomsLocal(), x, nullptr, masses, pbc, x_pbc, &comSumsTotal);
1075             }
1076             else
1077             {
1078 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
1079                 for (int t = 0; t < pull->nthreads; t++)
1080                 {
1081                     int ind_start = (pgrp->atomSet.numAtomsLocal() * (t + 0)) / pull->nthreads;
1082                     int ind_end   = (pgrp->atomSet.numAtomsLocal() * (t + 1)) / pull->nthreads;
1083                     sum_com_part(
1084                             pgrp, ind_start, ind_end, x, nullptr, masses, pbc, x_pbc, &pull->comSums[t]);
1085                 }
1086
1087                 /* Reduce the thread contributions to sum_com[0] */
1088                 for (int t = 1; t < pull->nthreads; t++)
1089                 {
1090                     comSumsTotal.sum_wm += pull->comSums[t].sum_wm;
1091                     comSumsTotal.sum_wwm += pull->comSums[t].sum_wwm;
1092                     dvec_inc(comSumsTotal.sum_wmx, pull->comSums[t].sum_wmx);
1093                     dvec_inc(comSumsTotal.sum_wmxp, pull->comSums[t].sum_wmxp);
1094                 }
1095             }
1096
1097             if (pgrp->localWeights.empty())
1098             {
1099                 comSumsTotal.sum_wwm = comSumsTotal.sum_wm;
1100             }
1101
1102             /* Copy local sums to a buffer for global summing */
1103             auto localSums = gmx::arrayRefFromArray(comm->comBuffer.data() + g * c_comBufferStride,
1104                                                     c_comBufferStride);
1105
1106             localSums[0]    = comSumsTotal.sum_wmx;
1107             localSums[1]    = comSumsTotal.sum_wmxp;
1108             localSums[2][0] = comSumsTotal.sum_wm;
1109             localSums[2][1] = comSumsTotal.sum_wwm;
1110             localSums[2][2] = 0;
1111         }
1112     }
1113
1114     pullAllReduce(cr, comm, ngroup * c_comBufferStride * DIM, static_cast<double*>(comm->comBuffer[0]));
1115
1116     for (size_t g = 0; g < ngroup; g++)
1117     {
1118         pull_group_work_t* pgrp;
1119
1120         pgrp = &pull->group[g];
1121         if (pgrp->needToCalcCom)
1122         {
1123             if (pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
1124             {
1125                 auto localSums = gmx::arrayRefFromArray(
1126                         comm->comBuffer.data() + g * c_comBufferStride, c_comBufferStride);
1127                 double wmass, wwmass;
1128
1129                 /* Determine the inverse mass */
1130                 wmass         = localSums[2][0];
1131                 wwmass        = localSums[2][1];
1132                 pgrp->mwscale = 1.0 / wmass;
1133                 /* invtm==0 signals a frozen group, so then we should keep it zero */
1134                 if (pgrp->invtm != 0)
1135                 {
1136                     pgrp->wscale = wmass / wwmass;
1137                     pgrp->invtm  = wwmass / (wmass * wmass);
1138                 }
1139                 /* Divide by the total mass */
1140                 for (int m = 0; m < DIM; m++)
1141                 {
1142                     pgrp->x[m] = localSums[0][m] * pgrp->mwscale;
1143                     pgrp->x[m] += comm->pbcAtomBuffer[g][m];
1144                 }
1145                 if (debug)
1146                 {
1147                     fprintf(debug, "Pull group %zu wmass %f invtm %f\n", g, 1.0 / pgrp->mwscale, pgrp->invtm);
1148                     fprintf(debug, "Initialising prev step COM of pull group %zu to", g);
1149                     for (int m = 0; m < DIM; m++)
1150                     {
1151                         fprintf(debug, " %f", pgrp->x[m]);
1152                     }
1153                     fprintf(debug, "\n");
1154                 }
1155                 copy_dvec(pgrp->x, pgrp->x_prev_step);
1156             }
1157         }
1158     }
1159 }