Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / pulling / pullutil.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <stdlib.h>
42
43 #include "gromacs/utility/futil.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "typedefs.h"
47 #include "types/commrec.h"
48 #include "names.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "macros.h"
51 #include "index.h"
52 #include "gromacs/fileio/confio.h"
53 #include "network.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "pull.h"
56 #include "gmx_ga2la.h"
57
58 static void pull_set_pbcatom(t_commrec *cr, t_pull_group *pgrp,
59                              t_mdatoms *md, rvec *x,
60                              rvec x_pbc)
61 {
62     int a, m;
63
64     if (cr && PAR(cr))
65     {
66         if (DOMAINDECOMP(cr))
67         {
68             if (!ga2la_get_home(cr->dd->ga2la, pgrp->pbcatom, &a))
69             {
70                 a = -1;
71             }
72         }
73         else
74         {
75             a = pgrp->pbcatom;
76         }
77
78         if (a >= 0 && a < md->homenr)
79         {
80             copy_rvec(x[a], x_pbc);
81         }
82         else
83         {
84             clear_rvec(x_pbc);
85         }
86     }
87     else
88     {
89         copy_rvec(x[pgrp->pbcatom], x_pbc);
90     }
91 }
92
93 static void pull_set_pbcatoms(t_commrec *cr, t_pull *pull,
94                               t_mdatoms *md, rvec *x,
95                               rvec *x_pbc)
96 {
97     int g, n, m;
98
99     n = 0;
100     for (g = 0; g < pull->ngroup; g++)
101     {
102         if ((g == 0 && PULL_CYL(pull)) || pull->group[g].pbcatom == -1)
103         {
104             clear_rvec(x_pbc[g]);
105         }
106         else
107         {
108             pull_set_pbcatom(cr, &pull->group[g], md, x, x_pbc[g]);
109             for (m = 0; m < DIM; m++)
110             {
111                 if (pull->dim[m] == 0)
112                 {
113                     x_pbc[g][m] = 0.0;
114                 }
115             }
116             n++;
117         }
118     }
119
120     if (cr && PAR(cr) && n > 0)
121     {
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);
124     }
125 }
126
127 /* switch function, x between r and w */
128 static real get_weight(real x, real r1, real r0)
129 {
130     real weight;
131
132     if (x >= r0)
133     {
134         weight = 0;
135     }
136     else if (x <= r1)
137     {
138         weight = 1;
139     }
140     else
141     {
142         weight = (r0 - x)/(r0 - r1);
143     }
144
145     return weight;
146 }
147
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)
150 {
151     int           c, i, ii, m, start, end;
152     rvec          g_x, dx, dir;
153     double        r0_2, sum_a, sum_ap, dr2, mass, weight, wmass, wwmass, inp;
154     t_pull_coord *pcrd;
155     t_pull_group *pref, *pgrp, *pdyna;
156     gmx_ga2la_t   ga2la = NULL;
157
158     if (pull->dbuf_cyl == NULL)
159     {
160         snew(pull->dbuf_cyl, pull->ncoord*4);
161     }
162
163     if (cr && DOMAINDECOMP(cr))
164     {
165         ga2la = cr->dd->ga2la;
166     }
167
168     start = 0;
169     end   = md->homenr;
170
171     r0_2 = dsqr(pull->cyl_r0);
172
173     /* loop over all groups to make a reference group for each*/
174     for (c = 0; c < pull->ncoord; c++)
175     {
176         pcrd  = &pull->coord[c];
177
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);
183         sum_a          = 0;
184         sum_ap         = 0;
185         wmass          = 0;
186         wwmass         = 0;
187         pdyna->nat_loc = 0;
188
189         for (m = 0; m < DIM; m++)
190         {
191             g_x[m] = pgrp->x[m] - pcrd->vec[m]*(pcrd->init + pcrd->rate*t);
192         }
193
194         /* loop over all atoms in the main ref group */
195         for (i = 0; i < pref->nat; i++)
196         {
197             ii = pref->ind[i];
198             if (ga2la)
199             {
200                 if (!ga2la_get_home(ga2la, pref->ind[i], &ii))
201                 {
202                     ii = -1;
203                 }
204             }
205             if (ii >= start && ii < end)
206             {
207                 pbc_dx_aiuc(pbc, x[ii], g_x, dx);
208                 inp = iprod(dir, dx);
209                 dr2 = 0;
210                 for (m = 0; m < DIM; m++)
211                 {
212                     dr2 += dsqr(dx[m] - inp*dir[m]);
213                 }
214
215                 if (dr2 < r0_2)
216                 {
217                     /* add to index, to sum of COM, to weight array */
218                     if (pdyna->nat_loc >= pdyna->nalloc_loc)
219                     {
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);
223                     }
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;
229                     if (xp)
230                     {
231                         pbc_dx_aiuc(pbc, xp[ii], g_x, dx);
232                         inp     = iprod(dir, dx);
233                         sum_ap += mass*weight*inp;
234                     }
235                     wmass  += mass*weight;
236                     wwmass += mass*sqr(weight);
237                     pdyna->nat_loc++;
238                 }
239             }
240         }
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;
245     }
246
247     if (cr && PAR(cr))
248     {
249         /* Sum the contributions over the nodes */
250         gmx_sumd(pull->ncoord*4, pull->dbuf_cyl, cr);
251     }
252
253     for (c = 0; c < pull->ncoord; c++)
254     {
255         pcrd  = &pull->coord[c];
256
257         pdyna = &pull->dyna[c];
258         pgrp  = &pull->group[pcrd->group[1]];
259
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);
264
265         for (m = 0; m < DIM; m++)
266         {
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;
269             if (xp)
270             {
271                 pdyna->xp[m] = g_x[m] + pcrd->vec[m]*pull->dbuf_cyl[c*4+3]/wmass;
272             }
273         }
274
275         if (debug)
276         {
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);
280         }
281     }
282 }
283
284 static double atan2_0_2pi(double y, double x)
285 {
286     double a;
287
288     a = atan2(y, x);
289     if (a < 0)
290     {
291         a += 2.0*M_PI;
292     }
293     return a;
294 }
295
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,
299                     rvec x[], rvec *xp)
300 {
301     int           g, i, ii, m;
302     real          mass, w, wm, twopi_box = 0;
303     double        wmass, wwmass, invwmass;
304     dvec          com, comp;
305     double        cm, sm, cmp, smp, ccm, csm, ssm, csw, snw;
306     rvec         *xx[2], x_pbc = {0, 0, 0}, dx;
307     t_pull_group *pgrp;
308
309     if (pull->rbuf == NULL)
310     {
311         snew(pull->rbuf, pull->ngroup);
312     }
313     if (pull->dbuf == NULL)
314     {
315         snew(pull->dbuf, 3*pull->ngroup);
316     }
317
318     if (pull->bRefAt)
319     {
320         pull_set_pbcatoms(cr, pull, md, x, pull->rbuf);
321     }
322
323     if (pull->cosdim >= 0)
324     {
325         for (m = pull->cosdim+1; m < pull->npbcdim; m++)
326         {
327             if (pbc->box[m][pull->cosdim] != 0)
328             {
329                 gmx_fatal(FARGS, "Can not do cosine weighting for trilinic dimensions");
330             }
331         }
332         twopi_box = 2.0*M_PI/pbc->box[pull->cosdim][pull->cosdim];
333     }
334
335     for (g = 0; g < pull->ngroup; g++)
336     {
337         pgrp = &pull->group[g];
338         clear_dvec(com);
339         clear_dvec(comp);
340         wmass  = 0;
341         wwmass = 0;
342         cm     = 0;
343         sm     = 0;
344         cmp    = 0;
345         smp    = 0;
346         ccm    = 0;
347         csm    = 0;
348         ssm    = 0;
349         if (!(g == 0 && PULL_CYL(pull)))
350         {
351             if (pgrp->epgrppbc == epgrppbcREFAT)
352             {
353                 /* Set the pbc atom */
354                 copy_rvec(pull->rbuf[g], x_pbc);
355             }
356             w = 1;
357             for (i = 0; i < pgrp->nat_loc; i++)
358             {
359                 ii   = pgrp->ind_loc[i];
360                 mass = md->massT[ii];
361                 if (pgrp->epgrppbc != epgrppbcCOS)
362                 {
363                     if (pgrp->weight_loc)
364                     {
365                         w = pgrp->weight_loc[i];
366                     }
367                     wm      = w*mass;
368                     wmass  += wm;
369                     wwmass += wm*w;
370                     if (pgrp->epgrppbc == epgrppbcNONE)
371                     {
372                         /* Plain COM: sum the coordinates */
373                         for (m = 0; m < DIM; m++)
374                         {
375                             com[m]    += wm*x[ii][m];
376                         }
377                         if (xp)
378                         {
379                             for (m = 0; m < DIM; m++)
380                             {
381                                 comp[m] += wm*xp[ii][m];
382                             }
383                         }
384                     }
385                     else
386                     {
387                         /* Sum the difference with the reference atom */
388                         pbc_dx(pbc, x[ii], x_pbc, dx);
389                         for (m = 0; m < DIM; m++)
390                         {
391                             com[m]    += wm*dx[m];
392                         }
393                         if (xp)
394                         {
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.
398                              */
399                             for (m = 0; m < DIM; m++)
400                             {
401                                 comp[m] += wm*(dx[m] + xp[ii][m] - x[ii][m]);
402                             }
403                         }
404                     }
405                 }
406                 else
407                 {
408                     /* Determine cos and sin sums */
409                     csw  = cos(x[ii][pull->cosdim]*twopi_box);
410                     snw  = sin(x[ii][pull->cosdim]*twopi_box);
411                     cm  += csw*mass;
412                     sm  += snw*mass;
413                     ccm += csw*csw*mass;
414                     csm += csw*snw*mass;
415                     ssm += snw*snw*mass;
416
417                     if (xp)
418                     {
419                         csw  = cos(xp[ii][pull->cosdim]*twopi_box);
420                         snw  = sin(xp[ii][pull->cosdim]*twopi_box);
421                         cmp += csw*mass;
422                         smp += snw*mass;
423                     }
424                 }
425             }
426         }
427
428         /* Copy local sums to a buffer for global summing */
429         switch (pgrp->epgrppbc)
430         {
431             case epgrppbcNONE:
432             case epgrppbcREFAT:
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;
438                 break;
439             case epgrppbcCOS:
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;
449                 break;
450         }
451     }
452
453     if (cr && PAR(cr))
454     {
455         /* Sum the contributions over the nodes */
456         gmx_sumd(pull->ngroup*3*DIM, pull->dbuf[0], cr);
457     }
458
459     for (g = 0; g < pull->ngroup; g++)
460     {
461         pgrp = &pull->group[g];
462         if (pgrp->nat > 0 && !(g == 0 && PULL_CYL(pull)))
463         {
464             if (pgrp->epgrppbc != epgrppbcCOS)
465             {
466                 /* Determine the inverse mass */
467                 wmass    = pull->dbuf[g*3+2][0];
468                 wwmass   = pull->dbuf[g*3+2][1];
469                 invwmass = 1/wmass;
470                 /* invtm==0 signals a frozen group, so then we should keep it zero */
471                 if (pgrp->invtm > 0)
472                 {
473                     pgrp->wscale = wmass/wwmass;
474                     pgrp->invtm  = 1.0/(pgrp->wscale*wmass);
475                 }
476                 /* Divide by the total mass */
477                 for (m = 0; m < DIM; m++)
478                 {
479                     pgrp->x[m]    = pull->dbuf[g*3  ][m]*invwmass;
480                     if (xp)
481                     {
482                         pgrp->xp[m] = pull->dbuf[g*3+1][m]*invwmass;
483                     }
484                     if (pgrp->epgrppbc == epgrppbcREFAT)
485                     {
486                         pgrp->x[m]    += pull->rbuf[g][m];
487                         if (xp)
488                         {
489                             pgrp->xp[m] += pull->rbuf[g][m];
490                         }
491                     }
492                 }
493             }
494             else
495             {
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 */
508                 csw *= pgrp->invtm;
509                 snw *= pgrp->invtm;
510                 for (i = 0; i < pgrp->nat_loc; i++)
511                 {
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]);
515                 }
516                 if (xp)
517                 {
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;
521                 }
522             }
523             if (debug)
524             {
525                 fprintf(debug, "Pull group %d wmass %f wwmass %f invtm %f\n",
526                         g, wmass, wwmass, pgrp->invtm);
527             }
528         }
529     }
530
531     if (PULL_CYL(pull))
532     {
533         /* Calculate the COMs for the cyclinder reference groups */
534         make_cyl_refgrps(cr, pull, md, pbc, t, x, xp);
535     }
536 }