Merge release-2019 into master
[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,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #include <cassert>
42 #include <cstdlib>
43
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/gmxlib/network.h"
46 #include "gromacs/math/functions.h"
47 #include "gromacs/math/utilities.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdtypes/commrec.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/mdtypes/mdatom.h"
53 #include "gromacs/mdtypes/state.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/pulling/pull.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/real.h"
60 #include "gromacs/utility/smalloc.h"
61
62 #include "pull_internal.h"
63
64 #if GMX_MPI
65
66 // Helper function to deduce MPI datatype from the type of data
67 gmx_unused static MPI_Datatype mpiDatatype(const float gmx_unused *data)
68 {
69     return MPI_FLOAT;
70 }
71
72 // Helper function to deduce MPI datatype from the type of data
73 gmx_unused static MPI_Datatype mpiDatatype(const double gmx_unused *data)
74 {
75     return MPI_DOUBLE;
76 }
77
78 #endif // GMX_MPI
79
80 #if !GMX_DOUBLE
81 // Helper function; note that gmx_sum(d) should actually be templated
82 gmx_unused static void gmxAllReduce(int n, real *data, const t_commrec *cr)
83 {
84     gmx_sum(n, data, cr);
85 }
86 #endif
87
88 // Helper function; note that gmx_sum(d) should actually be templated
89 gmx_unused static void gmxAllReduce(int n, double *data, const t_commrec *cr)
90 {
91     gmx_sumd(n, data, cr);
92 }
93
94 // Reduce data of n elements over all ranks currently participating in pull
95 template <typename T>
96 static void pullAllReduce(const t_commrec *cr,
97                           pull_comm_t     *comm,
98                           int              n,
99                           T               *data)
100 {
101     if (cr != nullptr && PAR(cr))
102     {
103         if (comm->bParticipateAll)
104         {
105             /* Sum the contributions over all DD ranks */
106             gmxAllReduce(n, data, cr);
107         }
108         else
109         {
110             /* Separate branch because gmx_sum uses cr->mpi_comm_mygroup */
111 #if GMX_MPI
112 #if MPI_IN_PLACE_EXISTS
113             MPI_Allreduce(MPI_IN_PLACE, data, n, mpiDatatype(data), MPI_SUM,
114                           comm->mpi_comm_com);
115 #else
116             std::vector<T> buf(n);
117
118             MPI_Allreduce(data, buf.data(), n, mpiDatatype(data), MPI_SUM,
119                           comm->mpi_comm_com);
120
121             /* Copy the result from the buffer to the input/output data */
122             for (int i = 0; i < n; i++)
123             {
124                 data[i] = buf[i];
125             }
126 #endif
127 #else
128             gmx_incons("comm->bParticipateAll=FALSE without GMX_MPI");
129 #endif
130         }
131     }
132 }
133
134 /* Copies the coordinates of the PBC atom of pgrp to x_pbc.
135  * When those coordinates are not available on this rank, clears x_pbc.
136  */
137 static void setPbcAtomCoords(const pull_group_work_t &pgrp,
138                              const rvec              *x,
139                              rvec                     x_pbc)
140 {
141     if (pgrp.pbcAtomSet != nullptr)
142     {
143         if (pgrp.pbcAtomSet->numAtomsLocal() > 0)
144         {
145             /* We have the atom locally, copy its coordinates */
146             copy_rvec(x[pgrp.pbcAtomSet->localIndex()[0]], x_pbc);
147         }
148         else
149         {
150             /* Another rank has it, clear the coordinates for MPI_Allreduce */
151             clear_rvec(x_pbc);
152         }
153     }
154     else
155     {
156         copy_rvec(x[pgrp.params.pbcatom], x_pbc);
157     }
158 }
159
160 static void pull_set_pbcatoms(const t_commrec *cr, struct pull_t *pull,
161                               const rvec *x,
162                               gmx::ArrayRef<gmx::RVec> x_pbc)
163 {
164     int numPbcAtoms = 0;
165     for (size_t g = 0; g < pull->group.size(); g++)
166     {
167         const pull_group_work_t &group = pull->group[g];
168         if (group.needToCalcCom && (group.epgrppbc == epgrppbcREFAT || group.epgrppbc == epgrppbcPREVSTEPCOM))
169         {
170             setPbcAtomCoords(pull->group[g], x, x_pbc[g]);
171             numPbcAtoms++;
172         }
173         else
174         {
175             clear_rvec(x_pbc[g]);
176         }
177     }
178
179     if (cr && PAR(cr) && numPbcAtoms > 0)
180     {
181         /* Sum over participating ranks to get x_pbc from the home ranks.
182          * This can be very expensive at high parallelization, so we only
183          * do this after each DD repartitioning.
184          */
185         pullAllReduce(cr, &pull->comm, pull->group.size()*DIM,
186                       static_cast<real *>(x_pbc[0]));
187     }
188 }
189
190 static void make_cyl_refgrps(const t_commrec *cr,
191                              pull_t          *pull,
192                              const t_mdatoms *md,
193                              t_pbc           *pbc,
194                              double           t,
195                              const rvec      *x)
196 {
197     pull_comm_t *comm = &pull->comm;
198
199     GMX_ASSERT(comm->cylinderBuffer.size() == pull->coord.size()*c_cylinderBufferStride, "cylinderBuffer should have the correct size");
200
201     double inv_cyl_r2 = 1.0/gmx::square(pull->params.cylinder_r);
202
203     /* loop over all groups to make a reference group for each*/
204     for (size_t c = 0; c < pull->coord.size(); c++)
205     {
206         pull_coord_work_t *pcrd;
207         double             sum_a, wmass, wwmass;
208         dvec               radf_fac0, radf_fac1;
209
210         pcrd   = &pull->coord[c];
211
212         sum_a  = 0;
213         wmass  = 0;
214         wwmass = 0;
215         clear_dvec(radf_fac0);
216         clear_dvec(radf_fac1);
217
218         if (pcrd->params.eGeom == epullgCYL)
219         {
220             /* pref will be the same group for all pull coordinates */
221             const pull_group_work_t &pref  = pull->group[pcrd->params.group[0]];
222             const pull_group_work_t &pgrp  = pull->group[pcrd->params.group[1]];
223             pull_group_work_t       &pdyna = pull->dyna[c];
224             rvec                     direction;
225             copy_dvec_to_rvec(pcrd->spatialData.vec, direction);
226
227             /* Since we have not calculated the COM of the cylinder group yet,
228              * we calculate distances with respect to location of the pull
229              * group minus the reference position along the vector.
230              * here we already have the COM of the pull group. This resolves
231              * any PBC issues and we don't need to use a PBC-atom here.
232              */
233             if (pcrd->params.rate != 0)
234             {
235                 /* With rate=0, value_ref is set initially */
236                 pcrd->value_ref = pcrd->params.init + pcrd->params.rate*t;
237             }
238             rvec reference;
239             for (int m = 0; m < DIM; m++)
240             {
241                 reference[m] = pgrp.x[m] - pcrd->spatialData.vec[m]*pcrd->value_ref;
242             }
243
244             auto localAtomIndices = pref.atomSet.localIndex();
245
246             /* This actually only needs to be done at init or DD time,
247              * but resizing with the same size does not cause much overhead.
248              */
249             pdyna.localWeights.resize(localAtomIndices.size());
250             pdyna.mdw.resize(localAtomIndices.size());
251             pdyna.dv.resize(localAtomIndices.size());
252
253             /* loop over all atoms in the main ref group */
254             for (gmx::index indexInSet = 0; indexInSet < localAtomIndices.ssize(); indexInSet++)
255             {
256                 int    atomIndex = localAtomIndices[indexInSet];
257                 rvec   dx;
258                 pbc_dx_aiuc(pbc, x[atomIndex], reference, dx);
259                 double axialLocation = iprod(direction, dx);
260                 dvec   radialLocation;
261                 double dr2 = 0;
262                 for (int m = 0; m < DIM; m++)
263                 {
264                     /* Determine the radial components */
265                     radialLocation[m]  = dx[m] - axialLocation*direction[m];
266                     dr2               += gmx::square(radialLocation[m]);
267                 }
268                 double dr2_rel = dr2*inv_cyl_r2;
269
270                 if (dr2_rel < 1)
271                 {
272                     /* add atom to sum of COM and to weight array */
273
274                     double mass                     = md->massT[atomIndex];
275                     /* The radial weight function is 1-2x^2+x^4,
276                      * where x=r/cylinder_r. Since this function depends
277                      * on the radial component, we also get radial forces
278                      * on both groups.
279                      */
280                     double weight                   =  1 + (-2 + dr2_rel)*dr2_rel;
281                     double dweight_r                = (-4 + 4*dr2_rel)*inv_cyl_r2;
282                     pdyna.localWeights[indexInSet]  = weight;
283                     sum_a                          += mass*weight*axialLocation;
284                     wmass                          += mass*weight;
285                     wwmass                         += mass*weight*weight;
286                     dvec mdw;
287                     dsvmul(mass*dweight_r, radialLocation, mdw);
288                     copy_dvec(mdw, pdyna.mdw[indexInSet]);
289                     /* Currently we only have the axial component of the
290                      * offset from the cylinder COM up to an unkown offset.
291                      * We add this offset after the reduction needed
292                      * for determining the COM of the cylinder group.
293                      */
294                     pdyna.dv[indexInSet] = axialLocation;
295                     for (int m = 0; m < DIM; m++)
296                     {
297                         radf_fac0[m] += mdw[m];
298                         radf_fac1[m] += mdw[m]*axialLocation;
299                     }
300                 }
301                 else
302                 {
303                     pdyna.localWeights[indexInSet] = 0;
304                 }
305             }
306         }
307
308         auto buffer = gmx::arrayRefFromArray(comm->cylinderBuffer.data() + c*c_cylinderBufferStride, c_cylinderBufferStride);
309
310         buffer[0] = wmass;
311         buffer[1] = wwmass;
312         buffer[2] = sum_a;
313
314         buffer[3] = radf_fac0[XX];
315         buffer[4] = radf_fac0[YY];
316         buffer[5] = radf_fac0[ZZ];
317
318         buffer[6] = radf_fac1[XX];
319         buffer[7] = radf_fac1[YY];
320         buffer[8] = radf_fac1[ZZ];
321     }
322
323     if (cr != nullptr && PAR(cr))
324     {
325         /* Sum the contributions over the ranks */
326         pullAllReduce(cr, comm, pull->coord.size()*c_cylinderBufferStride,
327                       comm->cylinderBuffer.data());
328     }
329
330     for (size_t c = 0; c < pull->coord.size(); c++)
331     {
332         pull_coord_work_t *pcrd;
333
334         pcrd  = &pull->coord[c];
335
336         if (pcrd->params.eGeom == epullgCYL)
337         {
338             pull_group_work_t    *pdyna       = &pull->dyna[c];
339             pull_group_work_t    *pgrp        = &pull->group[pcrd->params.group[1]];
340             PullCoordSpatialData &spatialData = pcrd->spatialData;
341
342             auto                  buffer      = gmx::constArrayRefFromArray(comm->cylinderBuffer.data() + c*c_cylinderBufferStride, c_cylinderBufferStride);
343             double                wmass       = buffer[0];
344             double                wwmass      = buffer[1];
345             pdyna->mwscale                    = 1.0/wmass;
346             /* Cylinder pulling can't be used with constraints, but we set
347              * wscale and invtm anyhow, in case someone would like to use them.
348              */
349             pdyna->wscale  = wmass/wwmass;
350             pdyna->invtm   = wwmass/(wmass*wmass);
351
352             /* We store the deviation of the COM from the reference location
353              * used above, since we need it when we apply the radial forces
354              * to the atoms in the cylinder group.
355              */
356             spatialData.cyl_dev = 0;
357             for (int m = 0; m < DIM; m++)
358             {
359                 double reference     = pgrp->x[m] - spatialData.vec[m]*pcrd->value_ref;
360                 double dist          = -spatialData.vec[m]*buffer[2]*pdyna->mwscale;
361                 pdyna->x[m]          = reference - dist;
362                 spatialData.cyl_dev += dist;
363             }
364             /* Now we know the exact COM of the cylinder reference group,
365              * we can determine the radial force factor (ffrad) that when
366              * multiplied with the axial pull force will give the radial
367              * force on the pulled (non-cylinder) group.
368              */
369             for (int m = 0; m < DIM; m++)
370             {
371                 spatialData.ffrad[m] = (buffer[6 + m] +
372                                         buffer[3 + m]*spatialData.cyl_dev)/wmass;
373             }
374
375             if (debug)
376             {
377                 fprintf(debug, "Pull cylinder group %zu:%8.3f%8.3f%8.3f m:%8.3f\n",
378                         c, pdyna->x[0], pdyna->x[1],
379                         pdyna->x[2], 1.0/pdyna->invtm);
380                 fprintf(debug, "ffrad %8.3f %8.3f %8.3f\n",
381                         spatialData.ffrad[XX], spatialData.ffrad[YY], spatialData.ffrad[ZZ]);
382             }
383         }
384     }
385 }
386
387 static double atan2_0_2pi(double y, double x)
388 {
389     double a;
390
391     a = atan2(y, x);
392     if (a < 0)
393     {
394         a += 2.0*M_PI;
395     }
396     return a;
397 }
398
399 static void sum_com_part(const pull_group_work_t *pgrp,
400                          int ind_start, int ind_end,
401                          const rvec *x, const rvec *xp,
402                          const real *mass,
403                          const t_pbc *pbc,
404                          const rvec x_pbc,
405                          ComSums *sum_com)
406 {
407     double sum_wm   = 0;
408     double sum_wwm  = 0;
409     dvec   sum_wmx  = { 0, 0, 0 };
410     dvec   sum_wmxp = { 0, 0, 0 };
411
412     auto   localAtomIndices = pgrp->atomSet.localIndex();
413     for (int i = ind_start; i < ind_end; i++)
414     {
415         int  ii = localAtomIndices[i];
416         real wm;
417         if (pgrp->localWeights.empty())
418         {
419             wm      = mass[ii];
420             sum_wm += wm;
421         }
422         else
423         {
424             real w;
425
426             w        = pgrp->localWeights[i];
427             wm       = w*mass[ii];
428             sum_wm  += wm;
429             sum_wwm += wm*w;
430         }
431         if (pgrp->epgrppbc == epgrppbcNONE)
432         {
433             /* Plain COM: sum the coordinates */
434             for (int d = 0; d < DIM; d++)
435             {
436                 sum_wmx[d]      += wm*x[ii][d];
437             }
438             if (xp)
439             {
440                 for (int d = 0; d < DIM; d++)
441                 {
442                     sum_wmxp[d] += wm*xp[ii][d];
443                 }
444             }
445         }
446         else
447         {
448             rvec dx;
449
450             /* Sum the difference with the reference atom */
451             pbc_dx(pbc, x[ii], x_pbc, dx);
452             for (int d = 0; d < DIM; d++)
453             {
454                 sum_wmx[d]     += wm*dx[d];
455             }
456             if (xp)
457             {
458                 /* For xp add the difference between xp and x to dx,
459                  * such that we use the same periodic image,
460                  * also when xp has a large displacement.
461                  */
462                 for (int d = 0; d < DIM; d++)
463                 {
464                     sum_wmxp[d] += wm*(dx[d] + xp[ii][d] - x[ii][d]);
465                 }
466             }
467         }
468     }
469
470     sum_com->sum_wm  = sum_wm;
471     sum_com->sum_wwm = sum_wwm;
472     copy_dvec(sum_wmx, sum_com->sum_wmx);
473     if (xp)
474     {
475         copy_dvec(sum_wmxp, sum_com->sum_wmxp);
476     }
477 }
478
479 static void sum_com_part_cosweight(const pull_group_work_t *pgrp,
480                                    int ind_start, int ind_end,
481                                    int cosdim, real twopi_box,
482                                    const rvec *x, const rvec *xp,
483                                    const real *mass,
484                                    ComSums *sum_com)
485 {
486     /* Cosine weighting geometry */
487     double sum_cm  = 0;
488     double sum_sm  = 0;
489     double sum_ccm = 0;
490     double sum_csm = 0;
491     double sum_ssm = 0;
492     double sum_cmp = 0;
493     double sum_smp = 0;
494
495     auto   localAtomIndices = pgrp->atomSet.localIndex();
496
497     for (int i = ind_start; i < ind_end; i++)
498     {
499         int  ii  = localAtomIndices[i];
500         real m   = mass[ii];
501         /* Determine cos and sin sums */
502         real cw  = std::cos(x[ii][cosdim]*twopi_box);
503         real sw  = std::sin(x[ii][cosdim]*twopi_box);
504         sum_cm  += static_cast<double>(cw*m);
505         sum_sm  += static_cast<double>(sw*m);
506         sum_ccm += static_cast<double>(cw*cw*m);
507         sum_csm += static_cast<double>(cw*sw*m);
508         sum_ssm += static_cast<double>(sw*sw*m);
509
510         if (xp != nullptr)
511         {
512             real cw  = std::cos(xp[ii][cosdim]*twopi_box);
513             real sw  = std::sin(xp[ii][cosdim]*twopi_box);
514             sum_cmp += static_cast<double>(cw*m);
515             sum_smp += static_cast<double>(sw*m);
516         }
517     }
518
519     sum_com->sum_cm  = sum_cm;
520     sum_com->sum_sm  = sum_sm;
521     sum_com->sum_ccm = sum_ccm;
522     sum_com->sum_csm = sum_csm;
523     sum_com->sum_ssm = sum_ssm;
524     sum_com->sum_cmp = sum_cmp;
525     sum_com->sum_smp = sum_smp;
526 }
527
528 /* calculates center of mass of selection index from all coordinates x */
529 void pull_calc_coms(const t_commrec *cr,
530                     pull_t *pull,
531                     const t_mdatoms *md,
532                     t_pbc *pbc,
533                     double t,
534                     const rvec x[], rvec *xp)
535 {
536     real         twopi_box = 0;
537     pull_comm_t *comm;
538
539     comm = &pull->comm;
540
541     GMX_ASSERT(comm->pbcAtomBuffer.size() == pull->group.size(), "pbcAtomBuffer should have size number of groups");
542     GMX_ASSERT(comm->comBuffer.size() == pull->group.size()*c_comBufferStride,
543                "comBuffer should have size #group*c_comBufferStride");
544
545     if (pull->bRefAt && pull->bSetPBCatoms)
546     {
547         pull_set_pbcatoms(cr, pull, x, comm->pbcAtomBuffer);
548
549         if (cr != nullptr && DOMAINDECOMP(cr))
550         {
551             /* We can keep these PBC reference coordinates fixed for nstlist
552              * steps, since atoms won't jump over PBC.
553              * This avoids a global reduction at the next nstlist-1 steps.
554              * Note that the exact values of the pbc reference coordinates
555              * are irrelevant, as long all atoms in the group are within
556              * half a box distance of the reference coordinate.
557              */
558             pull->bSetPBCatoms = FALSE;
559         }
560     }
561
562     if (pull->cosdim >= 0)
563     {
564         int m;
565
566         assert(pull->npbcdim <= DIM);
567
568         for (m = pull->cosdim+1; m < pull->npbcdim; m++)
569         {
570             if (pbc->box[m][pull->cosdim] != 0)
571             {
572                 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
573             }
574         }
575         twopi_box = 2.0*M_PI/pbc->box[pull->cosdim][pull->cosdim];
576     }
577
578     for (size_t g = 0; g < pull->group.size(); g++)
579     {
580         pull_group_work_t *pgrp      = &pull->group[g];
581
582         /* Cosine-weighted COMs behave different from all other weighted COMs
583          * in the sense that the weights depend on instantaneous coordinates,
584          * not on pre-set weights. Thus we resize the local weight buffer here.
585          */
586         if (pgrp->epgrppbc == epgrppbcCOS)
587         {
588             pgrp->localWeights.resize(pgrp->atomSet.localIndex().size());
589         }
590
591         auto               comBuffer =
592             gmx::arrayRefFromArray(comm->comBuffer.data() + g*c_comBufferStride, c_comBufferStride);
593
594         if (pgrp->needToCalcCom)
595         {
596             if (pgrp->epgrppbc != epgrppbcCOS)
597             {
598                 rvec x_pbc = { 0, 0, 0 };
599
600                 switch (pgrp->epgrppbc)
601                 {
602                     case epgrppbcREFAT:
603                         /* Set the pbc atom */
604                         copy_rvec(comm->pbcAtomBuffer[g], x_pbc);
605                         break;
606                     case epgrppbcPREVSTEPCOM:
607                         /* Set the pbc reference to the COM of the group of the last step */
608                         copy_dvec_to_rvec(pgrp->x_prev_step, comm->pbcAtomBuffer[g]);
609                         copy_dvec_to_rvec(pgrp->x_prev_step, x_pbc);
610                 }
611
612                 /* The final sums should end up in comSums[0] */
613                 ComSums &comSumsTotal = pull->comSums[0];
614
615                 /* If we have a single-atom group the mass is irrelevant, so
616                  * we can remove the mass factor to avoid division by zero.
617                  * Note that with constraint pulling the mass does matter, but
618                  * in that case a check group mass != 0 has been done before.
619                  */
620                 if (pgrp->params.nat == 1 &&
621                     pgrp->atomSet.numAtomsLocal() == 1 &&
622                     md->massT[pgrp->atomSet.localIndex()[0]] == 0)
623                 {
624                     GMX_ASSERT(xp == nullptr, "We should not have groups with zero mass with constraints, i.e. xp!=NULL");
625
626                     /* Copy the single atom coordinate */
627                     for (int d = 0; d < DIM; d++)
628                     {
629                         comSumsTotal.sum_wmx[d] = x[pgrp->atomSet.localIndex()[0]][d];
630                     }
631                     /* Set all mass factors to 1 to get the correct COM */
632                     comSumsTotal.sum_wm  = 1;
633                     comSumsTotal.sum_wwm = 1;
634                 }
635                 else if (pgrp->atomSet.numAtomsLocal() <= c_pullMaxNumLocalAtomsSingleThreaded)
636                 {
637                     sum_com_part(pgrp, 0, pgrp->atomSet.numAtomsLocal(),
638                                  x, xp, md->massT,
639                                  pbc, x_pbc,
640                                  &comSumsTotal);
641                 }
642                 else
643                 {
644 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
645                     for (int t = 0; t < pull->nthreads; t++)
646                     {
647                         int ind_start = (pgrp->atomSet.numAtomsLocal()*(t + 0))/pull->nthreads;
648                         int ind_end   = (pgrp->atomSet.numAtomsLocal()*(t + 1))/pull->nthreads;
649                         sum_com_part(pgrp, ind_start, ind_end,
650                                      x, xp, md->massT,
651                                      pbc, x_pbc,
652                                      &pull->comSums[t]);
653                     }
654
655                     /* Reduce the thread contributions to sum_com[0] */
656                     for (int t = 1; t < pull->nthreads; t++)
657                     {
658                         comSumsTotal.sum_wm  += pull->comSums[t].sum_wm;
659                         comSumsTotal.sum_wwm += pull->comSums[t].sum_wwm;
660                         dvec_inc(comSumsTotal.sum_wmx, pull->comSums[t].sum_wmx);
661                         dvec_inc(comSumsTotal.sum_wmxp, pull->comSums[t].sum_wmxp);
662                     }
663                 }
664
665                 if (pgrp->localWeights.empty())
666                 {
667                     comSumsTotal.sum_wwm = comSumsTotal.sum_wm;
668                 }
669
670                 /* Copy local sums to a buffer for global summing */
671                 copy_dvec(comSumsTotal.sum_wmx,  comBuffer[0]);
672
673                 copy_dvec(comSumsTotal.sum_wmxp, comBuffer[1]);
674
675                 comBuffer[2][0] = comSumsTotal.sum_wm;
676                 comBuffer[2][1] = comSumsTotal.sum_wwm;
677                 comBuffer[2][2] = 0;
678             }
679             else
680             {
681                 /* Cosine weighting geometry.
682                  * This uses a slab of the system, thus we always have many
683                  * atoms in the pull groups. Therefore, always use threads.
684                  */
685 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
686                 for (int t = 0; t < pull->nthreads; t++)
687                 {
688                     int ind_start = (pgrp->atomSet.numAtomsLocal()*(t + 0))/pull->nthreads;
689                     int ind_end   = (pgrp->atomSet.numAtomsLocal()*(t + 1))/pull->nthreads;
690                     sum_com_part_cosweight(pgrp, ind_start, ind_end,
691                                            pull->cosdim, twopi_box,
692                                            x, xp, md->massT,
693                                            &pull->comSums[t]);
694                 }
695
696                 /* Reduce the thread contributions to comSums[0] */
697                 ComSums &comSumsTotal = pull->comSums[0];
698                 for (int t = 1; t < pull->nthreads; t++)
699                 {
700                     comSumsTotal.sum_cm  += pull->comSums[t].sum_cm;
701                     comSumsTotal.sum_sm  += pull->comSums[t].sum_sm;
702                     comSumsTotal.sum_ccm += pull->comSums[t].sum_ccm;
703                     comSumsTotal.sum_csm += pull->comSums[t].sum_csm;
704                     comSumsTotal.sum_ssm += pull->comSums[t].sum_ssm;
705                     comSumsTotal.sum_cmp += pull->comSums[t].sum_cmp;
706                     comSumsTotal.sum_smp += pull->comSums[t].sum_smp;
707                 }
708
709                 /* Copy local sums to a buffer for global summing */
710                 comBuffer[0][0] = comSumsTotal.sum_cm;
711                 comBuffer[0][1] = comSumsTotal.sum_sm;
712                 comBuffer[0][2] = 0;
713                 comBuffer[1][0] = comSumsTotal.sum_ccm;
714                 comBuffer[1][1] = comSumsTotal.sum_csm;
715                 comBuffer[1][2] = comSumsTotal.sum_ssm;
716                 comBuffer[2][0] = comSumsTotal.sum_cmp;
717                 comBuffer[2][1] = comSumsTotal.sum_smp;
718                 comBuffer[2][2] = 0;
719             }
720         }
721         else
722         {
723             clear_dvec(comBuffer[0]);
724             clear_dvec(comBuffer[1]);
725             clear_dvec(comBuffer[2]);
726         }
727     }
728
729     pullAllReduce(cr, comm, pull->group.size()*c_comBufferStride*DIM,
730                   static_cast<double *>(comm->comBuffer[0]));
731
732     for (size_t g = 0; g < pull->group.size(); g++)
733     {
734         pull_group_work_t *pgrp;
735
736         pgrp = &pull->group[g];
737         if (pgrp->needToCalcCom)
738         {
739             GMX_ASSERT(pgrp->params.nat > 0, "Normal pull groups should have atoms, only group 0, which should have bCalcCom=FALSE has nat=0");
740
741             const auto comBuffer = gmx::constArrayRefFromArray(comm->comBuffer.data() + g*c_comBufferStride, c_comBufferStride);
742
743             if (pgrp->epgrppbc != epgrppbcCOS)
744             {
745                 double wmass, wwmass;
746                 int    m;
747
748                 /* Determine the inverse mass */
749                 wmass             = comBuffer[2][0];
750                 wwmass            = comBuffer[2][1];
751                 pgrp->mwscale     = 1.0/wmass;
752                 /* invtm==0 signals a frozen group, so then we should keep it zero */
753                 if (pgrp->invtm != 0)
754                 {
755                     pgrp->wscale  = wmass/wwmass;
756                     pgrp->invtm   = wwmass/(wmass*wmass);
757                 }
758                 /* Divide by the total mass */
759                 for (m = 0; m < DIM; m++)
760                 {
761                     pgrp->x[m]      = comBuffer[0][m]*pgrp->mwscale;
762                     if (xp)
763                     {
764                         pgrp->xp[m] = comBuffer[1][m]*pgrp->mwscale;
765                     }
766                     if (pgrp->epgrppbc == epgrppbcREFAT || pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
767                     {
768                         pgrp->x[m]      += comm->pbcAtomBuffer[g][m];
769                         if (xp)
770                         {
771                             pgrp->xp[m] += comm->pbcAtomBuffer[g][m];
772                         }
773                     }
774                 }
775             }
776             else
777             {
778                 /* Cosine weighting geometry */
779                 double csw, snw, wmass, wwmass;
780
781                 /* Determine the optimal location of the cosine weight */
782                 csw                   = comBuffer[0][0];
783                 snw                   = comBuffer[0][1];
784                 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
785                 /* Set the weights for the local atoms */
786                 wmass  = sqrt(csw*csw + snw*snw);
787                 wwmass = (comBuffer[1][0]*csw*csw +
788                           comBuffer[1][1]*csw*snw +
789                           comBuffer[1][2]*snw*snw)/(wmass*wmass);
790
791                 pgrp->mwscale = 1.0/wmass;
792                 pgrp->wscale  = wmass/wwmass;
793                 pgrp->invtm   = wwmass/(wmass*wmass);
794                 /* Set the weights for the local atoms */
795                 csw *= pgrp->invtm;
796                 snw *= pgrp->invtm;
797                 for (size_t i = 0; i < pgrp->atomSet.numAtomsLocal(); i++)
798                 {
799                     int ii = pgrp->atomSet.localIndex()[i];
800                     pgrp->localWeights[i] = csw*std::cos(twopi_box*x[ii][pull->cosdim]) +
801                         snw*std::sin(twopi_box*x[ii][pull->cosdim]);
802                 }
803                 if (xp)
804                 {
805                     csw                    = comBuffer[2][0];
806                     snw                    = comBuffer[2][1];
807                     pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
808                 }
809             }
810             if (debug)
811             {
812                 fprintf(debug, "Pull group %zu wmass %f invtm %f\n",
813                         g, 1.0/pgrp->mwscale, pgrp->invtm);
814             }
815         }
816     }
817
818     if (pull->bCylinder)
819     {
820         /* Calculate the COMs for the cyclinder reference groups */
821         make_cyl_refgrps(cr, pull, md, pbc, t, x);
822     }
823 }
824
825 using BoolVec = gmx::BasicVector<bool>;
826
827 /* Returns whether the pull group obeys the PBC restrictions */
828 static bool pullGroupObeysPbcRestrictions(const pull_group_work_t &group,
829                                           const BoolVec           &dimUsed,
830                                           const rvec              *x,
831                                           const t_pbc             &pbc,
832                                           const gmx::RVec         &x_pbc,
833                                           const real               pbcMargin)
834 {
835     /* Determine which dimensions are relevant for PBC */
836     BoolVec dimUsesPbc       = { false, false, false };
837     bool    pbcIsRectangular = true;
838     for (int d = 0; d < pbc.ndim_ePBC; d++)
839     {
840         if (dimUsed[d])
841         {
842             dimUsesPbc[d] = true;
843             /* All non-zero dimensions of vector v are involved in PBC */
844             for (int d2 = d + 1; d2 < pbc.ndim_ePBC; d2++)
845             {
846                 assert(d2 < DIM);
847                 if (pbc.box[d2][d] != 0)
848                 {
849                     dimUsesPbc[d2]   = true;
850                     pbcIsRectangular = false;
851                 }
852             }
853         }
854     }
855
856     rvec marginPerDim = {};
857     real marginDistance2 = 0;
858     if (pbcIsRectangular)
859     {
860         /* Use margins for dimensions independently */
861         for (int d = 0; d < pbc.ndim_ePBC; d++)
862         {
863             marginPerDim[d] = pbcMargin*pbc.hbox_diag[d];
864         }
865     }
866     else
867     {
868         /* Check the total distance along the relevant dimensions */
869         for (int d = 0; d < pbc.ndim_ePBC; d++)
870         {
871             if (dimUsesPbc[d])
872             {
873                 marginDistance2 += pbcMargin*gmx::square(0.5)*norm2(pbc.box[d]);
874             }
875         }
876     }
877
878     auto localAtomIndices = group.atomSet.localIndex();
879     for (gmx::index indexInSet = 0; indexInSet < localAtomIndices.ssize(); indexInSet++)
880     {
881         rvec dx;
882         pbc_dx(&pbc, x[localAtomIndices[indexInSet]], x_pbc, dx);
883
884         bool atomIsTooFar = false;
885         if (pbcIsRectangular)
886         {
887             for (int d = 0; d < pbc.ndim_ePBC; d++)
888             {
889                 if (dimUsesPbc[d] && (dx[d] < -marginPerDim[d] ||
890                                       dx[d] >  marginPerDim[d]))
891                 {
892                     atomIsTooFar = true;
893                 }
894             }
895         }
896         else
897         {
898             real pbcDistance2 = 0;
899             for (int d = 0; d < pbc.ndim_ePBC; d++)
900             {
901                 if (dimUsesPbc[d])
902                 {
903                     pbcDistance2 += gmx::square(dx[d]);
904                 }
905             }
906             atomIsTooFar = (pbcDistance2 > marginDistance2);
907         }
908         if (atomIsTooFar)
909         {
910             return false;
911         }
912     }
913
914     return true;
915 }
916
917 int pullCheckPbcWithinGroups(const pull_t &pull,
918                              const rvec   *x,
919                              const t_pbc  &pbc,
920                              real          pbcMargin)
921 {
922     if (pbc.ePBC == epbcNONE)
923     {
924         return -1;
925     }
926
927     /* Determine what dimensions are used for each group by pull coordinates */
928     std::vector<BoolVec> dimUsed(pull.group.size(), { false, false, false });
929     for (size_t c = 0; c < pull.coord.size(); c++)
930     {
931         const t_pull_coord &coordParams = pull.coord[c].params;
932         for (int groupIndex = 0; groupIndex < coordParams.ngroup; groupIndex++)
933         {
934             for (int d = 0; d < DIM; d++)
935             {
936                 if (coordParams.dim[d] &&
937                     !(coordParams.eGeom == epullgCYL && groupIndex == 0))
938                 {
939                     dimUsed[coordParams.group[groupIndex]][d] = true;
940                 }
941             }
942         }
943     }
944
945     /* Check PBC for every group that uses a PBC reference atom treatment */
946     for (size_t g = 0; g < pull.group.size(); g++)
947     {
948         const pull_group_work_t &group = pull.group[g];
949         if ((group.epgrppbc == epgrppbcREFAT || group.epgrppbc == epgrppbcPREVSTEPCOM) &&
950             !pullGroupObeysPbcRestrictions(group, dimUsed[g], x, pbc, pull.comm.pbcAtomBuffer[g], pbcMargin))
951         {
952             return g;
953         }
954     }
955
956     return -1;
957 }
958
959 bool pullCheckPbcWithinGroup(const pull_t                  &pull,
960                              gmx::ArrayRef<const gmx::RVec> x,
961                              const t_pbc                   &pbc,
962                              int                            groupNr,
963                              real                           pbcMargin)
964 {
965     if (pbc.ePBC == epbcNONE)
966     {
967         return true;
968     }
969     GMX_ASSERT(groupNr < gmx::ssize(pull.group), "groupNr is out of range");
970
971     /* Check PBC if the group uses a PBC reference atom treatment. */
972     const pull_group_work_t &group = pull.group[groupNr];
973     if (group.epgrppbc != epgrppbcREFAT && group.epgrppbc != epgrppbcPREVSTEPCOM)
974     {
975         return true;
976     }
977
978     /* Determine what dimensions are used for each group by pull coordinates */
979     BoolVec dimUsed = { false, false, false };
980     for (size_t c = 0; c < pull.coord.size(); c++)
981     {
982         const t_pull_coord &coordParams = pull.coord[c].params;
983         for (int groupIndex = 0; groupIndex < coordParams.ngroup; groupIndex++)
984         {
985             if (coordParams.group[groupIndex] == groupNr)
986             {
987                 for (int d = 0; d < DIM; d++)
988                 {
989                     if (coordParams.dim[d] &&
990                         !(coordParams.eGeom == epullgCYL && groupIndex == 0))
991                     {
992                         dimUsed[d] = true;
993                     }
994                 }
995             }
996         }
997     }
998
999     return (pullGroupObeysPbcRestrictions(group, dimUsed, as_rvec_array(x.data()), pbc, pull.comm.pbcAtomBuffer[groupNr], pbcMargin));
1000 }
1001
1002 void setPrevStepPullComFromState(struct pull_t *pull, const t_state *state)
1003 {
1004     for (size_t g = 0; g < pull->group.size(); g++)
1005     {
1006         for (int j = 0; j < DIM; j++)
1007         {
1008             pull->group[g].x_prev_step[j] = state->pull_com_prev_step[g*DIM+j];
1009         }
1010     }
1011 }
1012
1013 void updatePrevStepPullCom(struct pull_t *pull, t_state *state)
1014 {
1015     for (size_t g = 0; g < pull->group.size(); g++)
1016     {
1017         if (pull->group[g].needToCalcCom)
1018         {
1019             for (int j = 0; j < DIM; j++)
1020             {
1021                 pull->group[g].x_prev_step[j]      = pull->group[g].x[j];
1022                 state->pull_com_prev_step[g*DIM+j] = pull->group[g].x[j];
1023             }
1024         }
1025     }
1026 }
1027
1028 void allocStatePrevStepPullCom(t_state *state, const pull_t *pull)
1029 {
1030     if (!pull)
1031     {
1032         state->pull_com_prev_step.clear();
1033         return;
1034     }
1035     size_t ngroup = pull->group.size();
1036     if (state->pull_com_prev_step.size()/DIM != ngroup)
1037     {
1038         state->pull_com_prev_step.resize(ngroup * DIM, NAN);
1039     }
1040 }
1041
1042 void initPullComFromPrevStep(const t_commrec *cr,
1043                              pull_t          *pull,
1044                              const t_mdatoms *md,
1045                              t_pbc           *pbc,
1046                              const rvec       x[])
1047 {
1048     pull_comm_t *comm   = &pull->comm;
1049     size_t       ngroup = pull->group.size();
1050
1051     if (!comm->bParticipate)
1052     {
1053         return;
1054     }
1055
1056     GMX_ASSERT(comm->pbcAtomBuffer.size() == pull->group.size(), "pbcAtomBuffer should have size number of groups");
1057     GMX_ASSERT(comm->comBuffer.size() == pull->group.size()*c_comBufferStride,
1058                "comBuffer should have size #group*c_comBufferStride");
1059
1060     pull_set_pbcatoms(cr, pull, x, comm->pbcAtomBuffer);
1061
1062     for (size_t g = 0; g < ngroup; g++)
1063     {
1064         pull_group_work_t *pgrp;
1065
1066         pgrp = &pull->group[g];
1067
1068         if (pgrp->needToCalcCom && pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
1069         {
1070             GMX_ASSERT(pgrp->params.nat > 1, "Groups with no atoms, or only one atom, should not "
1071                        "use the COM from the previous step as reference.");
1072
1073             rvec x_pbc = { 0, 0, 0 };
1074             copy_rvec(comm->pbcAtomBuffer[g], x_pbc);
1075
1076             if (debug)
1077             {
1078                 fprintf(debug, "Initialising prev step COM of pull group %zu. x_pbc =", g);
1079                 for (int m = 0; m < DIM; m++)
1080                 {
1081                     fprintf(debug, " %f", x_pbc[m]);
1082                 }
1083                 fprintf(debug, "\n");
1084             }
1085
1086             /* The following is to a large extent similar to pull_calc_coms() */
1087
1088             /* The final sums should end up in sum_com[0] */
1089             ComSums &comSumsTotal = pull->comSums[0];
1090
1091             if (pgrp->atomSet.numAtomsLocal() <= c_pullMaxNumLocalAtomsSingleThreaded)
1092             {
1093                 sum_com_part(pgrp, 0, pgrp->atomSet.numAtomsLocal(),
1094                              x, nullptr, md->massT,
1095                              pbc, x_pbc,
1096                              &comSumsTotal);
1097             }
1098             else
1099             {
1100 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
1101                 for (int t = 0; t < pull->nthreads; t++)
1102                 {
1103                     int ind_start = (pgrp->atomSet.numAtomsLocal()*(t + 0))/pull->nthreads;
1104                     int ind_end   = (pgrp->atomSet.numAtomsLocal()*(t + 1))/pull->nthreads;
1105                     sum_com_part(pgrp, ind_start, ind_end,
1106                                  x, nullptr, md->massT,
1107                                  pbc, x_pbc,
1108                                  &pull->comSums[t]);
1109                 }
1110
1111                 /* Reduce the thread contributions to sum_com[0] */
1112                 for (int t = 1; t < pull->nthreads; t++)
1113                 {
1114                     comSumsTotal.sum_wm  += pull->comSums[t].sum_wm;
1115                     comSumsTotal.sum_wwm += pull->comSums[t].sum_wwm;
1116                     dvec_inc(comSumsTotal.sum_wmx, pull->comSums[t].sum_wmx);
1117                     dvec_inc(comSumsTotal.sum_wmxp, pull->comSums[t].sum_wmxp);
1118                 }
1119             }
1120
1121             if (pgrp->localWeights.empty())
1122             {
1123                 comSumsTotal.sum_wwm = comSumsTotal.sum_wm;
1124             }
1125
1126             /* Copy local sums to a buffer for global summing */
1127             auto localSums =
1128                 gmx::arrayRefFromArray(comm->comBuffer.data() + g*c_comBufferStride, c_comBufferStride);
1129
1130             localSums[0]    = comSumsTotal.sum_wmx;
1131             localSums[1]    = comSumsTotal.sum_wmxp;
1132             localSums[2][0] = comSumsTotal.sum_wm;
1133             localSums[2][1] = comSumsTotal.sum_wwm;
1134             localSums[2][2] = 0;
1135         }
1136     }
1137
1138     pullAllReduce(cr, comm, ngroup*c_comBufferStride*DIM, static_cast<double *>(comm->comBuffer[0]));
1139
1140     for (size_t g = 0; g < ngroup; g++)
1141     {
1142         pull_group_work_t *pgrp;
1143
1144         pgrp = &pull->group[g];
1145         if (pgrp->needToCalcCom)
1146         {
1147             if (pgrp->epgrppbc == epgrppbcPREVSTEPCOM)
1148             {
1149                 auto   localSums =
1150                     gmx::arrayRefFromArray(comm->comBuffer.data() + g*c_comBufferStride, c_comBufferStride);
1151                 double wmass, wwmass;
1152
1153                 /* Determine the inverse mass */
1154                 wmass             = localSums[2][0];
1155                 wwmass            = localSums[2][1];
1156                 pgrp->mwscale     = 1.0/wmass;
1157                 /* invtm==0 signals a frozen group, so then we should keep it zero */
1158                 if (pgrp->invtm != 0)
1159                 {
1160                     pgrp->wscale  = wmass/wwmass;
1161                     pgrp->invtm   = wwmass/(wmass*wmass);
1162                 }
1163                 /* Divide by the total mass */
1164                 for (int m = 0; m < DIM; m++)
1165                 {
1166                     pgrp->x[m]  = localSums[0][m]*pgrp->mwscale;
1167                     pgrp->x[m] += comm->pbcAtomBuffer[g][m];
1168                 }
1169                 if (debug)
1170                 {
1171                     fprintf(debug, "Pull group %zu wmass %f invtm %f\n",
1172                             g, 1.0/pgrp->mwscale, pgrp->invtm);
1173                     fprintf(debug, "Initialising prev step COM of pull group %zu to", g);
1174                     for (int m = 0; m < DIM; m++)
1175                     {
1176                         fprintf(debug, " %f", pgrp->x[m]);
1177                     }
1178                     fprintf(debug, "\n");
1179                 }
1180                 copy_dvec(pgrp->x, pgrp->x_prev_step);
1181             }
1182         }
1183     }
1184 }