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