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"
50 #include "gmx_fatal.h"
54 #include "gromacs/fileio/confio.h"
58 #include "gmx_ga2la.h"
60 static void pull_set_pbcatom(t_commrec *cr, t_pull_group *pgrp,
61 t_mdatoms *md, rvec *x,
70 if (!ga2la_get_home(cr->dd->ga2la, pgrp->pbcatom, &a))
80 if (a >= 0 && a < md->homenr)
82 copy_rvec(x[a], x_pbc);
91 copy_rvec(x[pgrp->pbcatom], x_pbc);
95 static void pull_set_pbcatoms(t_commrec *cr, t_pull *pull,
96 t_mdatoms *md, rvec *x,
102 for (g = 0; g < pull->ngroup; g++)
104 if ((g == 0 && PULL_CYL(pull)) || pull->group[g].pbcatom == -1)
106 clear_rvec(x_pbc[g]);
110 pull_set_pbcatom(cr, &pull->group[g], md, x, x_pbc[g]);
111 for (m = 0; m < DIM; m++)
113 if (pull->dim[m] == 0)
122 if (cr && PAR(cr) && n > 0)
124 /* Sum over the nodes to get x_pbc from the home node of pbcatom */
125 gmx_sum(pull->ngroup*DIM, x_pbc[0], cr);
129 /* switch function, x between r and w */
130 static real get_weight(real x, real r1, real r0)
144 weight = (r0 - x)/(r0 - r1);
150 static void make_cyl_refgrps(t_commrec *cr, t_pull *pull, t_mdatoms *md,
151 t_pbc *pbc, double t, rvec *x, rvec *xp)
153 int c, i, ii, m, start, end;
155 double r0_2, sum_a, sum_ap, dr2, mass, weight, wmass, wwmass, inp;
157 t_pull_group *pref, *pgrp, *pdyna;
158 gmx_ga2la_t ga2la = NULL;
160 if (pull->dbuf_cyl == NULL)
162 snew(pull->dbuf_cyl, pull->ncoord*4);
165 if (cr && DOMAINDECOMP(cr))
167 ga2la = cr->dd->ga2la;
173 r0_2 = dsqr(pull->cyl_r0);
175 /* loop over all groups to make a reference group for each*/
176 for (c = 0; c < pull->ncoord; c++)
178 pcrd = &pull->coord[c];
180 /* pref will be the same group for all pull coordinates */
181 pref = &pull->group[pcrd->group[0]];
182 pgrp = &pull->group[pcrd->group[1]];
183 pdyna = &pull->dyna[c];
184 copy_rvec(pcrd->vec, dir);
191 for (m = 0; m < DIM; m++)
193 g_x[m] = pgrp->x[m] - pcrd->vec[m]*(pcrd->init + pcrd->rate*t);
196 /* loop over all atoms in the main ref group */
197 for (i = 0; i < pref->nat; i++)
202 if (!ga2la_get_home(ga2la, pref->ind[i], &ii))
207 if (ii >= start && ii < end)
209 pbc_dx_aiuc(pbc, x[ii], g_x, dx);
210 inp = iprod(dir, dx);
212 for (m = 0; m < DIM; m++)
214 dr2 += dsqr(dx[m] - inp*dir[m]);
219 /* add to index, to sum of COM, to weight array */
220 if (pdyna->nat_loc >= pdyna->nalloc_loc)
222 pdyna->nalloc_loc = over_alloc_large(pdyna->nat_loc+1);
223 srenew(pdyna->ind_loc, pdyna->nalloc_loc);
224 srenew(pdyna->weight_loc, pdyna->nalloc_loc);
226 pdyna->ind_loc[pdyna->nat_loc] = ii;
227 mass = md->massT[ii];
228 weight = get_weight(sqrt(dr2), pull->cyl_r1, pull->cyl_r0);
229 pdyna->weight_loc[pdyna->nat_loc] = weight;
230 sum_a += mass*weight*inp;
233 pbc_dx_aiuc(pbc, xp[ii], g_x, dx);
234 inp = iprod(dir, dx);
235 sum_ap += mass*weight*inp;
237 wmass += mass*weight;
238 wwmass += mass*sqr(weight);
243 pull->dbuf_cyl[c*4+0] = wmass;
244 pull->dbuf_cyl[c*4+1] = wwmass;
245 pull->dbuf_cyl[c*4+2] = sum_a;
246 pull->dbuf_cyl[c*4+3] = sum_ap;
251 /* Sum the contributions over the nodes */
252 gmx_sumd(pull->ncoord*4, pull->dbuf_cyl, cr);
255 for (c = 0; c < pull->ncoord; c++)
257 pcrd = &pull->coord[c];
259 pdyna = &pull->dyna[c];
260 pgrp = &pull->group[pcrd->group[1]];
262 wmass = pull->dbuf_cyl[c*4+0];
263 wwmass = pull->dbuf_cyl[c*4+1];
264 pdyna->wscale = wmass/wwmass;
265 pdyna->invtm = 1.0/(pdyna->wscale*wmass);
267 for (m = 0; m < DIM; m++)
269 g_x[m] = pgrp->x[m] - pcrd->vec[m]*(pcrd->init + pcrd->rate*t);
270 pdyna->x[m] = g_x[m] + pcrd->vec[m]*pull->dbuf_cyl[c*4+2]/wmass;
273 pdyna->xp[m] = g_x[m] + pcrd->vec[m]*pull->dbuf_cyl[c*4+3]/wmass;
279 fprintf(debug, "Pull cylinder group %d:%8.3f%8.3f%8.3f m:%8.3f\n",
280 c, pdyna->x[0], pdyna->x[1],
281 pdyna->x[2], 1.0/pdyna->invtm);
286 static double atan2_0_2pi(double y, double x)
298 /* calculates center of mass of selection index from all coordinates x */
299 void pull_calc_coms(t_commrec *cr,
300 t_pull *pull, t_mdatoms *md, t_pbc *pbc, double t,
304 real mass, w, wm, twopi_box = 0;
305 double wmass, wwmass, invwmass;
307 double cm, sm, cmp, smp, ccm, csm, ssm, csw, snw;
308 rvec *xx[2], x_pbc = {0, 0, 0}, dx;
311 if (pull->rbuf == NULL)
313 snew(pull->rbuf, pull->ngroup);
315 if (pull->dbuf == NULL)
317 snew(pull->dbuf, 3*pull->ngroup);
322 pull_set_pbcatoms(cr, pull, md, x, pull->rbuf);
325 if (pull->cosdim >= 0)
327 for (m = pull->cosdim+1; m < pull->npbcdim; m++)
329 if (pbc->box[m][pull->cosdim] != 0)
331 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
334 twopi_box = 2.0*M_PI/pbc->box[pull->cosdim][pull->cosdim];
337 for (g = 0; g < pull->ngroup; g++)
339 pgrp = &pull->group[g];
351 if (!(g == 0 && PULL_CYL(pull)))
353 if (pgrp->epgrppbc == epgrppbcREFAT)
355 /* Set the pbc atom */
356 copy_rvec(pull->rbuf[g], x_pbc);
359 for (i = 0; i < pgrp->nat_loc; i++)
361 ii = pgrp->ind_loc[i];
362 mass = md->massT[ii];
363 if (pgrp->epgrppbc != epgrppbcCOS)
365 if (pgrp->weight_loc)
367 w = pgrp->weight_loc[i];
372 if (pgrp->epgrppbc == epgrppbcNONE)
374 /* Plain COM: sum the coordinates */
375 for (m = 0; m < DIM; m++)
377 com[m] += wm*x[ii][m];
381 for (m = 0; m < DIM; m++)
383 comp[m] += wm*xp[ii][m];
389 /* Sum the difference with the reference atom */
390 pbc_dx(pbc, x[ii], x_pbc, dx);
391 for (m = 0; m < DIM; m++)
397 /* For xp add the difference between xp and x to dx,
398 * such that we use the same periodic image,
399 * also when xp has a large displacement.
401 for (m = 0; m < DIM; m++)
403 comp[m] += wm*(dx[m] + xp[ii][m] - x[ii][m]);
410 /* Determine cos and sin sums */
411 csw = cos(x[ii][pull->cosdim]*twopi_box);
412 snw = sin(x[ii][pull->cosdim]*twopi_box);
421 csw = cos(xp[ii][pull->cosdim]*twopi_box);
422 snw = sin(xp[ii][pull->cosdim]*twopi_box);
430 /* Copy local sums to a buffer for global summing */
431 switch (pgrp->epgrppbc)
435 copy_dvec(com, pull->dbuf[g*3]);
436 copy_dvec(comp, pull->dbuf[g*3+1]);
437 pull->dbuf[g*3+2][0] = wmass;
438 pull->dbuf[g*3+2][1] = wwmass;
439 pull->dbuf[g*3+2][2] = 0;
442 pull->dbuf[g*3 ][0] = cm;
443 pull->dbuf[g*3 ][1] = sm;
444 pull->dbuf[g*3 ][2] = 0;
445 pull->dbuf[g*3+1][0] = ccm;
446 pull->dbuf[g*3+1][1] = csm;
447 pull->dbuf[g*3+1][2] = ssm;
448 pull->dbuf[g*3+2][0] = cmp;
449 pull->dbuf[g*3+2][1] = smp;
450 pull->dbuf[g*3+2][2] = 0;
457 /* Sum the contributions over the nodes */
458 gmx_sumd(pull->ngroup*3*DIM, pull->dbuf[0], cr);
461 for (g = 0; g < pull->ngroup; g++)
463 pgrp = &pull->group[g];
464 if (pgrp->nat > 0 && !(g == 0 && PULL_CYL(pull)))
466 if (pgrp->epgrppbc != epgrppbcCOS)
468 /* Determine the inverse mass */
469 wmass = pull->dbuf[g*3+2][0];
470 wwmass = pull->dbuf[g*3+2][1];
472 /* invtm==0 signals a frozen group, so then we should keep it zero */
475 pgrp->wscale = wmass/wwmass;
476 pgrp->invtm = 1.0/(pgrp->wscale*wmass);
478 /* Divide by the total mass */
479 for (m = 0; m < DIM; m++)
481 pgrp->x[m] = pull->dbuf[g*3 ][m]*invwmass;
484 pgrp->xp[m] = pull->dbuf[g*3+1][m]*invwmass;
486 if (pgrp->epgrppbc == epgrppbcREFAT)
488 pgrp->x[m] += pull->rbuf[g][m];
491 pgrp->xp[m] += pull->rbuf[g][m];
498 /* Determine the optimal location of the cosine weight */
499 csw = pull->dbuf[g*3][0];
500 snw = pull->dbuf[g*3][1];
501 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
502 /* Set the weights for the local atoms */
503 wmass = sqrt(csw*csw + snw*snw);
504 wwmass = (pull->dbuf[g*3+1][0]*csw*csw +
505 pull->dbuf[g*3+1][1]*csw*snw +
506 pull->dbuf[g*3+1][2]*snw*snw)/(wmass*wmass);
507 pgrp->wscale = wmass/wwmass;
508 pgrp->invtm = 1.0/(pgrp->wscale*wmass);
509 /* Set the weights for the local atoms */
512 for (i = 0; i < pgrp->nat_loc; i++)
514 ii = pgrp->ind_loc[i];
515 pgrp->weight_loc[i] = csw*cos(twopi_box*x[ii][pull->cosdim]) +
516 snw*sin(twopi_box*x[ii][pull->cosdim]);
520 csw = pull->dbuf[g*3+2][0];
521 snw = pull->dbuf[g*3+2][1];
522 pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
527 fprintf(debug, "Pull group %d wmass %f wwmass %f invtm %f\n",
528 g, wmass, wwmass, pgrp->invtm);
535 /* Calculate the COMs for the cyclinder reference groups */
536 make_cyl_refgrps(cr, pull, md, pbc, t, x, xp);