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