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.
163 /* With rate=0, value_ref is set initially */
164 pcrd->value_ref = pcrd->init + pcrd->rate*t;
166 for (m = 0; m < DIM; m++)
168 g_x[m] = pgrp->x[m] - pcrd->vec[m]*pcrd->value_ref;
171 /* loop over all atoms in the main ref group */
172 for (i = 0; i < pref->nat; i++)
177 if (!ga2la_get_home(ga2la, pref->ind[i], &ii))
182 if (ii >= start && ii < end)
184 double dr2, dr2_rel, inp;
187 pbc_dx_aiuc(pbc, x[ii], g_x, dx);
188 inp = iprod(dir, dx);
190 for (m = 0; m < DIM; m++)
192 /* Determine the radial components */
193 dr[m] = dx[m] - inp*dir[m];
196 dr2_rel = dr2*inv_cyl_r2;
200 double mass, weight, dweight_r;
203 /* add to index, to sum of COM, to weight array */
204 if (pdyna->nat_loc >= pdyna->nalloc_loc)
206 pdyna->nalloc_loc = over_alloc_large(pdyna->nat_loc+1);
207 srenew(pdyna->ind_loc, pdyna->nalloc_loc);
208 srenew(pdyna->weight_loc, pdyna->nalloc_loc);
209 srenew(pdyna->mdw, pdyna->nalloc_loc);
210 srenew(pdyna->dv, pdyna->nalloc_loc);
212 pdyna->ind_loc[pdyna->nat_loc] = ii;
214 mass = md->massT[ii];
215 /* The radial weight function is 1-2x^2+x^4,
216 * where x=r/cylinder_r. Since this function depends
217 * on the radial component, we also get radial forces
220 weight = 1 + (-2 + dr2_rel)*dr2_rel;
221 dweight_r = (-4 + 4*dr2_rel)*inv_cyl_r2;
222 pdyna->weight_loc[pdyna->nat_loc] = weight;
223 sum_a += mass*weight*inp;
224 wmass += mass*weight;
225 wwmass += mass*weight*weight;
226 dsvmul(mass*dweight_r, dr, mdw);
227 copy_dvec(mdw, pdyna->mdw[pdyna->nat_loc]);
228 /* Currently we only have the axial component of the
229 * distance (inp) up to an unkown offset. We add this
230 * offset after the reduction needs to determine the
231 * COM of the cylinder group.
233 pdyna->dv[pdyna->nat_loc] = inp;
234 for (m = 0; m < DIM; m++)
236 radf_fac0[m] += mdw[m];
237 radf_fac1[m] += mdw[m]*inp;
244 pull->dbuf_cyl[c*stride+0] = wmass;
245 pull->dbuf_cyl[c*stride+1] = wwmass;
246 pull->dbuf_cyl[c*stride+2] = sum_a;
247 pull->dbuf_cyl[c*stride+3] = radf_fac0[XX];
248 pull->dbuf_cyl[c*stride+4] = radf_fac0[YY];
249 pull->dbuf_cyl[c*stride+5] = radf_fac0[ZZ];
250 pull->dbuf_cyl[c*stride+6] = radf_fac1[XX];
251 pull->dbuf_cyl[c*stride+7] = radf_fac1[YY];
252 pull->dbuf_cyl[c*stride+8] = radf_fac1[ZZ];
255 if (cr != NULL && PAR(cr))
257 /* Sum the contributions over the ranks */
258 gmx_sumd(pull->ncoord*stride, pull->dbuf_cyl, cr);
261 for (c = 0; c < pull->ncoord; c++)
263 pcrd = &pull->coord[c];
265 if (pcrd->eGeom == epullgCYL)
267 double wmass, wwmass, dist;
269 pdyna = &pull->dyna[c];
270 pgrp = &pull->group[pcrd->group[1]];
272 wmass = pull->dbuf_cyl[c*stride+0];
273 wwmass = pull->dbuf_cyl[c*stride+1];
274 pdyna->mwscale = 1.0/wmass;
275 /* Cylinder pulling can't be used with constraints, but we set
276 * wscale and invtm anyhow, in case someone would like to use them.
278 pdyna->wscale = wmass/wwmass;
279 pdyna->invtm = wwmass/(wmass*wmass);
281 /* We store the deviation of the COM from the reference location
282 * used above, since we need it when we apply the radial forces
283 * to the atoms in the cylinder group.
286 for (m = 0; m < DIM; m++)
288 g_x[m] = pgrp->x[m] - pcrd->vec[m]*pcrd->value_ref;
289 dist = -pcrd->vec[m]*pull->dbuf_cyl[c*stride+2]*pdyna->mwscale;
290 pdyna->x[m] = g_x[m] - dist;
291 pcrd->cyl_dev += dist;
293 /* Now we know the exact COM of the cylinder reference group,
294 * we can determine the radial force factor (ffrad) that when
295 * multiplied with the axial pull force will give the radial
296 * force on the pulled (non-cylinder) group.
298 for (m = 0; m < DIM; m++)
300 pcrd->ffrad[m] = (pull->dbuf_cyl[c*stride+6+m] +
301 pull->dbuf_cyl[c*stride+3+m]*pcrd->cyl_dev)/wmass;
306 fprintf(debug, "Pull cylinder group %d:%8.3f%8.3f%8.3f m:%8.3f\n",
307 c, pdyna->x[0], pdyna->x[1],
308 pdyna->x[2], 1.0/pdyna->invtm);
309 fprintf(debug, "ffrad %8.3f %8.3f %8.3f\n",
310 pcrd->ffrad[XX], pcrd->ffrad[YY], pcrd->ffrad[ZZ]);
316 static double atan2_0_2pi(double y, double x)
328 /* calculates center of mass of selection index from all coordinates x */
329 void pull_calc_coms(t_commrec *cr,
330 t_pull *pull, t_mdatoms *md, t_pbc *pbc, double t,
336 if (pull->rbuf == NULL)
338 snew(pull->rbuf, pull->ngroup);
340 if (pull->dbuf == NULL)
342 snew(pull->dbuf, 3*pull->ngroup);
345 if (pull->bRefAt && pull->bSetPBCatoms)
347 pull_set_pbcatoms(cr, pull, x, pull->rbuf);
349 if (cr != NULL && DOMAINDECOMP(cr))
351 /* We can keep these PBC reference coordinates fixed for nstlist
352 * steps, since atoms won't jump over PBC.
353 * This avoids a global reduction at the next nstlist-1 steps.
354 * Note that the exact values of the pbc reference coordinates
355 * are irrelevant, as long all atoms in the group are within
356 * half a box distance of the reference coordinate.
358 pull->bSetPBCatoms = FALSE;
362 if (pull->cosdim >= 0)
366 assert(pull->npbcdim <= DIM);
368 for (m = pull->cosdim+1; m < pull->npbcdim; m++)
370 if (pbc->box[m][pull->cosdim] != 0)
372 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
375 twopi_box = 2.0*M_PI/pbc->box[pull->cosdim][pull->cosdim];
378 for (g = 0; g < pull->ngroup; g++)
382 pgrp = &pull->group[g];
386 if (pgrp->epgrppbc != epgrppbcCOS)
389 double wmass, wwmass;
390 rvec x_pbc = { 0, 0, 0 };
398 if (pgrp->epgrppbc == epgrppbcREFAT)
400 /* Set the pbc atom */
401 copy_rvec(pull->rbuf[g], x_pbc);
404 for (i = 0; i < pgrp->nat_loc; i++)
409 ii = pgrp->ind_loc[i];
410 mass = md->massT[ii];
411 if (pgrp->weight_loc == NULL)
420 w = pgrp->weight_loc[i];
425 if (pgrp->epgrppbc == epgrppbcNONE)
427 /* Plain COM: sum the coordinates */
428 for (m = 0; m < DIM; m++)
430 com[m] += wm*x[ii][m];
434 for (m = 0; m < DIM; m++)
436 comp[m] += wm*xp[ii][m];
444 /* Sum the difference with the reference atom */
445 pbc_dx(pbc, x[ii], x_pbc, dx);
446 for (m = 0; m < DIM; m++)
452 /* For xp add the difference between xp and x to dx,
453 * such that we use the same periodic image,
454 * also when xp has a large displacement.
456 for (m = 0; m < DIM; m++)
458 comp[m] += wm*(dx[m] + xp[ii][m] - x[ii][m]);
464 /* We do this check after the loop above to avoid more nesting.
465 * If we have a single-atom group the mass is irrelevant, so
466 * we can remove the mass factor to avoid division by zero.
467 * Note that with constraint pulling the mass does matter, but
468 * in that case a check group mass != 0 has been done before.
470 if (pgrp->nat == 1 && pgrp->nat_loc == 1 && wmass == 0)
474 /* Copy the single atom coordinate */
475 for (m = 0; m < DIM; m++)
477 com[m] = x[pgrp->ind_loc[0]][m];
479 /* Set all mass factors to 1 to get the correct COM */
484 if (pgrp->weight_loc == NULL)
489 /* Copy local sums to a buffer for global summing */
490 copy_dvec(com, pull->dbuf[g*3]);
491 copy_dvec(comp, pull->dbuf[g*3+1]);
492 pull->dbuf[g*3+2][0] = wmass;
493 pull->dbuf[g*3+2][1] = wwmass;
494 pull->dbuf[g*3+2][2] = 0;
498 /* Cosine weighting geometry */
499 double cm, sm, cmp, smp, ccm, csm, ssm, csw, snw;
510 for (i = 0; i < pgrp->nat_loc; i++)
515 ii = pgrp->ind_loc[i];
516 mass = md->massT[ii];
517 /* Determine cos and sin sums */
518 csw = cos(x[ii][pull->cosdim]*twopi_box);
519 snw = sin(x[ii][pull->cosdim]*twopi_box);
528 csw = cos(xp[ii][pull->cosdim]*twopi_box);
529 snw = sin(xp[ii][pull->cosdim]*twopi_box);
535 /* Copy local sums to a buffer for global summing */
536 pull->dbuf[g*3 ][0] = cm;
537 pull->dbuf[g*3 ][1] = sm;
538 pull->dbuf[g*3 ][2] = 0;
539 pull->dbuf[g*3+1][0] = ccm;
540 pull->dbuf[g*3+1][1] = csm;
541 pull->dbuf[g*3+1][2] = ssm;
542 pull->dbuf[g*3+2][0] = cmp;
543 pull->dbuf[g*3+2][1] = smp;
544 pull->dbuf[g*3+2][2] = 0;
551 /* Sum the contributions over the nodes */
552 gmx_sumd(pull->ngroup*3*DIM, pull->dbuf[0], cr);
555 for (g = 0; g < pull->ngroup; g++)
559 pgrp = &pull->group[g];
560 if (pgrp->nat > 0 && pgrp->bCalcCOM)
562 if (pgrp->epgrppbc != epgrppbcCOS)
564 double wmass, wwmass;
567 /* Determine the inverse mass */
568 wmass = pull->dbuf[g*3+2][0];
569 wwmass = pull->dbuf[g*3+2][1];
570 pgrp->mwscale = 1.0/wmass;
571 /* invtm==0 signals a frozen group, so then we should keep it zero */
572 if (pgrp->invtm != 0)
574 pgrp->wscale = wmass/wwmass;
575 pgrp->invtm = wwmass/(wmass*wmass);
577 /* Divide by the total mass */
578 for (m = 0; m < DIM; m++)
580 pgrp->x[m] = pull->dbuf[g*3 ][m]*pgrp->mwscale;
583 pgrp->xp[m] = pull->dbuf[g*3+1][m]*pgrp->mwscale;
585 if (pgrp->epgrppbc == epgrppbcREFAT)
587 pgrp->x[m] += pull->rbuf[g][m];
590 pgrp->xp[m] += pull->rbuf[g][m];
597 /* Cosine weighting geometry */
598 double csw, snw, wmass, wwmass;
601 /* Determine the optimal location of the cosine weight */
602 csw = pull->dbuf[g*3][0];
603 snw = pull->dbuf[g*3][1];
604 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
605 /* Set the weights for the local atoms */
606 wmass = sqrt(csw*csw + snw*snw);
607 wwmass = (pull->dbuf[g*3+1][0]*csw*csw +
608 pull->dbuf[g*3+1][1]*csw*snw +
609 pull->dbuf[g*3+1][2]*snw*snw)/(wmass*wmass);
611 pgrp->mwscale = 1.0/wmass;
612 pgrp->wscale = wmass/wwmass;
613 pgrp->invtm = wwmass/(wmass*wmass);
614 /* Set the weights for the local atoms */
617 for (i = 0; i < pgrp->nat_loc; i++)
619 ii = pgrp->ind_loc[i];
620 pgrp->weight_loc[i] = csw*cos(twopi_box*x[ii][pull->cosdim]) +
621 snw*sin(twopi_box*x[ii][pull->cosdim]);
625 csw = pull->dbuf[g*3+2][0];
626 snw = pull->dbuf[g*3+2][1];
627 pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
632 fprintf(debug, "Pull group %d wmass %f invtm %f\n",
633 g, 1.0/pgrp->mwscale, pgrp->invtm);
640 /* Calculate the COMs for the cyclinder reference groups */
641 make_cyl_refgrps(cr, pull, md, pbc, t, x);