2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018, 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.
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.
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.
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.
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.
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.
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/md_enums.h"
51 #include "gromacs/mdtypes/mdatom.h"
52 #include "gromacs/pbcutil/pbc.h"
53 #include "gromacs/pulling/pull.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/utility/gmxassert.h"
57 #include "gromacs/utility/real.h"
58 #include "gromacs/utility/smalloc.h"
60 #include "pull_internal.h"
64 // Helper function to deduce MPI datatype from the type of data
65 gmx_unused static MPI_Datatype mpiDatatype(const float gmx_unused *data)
70 // Helper function to deduce MPI datatype from the type of data
71 gmx_unused static MPI_Datatype mpiDatatype(const double gmx_unused *data)
79 // Helper function; note that gmx_sum(d) should actually be templated
80 gmx_unused static void gmxAllReduce(int n, real *data, const t_commrec *cr)
86 // Helper function; note that gmx_sum(d) should actually be templated
87 gmx_unused static void gmxAllReduce(int n, double *data, const t_commrec *cr)
89 gmx_sumd(n, data, cr);
92 // Reduce data of n elements over all ranks currently participating in pull
94 static void pullAllReduce(const t_commrec *cr,
99 if (cr != nullptr && PAR(cr))
101 if (comm->bParticipateAll)
103 /* Sum the contributions over all DD ranks */
104 gmxAllReduce(n, data, cr);
108 /* Separate branch because gmx_sum uses cr->mpi_comm_mygroup */
110 #if MPI_IN_PLACE_EXISTS
111 MPI_Allreduce(MPI_IN_PLACE, data, n, mpiDatatype(data), MPI_SUM,
114 std::vector<T> buf(n);
116 MPI_Allreduce(data, buf, n, mpiDatatype(data), MPI_SUM,
119 /* Copy the result from the buffer to the input/output data */
120 for (int i = 0; i < n; i++)
126 gmx_incons("comm->bParticipateAll=FALSE without GMX_MPI");
132 /* Copies the coordinates of the PBC atom of pgrp to x_pbc.
133 * When those coordinates are not available on this rank, clears x_pbc.
135 static void setPbcAtomCoords(const pull_group_work_t &pgrp,
139 if (pgrp.pbcAtomSet != nullptr)
141 if (pgrp.pbcAtomSet->numAtomsLocal() > 0)
143 /* We have the atom locally, copy its coordinates */
144 copy_rvec(x[pgrp.pbcAtomSet->localIndex()[0]], x_pbc);
148 /* Another rank has it, clear the coordinates for MPI_Allreduce */
154 copy_rvec(x[pgrp.params.pbcatom], x_pbc);
158 static void pull_set_pbcatoms(const t_commrec *cr, struct pull_t *pull,
163 for (size_t g = 0; g < pull->group.size(); g++)
165 const pull_group_work_t &group = pull->group[g];
166 if (group.needToCalcCom && group.epgrppbc == epgrppbcREFAT)
168 setPbcAtomCoords(pull->group[g], x, x_pbc[g]);
173 clear_rvec(x_pbc[g]);
177 if (cr && PAR(cr) && numPbcAtoms > 0)
179 /* Sum over participating ranks to get x_pbc from the home ranks.
180 * This can be very expensive at high parallelization, so we only
181 * do this after each DD repartitioning.
183 pullAllReduce(cr, &pull->comm, pull->group.size()*DIM, x_pbc[0]);
187 static void make_cyl_refgrps(const t_commrec *cr,
194 /* The size and stride per coord for the reduction buffer */
195 const int stride = 9;
201 if (comm->dbuf_cyl == nullptr)
203 snew(comm->dbuf_cyl, pull->coord.size()*stride);
206 inv_cyl_r2 = 1.0/gmx::square(pull->params.cylinder_r);
208 /* loop over all groups to make a reference group for each*/
209 for (size_t c = 0; c < pull->coord.size(); c++)
211 pull_coord_work_t *pcrd;
212 double sum_a, wmass, wwmass;
213 dvec radf_fac0, radf_fac1;
215 pcrd = &pull->coord[c];
220 clear_dvec(radf_fac0);
221 clear_dvec(radf_fac1);
223 if (pcrd->params.eGeom == epullgCYL)
225 /* pref will be the same group for all pull coordinates */
226 const pull_group_work_t &pref = pull->group[pcrd->params.group[0]];
227 const pull_group_work_t &pgrp = pull->group[pcrd->params.group[1]];
228 pull_group_work_t &pdyna = pull->dyna[c];
230 copy_dvec_to_rvec(pcrd->spatialData.vec, direction);
232 /* Since we have not calculated the COM of the cylinder group yet,
233 * we calculate distances with respect to location of the pull
234 * group minus the reference position along the vector.
235 * here we already have the COM of the pull group. This resolves
236 * any PBC issues and we don't need to use a PBC-atom here.
238 if (pcrd->params.rate != 0)
240 /* With rate=0, value_ref is set initially */
241 pcrd->value_ref = pcrd->params.init + pcrd->params.rate*t;
244 for (int m = 0; m < DIM; m++)
246 reference[m] = pgrp.x[m] - pcrd->spatialData.vec[m]*pcrd->value_ref;
249 auto localAtomIndices = pref.atomSet.localIndex();
251 /* This actually only needs to be done at init or DD time,
252 * but resizing with the same size does not cause much overhead.
254 pdyna.localWeights.resize(localAtomIndices.size());
255 pdyna.mdw.resize(localAtomIndices.size());
256 pdyna.dv.resize(localAtomIndices.size());
258 /* loop over all atoms in the main ref group */
259 for (gmx::index indexInSet = 0; indexInSet < localAtomIndices.size(); indexInSet++)
261 int atomIndex = localAtomIndices[indexInSet];
263 pbc_dx_aiuc(pbc, x[atomIndex], reference, dx);
264 double axialLocation = iprod(direction, dx);
267 for (int m = 0; m < DIM; m++)
269 /* Determine the radial components */
270 radialLocation[m] = dx[m] - axialLocation*direction[m];
271 dr2 += gmx::square(radialLocation[m]);
273 double dr2_rel = dr2*inv_cyl_r2;
277 /* add atom to sum of COM and to weight array */
279 double mass = md->massT[atomIndex];
280 /* The radial weight function is 1-2x^2+x^4,
281 * where x=r/cylinder_r. Since this function depends
282 * on the radial component, we also get radial forces
285 double weight = 1 + (-2 + dr2_rel)*dr2_rel;
286 double dweight_r = (-4 + 4*dr2_rel)*inv_cyl_r2;
287 pdyna.localWeights[indexInSet] = weight;
288 sum_a += mass*weight*axialLocation;
289 wmass += mass*weight;
290 wwmass += mass*weight*weight;
292 dsvmul(mass*dweight_r, radialLocation, mdw);
293 copy_dvec(mdw, pdyna.mdw[indexInSet]);
294 /* Currently we only have the axial component of the
295 * offset from the cylinder COM up to an unkown offset.
296 * We add this offset after the reduction needed
297 * for determining the COM of the cylinder group.
299 pdyna.dv[indexInSet] = axialLocation;
300 for (int m = 0; m < DIM; m++)
302 radf_fac0[m] += mdw[m];
303 radf_fac1[m] += mdw[m]*axialLocation;
308 pdyna.localWeights[indexInSet] = 0;
312 comm->dbuf_cyl[c*stride+0] = wmass;
313 comm->dbuf_cyl[c*stride+1] = wwmass;
314 comm->dbuf_cyl[c*stride+2] = sum_a;
315 comm->dbuf_cyl[c*stride+3] = radf_fac0[XX];
316 comm->dbuf_cyl[c*stride+4] = radf_fac0[YY];
317 comm->dbuf_cyl[c*stride+5] = radf_fac0[ZZ];
318 comm->dbuf_cyl[c*stride+6] = radf_fac1[XX];
319 comm->dbuf_cyl[c*stride+7] = radf_fac1[YY];
320 comm->dbuf_cyl[c*stride+8] = radf_fac1[ZZ];
323 if (cr != nullptr && PAR(cr))
325 /* Sum the contributions over the ranks */
326 pullAllReduce(cr, comm, pull->coord.size()*stride, comm->dbuf_cyl);
329 for (size_t c = 0; c < pull->coord.size(); c++)
331 pull_coord_work_t *pcrd;
333 pcrd = &pull->coord[c];
335 if (pcrd->params.eGeom == epullgCYL)
337 pull_group_work_t *pdyna = &pull->dyna[c];
338 pull_group_work_t *pgrp = &pull->group[pcrd->params.group[1]];
339 PullCoordSpatialData &spatialData = pcrd->spatialData;
341 double wmass = comm->dbuf_cyl[c*stride+0];
342 double wwmass = comm->dbuf_cyl[c*stride+1];
343 pdyna->mwscale = 1.0/wmass;
344 /* Cylinder pulling can't be used with constraints, but we set
345 * wscale and invtm anyhow, in case someone would like to use them.
347 pdyna->wscale = wmass/wwmass;
348 pdyna->invtm = wwmass/(wmass*wmass);
350 /* We store the deviation of the COM from the reference location
351 * used above, since we need it when we apply the radial forces
352 * to the atoms in the cylinder group.
354 spatialData.cyl_dev = 0;
355 for (int m = 0; m < DIM; m++)
357 double reference = pgrp->x[m] - spatialData.vec[m]*pcrd->value_ref;
358 double dist = -spatialData.vec[m]*comm->dbuf_cyl[c*stride+2]*pdyna->mwscale;
359 pdyna->x[m] = reference - dist;
360 spatialData.cyl_dev += dist;
362 /* Now we know the exact COM of the cylinder reference group,
363 * we can determine the radial force factor (ffrad) that when
364 * multiplied with the axial pull force will give the radial
365 * force on the pulled (non-cylinder) group.
367 for (int m = 0; m < DIM; m++)
369 spatialData.ffrad[m] = (comm->dbuf_cyl[c*stride+6+m] +
370 comm->dbuf_cyl[c*stride+3+m]*spatialData.cyl_dev)/wmass;
375 fprintf(debug, "Pull cylinder group %zu:%8.3f%8.3f%8.3f m:%8.3f\n",
376 c, pdyna->x[0], pdyna->x[1],
377 pdyna->x[2], 1.0/pdyna->invtm);
378 fprintf(debug, "ffrad %8.3f %8.3f %8.3f\n",
379 spatialData.ffrad[XX], spatialData.ffrad[YY], spatialData.ffrad[ZZ]);
385 static double atan2_0_2pi(double y, double x)
397 static void sum_com_part(const pull_group_work_t *pgrp,
398 int ind_start, int ind_end,
399 const rvec *x, const rvec *xp,
403 pull_sum_com_t *sum_com)
407 dvec sum_wmx = { 0, 0, 0 };
408 dvec sum_wmxp = { 0, 0, 0 };
410 auto localAtomIndices = pgrp->atomSet.localIndex();
411 for (int i = ind_start; i < ind_end; i++)
413 int ii = localAtomIndices[i];
415 if (pgrp->localWeights.empty())
424 w = pgrp->localWeights[i];
429 if (pgrp->epgrppbc == epgrppbcNONE)
431 /* Plain COM: sum the coordinates */
432 for (int d = 0; d < DIM; d++)
434 sum_wmx[d] += wm*x[ii][d];
438 for (int d = 0; d < DIM; d++)
440 sum_wmxp[d] += wm*xp[ii][d];
448 /* Sum the difference with the reference atom */
449 pbc_dx(pbc, x[ii], x_pbc, dx);
450 for (int d = 0; d < DIM; d++)
452 sum_wmx[d] += wm*dx[d];
456 /* For xp add the difference between xp and x to dx,
457 * such that we use the same periodic image,
458 * also when xp has a large displacement.
460 for (int d = 0; d < DIM; d++)
462 sum_wmxp[d] += wm*(dx[d] + xp[ii][d] - x[ii][d]);
468 sum_com->sum_wm = sum_wm;
469 sum_com->sum_wwm = sum_wwm;
470 copy_dvec(sum_wmx, sum_com->sum_wmx);
473 copy_dvec(sum_wmxp, sum_com->sum_wmxp);
477 static void sum_com_part_cosweight(const pull_group_work_t *pgrp,
478 int ind_start, int ind_end,
479 int cosdim, real twopi_box,
480 const rvec *x, const rvec *xp,
482 pull_sum_com_t *sum_com)
484 /* Cosine weighting geometry */
493 auto localAtomIndices = pgrp->atomSet.localIndex();
495 for (int i = ind_start; i < ind_end; i++)
497 int ii = localAtomIndices[i];
499 /* Determine cos and sin sums */
500 real cw = std::cos(x[ii][cosdim]*twopi_box);
501 real sw = std::sin(x[ii][cosdim]*twopi_box);
502 sum_cm += static_cast<double>(cw*m);
503 sum_sm += static_cast<double>(sw*m);
504 sum_ccm += static_cast<double>(cw*cw*m);
505 sum_csm += static_cast<double>(cw*sw*m);
506 sum_ssm += static_cast<double>(sw*sw*m);
510 real cw = std::cos(xp[ii][cosdim]*twopi_box);
511 real sw = std::sin(xp[ii][cosdim]*twopi_box);
512 sum_cmp += static_cast<double>(cw*m);
513 sum_smp += static_cast<double>(sw*m);
517 sum_com->sum_cm = sum_cm;
518 sum_com->sum_sm = sum_sm;
519 sum_com->sum_ccm = sum_ccm;
520 sum_com->sum_csm = sum_csm;
521 sum_com->sum_ssm = sum_ssm;
522 sum_com->sum_cmp = sum_cmp;
523 sum_com->sum_smp = sum_smp;
526 /* calculates center of mass of selection index from all coordinates x */
527 void pull_calc_coms(const t_commrec *cr,
532 const rvec x[], rvec *xp)
539 if (comm->rbuf == nullptr)
541 snew(comm->rbuf, pull->group.size());
543 if (comm->dbuf == nullptr)
545 snew(comm->dbuf, pull->group.size()*DIM);
548 if (pull->bRefAt && pull->bSetPBCatoms)
550 pull_set_pbcatoms(cr, pull, x, comm->rbuf);
552 if (cr != nullptr && DOMAINDECOMP(cr))
554 /* We can keep these PBC reference coordinates fixed for nstlist
555 * steps, since atoms won't jump over PBC.
556 * This avoids a global reduction at the next nstlist-1 steps.
557 * Note that the exact values of the pbc reference coordinates
558 * are irrelevant, as long all atoms in the group are within
559 * half a box distance of the reference coordinate.
561 pull->bSetPBCatoms = FALSE;
565 if (pull->cosdim >= 0)
569 assert(pull->npbcdim <= DIM);
571 for (m = pull->cosdim+1; m < pull->npbcdim; m++)
573 if (pbc->box[m][pull->cosdim] != 0)
575 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
578 twopi_box = 2.0*M_PI/pbc->box[pull->cosdim][pull->cosdim];
581 for (size_t g = 0; g < pull->group.size(); g++)
583 pull_group_work_t *pgrp;
585 pgrp = &pull->group[g];
587 if (pgrp->needToCalcCom)
589 if (pgrp->epgrppbc != epgrppbcCOS)
591 rvec x_pbc = { 0, 0, 0 };
593 if (pgrp->epgrppbc == epgrppbcREFAT)
595 /* Set the pbc atom */
596 copy_rvec(comm->rbuf[g], x_pbc);
599 /* The final sums should end up in sum_com[0] */
600 pull_sum_com_t *sum_com = &pull->sum_com[0];
602 /* If we have a single-atom group the mass is irrelevant, so
603 * we can remove the mass factor to avoid division by zero.
604 * Note that with constraint pulling the mass does matter, but
605 * in that case a check group mass != 0 has been done before.
607 if (pgrp->params.nat == 1 &&
608 pgrp->atomSet.numAtomsLocal() == 1 &&
609 md->massT[pgrp->atomSet.localIndex()[0]] == 0)
611 GMX_ASSERT(xp == nullptr, "We should not have groups with zero mass with constraints, i.e. xp!=NULL");
613 /* Copy the single atom coordinate */
614 for (int d = 0; d < DIM; d++)
616 sum_com->sum_wmx[d] = x[pgrp->atomSet.localIndex()[0]][d];
618 /* Set all mass factors to 1 to get the correct COM */
620 sum_com->sum_wwm = 1;
622 else if (pgrp->atomSet.numAtomsLocal() <= c_pullMaxNumLocalAtomsSingleThreaded)
624 sum_com_part(pgrp, 0, pgrp->atomSet.numAtomsLocal(),
631 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
632 for (int t = 0; t < pull->nthreads; t++)
634 int ind_start = (pgrp->atomSet.numAtomsLocal()*(t + 0))/pull->nthreads;
635 int ind_end = (pgrp->atomSet.numAtomsLocal()*(t + 1))/pull->nthreads;
636 sum_com_part(pgrp, ind_start, ind_end,
642 /* Reduce the thread contributions to sum_com[0] */
643 for (int t = 1; t < pull->nthreads; t++)
645 sum_com->sum_wm += pull->sum_com[t].sum_wm;
646 sum_com->sum_wwm += pull->sum_com[t].sum_wwm;
647 dvec_inc(sum_com->sum_wmx, pull->sum_com[t].sum_wmx);
648 dvec_inc(sum_com->sum_wmxp, pull->sum_com[t].sum_wmxp);
652 if (pgrp->localWeights.empty())
654 sum_com->sum_wwm = sum_com->sum_wm;
657 /* Copy local sums to a buffer for global summing */
658 copy_dvec(sum_com->sum_wmx, comm->dbuf[g*3]);
659 copy_dvec(sum_com->sum_wmxp, comm->dbuf[g*3 + 1]);
660 comm->dbuf[g*3 + 2][0] = sum_com->sum_wm;
661 comm->dbuf[g*3 + 2][1] = sum_com->sum_wwm;
662 comm->dbuf[g*3 + 2][2] = 0;
666 /* Cosine weighting geometry.
667 * This uses a slab of the system, thus we always have many
668 * atoms in the pull groups. Therefore, always use threads.
670 #pragma omp parallel for num_threads(pull->nthreads) schedule(static)
671 for (int t = 0; t < pull->nthreads; t++)
673 int ind_start = (pgrp->atomSet.numAtomsLocal()*(t + 0))/pull->nthreads;
674 int ind_end = (pgrp->atomSet.numAtomsLocal()*(t + 1))/pull->nthreads;
675 sum_com_part_cosweight(pgrp, ind_start, ind_end,
676 pull->cosdim, twopi_box,
681 /* Reduce the thread contributions to sum_com[0] */
682 pull_sum_com_t *sum_com = &pull->sum_com[0];
683 for (int t = 1; t < pull->nthreads; t++)
685 sum_com->sum_cm += pull->sum_com[t].sum_cm;
686 sum_com->sum_sm += pull->sum_com[t].sum_sm;
687 sum_com->sum_ccm += pull->sum_com[t].sum_ccm;
688 sum_com->sum_csm += pull->sum_com[t].sum_csm;
689 sum_com->sum_ssm += pull->sum_com[t].sum_ssm;
690 sum_com->sum_cmp += pull->sum_com[t].sum_cmp;
691 sum_com->sum_smp += pull->sum_com[t].sum_smp;
694 /* Copy local sums to a buffer for global summing */
695 comm->dbuf[g*3 ][0] = sum_com->sum_cm;
696 comm->dbuf[g*3 ][1] = sum_com->sum_sm;
697 comm->dbuf[g*3 ][2] = 0;
698 comm->dbuf[g*3 + 1][0] = sum_com->sum_ccm;
699 comm->dbuf[g*3 + 1][1] = sum_com->sum_csm;
700 comm->dbuf[g*3 + 1][2] = sum_com->sum_ssm;
701 comm->dbuf[g*3 + 2][0] = sum_com->sum_cmp;
702 comm->dbuf[g*3 + 2][1] = sum_com->sum_smp;
703 comm->dbuf[g*3 + 2][2] = 0;
708 pullAllReduce(cr, comm, pull->group.size()*3*DIM, comm->dbuf[0]);
710 for (size_t g = 0; g < pull->group.size(); g++)
712 pull_group_work_t *pgrp;
714 pgrp = &pull->group[g];
715 if (pgrp->needToCalcCom)
717 GMX_ASSERT(pgrp->params.nat > 0, "Normal pull groups should have atoms, only group 0, which should have bCalcCom=FALSE has nat=0");
719 if (pgrp->epgrppbc != epgrppbcCOS)
721 double wmass, wwmass;
724 /* Determine the inverse mass */
725 wmass = comm->dbuf[g*3+2][0];
726 wwmass = comm->dbuf[g*3+2][1];
727 pgrp->mwscale = 1.0/wmass;
728 /* invtm==0 signals a frozen group, so then we should keep it zero */
729 if (pgrp->invtm != 0)
731 pgrp->wscale = wmass/wwmass;
732 pgrp->invtm = wwmass/(wmass*wmass);
734 /* Divide by the total mass */
735 for (m = 0; m < DIM; m++)
737 pgrp->x[m] = comm->dbuf[g*3 ][m]*pgrp->mwscale;
740 pgrp->xp[m] = comm->dbuf[g*3+1][m]*pgrp->mwscale;
742 if (pgrp->epgrppbc == epgrppbcREFAT)
744 pgrp->x[m] += comm->rbuf[g][m];
747 pgrp->xp[m] += comm->rbuf[g][m];
754 /* Cosine weighting geometry */
755 double csw, snw, wmass, wwmass;
757 /* Determine the optimal location of the cosine weight */
758 csw = comm->dbuf[g*3][0];
759 snw = comm->dbuf[g*3][1];
760 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
761 /* Set the weights for the local atoms */
762 wmass = sqrt(csw*csw + snw*snw);
763 wwmass = (comm->dbuf[g*3+1][0]*csw*csw +
764 comm->dbuf[g*3+1][1]*csw*snw +
765 comm->dbuf[g*3+1][2]*snw*snw)/(wmass*wmass);
767 pgrp->mwscale = 1.0/wmass;
768 pgrp->wscale = wmass/wwmass;
769 pgrp->invtm = wwmass/(wmass*wmass);
770 /* Set the weights for the local atoms */
773 for (size_t i = 0; i < pgrp->atomSet.numAtomsLocal(); i++)
775 int ii = pgrp->atomSet.localIndex()[i];
776 pgrp->localWeights[i] = csw*std::cos(twopi_box*x[ii][pull->cosdim]) +
777 snw*std::sin(twopi_box*x[ii][pull->cosdim]);
781 csw = comm->dbuf[g*3+2][0];
782 snw = comm->dbuf[g*3+2][1];
783 pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
788 fprintf(debug, "Pull group %zu wmass %f invtm %f\n",
789 g, 1.0/pgrp->mwscale, pgrp->invtm);
796 /* Calculate the COMs for the cyclinder reference groups */
797 make_cyl_refgrps(cr, pull, md, pbc, t, x);