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/pulling/pull_internal.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
57 static void pull_set_pbcatom(t_commrec *cr, pull_group_work_t *pgrp,
63 if (cr != NULL && DOMAINDECOMP(cr))
65 if (ga2la_get_home(cr->dd->ga2la, pgrp->params.pbcatom, &a))
67 copy_rvec(x[a], x_pbc);
76 copy_rvec(x[pgrp->params.pbcatom], x_pbc);
80 static void pull_set_pbcatoms(t_commrec *cr, struct pull_t *pull,
87 for (g = 0; g < pull->ngroup; g++)
89 if (!pull->group[g].bCalcCOM || pull->group[g].params.pbcatom == -1)
95 pull_set_pbcatom(cr, &pull->group[g], x, x_pbc[g]);
100 if (cr && PAR(cr) && n > 0)
102 /* Sum over the nodes to get x_pbc from the home node of pbcatom */
103 gmx_sum(pull->ngroup*DIM, x_pbc[0], cr);
107 static void make_cyl_refgrps(t_commrec *cr, struct pull_t *pull, t_mdatoms *md,
108 t_pbc *pbc, double t, rvec *x)
110 /* The size and stride per coord for the reduction buffer */
111 const int stride = 9;
112 int c, i, ii, m, start, end;
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->params.cylinder_r);
132 /* loop over all groups to make a reference group for each*/
133 for (c = 0; c < pull->ncoord; c++)
135 pull_coord_work_t *pcrd;
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->params.eGeom == epullgCYL)
149 pull_group_work_t *pref, *pgrp, *pdyna;
151 /* pref will be the same group for all pull coordinates */
152 pref = &pull->group[pcrd->params.group[0]];
153 pgrp = &pull->group[pcrd->params.group[1]];
154 pdyna = &pull->dyna[c];
155 copy_rvec(pcrd->vec, dir);
158 /* We calculate distances with respect to the reference location
159 * of this cylinder group (g_x), which we already have now since
160 * we reduced the other group COM over the ranks. This resolves
161 * any PBC issues and we don't need to use a PBC-atom here.
163 if (pcrd->params.rate != 0)
165 /* With rate=0, value_ref is set initially */
166 pcrd->value_ref = pcrd->params.init + pcrd->params.rate*t;
168 for (m = 0; m < DIM; m++)
170 g_x[m] = pgrp->x[m] - pcrd->vec[m]*pcrd->value_ref;
173 /* loop over all atoms in the main ref group */
174 for (i = 0; i < pref->params.nat; i++)
176 ii = pref->params.ind[i];
179 if (!ga2la_get_home(ga2la, pref->params.ind[i], &ii))
184 if (ii >= start && ii < end)
186 double dr2, dr2_rel, inp;
189 pbc_dx_aiuc(pbc, x[ii], g_x, dx);
190 inp = iprod(dir, dx);
192 for (m = 0; m < DIM; m++)
194 /* Determine the radial components */
195 dr[m] = dx[m] - inp*dir[m];
198 dr2_rel = dr2*inv_cyl_r2;
202 double mass, weight, dweight_r;
205 /* add to index, to sum of COM, to weight array */
206 if (pdyna->nat_loc >= pdyna->nalloc_loc)
208 pdyna->nalloc_loc = over_alloc_large(pdyna->nat_loc+1);
209 srenew(pdyna->ind_loc, pdyna->nalloc_loc);
210 srenew(pdyna->weight_loc, pdyna->nalloc_loc);
211 srenew(pdyna->mdw, pdyna->nalloc_loc);
212 srenew(pdyna->dv, pdyna->nalloc_loc);
214 pdyna->ind_loc[pdyna->nat_loc] = ii;
216 mass = md->massT[ii];
217 /* The radial weight function is 1-2x^2+x^4,
218 * where x=r/cylinder_r. Since this function depends
219 * on the radial component, we also get radial forces
222 weight = 1 + (-2 + dr2_rel)*dr2_rel;
223 dweight_r = (-4 + 4*dr2_rel)*inv_cyl_r2;
224 pdyna->weight_loc[pdyna->nat_loc] = weight;
225 sum_a += mass*weight*inp;
226 wmass += mass*weight;
227 wwmass += mass*weight*weight;
228 dsvmul(mass*dweight_r, dr, mdw);
229 copy_dvec(mdw, pdyna->mdw[pdyna->nat_loc]);
230 /* Currently we only have the axial component of the
231 * distance (inp) up to an unkown offset. We add this
232 * offset after the reduction needs to determine the
233 * COM of the cylinder group.
235 pdyna->dv[pdyna->nat_loc] = inp;
236 for (m = 0; m < DIM; m++)
238 radf_fac0[m] += mdw[m];
239 radf_fac1[m] += mdw[m]*inp;
246 pull->dbuf_cyl[c*stride+0] = wmass;
247 pull->dbuf_cyl[c*stride+1] = wwmass;
248 pull->dbuf_cyl[c*stride+2] = sum_a;
249 pull->dbuf_cyl[c*stride+3] = radf_fac0[XX];
250 pull->dbuf_cyl[c*stride+4] = radf_fac0[YY];
251 pull->dbuf_cyl[c*stride+5] = radf_fac0[ZZ];
252 pull->dbuf_cyl[c*stride+6] = radf_fac1[XX];
253 pull->dbuf_cyl[c*stride+7] = radf_fac1[YY];
254 pull->dbuf_cyl[c*stride+8] = radf_fac1[ZZ];
257 if (cr != NULL && PAR(cr))
259 /* Sum the contributions over the ranks */
260 gmx_sumd(pull->ncoord*stride, pull->dbuf_cyl, cr);
263 for (c = 0; c < pull->ncoord; c++)
265 pull_coord_work_t *pcrd;
267 pcrd = &pull->coord[c];
269 if (pcrd->params.eGeom == epullgCYL)
271 pull_group_work_t *pdyna, *pgrp;
272 double wmass, wwmass, dist;
274 pdyna = &pull->dyna[c];
275 pgrp = &pull->group[pcrd->params.group[1]];
277 wmass = pull->dbuf_cyl[c*stride+0];
278 wwmass = pull->dbuf_cyl[c*stride+1];
279 pdyna->mwscale = 1.0/wmass;
280 /* Cylinder pulling can't be used with constraints, but we set
281 * wscale and invtm anyhow, in case someone would like to use them.
283 pdyna->wscale = wmass/wwmass;
284 pdyna->invtm = wwmass/(wmass*wmass);
286 /* We store the deviation of the COM from the reference location
287 * used above, since we need it when we apply the radial forces
288 * to the atoms in the cylinder group.
291 for (m = 0; m < DIM; m++)
293 g_x[m] = pgrp->x[m] - pcrd->vec[m]*pcrd->value_ref;
294 dist = -pcrd->vec[m]*pull->dbuf_cyl[c*stride+2]*pdyna->mwscale;
295 pdyna->x[m] = g_x[m] - dist;
296 pcrd->cyl_dev += dist;
298 /* Now we know the exact COM of the cylinder reference group,
299 * we can determine the radial force factor (ffrad) that when
300 * multiplied with the axial pull force will give the radial
301 * force on the pulled (non-cylinder) group.
303 for (m = 0; m < DIM; m++)
305 pcrd->ffrad[m] = (pull->dbuf_cyl[c*stride+6+m] +
306 pull->dbuf_cyl[c*stride+3+m]*pcrd->cyl_dev)/wmass;
311 fprintf(debug, "Pull cylinder group %d:%8.3f%8.3f%8.3f m:%8.3f\n",
312 c, pdyna->x[0], pdyna->x[1],
313 pdyna->x[2], 1.0/pdyna->invtm);
314 fprintf(debug, "ffrad %8.3f %8.3f %8.3f\n",
315 pcrd->ffrad[XX], pcrd->ffrad[YY], pcrd->ffrad[ZZ]);
321 static double atan2_0_2pi(double y, double x)
333 /* calculates center of mass of selection index from all coordinates x */
334 void pull_calc_coms(t_commrec *cr,
335 struct pull_t *pull, t_mdatoms *md, t_pbc *pbc, double t,
341 if (pull->rbuf == NULL)
343 snew(pull->rbuf, pull->ngroup);
345 if (pull->dbuf == NULL)
347 snew(pull->dbuf, 3*pull->ngroup);
350 if (pull->bRefAt && pull->bSetPBCatoms)
352 pull_set_pbcatoms(cr, pull, x, pull->rbuf);
354 if (cr != NULL && DOMAINDECOMP(cr))
356 /* We can keep these PBC reference coordinates fixed for nstlist
357 * steps, since atoms won't jump over PBC.
358 * This avoids a global reduction at the next nstlist-1 steps.
359 * Note that the exact values of the pbc reference coordinates
360 * are irrelevant, as long all atoms in the group are within
361 * half a box distance of the reference coordinate.
363 pull->bSetPBCatoms = FALSE;
367 if (pull->cosdim >= 0)
371 assert(pull->npbcdim <= DIM);
373 for (m = pull->cosdim+1; m < pull->npbcdim; m++)
375 if (pbc->box[m][pull->cosdim] != 0)
377 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
380 twopi_box = 2.0*M_PI/pbc->box[pull->cosdim][pull->cosdim];
383 for (g = 0; g < pull->ngroup; g++)
385 pull_group_work_t *pgrp;
387 pgrp = &pull->group[g];
391 if (pgrp->epgrppbc != epgrppbcCOS)
394 double wmass, wwmass;
395 rvec x_pbc = { 0, 0, 0 };
403 if (pgrp->epgrppbc == epgrppbcREFAT)
405 /* Set the pbc atom */
406 copy_rvec(pull->rbuf[g], x_pbc);
409 for (i = 0; i < pgrp->nat_loc; i++)
414 ii = pgrp->ind_loc[i];
415 mass = md->massT[ii];
416 if (pgrp->weight_loc == NULL)
425 w = pgrp->weight_loc[i];
430 if (pgrp->epgrppbc == epgrppbcNONE)
432 /* Plain COM: sum the coordinates */
433 for (m = 0; m < DIM; m++)
435 com[m] += wm*x[ii][m];
439 for (m = 0; m < DIM; m++)
441 comp[m] += wm*xp[ii][m];
449 /* Sum the difference with the reference atom */
450 pbc_dx(pbc, x[ii], x_pbc, dx);
451 for (m = 0; m < DIM; m++)
457 /* For xp add the difference between xp and x to dx,
458 * such that we use the same periodic image,
459 * also when xp has a large displacement.
461 for (m = 0; m < DIM; m++)
463 comp[m] += wm*(dx[m] + xp[ii][m] - x[ii][m]);
469 /* We do this check after the loop above to avoid more nesting.
470 * If we have a single-atom group the mass is irrelevant, so
471 * we can remove the mass factor to avoid division by zero.
472 * Note that with constraint pulling the mass does matter, but
473 * in that case a check group mass != 0 has been done before.
475 if (pgrp->params.nat == 1 && pgrp->nat_loc == 1 && wmass == 0)
479 /* Copy the single atom coordinate */
480 for (m = 0; m < DIM; m++)
482 com[m] = x[pgrp->ind_loc[0]][m];
484 /* Set all mass factors to 1 to get the correct COM */
489 if (pgrp->weight_loc == NULL)
494 /* Copy local sums to a buffer for global summing */
495 copy_dvec(com, pull->dbuf[g*3]);
496 copy_dvec(comp, pull->dbuf[g*3+1]);
497 pull->dbuf[g*3+2][0] = wmass;
498 pull->dbuf[g*3+2][1] = wwmass;
499 pull->dbuf[g*3+2][2] = 0;
503 /* Cosine weighting geometry */
504 double cm, sm, cmp, smp, ccm, csm, ssm, csw, snw;
515 for (i = 0; i < pgrp->nat_loc; i++)
520 ii = pgrp->ind_loc[i];
521 mass = md->massT[ii];
522 /* Determine cos and sin sums */
523 csw = cos(x[ii][pull->cosdim]*twopi_box);
524 snw = sin(x[ii][pull->cosdim]*twopi_box);
533 csw = cos(xp[ii][pull->cosdim]*twopi_box);
534 snw = sin(xp[ii][pull->cosdim]*twopi_box);
540 /* Copy local sums to a buffer for global summing */
541 pull->dbuf[g*3 ][0] = cm;
542 pull->dbuf[g*3 ][1] = sm;
543 pull->dbuf[g*3 ][2] = 0;
544 pull->dbuf[g*3+1][0] = ccm;
545 pull->dbuf[g*3+1][1] = csm;
546 pull->dbuf[g*3+1][2] = ssm;
547 pull->dbuf[g*3+2][0] = cmp;
548 pull->dbuf[g*3+2][1] = smp;
549 pull->dbuf[g*3+2][2] = 0;
556 /* Sum the contributions over the nodes */
557 gmx_sumd(pull->ngroup*3*DIM, pull->dbuf[0], cr);
560 for (g = 0; g < pull->ngroup; g++)
562 pull_group_work_t *pgrp;
564 pgrp = &pull->group[g];
565 if (pgrp->params.nat > 0 && pgrp->bCalcCOM)
567 if (pgrp->epgrppbc != epgrppbcCOS)
569 double wmass, wwmass;
572 /* Determine the inverse mass */
573 wmass = pull->dbuf[g*3+2][0];
574 wwmass = pull->dbuf[g*3+2][1];
575 pgrp->mwscale = 1.0/wmass;
576 /* invtm==0 signals a frozen group, so then we should keep it zero */
577 if (pgrp->invtm != 0)
579 pgrp->wscale = wmass/wwmass;
580 pgrp->invtm = wwmass/(wmass*wmass);
582 /* Divide by the total mass */
583 for (m = 0; m < DIM; m++)
585 pgrp->x[m] = pull->dbuf[g*3 ][m]*pgrp->mwscale;
588 pgrp->xp[m] = pull->dbuf[g*3+1][m]*pgrp->mwscale;
590 if (pgrp->epgrppbc == epgrppbcREFAT)
592 pgrp->x[m] += pull->rbuf[g][m];
595 pgrp->xp[m] += pull->rbuf[g][m];
602 /* Cosine weighting geometry */
603 double csw, snw, wmass, wwmass;
606 /* Determine the optimal location of the cosine weight */
607 csw = pull->dbuf[g*3][0];
608 snw = pull->dbuf[g*3][1];
609 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
610 /* Set the weights for the local atoms */
611 wmass = sqrt(csw*csw + snw*snw);
612 wwmass = (pull->dbuf[g*3+1][0]*csw*csw +
613 pull->dbuf[g*3+1][1]*csw*snw +
614 pull->dbuf[g*3+1][2]*snw*snw)/(wmass*wmass);
616 pgrp->mwscale = 1.0/wmass;
617 pgrp->wscale = wmass/wwmass;
618 pgrp->invtm = wwmass/(wmass*wmass);
619 /* Set the weights for the local atoms */
622 for (i = 0; i < pgrp->nat_loc; i++)
624 ii = pgrp->ind_loc[i];
625 pgrp->weight_loc[i] = csw*cos(twopi_box*x[ii][pull->cosdim]) +
626 snw*sin(twopi_box*x[ii][pull->cosdim]);
630 csw = pull->dbuf[g*3+2][0];
631 snw = pull->dbuf[g*3+2][1];
632 pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
637 fprintf(debug, "Pull group %d wmass %f invtm %f\n",
638 g, 1.0/pgrp->mwscale, pgrp->invtm);
645 /* Calculate the COMs for the cyclinder reference groups */
646 make_cyl_refgrps(cr, pull, md, pbc, t, x);