Tidy: modernize-use-nullptr
[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, 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 <assert.h>
42 #include <stdlib.h>
43
44 #include "gromacs/domdec/domdec_struct.h"
45 #include "gromacs/domdec/ga2la.h"
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/gmxlib/network.h"
48 #include "gromacs/math/functions.h"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdtypes/commrec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/mdtypes/mdatom.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/pulling/pull.h"
56 #include "gromacs/pulling/pull_internal.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 static void pull_reduce_real(t_commrec   *cr,
64                              pull_comm_t *comm,
65                              int          n,
66                              real        *data)
67 {
68     if (cr != nullptr && PAR(cr))
69     {
70         if (comm->bParticipateAll)
71         {
72             /* Sum the contributions over all DD ranks */
73             gmx_sum(n, data, cr);
74         }
75         else
76         {
77 #if GMX_MPI
78 #if MPI_IN_PLACE_EXISTS
79             MPI_Allreduce(MPI_IN_PLACE, data, n, GMX_MPI_REAL, MPI_SUM,
80                           comm->mpi_comm_com);
81 #else
82             real *buf;
83
84             snew(buf, n);
85
86             MPI_Allreduce(data, buf, n, GMX_MPI_REAL, MPI_SUM,
87                           comm->mpi_comm_com);
88
89             /* Copy the result from the buffer to the input/output data */
90             for (int i = 0; i < n; i++)
91             {
92                 data[i] = buf[i];
93             }
94             sfree(buf);
95 #endif
96 #else
97             gmx_incons("comm->bParticipateAll=FALSE without GMX_MPI");
98 #endif
99         }
100     }
101 }
102
103 static void pull_reduce_double(t_commrec   *cr,
104                                pull_comm_t *comm,
105                                int          n,
106                                double      *data)
107 {
108     if (cr != nullptr && PAR(cr))
109     {
110         if (comm->bParticipateAll)
111         {
112             /* Sum the contributions over all DD ranks */
113             gmx_sumd(n, data, cr);
114         }
115         else
116         {
117 #if GMX_MPI
118 #if MPI_IN_PLACE_EXISTS
119             MPI_Allreduce(MPI_IN_PLACE, data, n, MPI_DOUBLE, MPI_SUM,
120                           comm->mpi_comm_com);
121 #else
122             double *buf;
123
124             snew(buf, n);
125
126             MPI_Allreduce(data, buf, n, MPI_DOUBLE, MPI_SUM,
127                           comm->mpi_comm_com);
128
129             /* Copy the result from the buffer to the input/output data */
130             for (int i = 0; i < n; i++)
131             {
132                 data[i] = buf[i];
133             }
134             sfree(buf);
135 #endif
136 #else
137             gmx_incons("comm->bParticipateAll=FALSE without GMX_MPI");
138 #endif
139         }
140     }
141 }
142
143 static void pull_set_pbcatom(t_commrec *cr, pull_group_work_t *pgrp,
144                              rvec *x,
145                              rvec x_pbc)
146 {
147     int a;
148
149     if (cr != nullptr && DOMAINDECOMP(cr))
150     {
151         if (ga2la_get_home(cr->dd->ga2la, pgrp->params.pbcatom, &a))
152         {
153             copy_rvec(x[a], x_pbc);
154         }
155         else
156         {
157             clear_rvec(x_pbc);
158         }
159     }
160     else
161     {
162         copy_rvec(x[pgrp->params.pbcatom], x_pbc);
163     }
164 }
165
166 static void pull_set_pbcatoms(t_commrec *cr, struct pull_t *pull,
167                               rvec *x,
168                               rvec *x_pbc)
169 {
170     int g, n;
171
172     n = 0;
173     for (g = 0; g < pull->ngroup; g++)
174     {
175         if (!pull->group[g].bCalcCOM || pull->group[g].params.pbcatom == -1)
176         {
177             clear_rvec(x_pbc[g]);
178         }
179         else
180         {
181             pull_set_pbcatom(cr, &pull->group[g], x, x_pbc[g]);
182             n++;
183         }
184     }
185
186     if (cr && PAR(cr) && n > 0)
187     {
188         /* Sum over participating ranks to get x_pbc from the home ranks.
189          * This can be very expensive at high parallelization, so we only
190          * do this after each DD repartitioning.
191          */
192         pull_reduce_real(cr, &pull->comm, pull->ngroup*DIM, x_pbc[0]);
193     }
194 }
195
196 static void make_cyl_refgrps(t_commrec *cr, struct pull_t *pull, t_mdatoms *md,
197                              t_pbc *pbc, double t, rvec *x)
198 {
199     /* The size and stride per coord for the reduction buffer */
200     const int       stride = 9;
201     int             c, i, ii, m, start, end;
202     rvec            g_x, dx, dir;
203     double          inv_cyl_r2;
204     pull_comm_t    *comm;
205     gmx_ga2la_t    *ga2la = nullptr;
206
207     comm = &pull->comm;
208
209     if (comm->dbuf_cyl == nullptr)
210     {
211         snew(comm->dbuf_cyl, pull->ncoord*stride);
212     }
213
214     if (cr && DOMAINDECOMP(cr))
215     {
216         ga2la = cr->dd->ga2la;
217     }
218
219     start = 0;
220     end   = md->homenr;
221
222     inv_cyl_r2 = 1.0/gmx::square(pull->params.cylinder_r);
223
224     /* loop over all groups to make a reference group for each*/
225     for (c = 0; c < pull->ncoord; c++)
226     {
227         pull_coord_work_t *pcrd;
228         double             sum_a, wmass, wwmass;
229         dvec               radf_fac0, radf_fac1;
230
231         pcrd   = &pull->coord[c];
232
233         sum_a  = 0;
234         wmass  = 0;
235         wwmass = 0;
236         clear_dvec(radf_fac0);
237         clear_dvec(radf_fac1);
238
239         if (pcrd->params.eGeom == epullgCYL)
240         {
241             pull_group_work_t *pref, *pgrp, *pdyna;
242
243             /* pref will be the same group for all pull coordinates */
244             pref  = &pull->group[pcrd->params.group[0]];
245             pgrp  = &pull->group[pcrd->params.group[1]];
246             pdyna = &pull->dyna[c];
247             copy_dvec_to_rvec(pcrd->vec, dir);
248             pdyna->nat_loc = 0;
249
250             /* We calculate distances with respect to the reference location
251              * of this cylinder group (g_x), which we already have now since
252              * we reduced the other group COM over the ranks. This resolves
253              * any PBC issues and we don't need to use a PBC-atom here.
254              */
255             if (pcrd->params.rate != 0)
256             {
257                 /* With rate=0, value_ref is set initially */
258                 pcrd->value_ref = pcrd->params.init + pcrd->params.rate*t;
259             }
260             for (m = 0; m < DIM; m++)
261             {
262                 g_x[m] = pgrp->x[m] - pcrd->vec[m]*pcrd->value_ref;
263             }
264
265             /* loop over all atoms in the main ref group */
266             for (i = 0; i < pref->params.nat; i++)
267             {
268                 ii = pref->params.ind[i];
269                 if (ga2la)
270                 {
271                     if (!ga2la_get_home(ga2la, pref->params.ind[i], &ii))
272                     {
273                         ii = -1;
274                     }
275                 }
276                 if (ii >= start && ii < end)
277                 {
278                     double dr2, dr2_rel, inp;
279                     dvec   dr;
280
281                     pbc_dx_aiuc(pbc, x[ii], g_x, dx);
282                     inp = iprod(dir, dx);
283                     dr2 = 0;
284                     for (m = 0; m < DIM; m++)
285                     {
286                         /* Determine the radial components */
287                         dr[m] = dx[m] - inp*dir[m];
288                         dr2  += dr[m]*dr[m];
289                     }
290                     dr2_rel = dr2*inv_cyl_r2;
291
292                     if (dr2_rel < 1)
293                     {
294                         double mass, weight, dweight_r;
295                         dvec   mdw;
296
297                         /* add to index, to sum of COM, to weight array */
298                         if (pdyna->nat_loc >= pdyna->nalloc_loc)
299                         {
300                             pdyna->nalloc_loc = over_alloc_large(pdyna->nat_loc+1);
301                             srenew(pdyna->ind_loc,    pdyna->nalloc_loc);
302                             srenew(pdyna->weight_loc, pdyna->nalloc_loc);
303                             srenew(pdyna->mdw,        pdyna->nalloc_loc);
304                             srenew(pdyna->dv,         pdyna->nalloc_loc);
305                         }
306                         pdyna->ind_loc[pdyna->nat_loc] = ii;
307
308                         mass      = md->massT[ii];
309                         /* The radial weight function is 1-2x^2+x^4,
310                          * where x=r/cylinder_r. Since this function depends
311                          * on the radial component, we also get radial forces
312                          * on both groups.
313                          */
314                         weight    = 1 + (-2 + dr2_rel)*dr2_rel;
315                         dweight_r = (-4 + 4*dr2_rel)*inv_cyl_r2;
316                         pdyna->weight_loc[pdyna->nat_loc] = weight;
317                         sum_a    += mass*weight*inp;
318                         wmass    += mass*weight;
319                         wwmass   += mass*weight*weight;
320                         dsvmul(mass*dweight_r, dr, mdw);
321                         copy_dvec(mdw, pdyna->mdw[pdyna->nat_loc]);
322                         /* Currently we only have the axial component of the
323                          * distance (inp) up to an unkown offset. We add this
324                          * offset after the reduction needs to determine the
325                          * COM of the cylinder group.
326                          */
327                         pdyna->dv[pdyna->nat_loc] = inp;
328                         for (m = 0; m < DIM; m++)
329                         {
330                             radf_fac0[m] += mdw[m];
331                             radf_fac1[m] += mdw[m]*inp;
332                         }
333                         pdyna->nat_loc++;
334                     }
335                 }
336             }
337         }
338         comm->dbuf_cyl[c*stride+0] = wmass;
339         comm->dbuf_cyl[c*stride+1] = wwmass;
340         comm->dbuf_cyl[c*stride+2] = sum_a;
341         comm->dbuf_cyl[c*stride+3] = radf_fac0[XX];
342         comm->dbuf_cyl[c*stride+4] = radf_fac0[YY];
343         comm->dbuf_cyl[c*stride+5] = radf_fac0[ZZ];
344         comm->dbuf_cyl[c*stride+6] = radf_fac1[XX];
345         comm->dbuf_cyl[c*stride+7] = radf_fac1[YY];
346         comm->dbuf_cyl[c*stride+8] = radf_fac1[ZZ];
347     }
348
349     if (cr != nullptr && PAR(cr))
350     {
351         /* Sum the contributions over the ranks */
352         pull_reduce_double(cr, comm, pull->ncoord*stride, comm->dbuf_cyl);
353     }
354
355     for (c = 0; c < pull->ncoord; c++)
356     {
357         pull_coord_work_t *pcrd;
358
359         pcrd  = &pull->coord[c];
360
361         if (pcrd->params.eGeom == epullgCYL)
362         {
363             pull_group_work_t *pdyna, *pgrp;
364             double             wmass, wwmass, dist;
365
366             pdyna = &pull->dyna[c];
367             pgrp  = &pull->group[pcrd->params.group[1]];
368
369             wmass          = comm->dbuf_cyl[c*stride+0];
370             wwmass         = comm->dbuf_cyl[c*stride+1];
371             pdyna->mwscale = 1.0/wmass;
372             /* Cylinder pulling can't be used with constraints, but we set
373              * wscale and invtm anyhow, in case someone would like to use them.
374              */
375             pdyna->wscale  = wmass/wwmass;
376             pdyna->invtm   = wwmass/(wmass*wmass);
377
378             /* We store the deviation of the COM from the reference location
379              * used above, since we need it when we apply the radial forces
380              * to the atoms in the cylinder group.
381              */
382             pcrd->cyl_dev  = 0;
383             for (m = 0; m < DIM; m++)
384             {
385                 g_x[m]         = pgrp->x[m] - pcrd->vec[m]*pcrd->value_ref;
386                 dist           = -pcrd->vec[m]*comm->dbuf_cyl[c*stride+2]*pdyna->mwscale;
387                 pdyna->x[m]    = g_x[m] - dist;
388                 pcrd->cyl_dev += dist;
389             }
390             /* Now we know the exact COM of the cylinder reference group,
391              * we can determine the radial force factor (ffrad) that when
392              * multiplied with the axial pull force will give the radial
393              * force on the pulled (non-cylinder) group.
394              */
395             for (m = 0; m < DIM; m++)
396             {
397                 pcrd->ffrad[m] = (comm->dbuf_cyl[c*stride+6+m] +
398                                   comm->dbuf_cyl[c*stride+3+m]*pcrd->cyl_dev)/wmass;
399             }
400
401             if (debug)
402             {
403                 fprintf(debug, "Pull cylinder group %d:%8.3f%8.3f%8.3f m:%8.3f\n",
404                         c, pdyna->x[0], pdyna->x[1],
405                         pdyna->x[2], 1.0/pdyna->invtm);
406                 fprintf(debug, "ffrad %8.3f %8.3f %8.3f\n",
407                         pcrd->ffrad[XX], pcrd->ffrad[YY], pcrd->ffrad[ZZ]);
408             }
409         }
410     }
411 }
412
413 static double atan2_0_2pi(double y, double x)
414 {
415     double a;
416
417     a = atan2(y, x);
418     if (a < 0)
419     {
420         a += 2.0*M_PI;
421     }
422     return a;
423 }
424
425 static void sum_com_part(const pull_group_work_t *pgrp,
426                          int ind_start, int ind_end,
427                          const rvec *x, const rvec *xp,
428                          const real *mass,
429                          const t_pbc *pbc,
430                          const rvec x_pbc,
431                          pull_sum_com_t *sum_com)
432 {
433     double sum_wm   = 0;
434     double sum_wwm  = 0;
435     dvec   sum_wmx  = { 0, 0, 0 };
436     dvec   sum_wmxp = { 0, 0, 0 };
437
438     for (int i = ind_start; i < ind_end; i++)
439     {
440         int  ii = pgrp->ind_loc[i];
441         real wm;
442         if (pgrp->weight_loc == nullptr)
443         {
444             wm      = mass[ii];
445             sum_wm += wm;
446         }
447         else
448         {
449             real w;
450
451             w        = pgrp->weight_loc[i];
452             wm       = w*mass[ii];
453             sum_wm  += wm;
454             sum_wwm += wm*w;
455         }
456         if (pgrp->epgrppbc == epgrppbcNONE)
457         {
458             /* Plain COM: sum the coordinates */
459             for (int d = 0; d < DIM; d++)
460             {
461                 sum_wmx[d]      += wm*x[ii][d];
462             }
463             if (xp)
464             {
465                 for (int d = 0; d < DIM; d++)
466                 {
467                     sum_wmxp[d] += wm*xp[ii][d];
468                 }
469             }
470         }
471         else
472         {
473             rvec dx;
474
475             /* Sum the difference with the reference atom */
476             pbc_dx(pbc, x[ii], x_pbc, dx);
477             for (int d = 0; d < DIM; d++)
478             {
479                 sum_wmx[d]     += wm*dx[d];
480             }
481             if (xp)
482             {
483                 /* For xp add the difference between xp and x to dx,
484                  * such that we use the same periodic image,
485                  * also when xp has a large displacement.
486                  */
487                 for (int d = 0; d < DIM; d++)
488                 {
489                     sum_wmxp[d] += wm*(dx[d] + xp[ii][d] - x[ii][d]);
490                 }
491             }
492         }
493     }
494
495     sum_com->sum_wm  = sum_wm;
496     sum_com->sum_wwm = sum_wwm;
497     copy_dvec(sum_wmx, sum_com->sum_wmx);
498     if (xp)
499     {
500         copy_dvec(sum_wmxp, sum_com->sum_wmxp);
501     }
502 }
503
504 static void sum_com_part_cosweight(const pull_group_work_t *pgrp,
505                                    int ind_start, int ind_end,
506                                    int cosdim, real twopi_box,
507                                    const rvec *x, const rvec *xp,
508                                    const real *mass,
509                                    pull_sum_com_t *sum_com)
510 {
511     /* Cosine weighting geometry */
512     double sum_cm  = 0;
513     double sum_sm  = 0;
514     double sum_ccm = 0;
515     double sum_csm = 0;
516     double sum_ssm = 0;
517     double sum_cmp = 0;
518     double sum_smp = 0;
519
520     for (int i = ind_start; i < ind_end; i++)
521     {
522         int  ii  = pgrp->ind_loc[i];
523         real m   = mass[ii];
524         /* Determine cos and sin sums */
525         real cw  = std::cos(x[ii][cosdim]*twopi_box);
526         real sw  = std::sin(x[ii][cosdim]*twopi_box);
527         sum_cm  += static_cast<double>(cw*m);
528         sum_sm  += static_cast<double>(sw*m);
529         sum_ccm += static_cast<double>(cw*cw*m);
530         sum_csm += static_cast<double>(cw*sw*m);
531         sum_ssm += static_cast<double>(sw*sw*m);
532
533         if (xp != nullptr)
534         {
535             real cw  = std::cos(xp[ii][cosdim]*twopi_box);
536             real sw  = std::sin(xp[ii][cosdim]*twopi_box);
537             sum_cmp += static_cast<double>(cw*m);
538             sum_smp += static_cast<double>(sw*m);
539         }
540     }
541
542     sum_com->sum_cm  = sum_cm;
543     sum_com->sum_sm  = sum_sm;
544     sum_com->sum_ccm = sum_ccm;
545     sum_com->sum_csm = sum_csm;
546     sum_com->sum_ssm = sum_ssm;
547     sum_com->sum_cmp = sum_cmp;
548     sum_com->sum_smp = sum_smp;
549 }
550
551 /* calculates center of mass of selection index from all coordinates x */
552 void pull_calc_coms(t_commrec *cr,
553                     struct pull_t *pull, t_mdatoms *md, t_pbc *pbc, double t,
554                     rvec x[], rvec *xp)
555 {
556     int          g;
557     real         twopi_box = 0;
558     pull_comm_t *comm;
559
560     comm = &pull->comm;
561
562     if (comm->rbuf == nullptr)
563     {
564         snew(comm->rbuf, pull->ngroup);
565     }
566     if (comm->dbuf == nullptr)
567     {
568         snew(comm->dbuf, 3*pull->ngroup);
569     }
570
571     if (pull->bRefAt && pull->bSetPBCatoms)
572     {
573         pull_set_pbcatoms(cr, pull, x, comm->rbuf);
574
575         if (cr != nullptr && DOMAINDECOMP(cr))
576         {
577             /* We can keep these PBC reference coordinates fixed for nstlist
578              * steps, since atoms won't jump over PBC.
579              * This avoids a global reduction at the next nstlist-1 steps.
580              * Note that the exact values of the pbc reference coordinates
581              * are irrelevant, as long all atoms in the group are within
582              * half a box distance of the reference coordinate.
583              */
584             pull->bSetPBCatoms = FALSE;
585         }
586     }
587
588     if (pull->cosdim >= 0)
589     {
590         int m;
591
592         assert(pull->npbcdim <= DIM);
593
594         for (m = pull->cosdim+1; m < pull->npbcdim; m++)
595         {
596             if (pbc->box[m][pull->cosdim] != 0)
597             {
598                 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
599             }
600         }
601         twopi_box = 2.0*M_PI/pbc->box[pull->cosdim][pull->cosdim];
602     }
603
604     for (g = 0; g < pull->ngroup; g++)
605     {
606         pull_group_work_t *pgrp;
607
608         pgrp = &pull->group[g];
609
610         if (pgrp->bCalcCOM)
611         {
612             if (pgrp->epgrppbc != epgrppbcCOS)
613             {
614                 rvec   x_pbc = { 0, 0, 0 };
615
616                 if (pgrp->epgrppbc == epgrppbcREFAT)
617                 {
618                     /* Set the pbc atom */
619                     copy_rvec(comm->rbuf[g], x_pbc);
620                 }
621
622                 /* The final sums should end up in sum_com[0] */
623                 pull_sum_com_t *sum_com = &pull->sum_com[0];
624
625                 /* If we have a single-atom group the mass is irrelevant, so
626                  * we can remove the mass factor to avoid division by zero.
627                  * Note that with constraint pulling the mass does matter, but
628                  * in that case a check group mass != 0 has been done before.
629                  */
630                 if (pgrp->params.nat == 1 &&
631                     pgrp->nat_loc == 1 &&
632                     md->massT[pgrp->ind_loc[0]] == 0)
633                 {
634                     GMX_ASSERT(xp == NULL, "We should not have groups with zero mass with constraints, i.e. xp!=NULL");
635
636                     /* Copy the single atom coordinate */
637                     for (int d = 0; d < DIM; d++)
638                     {
639                         sum_com->sum_wmx[d] = x[pgrp->ind_loc[0]][d];
640                     }
641                     /* Set all mass factors to 1 to get the correct COM */
642                     sum_com->sum_wm  = 1;
643                     sum_com->sum_wwm = 1;
644                 }
645                 else if (pgrp->nat_loc <= c_pullMaxNumLocalAtomsSingleThreaded)
646                 {
647                     sum_com_part(pgrp, 0, pgrp->nat_loc,
648                                  x, xp, md->massT,
649                                  pbc, x_pbc,
650                                  sum_com);
651                 }
652                 else
653                 {
654 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
655                     for (int t = 0; t < pull->nthreads; t++)
656                     {
657                         int ind_start = (pgrp->nat_loc*(t + 0))/pull->nthreads;
658                         int ind_end   = (pgrp->nat_loc*(t + 1))/pull->nthreads;
659                         sum_com_part(pgrp, ind_start, ind_end,
660                                      x, xp, md->massT,
661                                      pbc, x_pbc,
662                                      &pull->sum_com[t]);
663                     }
664
665                     /* Reduce the thread contributions to sum_com[0] */
666                     for (int t = 1; t < pull->nthreads; t++)
667                     {
668                         sum_com->sum_wm  += pull->sum_com[t].sum_wm;
669                         sum_com->sum_wwm += pull->sum_com[t].sum_wwm;
670                         dvec_inc(sum_com->sum_wmx, pull->sum_com[t].sum_wmx);
671                         dvec_inc(sum_com->sum_wmxp, pull->sum_com[t].sum_wmxp);
672                     }
673                 }
674
675                 if (pgrp->weight_loc == nullptr)
676                 {
677                     sum_com->sum_wwm = sum_com->sum_wm;
678                 }
679
680                 /* Copy local sums to a buffer for global summing */
681                 copy_dvec(sum_com->sum_wmx,  comm->dbuf[g*3]);
682                 copy_dvec(sum_com->sum_wmxp, comm->dbuf[g*3 + 1]);
683                 comm->dbuf[g*3 + 2][0] = sum_com->sum_wm;
684                 comm->dbuf[g*3 + 2][1] = sum_com->sum_wwm;
685                 comm->dbuf[g*3 + 2][2] = 0;
686             }
687             else
688             {
689                 /* Cosine weighting geometry.
690                  * This uses a slab of the system, thus we always have many
691                  * atoms in the pull groups. Therefore, always use threads.
692                  */
693 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
694                 for (int t = 0; t < pull->nthreads; t++)
695                 {
696                     int ind_start = (pgrp->nat_loc*(t + 0))/pull->nthreads;
697                     int ind_end   = (pgrp->nat_loc*(t + 1))/pull->nthreads;
698                     sum_com_part_cosweight(pgrp, ind_start, ind_end,
699                                            pull->cosdim, twopi_box,
700                                            x, xp, md->massT,
701                                            &pull->sum_com[t]);
702                 }
703
704                 /* Reduce the thread contributions to sum_com[0] */
705                 pull_sum_com_t *sum_com = &pull->sum_com[0];
706                 for (int t = 1; t < pull->nthreads; t++)
707                 {
708                     sum_com->sum_cm  += pull->sum_com[t].sum_cm;
709                     sum_com->sum_sm  += pull->sum_com[t].sum_sm;
710                     sum_com->sum_ccm += pull->sum_com[t].sum_ccm;
711                     sum_com->sum_csm += pull->sum_com[t].sum_csm;
712                     sum_com->sum_ssm += pull->sum_com[t].sum_ssm;
713                     sum_com->sum_cmp += pull->sum_com[t].sum_cmp;
714                     sum_com->sum_smp += pull->sum_com[t].sum_smp;
715                 }
716
717                 /* Copy local sums to a buffer for global summing */
718                 comm->dbuf[g*3    ][0] = sum_com->sum_cm;
719                 comm->dbuf[g*3    ][1] = sum_com->sum_sm;
720                 comm->dbuf[g*3    ][2] = 0;
721                 comm->dbuf[g*3 + 1][0] = sum_com->sum_ccm;
722                 comm->dbuf[g*3 + 1][1] = sum_com->sum_csm;
723                 comm->dbuf[g*3 + 1][2] = sum_com->sum_ssm;
724                 comm->dbuf[g*3 + 2][0] = sum_com->sum_cmp;
725                 comm->dbuf[g*3 + 2][1] = sum_com->sum_smp;
726                 comm->dbuf[g*3 + 2][2] = 0;
727             }
728         }
729     }
730
731     pull_reduce_double(cr, comm, pull->ngroup*3*DIM, comm->dbuf[0]);
732
733     for (g = 0; g < pull->ngroup; g++)
734     {
735         pull_group_work_t *pgrp;
736
737         pgrp = &pull->group[g];
738         if (pgrp->params.nat > 0 && pgrp->bCalcCOM)
739         {
740             if (pgrp->epgrppbc != epgrppbcCOS)
741             {
742                 double wmass, wwmass;
743                 int    m;
744
745                 /* Determine the inverse mass */
746                 wmass             = comm->dbuf[g*3+2][0];
747                 wwmass            = comm->dbuf[g*3+2][1];
748                 pgrp->mwscale     = 1.0/wmass;
749                 /* invtm==0 signals a frozen group, so then we should keep it zero */
750                 if (pgrp->invtm != 0)
751                 {
752                     pgrp->wscale  = wmass/wwmass;
753                     pgrp->invtm   = wwmass/(wmass*wmass);
754                 }
755                 /* Divide by the total mass */
756                 for (m = 0; m < DIM; m++)
757                 {
758                     pgrp->x[m]      = comm->dbuf[g*3  ][m]*pgrp->mwscale;
759                     if (xp)
760                     {
761                         pgrp->xp[m] = comm->dbuf[g*3+1][m]*pgrp->mwscale;
762                     }
763                     if (pgrp->epgrppbc == epgrppbcREFAT)
764                     {
765                         pgrp->x[m]      += comm->rbuf[g][m];
766                         if (xp)
767                         {
768                             pgrp->xp[m] += comm->rbuf[g][m];
769                         }
770                     }
771                 }
772             }
773             else
774             {
775                 /* Cosine weighting geometry */
776                 double csw, snw, wmass, wwmass;
777                 int    i, ii;
778
779                 /* Determine the optimal location of the cosine weight */
780                 csw                   = comm->dbuf[g*3][0];
781                 snw                   = comm->dbuf[g*3][1];
782                 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
783                 /* Set the weights for the local atoms */
784                 wmass  = sqrt(csw*csw + snw*snw);
785                 wwmass = (comm->dbuf[g*3+1][0]*csw*csw +
786                           comm->dbuf[g*3+1][1]*csw*snw +
787                           comm->dbuf[g*3+1][2]*snw*snw)/(wmass*wmass);
788
789                 pgrp->mwscale = 1.0/wmass;
790                 pgrp->wscale  = wmass/wwmass;
791                 pgrp->invtm   = wwmass/(wmass*wmass);
792                 /* Set the weights for the local atoms */
793                 csw *= pgrp->invtm;
794                 snw *= pgrp->invtm;
795                 for (i = 0; i < pgrp->nat_loc; i++)
796                 {
797                     ii                  = pgrp->ind_loc[i];
798                     pgrp->weight_loc[i] = csw*cos(twopi_box*x[ii][pull->cosdim]) +
799                         snw*sin(twopi_box*x[ii][pull->cosdim]);
800                 }
801                 if (xp)
802                 {
803                     csw                    = comm->dbuf[g*3+2][0];
804                     snw                    = comm->dbuf[g*3+2][1];
805                     pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
806                 }
807             }
808             if (debug)
809             {
810                 fprintf(debug, "Pull group %d wmass %f invtm %f\n",
811                         g, 1.0/pgrp->mwscale, pgrp->invtm);
812             }
813         }
814     }
815
816     if (pull->bCylinder)
817     {
818         /* Calculate the COMs for the cyclinder reference groups */
819         make_cyl_refgrps(cr, pull, md, pbc, t, x);
820     }
821 }