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