64bdbfd3a5b3c96264a23a1a43b3b17a9cfb220b
[alexxy/gromacs.git] / src / gromacs / pulling / pull.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
42 #include <math.h>
43 #include <stdio.h>
44 #include <stdlib.h>
45 #include "gromacs/fileio/futil.h"
46 #include "index.h"
47 #include "gromacs/fileio/gmxfio.h"
48 #include "vec.h"
49 #include "typedefs.h"
50 #include "types/commrec.h"
51 #include "network.h"
52 #include "gromacs/fileio/filenm.h"
53 #include <string.h>
54 #include "gromacs/utility/smalloc.h"
55 #include "pull.h"
56 #include "xvgr.h"
57 #include "names.h"
58 #include "pbc.h"
59 #include "mtop_util.h"
60 #include "mdrun.h"
61 #include "gmx_ga2la.h"
62 #include "copyrite.h"
63 #include "macros.h"
64 #include "vec.h"
65
66 static void pull_print_group_x(FILE *out, ivec dim, const t_pull_group *pgrp)
67 {
68     int m;
69
70     for (m = 0; m < DIM; m++)
71     {
72         if (dim[m])
73         {
74             fprintf(out, "\t%g", pgrp->x[m]);
75         }
76     }
77 }
78
79 static void pull_print_coord_dr(FILE *out, ivec dim, const t_pull_coord *pcrd)
80 {
81     int m;
82
83     for (m = 0; m < DIM; m++)
84     {
85         if (dim[m])
86         {
87             fprintf(out, "\t%g", pcrd->dr[m]);
88         }
89     }
90 }
91
92 static void pull_print_x(FILE *out, t_pull *pull, double t)
93 {
94     int                 c;
95     const t_pull_coord *pcrd;
96
97     fprintf(out, "%.4f", t);
98
99     for (c = 0; c < pull->ncoord; c++)
100     {
101         pcrd = &pull->coord[c];
102
103         if (pull->bPrintRef)
104         {
105             if (PULL_CYL(pull))
106             {
107                 pull_print_group_x(out, pull->dim, &pull->dyna[c]);
108             }
109             else
110             {
111                 pull_print_group_x(out, pull->dim, &pull->group[pcrd->group[0]]);
112             }
113         }
114         pull_print_coord_dr(out, pull->dim, pcrd);
115     }
116     fprintf(out, "\n");
117 }
118
119 static void pull_print_f(FILE *out, t_pull *pull, double t)
120 {
121     int c, d;
122
123     fprintf(out, "%.4f", t);
124
125     for (c = 0; c < pull->ncoord; c++)
126     {
127         fprintf(out, "\t%g", pull->coord[c].f_scal);
128     }
129     fprintf(out, "\n");
130 }
131
132 void pull_print_output(t_pull *pull, gmx_int64_t step, double time)
133 {
134     if ((pull->nstxout != 0) && (step % pull->nstxout == 0))
135     {
136         pull_print_x(pull->out_x, pull, time);
137     }
138
139     if ((pull->nstfout != 0) && (step % pull->nstfout == 0))
140     {
141         pull_print_f(pull->out_f, pull, time);
142     }
143 }
144
145 static FILE *open_pull_out(const char *fn, t_pull *pull, const output_env_t oenv,
146                            gmx_bool bCoord, unsigned long Flags)
147 {
148     FILE  *fp;
149     int    nsets, c, m;
150     char **setname, buf[10];
151
152     if (Flags & MD_APPENDFILES)
153     {
154         fp = gmx_fio_fopen(fn, "a+");
155     }
156     else
157     {
158         fp = gmx_fio_fopen(fn, "w+");
159         if (bCoord)
160         {
161             xvgr_header(fp, "Pull COM",  "Time (ps)", "Position (nm)",
162                         exvggtXNY, oenv);
163         }
164         else
165         {
166             xvgr_header(fp, "Pull force", "Time (ps)", "Force (kJ/mol/nm)",
167                         exvggtXNY, oenv);
168         }
169
170         snew(setname, 2*pull->ncoord*DIM);
171         nsets = 0;
172         for (c = 0; c < pull->ncoord; c++)
173         {
174             if (bCoord)
175             {
176                 if (pull->bPrintRef)
177                 {
178                     for (m = 0; m < DIM; m++)
179                     {
180                         if (pull->dim[m])
181                         {
182                             sprintf(buf, "%d %s%c", c+1, "c", 'X'+m);
183                             setname[nsets] = strdup(buf);
184                             nsets++;
185                         }
186                     }
187                 }
188                 for (m = 0; m < DIM; m++)
189                 {
190                     if (pull->dim[m])
191                     {
192                         sprintf(buf, "%d %s%c", c+1, "d", 'X'+m);
193                         setname[nsets] = strdup(buf);
194                         nsets++;
195                     }
196                 }
197             }
198             else
199             {
200                 sprintf(buf, "%d", c+1);
201                 setname[nsets] = strdup(buf);
202                 nsets++;
203             }
204         }
205         if (nsets > 1)
206         {
207             xvgr_legend(fp, nsets, (const char**)setname, oenv);
208         }
209         for (c = 0; c < nsets; c++)
210         {
211             sfree(setname[c]);
212         }
213         sfree(setname);
214     }
215
216     return fp;
217 }
218
219 /* Apply forces in a mass weighted fashion */
220 static void apply_forces_grp(const t_pull_group *pgrp, const t_mdatoms *md,
221                              const dvec f_pull, int sign, rvec *f)
222 {
223     int    i, ii, m;
224     double wmass, inv_wm;
225
226     inv_wm = pgrp->wscale*pgrp->invtm;
227
228     for (i = 0; i < pgrp->nat_loc; i++)
229     {
230         ii    = pgrp->ind_loc[i];
231         wmass = md->massT[ii];
232         if (pgrp->weight_loc)
233         {
234             wmass *= pgrp->weight_loc[i];
235         }
236
237         for (m = 0; m < DIM; m++)
238         {
239             f[ii][m] += sign * wmass * f_pull[m] * inv_wm;
240         }
241     }
242 }
243
244 /* Apply forces in a mass weighted fashion */
245 static void apply_forces(t_pull * pull, t_mdatoms * md, rvec *f)
246 {
247     int                 c;
248     const t_pull_coord *pcrd;
249
250     for (c = 0; c < pull->ncoord; c++)
251     {
252         pcrd = &pull->coord[c];
253
254         if (PULL_CYL(pull))
255         {
256             apply_forces_grp(&pull->dyna[c], md, pcrd->f, -1, f);
257         }
258         else
259         {
260             if (pull->group[pcrd->group[0]].nat > 0)
261             {
262                 apply_forces_grp(&pull->group[pcrd->group[0]], md, pcrd->f, -1, f);
263             }
264         }
265         apply_forces_grp(&pull->group[pcrd->group[1]], md, pcrd->f, 1, f);
266     }
267 }
268
269 static double max_pull_distance2(const t_pull *pull, const t_pbc *pbc)
270 {
271     double max_d2;
272     int    m;
273
274     max_d2 = GMX_DOUBLE_MAX;
275
276     if (pull->eGeom != epullgDIRPBC)
277     {
278         for (m = 0; m < pbc->ndim_ePBC; m++)
279         {
280             if (pull->dim[m] != 0)
281             {
282                 max_d2 = min(max_d2, norm2(pbc->box[m]));
283             }
284         }
285     }
286
287     return 0.25*max_d2;
288 }
289
290 static void low_get_pull_coord_dr(const t_pull *pull,
291                                   const t_pull_coord *pcrd,
292                                   const t_pbc *pbc, double t,
293                                   dvec xg, dvec xref, double max_dist2,
294                                   dvec dr)
295 {
296     const t_pull_group *pgrp0, *pgrp1;
297     int                 m;
298     dvec                xrefr, dref = {0, 0, 0};
299     double              dr2;
300
301     pgrp0 = &pull->group[pcrd->group[0]];
302     pgrp1 = &pull->group[pcrd->group[1]];
303
304     /* Only the first group can be an absolute reference, in that case nat=0 */
305     if (pgrp0->nat == 0)
306     {
307         for (m = 0; m < DIM; m++)
308         {
309             xref[m] = pcrd->origin[m];
310         }
311     }
312
313     copy_dvec(xref, xrefr);
314
315     if (pull->eGeom == epullgDIRPBC)
316     {
317         for (m = 0; m < DIM; m++)
318         {
319             dref[m] = (pcrd->init + pcrd->rate*t)*pcrd->vec[m];
320         }
321         /* Add the reference position, so we use the correct periodic image */
322         dvec_inc(xrefr, dref);
323     }
324
325     pbc_dx_d(pbc, xg, xrefr, dr);
326     dr2 = 0;
327     for (m = 0; m < DIM; m++)
328     {
329         dr[m] *= pull->dim[m];
330         dr2   += dr[m]*dr[m];
331     }
332     if (max_dist2 >= 0 && dr2 > 0.98*0.98*max_dist2)
333     {
334         gmx_fatal(FARGS, "Distance between pull groups %d and %d (%f nm) is larger than 0.49 times the box size (%f).\nYou might want to consider using \"pull-geometry = direction-periodic\" instead.\n",
335                   pcrd->group[0], pcrd->group[1], sqrt(dr2), sqrt(max_dist2));
336     }
337
338     if (pull->eGeom == epullgDIRPBC)
339     {
340         dvec_inc(dr, dref);
341     }
342 }
343
344 static void get_pull_coord_dr(const t_pull *pull,
345                               int coord_ind,
346                               const t_pbc *pbc, double t,
347                               dvec dr)
348 {
349     double              md2;
350     const t_pull_coord *pcrd;
351
352     if (pull->eGeom == epullgDIRPBC)
353     {
354         md2 = -1;
355     }
356     else
357     {
358         md2 = max_pull_distance2(pull, pbc);
359     }
360
361     pcrd = &pull->coord[coord_ind];
362
363     low_get_pull_coord_dr(pull, pcrd, pbc, t,
364                           pull->group[pcrd->group[1]].x,
365                           PULL_CYL(pull) ? pull->dyna[coord_ind].x : pull->group[pcrd->group[0]].x,
366                           md2,
367                           dr);
368 }
369
370 void get_pull_coord_distance(const t_pull *pull,
371                              int coord_ind,
372                              const t_pbc *pbc, double t,
373                              dvec dr, double *dev)
374 {
375     static gmx_bool     bWarned = FALSE; /* TODO: this should be fixed for thread-safety,
376                                             but is fairly benign */
377     const t_pull_coord *pcrd;
378     int                 m;
379     double              ref, drs, inpr;
380
381     pcrd = &pull->coord[coord_ind];
382
383     get_pull_coord_dr(pull, coord_ind, pbc, t, dr);
384
385     ref = pcrd->init + pcrd->rate*t;
386
387     switch (pull->eGeom)
388     {
389         case epullgDIST:
390             /* Pull along the vector between the com's */
391             if (ref < 0 && !bWarned)
392             {
393                 fprintf(stderr, "\nPull reference distance for coordinate %d is negative (%f)\n", coord_ind+1, ref);
394                 bWarned = TRUE;
395             }
396             drs = dnorm(dr);
397             if (drs == 0)
398             {
399                 /* With no vector we can not determine the direction for the force,
400                  * so we set the force to zero.
401                  */
402                 *dev = 0;
403             }
404             else
405             {
406                 /* Determine the deviation */
407                 *dev = drs - ref;
408             }
409             break;
410         case epullgDIR:
411         case epullgDIRPBC:
412         case epullgCYL:
413             /* Pull along vec */
414             inpr = 0;
415             for (m = 0; m < DIM; m++)
416             {
417                 inpr += pcrd->vec[m]*dr[m];
418             }
419             *dev = inpr - ref;
420             break;
421     }
422 }
423
424 void clear_pull_forces(t_pull *pull)
425 {
426     int i;
427
428     /* Zeroing the forces is only required for constraint pulling.
429      * It can happen that multiple constraint steps need to be applied
430      * and therefore the constraint forces need to be accumulated.
431      */
432     for (i = 0; i < pull->ncoord; i++)
433     {
434         clear_dvec(pull->coord[i].f);
435         pull->coord[i].f_scal = 0;
436     }
437 }
438
439 /* Apply constraint using SHAKE */
440 static void do_constraint(t_pull *pull, t_pbc *pbc,
441                           rvec *x, rvec *v,
442                           gmx_bool bMaster, tensor vir,
443                           double dt, double t)
444 {
445
446     dvec         *r_ij;   /* x[i] com of i in prev. step. Obeys constr. -> r_ij[i] */
447     dvec          unc_ij; /* xp[i] com of i this step, before constr.   -> unc_ij  */
448     dvec         *rnew;   /* current 'new' positions of the groups */
449     double       *dr_tot; /* the total update of the coords */
450     double        ref;
451     dvec          vec;
452     double        d0, inpr;
453     double        lambda, rm, mass, invdt = 0;
454     gmx_bool      bConverged_all, bConverged = FALSE;
455     int           niter = 0, g, c, ii, j, m, max_iter = 100;
456     double        a;
457     dvec          f;       /* the pull force */
458     dvec          tmp, tmp3;
459     t_pull_group *pdyna, *pgrp0, *pgrp1;
460     t_pull_coord *pcrd;
461
462     snew(r_ij,   pull->ncoord);
463     snew(dr_tot, pull->ncoord);
464
465     snew(rnew, pull->ngroup);
466
467     /* copy the current unconstrained positions for use in iterations. We
468        iterate until rinew[i] and rjnew[j] obey the constraints. Then
469        rinew - pull.x_unc[i] is the correction dr to group i */
470     for (g = 0; g < pull->ngroup; g++)
471     {
472         copy_dvec(pull->group[g].xp, rnew[g]);
473     }
474     if (PULL_CYL(pull))
475     {
476         /* There is only one pull coordinate and reference group */
477         copy_dvec(pull->dyna[0].xp, rnew[pull->coord[0].group[0]]);
478     }
479
480     /* Determine the constraint directions from the old positions */
481     for (c = 0; c < pull->ncoord; c++)
482     {
483         get_pull_coord_dr(pull, c, pbc, t, r_ij[c]);
484         /* Store the difference vector at time t for printing */
485         copy_dvec(r_ij[c], pull->coord[c].dr);
486         if (debug)
487         {
488             fprintf(debug, "Pull coord %d dr %f %f %f\n",
489                     c, r_ij[c][XX], r_ij[c][YY], r_ij[c][ZZ]);
490         }
491
492         if (pull->eGeom == epullgDIR || pull->eGeom == epullgDIRPBC)
493         {
494             /* Select the component along vec */
495             a = 0;
496             for (m = 0; m < DIM; m++)
497             {
498                 a += pull->coord[c].vec[m]*r_ij[c][m];
499             }
500             for (m = 0; m < DIM; m++)
501             {
502                 r_ij[c][m] = a*pull->coord[c].vec[m];
503             }
504         }
505     }
506
507     bConverged_all = FALSE;
508     while (!bConverged_all && niter < max_iter)
509     {
510         bConverged_all = TRUE;
511
512         /* loop over all constraints */
513         for (c = 0; c < pull->ncoord; c++)
514         {
515             dvec dr0, dr1;
516
517             pcrd  = &pull->coord[c];
518             pgrp0 = &pull->group[pcrd->group[0]];
519             pgrp1 = &pull->group[pcrd->group[1]];
520
521             /* Get the current difference vector */
522             low_get_pull_coord_dr(pull, pcrd, pbc, t,
523                                   rnew[pcrd->group[1]],
524                                   rnew[pcrd->group[0]],
525                                   -1, unc_ij);
526
527             ref = pcrd->init + pcrd->rate*t;
528
529             if (debug)
530             {
531                 fprintf(debug, "Pull coord %d, iteration %d\n", c, niter);
532             }
533
534             rm = 1.0/(pgrp0->invtm + pgrp1->invtm);
535
536             switch (pull->eGeom)
537             {
538                 case epullgDIST:
539                     if (ref <= 0)
540                     {
541                         gmx_fatal(FARGS, "The pull constraint reference distance for group %d is <= 0 (%f)", c, ref);
542                     }
543
544                     {
545                         double q, c_a, c_b, c_c;
546
547                         c_a = diprod(r_ij[c], r_ij[c]);
548                         c_b = diprod(unc_ij, r_ij[c])*2;
549                         c_c = diprod(unc_ij, unc_ij) - dsqr(ref);
550
551                         if (c_b < 0)
552                         {
553                             q      = -0.5*(c_b - sqrt(c_b*c_b - 4*c_a*c_c));
554                             lambda = -q/c_a;
555                         }
556                         else
557                         {
558                             q      = -0.5*(c_b + sqrt(c_b*c_b - 4*c_a*c_c));
559                             lambda = -c_c/q;
560                         }
561
562                         if (debug)
563                         {
564                             fprintf(debug,
565                                     "Pull ax^2+bx+c=0: a=%e b=%e c=%e lambda=%e\n",
566                                     c_a, c_b, c_c, lambda);
567                         }
568                     }
569
570                     /* The position corrections dr due to the constraints */
571                     dsvmul(-lambda*rm*pgrp1->invtm, r_ij[c], dr1);
572                     dsvmul( lambda*rm*pgrp0->invtm, r_ij[c], dr0);
573                     dr_tot[c] += -lambda*dnorm(r_ij[c]);
574                     break;
575                 case epullgDIR:
576                 case epullgDIRPBC:
577                 case epullgCYL:
578                     /* A 1-dimensional constraint along a vector */
579                     a = 0;
580                     for (m = 0; m < DIM; m++)
581                     {
582                         vec[m] = pcrd->vec[m];
583                         a     += unc_ij[m]*vec[m];
584                     }
585                     /* Select only the component along the vector */
586                     dsvmul(a, vec, unc_ij);
587                     lambda = a - ref;
588                     if (debug)
589                     {
590                         fprintf(debug, "Pull inpr %e lambda: %e\n", a, lambda);
591                     }
592
593                     /* The position corrections dr due to the constraints */
594                     dsvmul(-lambda*rm*pgrp1->invtm, vec, dr1);
595                     dsvmul( lambda*rm*pgrp0->invtm, vec, dr0);
596                     dr_tot[c] += -lambda;
597                     break;
598             }
599
600             /* DEBUG */
601             if (debug)
602             {
603                 int g0, g1;
604
605                 g0 = pcrd->group[0];
606                 g1 = pcrd->group[1];
607                 low_get_pull_coord_dr(pull, pcrd, pbc, t, rnew[g1], rnew[g0], -1, tmp);
608                 low_get_pull_coord_dr(pull, pcrd, pbc, t, dr1, dr0, -1, tmp3);
609                 fprintf(debug,
610                         "Pull cur %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
611                         rnew[g0][0], rnew[g0][1], rnew[g0][2],
612                         rnew[g1][0], rnew[g1][1], rnew[g1][2], dnorm(tmp));
613                 fprintf(debug,
614                         "Pull ref %8s %8s %8s   %8s %8s %8s d: %8.5f\n",
615                         "", "", "", "", "", "", ref);
616                 fprintf(debug,
617                         "Pull cor %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
618                         dr0[0], dr0[1], dr0[2],
619                         dr1[0], dr1[1], dr1[2],
620                         dnorm(tmp3));
621             } /* END DEBUG */
622
623             /* Update the COMs with dr */
624             dvec_inc(rnew[pcrd->group[1]], dr1);
625             dvec_inc(rnew[pcrd->group[0]], dr0);
626         }
627
628         /* Check if all constraints are fullfilled now */
629         for (c = 0; c < pull->ncoord; c++)
630         {
631             pcrd = &pull->coord[c];
632             ref  = pcrd->init + pcrd->rate*t;
633
634             low_get_pull_coord_dr(pull, pcrd, pbc, t,
635                                   rnew[pcrd->group[1]],
636                                   rnew[pcrd->group[0]],
637                                   -1, unc_ij);
638
639             switch (pull->eGeom)
640             {
641                 case epullgDIST:
642                     bConverged = fabs(dnorm(unc_ij) - ref) < pull->constr_tol;
643                     break;
644                 case epullgDIR:
645                 case epullgDIRPBC:
646                 case epullgCYL:
647                     for (m = 0; m < DIM; m++)
648                     {
649                         vec[m] = pcrd->vec[m];
650                     }
651                     inpr = diprod(unc_ij, vec);
652                     dsvmul(inpr, vec, unc_ij);
653                     bConverged =
654                         fabs(diprod(unc_ij, vec) - ref) < pull->constr_tol;
655                     break;
656             }
657
658             if (!bConverged)
659             {
660                 if (debug)
661                 {
662                     fprintf(debug, "NOT CONVERGED YET: Group %d:"
663                             "d_ref = %f, current d = %f\n",
664                             g, ref, dnorm(unc_ij));
665                 }
666
667                 bConverged_all = FALSE;
668             }
669         }
670
671         niter++;
672         /* if after all constraints are dealt with and bConverged is still TRUE
673            we're finished, if not we do another iteration */
674     }
675     if (niter > max_iter)
676     {
677         gmx_fatal(FARGS, "Too many iterations for constraint run: %d", niter);
678     }
679
680     /* DONE ITERATING, NOW UPDATE COORDINATES AND CALC. CONSTRAINT FORCES */
681
682     if (v)
683     {
684         invdt = 1/dt;
685     }
686
687     /* update atoms in the groups */
688     for (g = 0; g < pull->ngroup; g++)
689     {
690         const t_pull_group *pgrp;
691         dvec                dr;
692
693         if (PULL_CYL(pull) && g == pull->coord[0].group[0])
694         {
695             pgrp = &pull->dyna[0];
696         }
697         else
698         {
699             pgrp = &pull->group[g];
700         }
701
702         /* get the final constraint displacement dr for group g */
703         dvec_sub(rnew[g], pgrp->xp, dr);
704         /* select components of dr */
705         for (m = 0; m < DIM; m++)
706         {
707             dr[m] *= pull->dim[m];
708         }
709
710         /* update the atom positions */
711         copy_dvec(dr, tmp);
712         for (j = 0; j < pgrp->nat_loc; j++)
713         {
714             ii = pgrp->ind_loc[j];
715             if (pgrp->weight_loc)
716             {
717                 dsvmul(pgrp->wscale*pgrp->weight_loc[j], dr, tmp);
718             }
719             for (m = 0; m < DIM; m++)
720             {
721                 x[ii][m] += tmp[m];
722             }
723             if (v)
724             {
725                 for (m = 0; m < DIM; m++)
726                 {
727                     v[ii][m] += invdt*tmp[m];
728                 }
729             }
730         }
731     }
732
733     /* calculate the constraint forces, used for output and virial only */
734     for (c = 0; c < pull->ncoord; c++)
735     {
736         pcrd         = &pull->coord[c];
737         pcrd->f_scal = dr_tot[c]/((pull->group[pcrd->group[0]].invtm + pull->group[pcrd->group[1]].invtm)*dt*dt);
738
739         if (vir && bMaster)
740         {
741             double f_invr;
742
743             /* Add the pull contribution to the virial */
744             f_invr = pcrd->f_scal/dnorm(r_ij[c]);
745
746             for (j = 0; j < DIM; j++)
747             {
748                 for (m = 0; m < DIM; m++)
749                 {
750                     vir[j][m] -= 0.5*f_invr*r_ij[c][j]*r_ij[c][m];
751                 }
752             }
753         }
754     }
755
756     /* finished! I hope. Give back some memory */
757     sfree(r_ij);
758     sfree(dr_tot);
759     sfree(rnew);
760 }
761
762 /* Pulling with a harmonic umbrella potential or constant force */
763 static void do_pull_pot(int ePull,
764                         t_pull *pull, t_pbc *pbc, double t, real lambda,
765                         real *V, tensor vir, real *dVdl)
766 {
767     int           c, j, m;
768     double        dev, ndr, invdr;
769     real          k, dkdl;
770     t_pull_coord *pcrd;
771
772     /* loop over the pull coordinates */
773     *V    = 0;
774     *dVdl = 0;
775     for (c = 0; c < pull->ncoord; c++)
776     {
777         pcrd = &pull->coord[c];
778
779         get_pull_coord_distance(pull, c, pbc, t, pcrd->dr, &dev);
780
781         k    = (1.0 - lambda)*pcrd->k + lambda*pcrd->kB;
782         dkdl = pcrd->kB - pcrd->k;
783
784         switch (pull->eGeom)
785         {
786             case epullgDIST:
787                 ndr   = dnorm(pcrd->dr);
788                 invdr = 1/ndr;
789                 if (ePull == epullUMBRELLA)
790                 {
791                     pcrd->f_scal  =       -k*dev;
792                     *V           += 0.5*   k*dsqr(dev);
793                     *dVdl        += 0.5*dkdl*dsqr(dev);
794                 }
795                 else
796                 {
797                     pcrd->f_scal  =   -k;
798                     *V           +=    k*ndr;
799                     *dVdl        += dkdl*ndr;
800                 }
801                 for (m = 0; m < DIM; m++)
802                 {
803                     pcrd->f[m]    = pcrd->f_scal*pcrd->dr[m]*invdr;
804                 }
805                 break;
806             case epullgDIR:
807             case epullgDIRPBC:
808             case epullgCYL:
809                 if (ePull == epullUMBRELLA)
810                 {
811                     pcrd->f_scal  =       -k*dev;
812                     *V           += 0.5*   k*dsqr(dev);
813                     *dVdl        += 0.5*dkdl*dsqr(dev);
814                 }
815                 else
816                 {
817                     ndr = 0;
818                     for (m = 0; m < DIM; m++)
819                     {
820                         ndr += pcrd->vec[m]*pcrd->dr[m];
821                     }
822                     pcrd->f_scal  =   -k;
823                     *V           +=    k*ndr;
824                     *dVdl        += dkdl*ndr;
825                 }
826                 for (m = 0; m < DIM; m++)
827                 {
828                     pcrd->f[m]    = pcrd->f_scal*pcrd->vec[m];
829                 }
830                 break;
831         }
832
833         if (vir)
834         {
835             /* Add the pull contribution to the virial */
836             for (j = 0; j < DIM; j++)
837             {
838                 for (m = 0; m < DIM; m++)
839                 {
840                     vir[j][m] -= 0.5*pcrd->f[j]*pcrd->dr[m];
841                 }
842             }
843         }
844     }
845 }
846
847 real pull_potential(int ePull, t_pull *pull, t_mdatoms *md, t_pbc *pbc,
848                     t_commrec *cr, double t, real lambda,
849                     rvec *x, rvec *f, tensor vir, real *dvdlambda)
850 {
851     real V, dVdl;
852
853     pull_calc_coms(cr, pull, md, pbc, t, x, NULL);
854
855     do_pull_pot(ePull, pull, pbc, t, lambda,
856                 &V, pull->bVirial && MASTER(cr) ? vir : NULL, &dVdl);
857
858     /* Distribute forces over pulled groups */
859     apply_forces(pull, md, f);
860
861     if (MASTER(cr))
862     {
863         *dvdlambda += dVdl;
864     }
865
866     return (MASTER(cr) ? V : 0.0);
867 }
868
869 void pull_constraint(t_pull *pull, t_mdatoms *md, t_pbc *pbc,
870                      t_commrec *cr, double dt, double t,
871                      rvec *x, rvec *xp, rvec *v, tensor vir)
872 {
873     pull_calc_coms(cr, pull, md, pbc, t, x, xp);
874
875     do_constraint(pull, pbc, xp, v, pull->bVirial && MASTER(cr), vir, dt, t);
876 }
877
878 static void make_local_pull_group(gmx_ga2la_t ga2la,
879                                   t_pull_group *pg, int start, int end)
880 {
881     int i, ii;
882
883     pg->nat_loc = 0;
884     for (i = 0; i < pg->nat; i++)
885     {
886         ii = pg->ind[i];
887         if (ga2la)
888         {
889             if (!ga2la_get_home(ga2la, ii, &ii))
890             {
891                 ii = -1;
892             }
893         }
894         if (ii >= start && ii < end)
895         {
896             /* This is a home atom, add it to the local pull group */
897             if (pg->nat_loc >= pg->nalloc_loc)
898             {
899                 pg->nalloc_loc = over_alloc_dd(pg->nat_loc+1);
900                 srenew(pg->ind_loc, pg->nalloc_loc);
901                 if (pg->epgrppbc == epgrppbcCOS || pg->weight)
902                 {
903                     srenew(pg->weight_loc, pg->nalloc_loc);
904                 }
905             }
906             pg->ind_loc[pg->nat_loc] = ii;
907             if (pg->weight)
908             {
909                 pg->weight_loc[pg->nat_loc] = pg->weight[i];
910             }
911             pg->nat_loc++;
912         }
913     }
914 }
915
916 void dd_make_local_pull_groups(gmx_domdec_t *dd, t_pull *pull, t_mdatoms *md)
917 {
918     gmx_ga2la_t ga2la;
919     int         g;
920
921     if (dd)
922     {
923         ga2la = dd->ga2la;
924     }
925     else
926     {
927         ga2la = NULL;
928     }
929
930     for (g = 0; g < pull->ngroup; g++)
931     {
932         make_local_pull_group(ga2la, &pull->group[g],
933                               0, md->homenr);
934     }
935 }
936
937 static void init_pull_group_index(FILE *fplog, t_commrec *cr,
938                                   int g, t_pull_group *pg, ivec pulldims,
939                                   gmx_mtop_t *mtop, t_inputrec *ir, real lambda)
940 {
941     int                   i, ii, d, nfrozen, ndim;
942     real                  m, w, mbd;
943     double                tmass, wmass, wwmass;
944     gmx_groups_t         *groups;
945     gmx_mtop_atomlookup_t alook;
946     t_atom               *atom;
947
948     if (EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD)
949     {
950         /* There are no masses in the integrator.
951          * But we still want to have the correct mass-weighted COMs.
952          * So we store the real masses in the weights.
953          * We do not set nweight, so these weights do not end up in the tpx file.
954          */
955         if (pg->nweight == 0)
956         {
957             snew(pg->weight, pg->nat);
958         }
959     }
960
961     if (cr && PAR(cr))
962     {
963         pg->nat_loc    = 0;
964         pg->nalloc_loc = 0;
965         pg->ind_loc    = NULL;
966         pg->weight_loc = NULL;
967     }
968     else
969     {
970         pg->nat_loc = pg->nat;
971         pg->ind_loc = pg->ind;
972         if (pg->epgrppbc == epgrppbcCOS)
973         {
974             snew(pg->weight_loc, pg->nat);
975         }
976         else
977         {
978             pg->weight_loc = pg->weight;
979         }
980     }
981
982     groups = &mtop->groups;
983
984     alook = gmx_mtop_atomlookup_init(mtop);
985
986     nfrozen = 0;
987     tmass   = 0;
988     wmass   = 0;
989     wwmass  = 0;
990     for (i = 0; i < pg->nat; i++)
991     {
992         ii = pg->ind[i];
993         gmx_mtop_atomnr_to_atom(alook, ii, &atom);
994         if (ir->opts.nFreeze)
995         {
996             for (d = 0; d < DIM; d++)
997             {
998                 if (pulldims[d] && ir->opts.nFreeze[ggrpnr(groups, egcFREEZE, ii)][d])
999                 {
1000                     nfrozen++;
1001                 }
1002             }
1003         }
1004         if (ir->efep == efepNO)
1005         {
1006             m = atom->m;
1007         }
1008         else
1009         {
1010             m = (1 - lambda)*atom->m + lambda*atom->mB;
1011         }
1012         if (pg->nweight > 0)
1013         {
1014             w = pg->weight[i];
1015         }
1016         else
1017         {
1018             w = 1;
1019         }
1020         if (EI_ENERGY_MINIMIZATION(ir->eI))
1021         {
1022             /* Move the mass to the weight */
1023             w            *= m;
1024             m             = 1;
1025             pg->weight[i] = w;
1026         }
1027         else if (ir->eI == eiBD)
1028         {
1029             if (ir->bd_fric)
1030             {
1031                 mbd = ir->bd_fric*ir->delta_t;
1032             }
1033             else
1034             {
1035                 if (groups->grpnr[egcTC] == NULL)
1036                 {
1037                     mbd = ir->delta_t/ir->opts.tau_t[0];
1038                 }
1039                 else
1040                 {
1041                     mbd = ir->delta_t/ir->opts.tau_t[groups->grpnr[egcTC][ii]];
1042                 }
1043             }
1044             w            *= m/mbd;
1045             m             = mbd;
1046             pg->weight[i] = w;
1047         }
1048         tmass  += m;
1049         wmass  += m*w;
1050         wwmass += m*w*w;
1051     }
1052
1053     gmx_mtop_atomlookup_destroy(alook);
1054
1055     if (wmass == 0)
1056     {
1057         gmx_fatal(FARGS, "The total%s mass of pull group %d is zero",
1058                   pg->weight ? " weighted" : "", g);
1059     }
1060     if (fplog)
1061     {
1062         fprintf(fplog,
1063                 "Pull group %d: %5d atoms, mass %9.3f", g, pg->nat, tmass);
1064         if (pg->weight || EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD)
1065         {
1066             fprintf(fplog, ", weighted mass %9.3f", wmass*wmass/wwmass);
1067         }
1068         if (pg->epgrppbc == epgrppbcCOS)
1069         {
1070             fprintf(fplog, ", cosine weighting will be used");
1071         }
1072         fprintf(fplog, "\n");
1073     }
1074
1075     if (nfrozen == 0)
1076     {
1077         /* A value > 0 signals not frozen, it is updated later */
1078         pg->invtm  = 1.0;
1079     }
1080     else
1081     {
1082         ndim = 0;
1083         for (d = 0; d < DIM; d++)
1084         {
1085             ndim += pulldims[d]*pg->nat;
1086         }
1087         if (fplog && nfrozen > 0 && nfrozen < ndim)
1088         {
1089             fprintf(fplog,
1090                     "\nWARNING: In pull group %d some, but not all of the degrees of freedom\n"
1091                     "         that are subject to pulling are frozen.\n"
1092                     "         For pulling the whole group will be frozen.\n\n",
1093                     g);
1094         }
1095         pg->invtm  = 0.0;
1096         pg->wscale = 1.0;
1097     }
1098 }
1099
1100 void init_pull(FILE *fplog, t_inputrec *ir, int nfile, const t_filenm fnm[],
1101                gmx_mtop_t *mtop, t_commrec *cr, const output_env_t oenv, real lambda,
1102                gmx_bool bOutFile, unsigned long Flags)
1103 {
1104     t_pull       *pull;
1105     t_pull_group *pgrp;
1106     int           c, g, start = 0, end = 0, m;
1107
1108     pull = ir->pull;
1109
1110     pull->ePBC = ir->ePBC;
1111     switch (pull->ePBC)
1112     {
1113         case epbcNONE: pull->npbcdim = 0; break;
1114         case epbcXY:   pull->npbcdim = 2; break;
1115         default:       pull->npbcdim = 3; break;
1116     }
1117
1118     if (fplog)
1119     {
1120         gmx_bool bAbs, bCos;
1121
1122         bAbs = FALSE;
1123         for (c = 0; c < pull->ncoord; c++)
1124         {
1125             if (pull->group[pull->coord[c].group[0]].nat == 0 ||
1126                 pull->group[pull->coord[c].group[1]].nat == 0)
1127             {
1128                 bAbs = TRUE;
1129             }
1130         }
1131
1132         fprintf(fplog, "\nWill apply %s COM pulling in geometry '%s'\n",
1133                 EPULLTYPE(ir->ePull), EPULLGEOM(pull->eGeom));
1134         fprintf(fplog, "with %d pull coordinate%s and %d group%s\n",
1135                 pull->ncoord, pull->ncoord == 1 ? "" : "s",
1136                 pull->ngroup, pull->ngroup == 1 ? "" : "s");
1137         if (bAbs)
1138         {
1139             fprintf(fplog, "with an absolute reference\n");
1140         }
1141         bCos = FALSE;
1142         for (g = 0; g < pull->ngroup; g++)
1143         {
1144             if (pull->group[g].nat > 1 &&
1145                 pull->group[g].pbcatom < 0)
1146             {
1147                 /* We are using cosine weighting */
1148                 fprintf(fplog, "Cosine weighting is used for group %d\n", g);
1149                 bCos = TRUE;
1150             }
1151         }
1152         if (bCos)
1153         {
1154             please_cite(fplog, "Engin2010");
1155         }
1156     }
1157
1158     /* We always add the virial contribution,
1159      * except for geometry = direction_periodic where this is impossible.
1160      */
1161     pull->bVirial = (pull->eGeom != epullgDIRPBC);
1162     if (getenv("GMX_NO_PULLVIR") != NULL)
1163     {
1164         if (fplog)
1165         {
1166             fprintf(fplog, "Found env. var., will not add the virial contribution of the COM pull forces\n");
1167         }
1168         pull->bVirial = FALSE;
1169     }
1170
1171     pull->rbuf     = NULL;
1172     pull->dbuf     = NULL;
1173     pull->dbuf_cyl = NULL;
1174     pull->bRefAt   = FALSE;
1175     pull->cosdim   = -1;
1176     for (g = 0; g < pull->ngroup; g++)
1177     {
1178         pgrp           = &pull->group[g];
1179         pgrp->epgrppbc = epgrppbcNONE;
1180         if (pgrp->nat > 0)
1181         {
1182             /* Determine if we need to take PBC into account for calculating
1183              * the COM's of the pull groups.
1184              */
1185             for (m = 0; m < pull->npbcdim; m++)
1186             {
1187                 if (pull->dim[m] && pgrp->nat > 1)
1188                 {
1189                     if (pgrp->pbcatom >= 0)
1190                     {
1191                         pgrp->epgrppbc = epgrppbcREFAT;
1192                         pull->bRefAt   = TRUE;
1193                     }
1194                     else
1195                     {
1196                         if (pgrp->weight)
1197                         {
1198                             gmx_fatal(FARGS, "Pull groups can not have relative weights and cosine weighting at same time");
1199                         }
1200                         pgrp->epgrppbc = epgrppbcCOS;
1201                         if (pull->cosdim >= 0 && pull->cosdim != m)
1202                         {
1203                             gmx_fatal(FARGS, "Can only use cosine weighting with pulling in one dimension (use mdp option pull_dim)");
1204                         }
1205                         pull->cosdim = m;
1206                     }
1207                 }
1208             }
1209             /* Set the indices */
1210             init_pull_group_index(fplog, cr, g, pgrp, pull->dim, mtop, ir, lambda);
1211             if (PULL_CYL(pull) && pgrp->invtm == 0)
1212             {
1213                 gmx_fatal(FARGS, "Can not have frozen atoms in a cylinder pull group");
1214             }
1215         }
1216         else
1217         {
1218             /* Absolute reference, set the inverse mass to zero */
1219             pgrp->invtm  = 0;
1220             pgrp->wscale = 1;
1221         }
1222     }
1223
1224     /* if we use dynamic reference groups, do some initialising for them */
1225     if (PULL_CYL(pull))
1226     {
1227         if (ir->ePull == epullCONSTRAINT && pull->ncoord > 1)
1228         {
1229             /* We can't easily update the single reference group with multiple
1230              * constraints. This would require recalculating COMs.
1231              */
1232             gmx_fatal(FARGS, "Constraint COM pulling supports only one coordinate with geometry=cylinder, you can use umbrella pulling with multiple coordinates");
1233         }
1234
1235         for (c = 0; c < pull->ncoord; c++)
1236         {
1237             if (pull->group[pull->coord[c].group[0]].nat == 0)
1238             {
1239                 gmx_fatal(FARGS, "Dynamic reference groups are not supported when using absolute reference!\n");
1240             }
1241         }
1242
1243         snew(pull->dyna, pull->ncoord);
1244     }
1245
1246     /* Only do I/O when we are doing dynamics and if we are the MASTER */
1247     pull->out_x = NULL;
1248     pull->out_f = NULL;
1249     if (bOutFile)
1250     {
1251         if (pull->nstxout > 0)
1252         {
1253             pull->out_x = open_pull_out(opt2fn("-px", nfile, fnm), pull, oenv, TRUE, Flags);
1254         }
1255         if (pull->nstfout > 0)
1256         {
1257             pull->out_f = open_pull_out(opt2fn("-pf", nfile, fnm), pull, oenv,
1258                                         FALSE, Flags);
1259         }
1260     }
1261 }
1262
1263 void finish_pull(t_pull *pull)
1264 {
1265     if (pull->out_x)
1266     {
1267         gmx_fio_fclose(pull->out_x);
1268     }
1269     if (pull->out_f)
1270     {
1271         gmx_fio_fclose(pull->out_f);
1272     }
1273 }