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