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