Fix recent compilation issue in pullcode
[gromacs.git] / src / gromacs / pulling / pullutil.cpp
blobd4bd740ace781b2ea6e2a27bf97e0c731cf3bf4e
1 /*
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.
37 #include "gmxpre.h"
39 #include "config.h"
41 #include <assert.h>
42 #include <stdlib.h>
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/legacyheaders/gmx_ga2la.h"
46 #include "gromacs/legacyheaders/macros.h"
47 #include "gromacs/legacyheaders/names.h"
48 #include "gromacs/legacyheaders/network.h"
49 #include "gromacs/legacyheaders/typedefs.h"
50 #include "gromacs/legacyheaders/types/commrec.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/pbcutil/pbc.h"
53 #include "gromacs/pulling/pull.h"
54 #include "gromacs/pulling/pull_internal.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/futil.h"
57 #include "gromacs/utility/real.h"
58 #include "gromacs/utility/smalloc.h"
60 static void pull_reduce_real(t_commrec *cr,
61 pull_comm_t *comm,
62 int n,
63 real *data)
65 if (cr != NULL && PAR(cr))
67 if (comm->bParticipateAll)
69 /* Sum the contributions over all DD ranks */
70 gmx_sum(n, data, cr);
72 else
74 #ifdef GMX_MPI
75 #ifdef MPI_IN_PLACE_EXISTS
76 MPI_Allreduce(MPI_IN_PLACE, data, n, GMX_MPI_REAL, MPI_SUM,
77 comm->mpi_comm_com);
78 #else
79 real *buf;
81 snew(buf, n);
83 MPI_Allreduce(data, buf, n, GMX_MPI_REAL, MPI_SUM,
84 comm->mpi_comm_com);
86 /* Copy the result from the buffer to the input/output data */
87 for (int i = 0; i < n; i++)
89 data[i] = buf[i];
91 sfree(buf);
92 #endif
93 #else
94 gmx_incons("comm->bParticipateAll=FALSE without GMX_MPI");
95 #endif
100 static void pull_reduce_double(t_commrec *cr,
101 pull_comm_t *comm,
102 int n,
103 double *data)
105 if (cr != NULL && PAR(cr))
107 if (comm->bParticipateAll)
109 /* Sum the contributions over all DD ranks */
110 gmx_sumd(n, data, cr);
112 else
114 #ifdef GMX_MPI
115 #ifdef MPI_IN_PLACE_EXISTS
116 MPI_Allreduce(MPI_IN_PLACE, data, n, MPI_DOUBLE, MPI_SUM,
117 comm->mpi_comm_com);
118 #else
119 double *buf;
121 snew(buf, n);
123 MPI_Allreduce(data, buf, n, MPI_DOUBLE, MPI_SUM,
124 comm->mpi_comm_com);
126 /* Copy the result from the buffer to the input/output data */
127 for (int i = 0; i < n; i++)
129 data[i] = buf[i];
131 sfree(buf);
132 #endif
133 #else
134 gmx_incons("comm->bParticipateAll=FALSE without GMX_MPI");
135 #endif
140 static void pull_set_pbcatom(t_commrec *cr, pull_group_work_t *pgrp,
141 rvec *x,
142 rvec x_pbc)
144 int a;
146 if (cr != NULL && DOMAINDECOMP(cr))
148 if (ga2la_get_home(cr->dd->ga2la, pgrp->params.pbcatom, &a))
150 copy_rvec(x[a], x_pbc);
152 else
154 clear_rvec(x_pbc);
157 else
159 copy_rvec(x[pgrp->params.pbcatom], x_pbc);
163 static void pull_set_pbcatoms(t_commrec *cr, struct pull_t *pull,
164 rvec *x,
165 rvec *x_pbc)
167 int g, n;
169 n = 0;
170 for (g = 0; g < pull->ngroup; g++)
172 if (!pull->group[g].bCalcCOM || pull->group[g].params.pbcatom == -1)
174 clear_rvec(x_pbc[g]);
176 else
178 pull_set_pbcatom(cr, &pull->group[g], x, x_pbc[g]);
179 n++;
183 if (cr && PAR(cr) && n > 0)
185 /* Sum over participating ranks to get x_pbc from the home ranks.
186 * This can be very expensive at high parallelization, so we only
187 * do this after each DD repartitioning.
189 pull_reduce_real(cr, &pull->comm, pull->ngroup*DIM, x_pbc[0]);
193 static void make_cyl_refgrps(t_commrec *cr, struct pull_t *pull, t_mdatoms *md,
194 t_pbc *pbc, double t, rvec *x)
196 /* The size and stride per coord for the reduction buffer */
197 const int stride = 9;
198 int c, i, ii, m, start, end;
199 rvec g_x, dx, dir;
200 double inv_cyl_r2;
201 pull_comm_t *comm;
202 gmx_ga2la_t ga2la = NULL;
204 comm = &pull->comm;
206 if (comm->dbuf_cyl == NULL)
208 snew(comm->dbuf_cyl, pull->ncoord*stride);
211 if (cr && DOMAINDECOMP(cr))
213 ga2la = cr->dd->ga2la;
216 start = 0;
217 end = md->homenr;
219 inv_cyl_r2 = 1/dsqr(pull->params.cylinder_r);
221 /* loop over all groups to make a reference group for each*/
222 for (c = 0; c < pull->ncoord; c++)
224 pull_coord_work_t *pcrd;
225 double sum_a, wmass, wwmass;
226 dvec radf_fac0, radf_fac1;
228 pcrd = &pull->coord[c];
230 sum_a = 0;
231 wmass = 0;
232 wwmass = 0;
233 clear_dvec(radf_fac0);
234 clear_dvec(radf_fac1);
236 if (pcrd->params.eGeom == epullgCYL)
238 pull_group_work_t *pref, *pgrp, *pdyna;
240 /* pref will be the same group for all pull coordinates */
241 pref = &pull->group[pcrd->params.group[0]];
242 pgrp = &pull->group[pcrd->params.group[1]];
243 pdyna = &pull->dyna[c];
244 copy_rvec(pcrd->vec, dir);
245 pdyna->nat_loc = 0;
247 /* We calculate distances with respect to the reference location
248 * of this cylinder group (g_x), which we already have now since
249 * we reduced the other group COM over the ranks. This resolves
250 * any PBC issues and we don't need to use a PBC-atom here.
252 if (pcrd->params.rate != 0)
254 /* With rate=0, value_ref is set initially */
255 pcrd->value_ref = pcrd->params.init + pcrd->params.rate*t;
257 for (m = 0; m < DIM; m++)
259 g_x[m] = pgrp->x[m] - pcrd->vec[m]*pcrd->value_ref;
262 /* loop over all atoms in the main ref group */
263 for (i = 0; i < pref->params.nat; i++)
265 ii = pref->params.ind[i];
266 if (ga2la)
268 if (!ga2la_get_home(ga2la, pref->params.ind[i], &ii))
270 ii = -1;
273 if (ii >= start && ii < end)
275 double dr2, dr2_rel, inp;
276 dvec dr;
278 pbc_dx_aiuc(pbc, x[ii], g_x, dx);
279 inp = iprod(dir, dx);
280 dr2 = 0;
281 for (m = 0; m < DIM; m++)
283 /* Determine the radial components */
284 dr[m] = dx[m] - inp*dir[m];
285 dr2 += dr[m]*dr[m];
287 dr2_rel = dr2*inv_cyl_r2;
289 if (dr2_rel < 1)
291 double mass, weight, dweight_r;
292 dvec mdw;
294 /* add to index, to sum of COM, to weight array */
295 if (pdyna->nat_loc >= pdyna->nalloc_loc)
297 pdyna->nalloc_loc = over_alloc_large(pdyna->nat_loc+1);
298 srenew(pdyna->ind_loc, pdyna->nalloc_loc);
299 srenew(pdyna->weight_loc, pdyna->nalloc_loc);
300 srenew(pdyna->mdw, pdyna->nalloc_loc);
301 srenew(pdyna->dv, pdyna->nalloc_loc);
303 pdyna->ind_loc[pdyna->nat_loc] = ii;
305 mass = md->massT[ii];
306 /* The radial weight function is 1-2x^2+x^4,
307 * where x=r/cylinder_r. Since this function depends
308 * on the radial component, we also get radial forces
309 * on both groups.
311 weight = 1 + (-2 + dr2_rel)*dr2_rel;
312 dweight_r = (-4 + 4*dr2_rel)*inv_cyl_r2;
313 pdyna->weight_loc[pdyna->nat_loc] = weight;
314 sum_a += mass*weight*inp;
315 wmass += mass*weight;
316 wwmass += mass*weight*weight;
317 dsvmul(mass*dweight_r, dr, mdw);
318 copy_dvec(mdw, pdyna->mdw[pdyna->nat_loc]);
319 /* Currently we only have the axial component of the
320 * distance (inp) up to an unkown offset. We add this
321 * offset after the reduction needs to determine the
322 * COM of the cylinder group.
324 pdyna->dv[pdyna->nat_loc] = inp;
325 for (m = 0; m < DIM; m++)
327 radf_fac0[m] += mdw[m];
328 radf_fac1[m] += mdw[m]*inp;
330 pdyna->nat_loc++;
335 comm->dbuf_cyl[c*stride+0] = wmass;
336 comm->dbuf_cyl[c*stride+1] = wwmass;
337 comm->dbuf_cyl[c*stride+2] = sum_a;
338 comm->dbuf_cyl[c*stride+3] = radf_fac0[XX];
339 comm->dbuf_cyl[c*stride+4] = radf_fac0[YY];
340 comm->dbuf_cyl[c*stride+5] = radf_fac0[ZZ];
341 comm->dbuf_cyl[c*stride+6] = radf_fac1[XX];
342 comm->dbuf_cyl[c*stride+7] = radf_fac1[YY];
343 comm->dbuf_cyl[c*stride+8] = radf_fac1[ZZ];
346 if (cr != NULL && PAR(cr))
348 /* Sum the contributions over the ranks */
349 pull_reduce_double(cr, comm, pull->ncoord*stride, comm->dbuf_cyl);
352 for (c = 0; c < pull->ncoord; c++)
354 pull_coord_work_t *pcrd;
356 pcrd = &pull->coord[c];
358 if (pcrd->params.eGeom == epullgCYL)
360 pull_group_work_t *pdyna, *pgrp;
361 double wmass, wwmass, dist;
363 pdyna = &pull->dyna[c];
364 pgrp = &pull->group[pcrd->params.group[1]];
366 wmass = comm->dbuf_cyl[c*stride+0];
367 wwmass = comm->dbuf_cyl[c*stride+1];
368 pdyna->mwscale = 1.0/wmass;
369 /* Cylinder pulling can't be used with constraints, but we set
370 * wscale and invtm anyhow, in case someone would like to use them.
372 pdyna->wscale = wmass/wwmass;
373 pdyna->invtm = wwmass/(wmass*wmass);
375 /* We store the deviation of the COM from the reference location
376 * used above, since we need it when we apply the radial forces
377 * to the atoms in the cylinder group.
379 pcrd->cyl_dev = 0;
380 for (m = 0; m < DIM; m++)
382 g_x[m] = pgrp->x[m] - pcrd->vec[m]*pcrd->value_ref;
383 dist = -pcrd->vec[m]*comm->dbuf_cyl[c*stride+2]*pdyna->mwscale;
384 pdyna->x[m] = g_x[m] - dist;
385 pcrd->cyl_dev += dist;
387 /* Now we know the exact COM of the cylinder reference group,
388 * we can determine the radial force factor (ffrad) that when
389 * multiplied with the axial pull force will give the radial
390 * force on the pulled (non-cylinder) group.
392 for (m = 0; m < DIM; m++)
394 pcrd->ffrad[m] = (comm->dbuf_cyl[c*stride+6+m] +
395 comm->dbuf_cyl[c*stride+3+m]*pcrd->cyl_dev)/wmass;
398 if (debug)
400 fprintf(debug, "Pull cylinder group %d:%8.3f%8.3f%8.3f m:%8.3f\n",
401 c, pdyna->x[0], pdyna->x[1],
402 pdyna->x[2], 1.0/pdyna->invtm);
403 fprintf(debug, "ffrad %8.3f %8.3f %8.3f\n",
404 pcrd->ffrad[XX], pcrd->ffrad[YY], pcrd->ffrad[ZZ]);
410 static double atan2_0_2pi(double y, double x)
412 double a;
414 a = atan2(y, x);
415 if (a < 0)
417 a += 2.0*M_PI;
419 return a;
422 /* calculates center of mass of selection index from all coordinates x */
423 void pull_calc_coms(t_commrec *cr,
424 struct pull_t *pull, t_mdatoms *md, t_pbc *pbc, double t,
425 rvec x[], rvec *xp)
427 int g;
428 real twopi_box = 0;
429 pull_comm_t *comm;
431 comm = &pull->comm;
433 if (comm->rbuf == NULL)
435 snew(comm->rbuf, pull->ngroup);
437 if (comm->dbuf == NULL)
439 snew(comm->dbuf, 3*pull->ngroup);
442 if (pull->bRefAt && pull->bSetPBCatoms)
444 pull_set_pbcatoms(cr, pull, x, comm->rbuf);
446 if (cr != NULL && DOMAINDECOMP(cr))
448 /* We can keep these PBC reference coordinates fixed for nstlist
449 * steps, since atoms won't jump over PBC.
450 * This avoids a global reduction at the next nstlist-1 steps.
451 * Note that the exact values of the pbc reference coordinates
452 * are irrelevant, as long all atoms in the group are within
453 * half a box distance of the reference coordinate.
455 pull->bSetPBCatoms = FALSE;
459 if (pull->cosdim >= 0)
461 int m;
463 assert(pull->npbcdim <= DIM);
465 for (m = pull->cosdim+1; m < pull->npbcdim; m++)
467 if (pbc->box[m][pull->cosdim] != 0)
469 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
472 twopi_box = 2.0*M_PI/pbc->box[pull->cosdim][pull->cosdim];
475 for (g = 0; g < pull->ngroup; g++)
477 pull_group_work_t *pgrp;
479 pgrp = &pull->group[g];
481 if (pgrp->bCalcCOM)
483 if (pgrp->epgrppbc != epgrppbcCOS)
485 dvec com, comp;
486 double wmass, wwmass;
487 rvec x_pbc = { 0, 0, 0 };
488 int i;
490 clear_dvec(com);
491 clear_dvec(comp);
492 wmass = 0;
493 wwmass = 0;
495 if (pgrp->epgrppbc == epgrppbcREFAT)
497 /* Set the pbc atom */
498 copy_rvec(comm->rbuf[g], x_pbc);
501 for (i = 0; i < pgrp->nat_loc; i++)
503 int ii, m;
504 real mass, wm;
506 ii = pgrp->ind_loc[i];
507 mass = md->massT[ii];
508 if (pgrp->weight_loc == NULL)
510 wm = mass;
511 wmass += wm;
513 else
515 real w;
517 w = pgrp->weight_loc[i];
518 wm = w*mass;
519 wmass += wm;
520 wwmass += wm*w;
522 if (pgrp->epgrppbc == epgrppbcNONE)
524 /* Plain COM: sum the coordinates */
525 for (m = 0; m < DIM; m++)
527 com[m] += wm*x[ii][m];
529 if (xp)
531 for (m = 0; m < DIM; m++)
533 comp[m] += wm*xp[ii][m];
537 else
539 rvec dx;
541 /* Sum the difference with the reference atom */
542 pbc_dx(pbc, x[ii], x_pbc, dx);
543 for (m = 0; m < DIM; m++)
545 com[m] += wm*dx[m];
547 if (xp)
549 /* For xp add the difference between xp and x to dx,
550 * such that we use the same periodic image,
551 * also when xp has a large displacement.
553 for (m = 0; m < DIM; m++)
555 comp[m] += wm*(dx[m] + xp[ii][m] - x[ii][m]);
561 /* We do this check after the loop above to avoid more nesting.
562 * If we have a single-atom group the mass is irrelevant, so
563 * we can remove the mass factor to avoid division by zero.
564 * Note that with constraint pulling the mass does matter, but
565 * in that case a check group mass != 0 has been done before.
567 if (pgrp->params.nat == 1 && pgrp->nat_loc == 1 && wmass == 0)
569 int m;
571 /* Copy the single atom coordinate */
572 for (m = 0; m < DIM; m++)
574 com[m] = x[pgrp->ind_loc[0]][m];
576 /* Set all mass factors to 1 to get the correct COM */
577 wmass = 1;
578 wwmass = 1;
581 if (pgrp->weight_loc == NULL)
583 wwmass = wmass;
586 /* Copy local sums to a buffer for global summing */
587 copy_dvec(com, comm->dbuf[g*3]);
588 copy_dvec(comp, comm->dbuf[g*3 + 1]);
589 comm->dbuf[g*3 + 2][0] = wmass;
590 comm->dbuf[g*3 + 2][1] = wwmass;
591 comm->dbuf[g*3 + 2][2] = 0;
593 else
595 /* Cosine weighting geometry */
596 double cm, sm, cmp, smp, ccm, csm, ssm, csw, snw;
597 int i;
599 cm = 0;
600 sm = 0;
601 cmp = 0;
602 smp = 0;
603 ccm = 0;
604 csm = 0;
605 ssm = 0;
607 for (i = 0; i < pgrp->nat_loc; i++)
609 int ii;
610 real mass;
612 ii = pgrp->ind_loc[i];
613 mass = md->massT[ii];
614 /* Determine cos and sin sums */
615 csw = cos(x[ii][pull->cosdim]*twopi_box);
616 snw = sin(x[ii][pull->cosdim]*twopi_box);
617 cm += csw*mass;
618 sm += snw*mass;
619 ccm += csw*csw*mass;
620 csm += csw*snw*mass;
621 ssm += snw*snw*mass;
623 if (xp)
625 csw = cos(xp[ii][pull->cosdim]*twopi_box);
626 snw = sin(xp[ii][pull->cosdim]*twopi_box);
627 cmp += csw*mass;
628 smp += snw*mass;
632 /* Copy local sums to a buffer for global summing */
633 comm->dbuf[g*3 ][0] = cm;
634 comm->dbuf[g*3 ][1] = sm;
635 comm->dbuf[g*3 ][2] = 0;
636 comm->dbuf[g*3+1][0] = ccm;
637 comm->dbuf[g*3+1][1] = csm;
638 comm->dbuf[g*3+1][2] = ssm;
639 comm->dbuf[g*3+2][0] = cmp;
640 comm->dbuf[g*3+2][1] = smp;
641 comm->dbuf[g*3+2][2] = 0;
646 pull_reduce_double(cr, comm, pull->ngroup*3*DIM, comm->dbuf[0]);
648 for (g = 0; g < pull->ngroup; g++)
650 pull_group_work_t *pgrp;
652 pgrp = &pull->group[g];
653 if (pgrp->params.nat > 0 && pgrp->bCalcCOM)
655 if (pgrp->epgrppbc != epgrppbcCOS)
657 double wmass, wwmass;
658 int m;
660 /* Determine the inverse mass */
661 wmass = comm->dbuf[g*3+2][0];
662 wwmass = comm->dbuf[g*3+2][1];
663 pgrp->mwscale = 1.0/wmass;
664 /* invtm==0 signals a frozen group, so then we should keep it zero */
665 if (pgrp->invtm != 0)
667 pgrp->wscale = wmass/wwmass;
668 pgrp->invtm = wwmass/(wmass*wmass);
670 /* Divide by the total mass */
671 for (m = 0; m < DIM; m++)
673 pgrp->x[m] = comm->dbuf[g*3 ][m]*pgrp->mwscale;
674 if (xp)
676 pgrp->xp[m] = comm->dbuf[g*3+1][m]*pgrp->mwscale;
678 if (pgrp->epgrppbc == epgrppbcREFAT)
680 pgrp->x[m] += comm->rbuf[g][m];
681 if (xp)
683 pgrp->xp[m] += comm->rbuf[g][m];
688 else
690 /* Cosine weighting geometry */
691 double csw, snw, wmass, wwmass;
692 int i, ii;
694 /* Determine the optimal location of the cosine weight */
695 csw = comm->dbuf[g*3][0];
696 snw = comm->dbuf[g*3][1];
697 pgrp->x[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
698 /* Set the weights for the local atoms */
699 wmass = sqrt(csw*csw + snw*snw);
700 wwmass = (comm->dbuf[g*3+1][0]*csw*csw +
701 comm->dbuf[g*3+1][1]*csw*snw +
702 comm->dbuf[g*3+1][2]*snw*snw)/(wmass*wmass);
704 pgrp->mwscale = 1.0/wmass;
705 pgrp->wscale = wmass/wwmass;
706 pgrp->invtm = wwmass/(wmass*wmass);
707 /* Set the weights for the local atoms */
708 csw *= pgrp->invtm;
709 snw *= pgrp->invtm;
710 for (i = 0; i < pgrp->nat_loc; i++)
712 ii = pgrp->ind_loc[i];
713 pgrp->weight_loc[i] = csw*cos(twopi_box*x[ii][pull->cosdim]) +
714 snw*sin(twopi_box*x[ii][pull->cosdim]);
716 if (xp)
718 csw = comm->dbuf[g*3+2][0];
719 snw = comm->dbuf[g*3+2][1];
720 pgrp->xp[pull->cosdim] = atan2_0_2pi(snw, csw)/twopi_box;
723 if (debug)
725 fprintf(debug, "Pull group %d wmass %f invtm %f\n",
726 g, 1.0/pgrp->mwscale, pgrp->invtm);
731 if (pull->bCylinder)
733 /* Calculate the COMs for the cyclinder reference groups */
734 make_cyl_refgrps(cr, pull, md, pbc, t, x);