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.
43 #include "gromacs/utility/futil.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/utility/smalloc.h"
47 #include "types/commrec.h"
49 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/fileio/confio.h"
54 #include "gromacs/pbcutil/pbc.h"
56 #include "gmx_ga2la.h"
58 static void pull_set_pbcatom(t_commrec *cr, t_pull_group *pgrp,
59 t_mdatoms *md, rvec *x,
68 if (!ga2la_get_home(cr->dd->ga2la, pgrp->pbcatom, &a))
78 if (a >= 0 && a < md->homenr)
80 copy_rvec(x[a], x_pbc);
89 copy_rvec(x[pgrp->pbcatom], x_pbc);
93 static void pull_set_pbcatoms(t_commrec *cr, t_pull *pull,
94 t_mdatoms *md, rvec *x,
100 for (g = 0; g < pull->ngroup; g++)
102 if ((g == 0 && PULL_CYL(pull)) || pull->group[g].pbcatom == -1)
104 clear_rvec(x_pbc[g]);
108 pull_set_pbcatom(cr, &pull->group[g], md, x, x_pbc[g]);
109 for (m = 0; m < DIM; m++)
111 if (pull->dim[m] == 0)
120 if (cr && PAR(cr) && n > 0)
122 /* Sum over the nodes to get x_pbc from the home node of pbcatom */
123 gmx_sum(pull->ngroup*DIM, x_pbc[0], cr);
127 /* switch function, x between r and w */
128 static real get_weight(real x, real r1, real r0)
142 weight = (r0 - x)/(r0 - r1);
148 static void make_cyl_refgrps(t_commrec *cr, t_pull *pull, t_mdatoms *md,
149 t_pbc *pbc, double t, rvec *x, rvec *xp)
151 int c, i, ii, m, start, end;
153 double r0_2, sum_a, sum_ap, dr2, mass, weight, wmass, wwmass, inp;
155 t_pull_group *pref, *pgrp, *pdyna;
156 gmx_ga2la_t ga2la = NULL;
158 if (pull->dbuf_cyl == NULL)
160 snew(pull->dbuf_cyl, pull->ncoord*4);
163 if (cr && DOMAINDECOMP(cr))
165 ga2la = cr->dd->ga2la;
171 r0_2 = dsqr(pull->cyl_r0);
173 /* loop over all groups to make a reference group for each*/
174 for (c = 0; c < pull->ncoord; c++)
176 pcrd = &pull->coord[c];
178 /* pref will be the same group for all pull coordinates */
179 pref = &pull->group[pcrd->group[0]];
180 pgrp = &pull->group[pcrd->group[1]];
181 pdyna = &pull->dyna[c];
182 copy_rvec(pcrd->vec, dir);
189 for (m = 0; m < DIM; m++)
191 g_x[m] = pgrp->x[m] - pcrd->vec[m]*(pcrd->init + pcrd->rate*t);
194 /* loop over all atoms in the main ref group */
195 for (i = 0; i < pref->nat; i++)
200 if (!ga2la_get_home(ga2la, pref->ind[i], &ii))
205 if (ii >= start && ii < end)
207 pbc_dx_aiuc(pbc, x[ii], g_x, dx);
208 inp = iprod(dir, dx);
210 for (m = 0; m < DIM; m++)
212 dr2 += dsqr(dx[m] - inp*dir[m]);
217 /* add to index, to sum of COM, to weight array */
218 if (pdyna->nat_loc >= pdyna->nalloc_loc)
220 pdyna->nalloc_loc = over_alloc_large(pdyna->nat_loc+1);
221 srenew(pdyna->ind_loc, pdyna->nalloc_loc);
222 srenew(pdyna->weight_loc, pdyna->nalloc_loc);
224 pdyna->ind_loc[pdyna->nat_loc] = ii;
225 mass = md->massT[ii];
226 weight = get_weight(sqrt(dr2), pull->cyl_r1, pull->cyl_r0);
227 pdyna->weight_loc[pdyna->nat_loc] = weight;
228 sum_a += mass*weight*inp;
231 pbc_dx_aiuc(pbc, xp[ii], g_x, dx);
232 inp = iprod(dir, dx);
233 sum_ap += mass*weight*inp;
235 wmass += mass*weight;
236 wwmass += mass*sqr(weight);
241 pull->dbuf_cyl[c*4+0] = wmass;
242 pull->dbuf_cyl[c*4+1] = wwmass;
243 pull->dbuf_cyl[c*4+2] = sum_a;
244 pull->dbuf_cyl[c*4+3] = sum_ap;
249 /* Sum the contributions over the nodes */
250 gmx_sumd(pull->ncoord*4, pull->dbuf_cyl, cr);
253 for (c = 0; c < pull->ncoord; c++)
255 pcrd = &pull->coord[c];
257 pdyna = &pull->dyna[c];
258 pgrp = &pull->group[pcrd->group[1]];
260 wmass = pull->dbuf_cyl[c*4+0];
261 wwmass = pull->dbuf_cyl[c*4+1];
262 pdyna->wscale = wmass/wwmass;
263 pdyna->invtm = 1.0/(pdyna->wscale*wmass);
265 for (m = 0; m < DIM; m++)
267 g_x[m] = pgrp->x[m] - pcrd->vec[m]*(pcrd->init + pcrd->rate*t);
268 pdyna->x[m] = g_x[m] + pcrd->vec[m]*pull->dbuf_cyl[c*4+2]/wmass;
271 pdyna->xp[m] = g_x[m] + pcrd->vec[m]*pull->dbuf_cyl[c*4+3]/wmass;
277 fprintf(debug, "Pull cylinder group %d:%8.3f%8.3f%8.3f m:%8.3f\n",
278 c, pdyna->x[0], pdyna->x[1],
279 pdyna->x[2], 1.0/pdyna->invtm);
284 static double atan2_0_2pi(double y, double x)
296 /* calculates center of mass of selection index from all coordinates x */
297 void pull_calc_coms(t_commrec *cr,
298 t_pull *pull, t_mdatoms *md, t_pbc *pbc, double t,
302 real mass, w, wm, twopi_box = 0;
303 double wmass, wwmass, invwmass;
305 double cm, sm, cmp, smp, ccm, csm, ssm, csw, snw;
306 rvec *xx[2], x_pbc = {0, 0, 0}, dx;
309 if (pull->rbuf == NULL)
311 snew(pull->rbuf, pull->ngroup);
313 if (pull->dbuf == NULL)
315 snew(pull->dbuf, 3*pull->ngroup);
320 pull_set_pbcatoms(cr, pull, md, x, pull->rbuf);
323 if (pull->cosdim >= 0)
325 for (m = pull->cosdim+1; m < pull->npbcdim; m++)
327 if (pbc->box[m][pull->cosdim] != 0)
329 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
332 twopi_box = 2.0*M_PI/pbc->box[pull->cosdim][pull->cosdim];
335 for (g = 0; g < pull->ngroup; g++)
337 pgrp = &pull->group[g];
349 if (!(g == 0 && PULL_CYL(pull)))
351 if (pgrp->epgrppbc == epgrppbcREFAT)
353 /* Set the pbc atom */
354 copy_rvec(pull->rbuf[g], x_pbc);
357 for (i = 0; i < pgrp->nat_loc; i++)
359 ii = pgrp->ind_loc[i];
360 mass = md->massT[ii];
361 if (pgrp->epgrppbc != epgrppbcCOS)
363 if (pgrp->weight_loc)
365 w = pgrp->weight_loc[i];
370 if (pgrp->epgrppbc == epgrppbcNONE)
372 /* Plain COM: sum the coordinates */
373 for (m = 0; m < DIM; m++)
375 com[m] += wm*x[ii][m];
379 for (m = 0; m < DIM; m++)
381 comp[m] += wm*xp[ii][m];
387 /* Sum the difference with the reference atom */
388 pbc_dx(pbc, x[ii], x_pbc, dx);
389 for (m = 0; m < DIM; m++)
395 /* For xp add the difference between xp and x to dx,
396 * such that we use the same periodic image,
397 * also when xp has a large displacement.
399 for (m = 0; m < DIM; m++)
401 comp[m] += wm*(dx[m] + xp[ii][m] - x[ii][m]);
408 /* Determine cos and sin sums */
409 csw = cos(x[ii][pull->cosdim]*twopi_box);
410 snw = sin(x[ii][pull->cosdim]*twopi_box);
419 csw = cos(xp[ii][pull->cosdim]*twopi_box);
420 snw = sin(xp[ii][pull->cosdim]*twopi_box);
428 /* Copy local sums to a buffer for global summing */
429 switch (pgrp->epgrppbc)
433 copy_dvec(com, pull->dbuf[g*3]);
434 copy_dvec(comp, pull->dbuf[g*3+1]);
435 pull->dbuf[g*3+2][0] = wmass;
436 pull->dbuf[g*3+2][1] = wwmass;
437 pull->dbuf[g*3+2][2] = 0;
440 pull->dbuf[g*3 ][0] = cm;
441 pull->dbuf[g*3 ][1] = sm;
442 pull->dbuf[g*3 ][2] = 0;
443 pull->dbuf[g*3+1][0] = ccm;
444 pull->dbuf[g*3+1][1] = csm;
445 pull->dbuf[g*3+1][2] = ssm;
446 pull->dbuf[g*3+2][0] = cmp;
447 pull->dbuf[g*3+2][1] = smp;
448 pull->dbuf[g*3+2][2] = 0;
455 /* Sum the contributions over the nodes */
456 gmx_sumd(pull->ngroup*3*DIM, pull->dbuf[0], cr);
459 for (g = 0; g < pull->ngroup; g++)
461 pgrp = &pull->group[g];
462 if (pgrp->nat > 0 && !(g == 0 && PULL_CYL(pull)))
464 if (pgrp->epgrppbc != epgrppbcCOS)
466 /* Determine the inverse mass */
467 wmass = pull->dbuf[g*3+2][0];
468 wwmass = pull->dbuf[g*3+2][1];
470 /* invtm==0 signals a frozen group, so then we should keep it zero */
473 pgrp->wscale = wmass/wwmass;
474 pgrp->invtm = 1.0/(pgrp->wscale*wmass);
476 /* Divide by the total mass */
477 for (m = 0; m < DIM; m++)
479 pgrp->x[m] = pull->dbuf[g*3 ][m]*invwmass;
482 pgrp->xp[m] = pull->dbuf[g*3+1][m]*invwmass;
484 if (pgrp->epgrppbc == epgrppbcREFAT)
486 pgrp->x[m] += pull->rbuf[g][m];
489 pgrp->xp[m] += pull->rbuf[g][m];
496 /* Determine the optimal location of the cosine weight */
497 csw = pull->dbuf[g*3][0];
498 snw = pull->dbuf[g*3][1];
499 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
500 /* Set the weights for the local atoms */
501 wmass = sqrt(csw*csw + snw*snw);
502 wwmass = (pull->dbuf[g*3+1][0]*csw*csw +
503 pull->dbuf[g*3+1][1]*csw*snw +
504 pull->dbuf[g*3+1][2]*snw*snw)/(wmass*wmass);
505 pgrp->wscale = wmass/wwmass;
506 pgrp->invtm = 1.0/(pgrp->wscale*wmass);
507 /* Set the weights for the local atoms */
510 for (i = 0; i < pgrp->nat_loc; i++)
512 ii = pgrp->ind_loc[i];
513 pgrp->weight_loc[i] = csw*cos(twopi_box*x[ii][pull->cosdim]) +
514 snw*sin(twopi_box*x[ii][pull->cosdim]);
518 csw = pull->dbuf[g*3+2][0];
519 snw = pull->dbuf[g*3+2][1];
520 pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
525 fprintf(debug, "Pull group %d wmass %f wwmass %f invtm %f\n",
526 g, wmass, wwmass, pgrp->invtm);
533 /* Calculate the COMs for the cyclinder reference groups */
534 make_cyl_refgrps(cr, pull, md, pbc, t, x, xp);