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