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.
41 #include "gromacs/fileio/confio.h"
42 #include "gromacs/legacyheaders/gmx_ga2la.h"
43 #include "gromacs/legacyheaders/macros.h"
44 #include "gromacs/legacyheaders/names.h"
45 #include "gromacs/legacyheaders/network.h"
46 #include "gromacs/legacyheaders/typedefs.h"
47 #include "gromacs/legacyheaders/types/commrec.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/pbcutil/pbc.h"
50 #include "gromacs/pulling/pull.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/futil.h"
53 #include "gromacs/utility/smalloc.h"
55 static void pull_set_pbcatom(t_commrec *cr, t_pull_group *pgrp,
61 if (cr != NULL && DOMAINDECOMP(cr))
63 if (ga2la_get_home(cr->dd->ga2la, pgrp->pbcatom, &a))
65 copy_rvec(x[a], x_pbc);
74 copy_rvec(x[pgrp->pbcatom], x_pbc);
78 static void pull_set_pbcatoms(t_commrec *cr, t_pull *pull,
85 for (g = 0; g < pull->ngroup; g++)
87 if (!pull->group[g].bCalcCOM || pull->group[g].pbcatom == -1)
93 pull_set_pbcatom(cr, &pull->group[g], x, x_pbc[g]);
98 if (cr && PAR(cr) && n > 0)
100 /* Sum over the nodes to get x_pbc from the home node of pbcatom */
101 gmx_sum(pull->ngroup*DIM, x_pbc[0], cr);
105 static void make_cyl_refgrps(t_commrec *cr, t_pull *pull, t_mdatoms *md,
106 t_pbc *pbc, double t, rvec *x)
108 /* The size and stride per coord for the reduction buffer */
109 const int stride = 9;
110 int c, i, ii, m, start, end;
114 t_pull_group *pref, *pgrp, *pdyna;
115 gmx_ga2la_t ga2la = NULL;
117 if (pull->dbuf_cyl == NULL)
119 snew(pull->dbuf_cyl, pull->ncoord*stride);
122 if (cr && DOMAINDECOMP(cr))
124 ga2la = cr->dd->ga2la;
130 inv_cyl_r2 = 1/dsqr(pull->cylinder_r);
132 /* loop over all groups to make a reference group for each*/
133 for (c = 0; c < pull->ncoord; c++)
135 double sum_a, wmass, wwmass;
136 dvec radf_fac0, radf_fac1;
138 pcrd = &pull->coord[c];
143 clear_dvec(radf_fac0);
144 clear_dvec(radf_fac1);
146 if (pcrd->eGeom == epullgCYL)
148 /* pref will be the same group for all pull coordinates */
149 pref = &pull->group[pcrd->group[0]];
150 pgrp = &pull->group[pcrd->group[1]];
151 pdyna = &pull->dyna[c];
152 copy_rvec(pcrd->vec, dir);
155 /* We calculate distances with respect to the reference location
156 * of this cylinder group (g_x), which we already have now since
157 * we reduced the other group COM over the ranks. This resolves
158 * any PBC issues and we don't need to use a PBC-atom here.
160 for (m = 0; m < DIM; m++)
162 g_x[m] = pgrp->x[m] - pcrd->vec[m]*(pcrd->init + pcrd->rate*t);
165 /* loop over all atoms in the main ref group */
166 for (i = 0; i < pref->nat; i++)
171 if (!ga2la_get_home(ga2la, pref->ind[i], &ii))
176 if (ii >= start && ii < end)
178 double dr2, dr2_rel, inp;
181 pbc_dx_aiuc(pbc, x[ii], g_x, dx);
182 inp = iprod(dir, dx);
184 for (m = 0; m < DIM; m++)
186 /* Determine the radial components */
187 dr[m] = dx[m] - inp*dir[m];
190 dr2_rel = dr2*inv_cyl_r2;
194 double mass, weight, dweight_r;
197 /* add to index, to sum of COM, to weight array */
198 if (pdyna->nat_loc >= pdyna->nalloc_loc)
200 pdyna->nalloc_loc = over_alloc_large(pdyna->nat_loc+1);
201 srenew(pdyna->ind_loc, pdyna->nalloc_loc);
202 srenew(pdyna->weight_loc, pdyna->nalloc_loc);
203 srenew(pdyna->mdw, pdyna->nalloc_loc);
204 srenew(pdyna->dv, pdyna->nalloc_loc);
206 pdyna->ind_loc[pdyna->nat_loc] = ii;
208 mass = md->massT[ii];
209 /* The radial weight function is 1-2x^2+x^4,
210 * where x=r/cylinder_r. Since this function depends
211 * on the radial component, we also get radial forces
214 weight = 1 + (-2 + dr2_rel)*dr2_rel;
215 dweight_r = (-4 + 4*dr2_rel)*inv_cyl_r2;
216 pdyna->weight_loc[pdyna->nat_loc] = weight;
217 sum_a += mass*weight*inp;
218 wmass += mass*weight;
219 wwmass += mass*weight*weight;
220 dsvmul(mass*dweight_r, dr, mdw);
221 copy_dvec(mdw, pdyna->mdw[pdyna->nat_loc]);
222 /* Currently we only have the axial component of the
223 * distance (inp) up to an unkown offset. We add this
224 * offset after the reduction needs to determine the
225 * COM of the cylinder group.
227 pdyna->dv[pdyna->nat_loc] = inp;
228 for (m = 0; m < DIM; m++)
230 radf_fac0[m] += mdw[m];
231 radf_fac1[m] += mdw[m]*inp;
238 pull->dbuf_cyl[c*stride+0] = wmass;
239 pull->dbuf_cyl[c*stride+1] = wwmass;
240 pull->dbuf_cyl[c*stride+2] = sum_a;
241 pull->dbuf_cyl[c*stride+3] = radf_fac0[XX];
242 pull->dbuf_cyl[c*stride+4] = radf_fac0[YY];
243 pull->dbuf_cyl[c*stride+5] = radf_fac0[ZZ];
244 pull->dbuf_cyl[c*stride+6] = radf_fac1[XX];
245 pull->dbuf_cyl[c*stride+7] = radf_fac1[YY];
246 pull->dbuf_cyl[c*stride+8] = radf_fac1[ZZ];
249 if (cr != NULL && PAR(cr))
251 /* Sum the contributions over the ranks */
252 gmx_sumd(pull->ncoord*stride, pull->dbuf_cyl, cr);
255 for (c = 0; c < pull->ncoord; c++)
257 pcrd = &pull->coord[c];
259 if (pcrd->eGeom == epullgCYL)
261 double wmass, wwmass, inp, dist;
263 pdyna = &pull->dyna[c];
264 pgrp = &pull->group[pcrd->group[1]];
266 wmass = pull->dbuf_cyl[c*stride+0];
267 wwmass = pull->dbuf_cyl[c*stride+1];
268 pdyna->mwscale = 1.0/wmass;
269 /* Cylinder pulling can't be used with constraints, but we set
270 * wscale and invtm anyhow, in case someone would like to use them.
272 pdyna->wscale = wmass/wwmass;
273 pdyna->invtm = wwmass/(wmass*wmass);
275 /* We store the deviation of the COM from the reference location
276 * used above, since we need it when we apply the radial forces
277 * to the atoms in the cylinder group.
280 for (m = 0; m < DIM; m++)
282 g_x[m] = pgrp->x[m] - pcrd->vec[m]*(pcrd->init + pcrd->rate*t);
283 dist = -pcrd->vec[m]*pull->dbuf_cyl[c*stride+2]*pdyna->mwscale;
284 pdyna->x[m] = g_x[m] - dist;
285 pcrd->cyl_dev += dist;
287 /* Now we know the exact COM of the cylinder reference group,
288 * we can determine the radial force factor (ffrad) that when
289 * multiplied with the axial pull force will give the radial
290 * force on the pulled (non-cylinder) group.
292 for (m = 0; m < DIM; m++)
294 pcrd->ffrad[m] = (pull->dbuf_cyl[c*stride+6+m] +
295 pull->dbuf_cyl[c*stride+3+m]*pcrd->cyl_dev)/wmass;
300 fprintf(debug, "Pull cylinder group %d:%8.3f%8.3f%8.3f m:%8.3f\n",
301 c, pdyna->x[0], pdyna->x[1],
302 pdyna->x[2], 1.0/pdyna->invtm);
303 fprintf(debug, "ffrad %8.3f %8.3f %8.3f\n",
304 pcrd->ffrad[XX], pcrd->ffrad[YY], pcrd->ffrad[ZZ]);
310 static double atan2_0_2pi(double y, double x)
322 /* calculates center of mass of selection index from all coordinates x */
323 void pull_calc_coms(t_commrec *cr,
324 t_pull *pull, t_mdatoms *md, t_pbc *pbc, double t,
328 real mass, w, wm, twopi_box = 0;
329 double wmass, wwmass;
331 double cm, sm, cmp, smp, ccm, csm, ssm, csw, snw;
332 rvec x_pbc = {0, 0, 0}, dx;
335 if (pull->rbuf == NULL)
337 snew(pull->rbuf, pull->ngroup);
339 if (pull->dbuf == NULL)
341 snew(pull->dbuf, 3*pull->ngroup);
344 if (pull->bRefAt && pull->bSetPBCatoms)
346 pull_set_pbcatoms(cr, pull, x, pull->rbuf);
348 if (cr != NULL && DOMAINDECOMP(cr))
350 /* We can keep these PBC reference coordinates fixed for nstlist
351 * steps, since atoms won't jump over PBC.
352 * This avoids a global reduction at the next nstlist-1 steps.
353 * Note that the exact values of the pbc reference coordinates
354 * are irrelevant, as long all atoms in the group are within
355 * half a box distance of the reference coordinate.
357 pull->bSetPBCatoms = FALSE;
361 if (pull->cosdim >= 0)
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++)
375 pgrp = &pull->group[g];
389 if (pgrp->epgrppbc == epgrppbcREFAT)
391 /* Set the pbc atom */
392 copy_rvec(pull->rbuf[g], x_pbc);
395 for (i = 0; i < pgrp->nat_loc; i++)
397 ii = pgrp->ind_loc[i];
398 mass = md->massT[ii];
399 if (pgrp->epgrppbc != epgrppbcCOS)
401 if (pgrp->weight_loc)
403 w = pgrp->weight_loc[i];
408 if (pgrp->epgrppbc == epgrppbcNONE)
410 /* Plain COM: sum the coordinates */
411 for (m = 0; m < DIM; m++)
413 com[m] += wm*x[ii][m];
417 for (m = 0; m < DIM; m++)
419 comp[m] += wm*xp[ii][m];
425 /* Sum the difference with the reference atom */
426 pbc_dx(pbc, x[ii], x_pbc, dx);
427 for (m = 0; m < DIM; m++)
433 /* For xp add the difference between xp and x to dx,
434 * such that we use the same periodic image,
435 * also when xp has a large displacement.
437 for (m = 0; m < DIM; m++)
439 comp[m] += wm*(dx[m] + xp[ii][m] - x[ii][m]);
446 /* Determine cos and sin sums */
447 csw = cos(x[ii][pull->cosdim]*twopi_box);
448 snw = sin(x[ii][pull->cosdim]*twopi_box);
457 csw = cos(xp[ii][pull->cosdim]*twopi_box);
458 snw = sin(xp[ii][pull->cosdim]*twopi_box);
466 /* Copy local sums to a buffer for global summing */
467 switch (pgrp->epgrppbc)
471 copy_dvec(com, pull->dbuf[g*3]);
472 copy_dvec(comp, pull->dbuf[g*3+1]);
473 pull->dbuf[g*3+2][0] = wmass;
474 pull->dbuf[g*3+2][1] = wwmass;
475 pull->dbuf[g*3+2][2] = 0;
478 pull->dbuf[g*3 ][0] = cm;
479 pull->dbuf[g*3 ][1] = sm;
480 pull->dbuf[g*3 ][2] = 0;
481 pull->dbuf[g*3+1][0] = ccm;
482 pull->dbuf[g*3+1][1] = csm;
483 pull->dbuf[g*3+1][2] = ssm;
484 pull->dbuf[g*3+2][0] = cmp;
485 pull->dbuf[g*3+2][1] = smp;
486 pull->dbuf[g*3+2][2] = 0;
493 /* Sum the contributions over the nodes */
494 gmx_sumd(pull->ngroup*3*DIM, pull->dbuf[0], cr);
497 for (g = 0; g < pull->ngroup; g++)
499 pgrp = &pull->group[g];
500 if (pgrp->nat > 0 && pgrp->bCalcCOM)
502 if (pgrp->epgrppbc != epgrppbcCOS)
504 /* Determine the inverse mass */
505 wmass = pull->dbuf[g*3+2][0];
506 wwmass = pull->dbuf[g*3+2][1];
507 pgrp->mwscale = 1.0/wmass;
508 /* invtm==0 signals a frozen group, so then we should keep it zero */
509 if (pgrp->invtm != 0)
511 pgrp->wscale = wmass/wwmass;
512 pgrp->invtm = wwmass/(wmass*wmass);
514 /* Divide by the total mass */
515 for (m = 0; m < DIM; m++)
517 pgrp->x[m] = pull->dbuf[g*3 ][m]*pgrp->mwscale;
520 pgrp->xp[m] = pull->dbuf[g*3+1][m]*pgrp->mwscale;
522 if (pgrp->epgrppbc == epgrppbcREFAT)
524 pgrp->x[m] += pull->rbuf[g][m];
527 pgrp->xp[m] += pull->rbuf[g][m];
534 /* Determine the optimal location of the cosine weight */
535 csw = pull->dbuf[g*3][0];
536 snw = pull->dbuf[g*3][1];
537 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
538 /* Set the weights for the local atoms */
539 wmass = sqrt(csw*csw + snw*snw);
540 wwmass = (pull->dbuf[g*3+1][0]*csw*csw +
541 pull->dbuf[g*3+1][1]*csw*snw +
542 pull->dbuf[g*3+1][2]*snw*snw)/(wmass*wmass);
544 pgrp->mwscale = 1.0/wmass;
545 pgrp->wscale = wmass/wwmass;
546 pgrp->invtm = wwmass/(wmass*wmass);
547 /* Set the weights for the local atoms */
550 for (i = 0; i < pgrp->nat_loc; i++)
552 ii = pgrp->ind_loc[i];
553 pgrp->weight_loc[i] = csw*cos(twopi_box*x[ii][pull->cosdim]) +
554 snw*sin(twopi_box*x[ii][pull->cosdim]);
558 csw = pull->dbuf[g*3+2][0];
559 snw = pull->dbuf[g*3+2][1];
560 pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
565 fprintf(debug, "Pull group %d wmass %f wwmass %f invtm %f\n",
566 g, wmass, wwmass, pgrp->invtm);
573 /* Calculate the COMs for the cyclinder reference groups */
574 make_cyl_refgrps(cr, pull, md, pbc, t, x);