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