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, 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.
42 #include "gromacs/fileio/confio.h"
43 #include "gromacs/legacyheaders/gmx_ga2la.h"
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/legacyheaders/names.h"
46 #include "gromacs/legacyheaders/network.h"
47 #include "gromacs/legacyheaders/typedefs.h"
48 #include "gromacs/legacyheaders/types/commrec.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include "gromacs/pulling/pull.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/smalloc.h"
56 static void pull_set_pbcatom(t_commrec *cr, t_pull_group *pgrp,
62 if (cr != NULL && DOMAINDECOMP(cr))
64 if (ga2la_get_home(cr->dd->ga2la, pgrp->pbcatom, &a))
66 copy_rvec(x[a], x_pbc);
75 copy_rvec(x[pgrp->pbcatom], x_pbc);
79 static void pull_set_pbcatoms(t_commrec *cr, t_pull *pull,
86 for (g = 0; g < pull->ngroup; g++)
88 if (!pull->group[g].bCalcCOM || pull->group[g].pbcatom == -1)
94 pull_set_pbcatom(cr, &pull->group[g], x, x_pbc[g]);
99 if (cr && PAR(cr) && n > 0)
101 /* Sum over the nodes to get x_pbc from the home node of pbcatom */
102 gmx_sum(pull->ngroup*DIM, x_pbc[0], cr);
106 static void make_cyl_refgrps(t_commrec *cr, t_pull *pull, t_mdatoms *md,
107 t_pbc *pbc, double t, rvec *x)
109 /* The size and stride per coord for the reduction buffer */
110 const int stride = 9;
111 int c, i, ii, m, start, end;
115 t_pull_group *pref, *pgrp, *pdyna;
116 gmx_ga2la_t ga2la = NULL;
118 if (pull->dbuf_cyl == NULL)
120 snew(pull->dbuf_cyl, pull->ncoord*stride);
123 if (cr && DOMAINDECOMP(cr))
125 ga2la = cr->dd->ga2la;
131 inv_cyl_r2 = 1/dsqr(pull->cylinder_r);
133 /* loop over all groups to make a reference group for each*/
134 for (c = 0; c < pull->ncoord; c++)
136 double sum_a, wmass, wwmass;
137 dvec radf_fac0, radf_fac1;
139 pcrd = &pull->coord[c];
144 clear_dvec(radf_fac0);
145 clear_dvec(radf_fac1);
147 if (pcrd->eGeom == epullgCYL)
149 /* pref will be the same group for all pull coordinates */
150 pref = &pull->group[pcrd->group[0]];
151 pgrp = &pull->group[pcrd->group[1]];
152 pdyna = &pull->dyna[c];
153 copy_rvec(pcrd->vec, dir);
156 /* We calculate distances with respect to the reference location
157 * of this cylinder group (g_x), which we already have now since
158 * we reduced the other group COM over the ranks. This resolves
159 * any PBC issues and we don't need to use a PBC-atom here.
161 for (m = 0; m < DIM; m++)
163 g_x[m] = pgrp->x[m] - pcrd->vec[m]*(pcrd->init + pcrd->rate*t);
166 /* loop over all atoms in the main ref group */
167 for (i = 0; i < pref->nat; i++)
172 if (!ga2la_get_home(ga2la, pref->ind[i], &ii))
177 if (ii >= start && ii < end)
179 double dr2, dr2_rel, inp;
182 pbc_dx_aiuc(pbc, x[ii], g_x, dx);
183 inp = iprod(dir, dx);
185 for (m = 0; m < DIM; m++)
187 /* Determine the radial components */
188 dr[m] = dx[m] - inp*dir[m];
191 dr2_rel = dr2*inv_cyl_r2;
195 double mass, weight, dweight_r;
198 /* add to index, to sum of COM, to weight array */
199 if (pdyna->nat_loc >= pdyna->nalloc_loc)
201 pdyna->nalloc_loc = over_alloc_large(pdyna->nat_loc+1);
202 srenew(pdyna->ind_loc, pdyna->nalloc_loc);
203 srenew(pdyna->weight_loc, pdyna->nalloc_loc);
204 srenew(pdyna->mdw, pdyna->nalloc_loc);
205 srenew(pdyna->dv, pdyna->nalloc_loc);
207 pdyna->ind_loc[pdyna->nat_loc] = ii;
209 mass = md->massT[ii];
210 /* The radial weight function is 1-2x^2+x^4,
211 * where x=r/cylinder_r. Since this function depends
212 * on the radial component, we also get radial forces
215 weight = 1 + (-2 + dr2_rel)*dr2_rel;
216 dweight_r = (-4 + 4*dr2_rel)*inv_cyl_r2;
217 pdyna->weight_loc[pdyna->nat_loc] = weight;
218 sum_a += mass*weight*inp;
219 wmass += mass*weight;
220 wwmass += mass*weight*weight;
221 dsvmul(mass*dweight_r, dr, mdw);
222 copy_dvec(mdw, pdyna->mdw[pdyna->nat_loc]);
223 /* Currently we only have the axial component of the
224 * distance (inp) up to an unkown offset. We add this
225 * offset after the reduction needs to determine the
226 * COM of the cylinder group.
228 pdyna->dv[pdyna->nat_loc] = inp;
229 for (m = 0; m < DIM; m++)
231 radf_fac0[m] += mdw[m];
232 radf_fac1[m] += mdw[m]*inp;
239 pull->dbuf_cyl[c*stride+0] = wmass;
240 pull->dbuf_cyl[c*stride+1] = wwmass;
241 pull->dbuf_cyl[c*stride+2] = sum_a;
242 pull->dbuf_cyl[c*stride+3] = radf_fac0[XX];
243 pull->dbuf_cyl[c*stride+4] = radf_fac0[YY];
244 pull->dbuf_cyl[c*stride+5] = radf_fac0[ZZ];
245 pull->dbuf_cyl[c*stride+6] = radf_fac1[XX];
246 pull->dbuf_cyl[c*stride+7] = radf_fac1[YY];
247 pull->dbuf_cyl[c*stride+8] = radf_fac1[ZZ];
250 if (cr != NULL && PAR(cr))
252 /* Sum the contributions over the ranks */
253 gmx_sumd(pull->ncoord*stride, pull->dbuf_cyl, cr);
256 for (c = 0; c < pull->ncoord; c++)
258 pcrd = &pull->coord[c];
260 if (pcrd->eGeom == epullgCYL)
262 double wmass, wwmass, inp, dist;
264 pdyna = &pull->dyna[c];
265 pgrp = &pull->group[pcrd->group[1]];
267 wmass = pull->dbuf_cyl[c*stride+0];
268 wwmass = pull->dbuf_cyl[c*stride+1];
269 pdyna->mwscale = 1.0/wmass;
270 /* Cylinder pulling can't be used with constraints, but we set
271 * wscale and invtm anyhow, in case someone would like to use them.
273 pdyna->wscale = wmass/wwmass;
274 pdyna->invtm = wwmass/(wmass*wmass);
276 /* We store the deviation of the COM from the reference location
277 * used above, since we need it when we apply the radial forces
278 * to the atoms in the cylinder group.
281 for (m = 0; m < DIM; m++)
283 g_x[m] = pgrp->x[m] - pcrd->vec[m]*(pcrd->init + pcrd->rate*t);
284 dist = -pcrd->vec[m]*pull->dbuf_cyl[c*stride+2]*pdyna->mwscale;
285 pdyna->x[m] = g_x[m] - dist;
286 pcrd->cyl_dev += dist;
288 /* Now we know the exact COM of the cylinder reference group,
289 * we can determine the radial force factor (ffrad) that when
290 * multiplied with the axial pull force will give the radial
291 * force on the pulled (non-cylinder) group.
293 for (m = 0; m < DIM; m++)
295 pcrd->ffrad[m] = (pull->dbuf_cyl[c*stride+6+m] +
296 pull->dbuf_cyl[c*stride+3+m]*pcrd->cyl_dev)/wmass;
301 fprintf(debug, "Pull cylinder group %d:%8.3f%8.3f%8.3f m:%8.3f\n",
302 c, pdyna->x[0], pdyna->x[1],
303 pdyna->x[2], 1.0/pdyna->invtm);
304 fprintf(debug, "ffrad %8.3f %8.3f %8.3f\n",
305 pcrd->ffrad[XX], pcrd->ffrad[YY], pcrd->ffrad[ZZ]);
311 static double atan2_0_2pi(double y, double x)
323 /* calculates center of mass of selection index from all coordinates x */
324 void pull_calc_coms(t_commrec *cr,
325 t_pull *pull, t_mdatoms *md, t_pbc *pbc, double t,
331 if (pull->rbuf == NULL)
333 snew(pull->rbuf, pull->ngroup);
335 if (pull->dbuf == NULL)
337 snew(pull->dbuf, 3*pull->ngroup);
340 if (pull->bRefAt && pull->bSetPBCatoms)
342 pull_set_pbcatoms(cr, pull, x, pull->rbuf);
344 if (cr != NULL && DOMAINDECOMP(cr))
346 /* We can keep these PBC reference coordinates fixed for nstlist
347 * steps, since atoms won't jump over PBC.
348 * This avoids a global reduction at the next nstlist-1 steps.
349 * Note that the exact values of the pbc reference coordinates
350 * are irrelevant, as long all atoms in the group are within
351 * half a box distance of the reference coordinate.
353 pull->bSetPBCatoms = FALSE;
357 if (pull->cosdim >= 0)
361 assert(pull->npbcdim <= DIM);
363 for (m = pull->cosdim+1; m < pull->npbcdim; m++)
365 if (pbc->box[m][pull->cosdim] != 0)
367 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
370 twopi_box = 2.0*M_PI/pbc->box[pull->cosdim][pull->cosdim];
373 for (g = 0; g < pull->ngroup; g++)
377 pgrp = &pull->group[g];
381 if (pgrp->epgrppbc != epgrppbcCOS)
384 double wmass, wwmass;
385 rvec x_pbc = { 0, 0, 0 };
393 if (pgrp->epgrppbc == epgrppbcREFAT)
395 /* Set the pbc atom */
396 copy_rvec(pull->rbuf[g], x_pbc);
399 for (i = 0; i < pgrp->nat_loc; i++)
404 ii = pgrp->ind_loc[i];
405 mass = md->massT[ii];
406 if (pgrp->weight_loc == NULL)
415 w = pgrp->weight_loc[i];
420 if (pgrp->epgrppbc == epgrppbcNONE)
422 /* Plain COM: sum the coordinates */
423 for (m = 0; m < DIM; m++)
425 com[m] += wm*x[ii][m];
429 for (m = 0; m < DIM; m++)
431 comp[m] += wm*xp[ii][m];
439 /* Sum the difference with the reference atom */
440 pbc_dx(pbc, x[ii], x_pbc, dx);
441 for (m = 0; m < DIM; m++)
447 /* For xp add the difference between xp and x to dx,
448 * such that we use the same periodic image,
449 * also when xp has a large displacement.
451 for (m = 0; m < DIM; m++)
453 comp[m] += wm*(dx[m] + xp[ii][m] - x[ii][m]);
459 /* We do this check after the loop above to avoid more nesting.
460 * If we have a single-atom group the mass is irrelevant, so
461 * we can remove the mass factor to avoid division by zero.
462 * Note that with constraint pulling the mass does matter, but
463 * in that case a check group mass != 0 has been done before.
465 if (pgrp->nat == 1 && pgrp->nat_loc == 1 && wmass == 0)
469 /* Copy the single atom coordinate */
470 for (m = 0; m < DIM; m++)
472 com[m] = x[pgrp->ind_loc[0]][m];
474 /* Set all mass factors to 1 to get the correct COM */
479 if (pgrp->weight_loc == NULL)
484 /* Copy local sums to a buffer for global summing */
485 copy_dvec(com, pull->dbuf[g*3]);
486 copy_dvec(comp, pull->dbuf[g*3+1]);
487 pull->dbuf[g*3+2][0] = wmass;
488 pull->dbuf[g*3+2][1] = wwmass;
489 pull->dbuf[g*3+2][2] = 0;
493 /* Cosine weighting geometry */
494 double cm, sm, cmp, smp, ccm, csm, ssm, csw, snw;
505 for (i = 0; i < pgrp->nat_loc; i++)
510 ii = pgrp->ind_loc[i];
511 mass = md->massT[ii];
512 /* Determine cos and sin sums */
513 csw = cos(x[ii][pull->cosdim]*twopi_box);
514 snw = sin(x[ii][pull->cosdim]*twopi_box);
523 csw = cos(xp[ii][pull->cosdim]*twopi_box);
524 snw = sin(xp[ii][pull->cosdim]*twopi_box);
530 /* Copy local sums to a buffer for global summing */
531 pull->dbuf[g*3 ][0] = cm;
532 pull->dbuf[g*3 ][1] = sm;
533 pull->dbuf[g*3 ][2] = 0;
534 pull->dbuf[g*3+1][0] = ccm;
535 pull->dbuf[g*3+1][1] = csm;
536 pull->dbuf[g*3+1][2] = ssm;
537 pull->dbuf[g*3+2][0] = cmp;
538 pull->dbuf[g*3+2][1] = smp;
539 pull->dbuf[g*3+2][2] = 0;
546 /* Sum the contributions over the nodes */
547 gmx_sumd(pull->ngroup*3*DIM, pull->dbuf[0], cr);
550 for (g = 0; g < pull->ngroup; g++)
554 pgrp = &pull->group[g];
555 if (pgrp->nat > 0 && pgrp->bCalcCOM)
557 if (pgrp->epgrppbc != epgrppbcCOS)
559 double wmass, wwmass;
562 /* Determine the inverse mass */
563 wmass = pull->dbuf[g*3+2][0];
564 wwmass = pull->dbuf[g*3+2][1];
565 pgrp->mwscale = 1.0/wmass;
566 /* invtm==0 signals a frozen group, so then we should keep it zero */
567 if (pgrp->invtm != 0)
569 pgrp->wscale = wmass/wwmass;
570 pgrp->invtm = wwmass/(wmass*wmass);
572 /* Divide by the total mass */
573 for (m = 0; m < DIM; m++)
575 pgrp->x[m] = pull->dbuf[g*3 ][m]*pgrp->mwscale;
578 pgrp->xp[m] = pull->dbuf[g*3+1][m]*pgrp->mwscale;
580 if (pgrp->epgrppbc == epgrppbcREFAT)
582 pgrp->x[m] += pull->rbuf[g][m];
585 pgrp->xp[m] += pull->rbuf[g][m];
592 /* Cosine weighting geometry */
593 double csw, snw, wmass, wwmass;
596 /* Determine the optimal location of the cosine weight */
597 csw = pull->dbuf[g*3][0];
598 snw = pull->dbuf[g*3][1];
599 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
600 /* Set the weights for the local atoms */
601 wmass = sqrt(csw*csw + snw*snw);
602 wwmass = (pull->dbuf[g*3+1][0]*csw*csw +
603 pull->dbuf[g*3+1][1]*csw*snw +
604 pull->dbuf[g*3+1][2]*snw*snw)/(wmass*wmass);
606 pgrp->mwscale = 1.0/wmass;
607 pgrp->wscale = wmass/wwmass;
608 pgrp->invtm = wwmass/(wmass*wmass);
609 /* Set the weights for the local atoms */
612 for (i = 0; i < pgrp->nat_loc; i++)
614 ii = pgrp->ind_loc[i];
615 pgrp->weight_loc[i] = csw*cos(twopi_box*x[ii][pull->cosdim]) +
616 snw*sin(twopi_box*x[ii][pull->cosdim]);
620 csw = pull->dbuf[g*3+2][0];
621 snw = pull->dbuf[g*3+2][1];
622 pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
627 fprintf(debug, "Pull group %d wmass %f invtm %f\n",
628 g, 1.0/pgrp->mwscale, pgrp->invtm);
635 /* Calculate the COMs for the cyclinder reference groups */
636 make_cyl_refgrps(cr, pull, md, pbc, t, x);