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, 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.
45 #include "gromacs/fileio/futil.h"
47 #include "gromacs/utility/smalloc.h"
49 #include "types/commrec.h"
51 #include "gmx_fatal.h"
55 #include "gromacs/fileio/confio.h"
59 #include "gmx_ga2la.h"
61 static void pull_set_pbcatom(t_commrec *cr, t_pull_group *pgrp,
67 if (cr != NULL && DOMAINDECOMP(cr))
69 if (ga2la_get_home(cr->dd->ga2la, pgrp->pbcatom, &a))
71 copy_rvec(x[a], x_pbc);
80 copy_rvec(x[pgrp->pbcatom], x_pbc);
84 static void pull_set_pbcatoms(t_commrec *cr, t_pull *pull,
91 for (g = 0; g < pull->ngroup; g++)
93 if ((g == 0 && PULL_CYL(pull)) || pull->group[g].pbcatom == -1)
99 pull_set_pbcatom(cr, &pull->group[g], x, x_pbc[g]);
100 for (m = 0; m < DIM; m++)
102 if (pull->dim[m] == 0)
111 if (cr && PAR(cr) && n > 0)
113 /* Sum over the nodes to get x_pbc from the home node of pbcatom */
114 gmx_sum(pull->ngroup*DIM, x_pbc[0], cr);
118 /* switch function, x between r and w */
119 static real get_weight(real x, real r1, real r0)
133 weight = (r0 - x)/(r0 - r1);
139 static void make_cyl_refgrps(t_commrec *cr, t_pull *pull, t_mdatoms *md,
140 t_pbc *pbc, double t, rvec *x, rvec *xp)
142 int c, i, ii, m, start, end;
144 double r0_2, sum_a, sum_ap, dr2, mass, weight, wmass, wwmass, inp;
146 t_pull_group *pref, *pgrp, *pdyna;
147 gmx_ga2la_t ga2la = NULL;
149 if (pull->dbuf_cyl == NULL)
151 snew(pull->dbuf_cyl, pull->ncoord*4);
154 if (cr && DOMAINDECOMP(cr))
156 ga2la = cr->dd->ga2la;
162 r0_2 = dsqr(pull->cyl_r0);
164 /* loop over all groups to make a reference group for each*/
165 for (c = 0; c < pull->ncoord; c++)
167 pcrd = &pull->coord[c];
169 /* pref will be the same group for all pull coordinates */
170 pref = &pull->group[pcrd->group[0]];
171 pgrp = &pull->group[pcrd->group[1]];
172 pdyna = &pull->dyna[c];
173 copy_rvec(pcrd->vec, dir);
180 for (m = 0; m < DIM; m++)
182 g_x[m] = pgrp->x[m] - pcrd->vec[m]*(pcrd->init + pcrd->rate*t);
185 /* loop over all atoms in the main ref group */
186 for (i = 0; i < pref->nat; i++)
191 if (!ga2la_get_home(ga2la, pref->ind[i], &ii))
196 if (ii >= start && ii < end)
198 pbc_dx_aiuc(pbc, x[ii], g_x, dx);
199 inp = iprod(dir, dx);
201 for (m = 0; m < DIM; m++)
203 dr2 += dsqr(dx[m] - inp*dir[m]);
208 /* add to index, to sum of COM, to weight array */
209 if (pdyna->nat_loc >= pdyna->nalloc_loc)
211 pdyna->nalloc_loc = over_alloc_large(pdyna->nat_loc+1);
212 srenew(pdyna->ind_loc, pdyna->nalloc_loc);
213 srenew(pdyna->weight_loc, pdyna->nalloc_loc);
215 pdyna->ind_loc[pdyna->nat_loc] = ii;
216 mass = md->massT[ii];
217 weight = get_weight(sqrt(dr2), pull->cyl_r1, pull->cyl_r0);
218 pdyna->weight_loc[pdyna->nat_loc] = weight;
219 sum_a += mass*weight*inp;
222 pbc_dx_aiuc(pbc, xp[ii], g_x, dx);
223 inp = iprod(dir, dx);
224 sum_ap += mass*weight*inp;
226 wmass += mass*weight;
227 wwmass += mass*sqr(weight);
232 pull->dbuf_cyl[c*4+0] = wmass;
233 pull->dbuf_cyl[c*4+1] = wwmass;
234 pull->dbuf_cyl[c*4+2] = sum_a;
235 pull->dbuf_cyl[c*4+3] = sum_ap;
240 /* Sum the contributions over the nodes */
241 gmx_sumd(pull->ncoord*4, pull->dbuf_cyl, cr);
244 for (c = 0; c < pull->ncoord; c++)
246 pcrd = &pull->coord[c];
248 pdyna = &pull->dyna[c];
249 pgrp = &pull->group[pcrd->group[1]];
251 wmass = pull->dbuf_cyl[c*4+0];
252 wwmass = pull->dbuf_cyl[c*4+1];
253 pdyna->wscale = wmass/wwmass;
254 pdyna->invtm = 1.0/(pdyna->wscale*wmass);
256 for (m = 0; m < DIM; m++)
258 g_x[m] = pgrp->x[m] - pcrd->vec[m]*(pcrd->init + pcrd->rate*t);
259 pdyna->x[m] = g_x[m] + pcrd->vec[m]*pull->dbuf_cyl[c*4+2]/wmass;
262 pdyna->xp[m] = g_x[m] + pcrd->vec[m]*pull->dbuf_cyl[c*4+3]/wmass;
268 fprintf(debug, "Pull cylinder group %d:%8.3f%8.3f%8.3f m:%8.3f\n",
269 c, pdyna->x[0], pdyna->x[1],
270 pdyna->x[2], 1.0/pdyna->invtm);
275 static double atan2_0_2pi(double y, double x)
287 /* calculates center of mass of selection index from all coordinates x */
288 void pull_calc_coms(t_commrec *cr,
289 t_pull *pull, t_mdatoms *md, t_pbc *pbc, double t,
293 real mass, w, wm, twopi_box = 0;
294 double wmass, wwmass, invwmass;
296 double cm, sm, cmp, smp, ccm, csm, ssm, csw, snw;
297 rvec *xx[2], x_pbc = {0, 0, 0}, dx;
300 if (pull->rbuf == NULL)
302 snew(pull->rbuf, pull->ngroup);
304 if (pull->dbuf == NULL)
306 snew(pull->dbuf, 3*pull->ngroup);
309 if (pull->bRefAt && pull->bSetPBCatoms)
311 pull_set_pbcatoms(cr, pull, x, pull->rbuf);
313 if (cr != NULL && DOMAINDECOMP(cr))
315 /* We can keep these PBC reference coordinates fixed for nstlist
316 * steps, since atoms won't jump over PBC.
317 * This avoids a global reduction at the next nstlist-1 steps.
318 * Note that the exact values of the pbc reference coordinates
319 * are irrelevant, as long all atoms in the group are within
320 * half a box distance of the reference coordinate.
322 pull->bSetPBCatoms = FALSE;
326 if (pull->cosdim >= 0)
328 for (m = pull->cosdim+1; m < pull->npbcdim; m++)
330 if (pbc->box[m][pull->cosdim] != 0)
332 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
335 twopi_box = 2.0*M_PI/pbc->box[pull->cosdim][pull->cosdim];
338 for (g = 0; g < pull->ngroup; g++)
340 pgrp = &pull->group[g];
352 if (!(g == 0 && PULL_CYL(pull)))
354 if (pgrp->epgrppbc == epgrppbcREFAT)
356 /* Set the pbc atom */
357 copy_rvec(pull->rbuf[g], x_pbc);
360 for (i = 0; i < pgrp->nat_loc; i++)
362 ii = pgrp->ind_loc[i];
363 mass = md->massT[ii];
364 if (pgrp->epgrppbc != epgrppbcCOS)
366 if (pgrp->weight_loc)
368 w = pgrp->weight_loc[i];
373 if (pgrp->epgrppbc == epgrppbcNONE)
375 /* Plain COM: sum the coordinates */
376 for (m = 0; m < DIM; m++)
378 com[m] += wm*x[ii][m];
382 for (m = 0; m < DIM; m++)
384 comp[m] += wm*xp[ii][m];
390 /* Sum the difference with the reference atom */
391 pbc_dx(pbc, x[ii], x_pbc, dx);
392 for (m = 0; m < DIM; m++)
398 /* For xp add the difference between xp and x to dx,
399 * such that we use the same periodic image,
400 * also when xp has a large displacement.
402 for (m = 0; m < DIM; m++)
404 comp[m] += wm*(dx[m] + xp[ii][m] - x[ii][m]);
411 /* Determine cos and sin sums */
412 csw = cos(x[ii][pull->cosdim]*twopi_box);
413 snw = sin(x[ii][pull->cosdim]*twopi_box);
422 csw = cos(xp[ii][pull->cosdim]*twopi_box);
423 snw = sin(xp[ii][pull->cosdim]*twopi_box);
431 /* Copy local sums to a buffer for global summing */
432 switch (pgrp->epgrppbc)
436 copy_dvec(com, pull->dbuf[g*3]);
437 copy_dvec(comp, pull->dbuf[g*3+1]);
438 pull->dbuf[g*3+2][0] = wmass;
439 pull->dbuf[g*3+2][1] = wwmass;
440 pull->dbuf[g*3+2][2] = 0;
443 pull->dbuf[g*3 ][0] = cm;
444 pull->dbuf[g*3 ][1] = sm;
445 pull->dbuf[g*3 ][2] = 0;
446 pull->dbuf[g*3+1][0] = ccm;
447 pull->dbuf[g*3+1][1] = csm;
448 pull->dbuf[g*3+1][2] = ssm;
449 pull->dbuf[g*3+2][0] = cmp;
450 pull->dbuf[g*3+2][1] = smp;
451 pull->dbuf[g*3+2][2] = 0;
458 /* Sum the contributions over the nodes */
459 gmx_sumd(pull->ngroup*3*DIM, pull->dbuf[0], cr);
462 for (g = 0; g < pull->ngroup; g++)
464 pgrp = &pull->group[g];
465 if (pgrp->nat > 0 && !(g == 0 && PULL_CYL(pull)))
467 if (pgrp->epgrppbc != epgrppbcCOS)
469 /* Determine the inverse mass */
470 wmass = pull->dbuf[g*3+2][0];
471 wwmass = pull->dbuf[g*3+2][1];
473 /* invtm==0 signals a frozen group, so then we should keep it zero */
476 pgrp->wscale = wmass/wwmass;
477 pgrp->invtm = 1.0/(pgrp->wscale*wmass);
479 /* Divide by the total mass */
480 for (m = 0; m < DIM; m++)
482 pgrp->x[m] = pull->dbuf[g*3 ][m]*invwmass;
485 pgrp->xp[m] = pull->dbuf[g*3+1][m]*invwmass;
487 if (pgrp->epgrppbc == epgrppbcREFAT)
489 pgrp->x[m] += pull->rbuf[g][m];
492 pgrp->xp[m] += pull->rbuf[g][m];
499 /* Determine the optimal location of the cosine weight */
500 csw = pull->dbuf[g*3][0];
501 snw = pull->dbuf[g*3][1];
502 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
503 /* Set the weights for the local atoms */
504 wmass = sqrt(csw*csw + snw*snw);
505 wwmass = (pull->dbuf[g*3+1][0]*csw*csw +
506 pull->dbuf[g*3+1][1]*csw*snw +
507 pull->dbuf[g*3+1][2]*snw*snw)/(wmass*wmass);
508 pgrp->wscale = wmass/wwmass;
509 pgrp->invtm = 1.0/(pgrp->wscale*wmass);
510 /* Set the weights for the local atoms */
513 for (i = 0; i < pgrp->nat_loc; i++)
515 ii = pgrp->ind_loc[i];
516 pgrp->weight_loc[i] = csw*cos(twopi_box*x[ii][pull->cosdim]) +
517 snw*sin(twopi_box*x[ii][pull->cosdim]);
521 csw = pull->dbuf[g*3+2][0];
522 snw = pull->dbuf[g*3+2][1];
523 pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
528 fprintf(debug, "Pull group %d wmass %f wwmass %f invtm %f\n",
529 g, wmass, wwmass, pgrp->invtm);
536 /* Calculate the COMs for the cyclinder reference groups */
537 make_cyl_refgrps(cr, pull, md, pbc, t, x, xp);