Merging in free energy, exp. ensemble, & andersen t-control to 4.6
[gromacs.git] / src / mdlib / pull.c
blobbb85739bf923c303aa0f83fba349199861ca5993
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * GROwing Monsters And Cloning Shrimps
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
41 #include <math.h>
42 #include <stdio.h>
43 #include <stdlib.h>
44 #include "futil.h"
45 #include "rdgroup.h"
46 #include "statutil.h"
47 #include "gmxfio.h"
48 #include "vec.h"
49 #include "typedefs.h"
50 #include "network.h"
51 #include "filenm.h"
52 #include <string.h>
53 #include "smalloc.h"
54 #include "pull.h"
55 #include "xvgr.h"
56 #include "names.h"
57 #include "partdec.h"
58 #include "pbc.h"
59 #include "mtop_util.h"
60 #include "mdrun.h"
61 #include "gmx_ga2la.h"
62 #include "copyrite.h"
64 static void pull_print_x_grp(FILE *out,gmx_bool bRef,ivec dim,t_pullgrp *pgrp)
66 int m;
68 for(m=0; m<DIM; m++)
70 if (dim[m])
72 fprintf(out,"\t%g",bRef ? pgrp->x[m] : pgrp->dr[m]);
77 static void pull_print_x(FILE *out,t_pull *pull,double t)
79 int g;
81 fprintf(out, "%.4f", t);
83 if (PULL_CYL(pull))
85 for (g=1; g<1+pull->ngrp; g++)
87 pull_print_x_grp(out,TRUE ,pull->dim,&pull->dyna[g]);
88 pull_print_x_grp(out,FALSE,pull->dim,&pull->grp[g]);
91 else
93 for (g=0; g<1+pull->ngrp; g++)
95 if (pull->grp[g].nat > 0)
97 pull_print_x_grp(out,g==0,pull->dim,&pull->grp[g]);
101 fprintf(out,"\n");
104 static void pull_print_f(FILE *out,t_pull *pull,double t)
106 int g,d;
108 fprintf(out, "%.4f", t);
110 for(g=1; g<1+pull->ngrp; g++)
112 if (pull->eGeom == epullgPOS)
114 for(d=0; d<DIM; d++)
116 if (pull->dim[d])
118 fprintf(out,"\t%g",pull->grp[g].f[d]);
122 else
124 fprintf(out,"\t%g",pull->grp[g].f_scal);
127 fprintf(out,"\n");
130 void pull_print_output(t_pull *pull, gmx_large_int_t step, double time)
132 if ((pull->nstxout != 0) && (step % pull->nstxout == 0))
134 pull_print_x(pull->out_x,pull,time);
137 if ((pull->nstfout != 0) && (step % pull->nstfout == 0))
139 pull_print_f(pull->out_f,pull,time);
143 static FILE *open_pull_out(const char *fn,t_pull *pull,const output_env_t oenv,
144 gmx_bool bCoord, unsigned long Flags)
146 FILE *fp;
147 int nsets,g,m;
148 char **setname,buf[10];
150 if(Flags & MD_APPENDFILES)
152 fp = gmx_fio_fopen(fn,"a+");
154 else
156 fp = gmx_fio_fopen(fn,"w+");
157 if (bCoord)
159 xvgr_header(fp,"Pull COM", "Time (ps)","Position (nm)",
160 exvggtXNY,oenv);
162 else
164 xvgr_header(fp,"Pull force","Time (ps)","Force (kJ/mol/nm)",
165 exvggtXNY,oenv);
168 snew(setname,(1+pull->ngrp)*DIM);
169 nsets = 0;
170 for(g=0; g<1+pull->ngrp; g++)
172 if (pull->grp[g].nat > 0 &&
173 (g > 0 || (bCoord && !PULL_CYL(pull))))
175 if (bCoord || pull->eGeom == epullgPOS)
177 if (PULL_CYL(pull))
179 for(m=0; m<DIM; m++)
181 if (pull->dim[m])
183 sprintf(buf,"%d %s%c",g,"c",'X'+m);
184 setname[nsets] = strdup(buf);
185 nsets++;
189 for(m=0; m<DIM; m++)
191 if (pull->dim[m])
193 sprintf(buf,"%d %s%c",
194 g,(bCoord && g > 0)?"d":"",'X'+m);
195 setname[nsets] = strdup(buf);
196 nsets++;
200 else
202 sprintf(buf,"%d",g);
203 setname[nsets] = strdup(buf);
204 nsets++;
208 if (bCoord || nsets > 1)
210 xvgr_legend(fp,nsets,(const char**)setname,oenv);
212 for(g=0; g<nsets; g++)
214 sfree(setname[g]);
216 sfree(setname);
219 return fp;
222 /* Apply forces in a mass weighted fashion */
223 static void apply_forces_grp(t_pullgrp *pgrp, t_mdatoms * md,
224 gmx_ga2la_t ga2la,
225 dvec f_pull, int sign, rvec *f)
227 int i,ii,m,start,end;
228 double wmass,inv_wm;
230 start = md->start;
231 end = md->homenr + start;
233 inv_wm = pgrp->wscale*pgrp->invtm;
235 for(i=0; i<pgrp->nat_loc; i++)
237 ii = pgrp->ind_loc[i];
238 wmass = md->massT[ii];
239 if (pgrp->weight_loc)
241 wmass *= pgrp->weight_loc[i];
244 for(m=0; m<DIM; m++)
246 f[ii][m] += sign * wmass * f_pull[m] * inv_wm;
251 /* Apply forces in a mass weighted fashion */
252 static void apply_forces(t_pull * pull, t_mdatoms * md, gmx_ga2la_t ga2la,
253 rvec *f)
255 int i;
256 t_pullgrp *pgrp;
258 for(i=1; i<pull->ngrp+1; i++)
260 pgrp = &(pull->grp[i]);
261 apply_forces_grp(pgrp,md,ga2la,pgrp->f,1,f);
262 if (pull->grp[0].nat)
264 if (PULL_CYL(pull))
266 apply_forces_grp(&(pull->dyna[i]),md,ga2la,pgrp->f,-1,f);
268 else
270 apply_forces_grp(&(pull->grp[0]),md,ga2la,pgrp->f,-1,f);
276 static double max_pull_distance2(const t_pull *pull,const t_pbc *pbc)
278 double max_d2;
279 int m;
281 max_d2 = GMX_DOUBLE_MAX;
283 if (pull->eGeom != epullgDIRPBC)
285 for(m=0; m<pbc->ndim_ePBC; m++)
287 if (pull->dim[m] != 0)
289 max_d2 = min(max_d2,norm2(pbc->box[m]));
294 return 0.25*max_d2;
297 static void get_pullgrps_dr(const t_pull *pull,const t_pbc *pbc,int g,double t,
298 dvec xg,dvec xref,double max_dist2,
299 dvec dr)
301 t_pullgrp *pref,*pgrp;
302 int m;
303 dvec xrefr,dref={0,0,0};
304 double dr2;
306 pgrp = &pull->grp[g];
308 copy_dvec(xref,xrefr);
310 if (pull->eGeom == epullgDIRPBC)
312 for(m=0; m<DIM; m++)
314 dref[m] = (pgrp->init[0] + pgrp->rate*t)*pull->grp[g].vec[m];
316 /* Add the reference position, so we use the correct periodic image */
317 dvec_inc(xrefr,dref);
320 pbc_dx_d(pbc, xg, xrefr, dr);
321 dr2 = 0;
322 for(m=0; m<DIM; m++)
324 dr[m] *= pull->dim[m];
325 dr2 += dr[m]*dr[m];
327 if (max_dist2 >= 0 && dr2 > 0.98*0.98*max_dist2)
329 gmx_fatal(FARGS,"Distance of pull group %d (%f nm) is larger than 0.49 times the box size (%f)",g,sqrt(dr2),sqrt(max_dist2));
332 if (pull->eGeom == epullgDIRPBC)
334 dvec_inc(dr,dref);
338 static void get_pullgrp_dr(const t_pull *pull,const t_pbc *pbc,int g,double t,
339 dvec dr)
341 double md2;
343 if (pull->eGeom == epullgDIRPBC)
345 md2 = -1;
347 else
349 md2 = max_pull_distance2(pull,pbc);
352 get_pullgrps_dr(pull,pbc,g,t,
353 pull->grp[g].x,
354 PULL_CYL(pull) ? pull->dyna[g].x : pull->grp[0].x,
355 md2,
356 dr);
359 void get_pullgrp_distance(t_pull *pull,t_pbc *pbc,int g,double t,
360 dvec dr,dvec dev)
362 static gmx_bool bWarned=FALSE; /* TODO: this should be fixed for thread-safety,
363 but is fairly benign */
364 t_pullgrp *pgrp;
365 int m;
366 dvec ref;
367 double drs,inpr;
369 pgrp = &pull->grp[g];
371 get_pullgrp_dr(pull,pbc,g,t,dr);
373 if (pull->eGeom == epullgPOS)
375 for(m=0; m<DIM; m++)
377 ref[m] = pgrp->init[m] + pgrp->rate*t*pgrp->vec[m];
380 else
382 ref[0] = pgrp->init[0] + pgrp->rate*t;
385 switch (pull->eGeom)
387 case epullgDIST:
388 /* Pull along the vector between the com's */
389 if (ref[0] < 0 && !bWarned)
391 fprintf(stderr,"\nPull reference distance for group %d is negative (%f)\n",g,ref[0]);
392 bWarned = TRUE;
394 drs = dnorm(dr);
395 if (drs == 0)
397 /* With no vector we can not determine the direction for the force,
398 * so we set the force to zero.
400 dev[0] = 0;
402 else
404 /* Determine the deviation */
405 dev[0] = drs - ref[0];
407 break;
408 case epullgDIR:
409 case epullgDIRPBC:
410 case epullgCYL:
411 /* Pull along vec */
412 inpr = 0;
413 for(m=0; m<DIM; m++)
415 inpr += pgrp->vec[m]*dr[m];
417 dev[0] = inpr - ref[0];
418 break;
419 case epullgPOS:
420 /* Determine the difference of dr and ref along each dimension */
421 for(m=0; m<DIM; m++)
423 dev[m] = (dr[m] - ref[m])*pull->dim[m];
425 break;
429 void clear_pull_forces(t_pull *pull)
431 int i;
433 /* Zeroing the forces is only required for constraint pulling.
434 * It can happen that multiple constraint steps need to be applied
435 * and therefore the constraint forces need to be accumulated.
437 for(i=0; i<1+pull->ngrp; i++)
439 clear_dvec(pull->grp[i].f);
440 pull->grp[i].f_scal = 0;
444 /* Apply constraint using SHAKE */
445 static void do_constraint(t_pull *pull, t_mdatoms *md, t_pbc *pbc,
446 rvec *x, rvec *v,
447 gmx_bool bMaster, tensor vir,
448 double dt, double t)
451 dvec *r_ij; /* x[i] com of i in prev. step. Obeys constr. -> r_ij[i] */
452 dvec unc_ij; /* xp[i] com of i this step, before constr. -> unc_ij */
454 dvec *rinew; /* current 'new' position of group i */
455 dvec *rjnew; /* current 'new' position of group j */
456 dvec ref,vec;
457 double d0,inpr;
458 double lambda, rm, mass, invdt=0;
459 gmx_bool bConverged_all,bConverged=FALSE;
460 int niter=0,g,ii,j,m,max_iter=100;
461 double q,a,b,c; /* for solving the quadratic equation,
462 see Num. Recipes in C ed 2 p. 184 */
463 dvec *dr; /* correction for group i */
464 dvec ref_dr; /* correction for group j */
465 dvec f; /* the pull force */
466 dvec tmp,tmp3;
467 t_pullgrp *pdyna,*pgrp,*pref;
469 snew(r_ij,pull->ngrp+1);
470 if (PULL_CYL(pull))
472 snew(rjnew,pull->ngrp+1);
474 else
476 snew(rjnew,1);
478 snew(dr,pull->ngrp+1);
479 snew(rinew,pull->ngrp+1);
481 /* copy the current unconstrained positions for use in iterations. We
482 iterate until rinew[i] and rjnew[j] obey the constraints. Then
483 rinew - pull.x_unc[i] is the correction dr to group i */
484 for(g=1; g<1+pull->ngrp; g++)
486 copy_dvec(pull->grp[g].xp,rinew[g]);
488 if (PULL_CYL(pull))
490 for(g=1; g<1+pull->ngrp; g++)
492 copy_dvec(pull->dyna[g].xp,rjnew[g]);
495 else
497 copy_dvec(pull->grp[0].xp,rjnew[0]);
500 /* Determine the constraint directions from the old positions */
501 for(g=1; g<1+pull->ngrp; g++)
503 get_pullgrp_dr(pull,pbc,g,t,r_ij[g]);
504 /* Store the difference vector at time t for printing */
505 copy_dvec(r_ij[g],pull->grp[g].dr);
506 if (debug)
508 fprintf(debug,"Pull group %d dr %f %f %f\n",
509 g,r_ij[g][XX],r_ij[g][YY],r_ij[g][ZZ]);
512 if (pull->eGeom == epullgDIR || pull->eGeom == epullgDIRPBC)
514 /* Select the component along vec */
515 a = 0;
516 for(m=0; m<DIM; m++)
518 a += pull->grp[g].vec[m]*r_ij[g][m];
520 for(m=0; m<DIM; m++)
522 r_ij[g][m] = a*pull->grp[g].vec[m];
527 bConverged_all = FALSE;
528 while (!bConverged_all && niter < max_iter)
530 bConverged_all = TRUE;
532 /* loop over all constraints */
533 for(g=1; g<1+pull->ngrp; g++)
535 pgrp = &pull->grp[g];
536 if (PULL_CYL(pull))
537 pref = &pull->dyna[g];
538 else
539 pref = &pull->grp[0];
541 /* Get the current difference vector */
542 get_pullgrps_dr(pull,pbc,g,t,rinew[g],rjnew[PULL_CYL(pull) ? g : 0],
543 -1,unc_ij);
545 if (pull->eGeom == epullgPOS)
547 for(m=0; m<DIM; m++)
549 ref[m] = pgrp->init[m] + pgrp->rate*t*pgrp->vec[m];
552 else
554 ref[0] = pgrp->init[0] + pgrp->rate*t;
555 /* Keep the compiler happy */
556 ref[1] = 0;
557 ref[2] = 0;
560 if (debug)
562 fprintf(debug,"Pull group %d, iteration %d\n",g,niter);
565 rm = 1.0/(pull->grp[g].invtm + pref->invtm);
567 switch (pull->eGeom)
569 case epullgDIST:
570 if (ref[0] <= 0)
572 gmx_fatal(FARGS,"The pull constraint reference distance for group %d is <= 0 (%f)",g,ref[0]);
575 a = diprod(r_ij[g],r_ij[g]);
576 b = diprod(unc_ij,r_ij[g])*2;
577 c = diprod(unc_ij,unc_ij) - dsqr(ref[0]);
579 if (b < 0)
581 q = -0.5*(b - sqrt(b*b - 4*a*c));
582 lambda = -q/a;
584 else
586 q = -0.5*(b + sqrt(b*b - 4*a*c));
587 lambda = -c/q;
590 if (debug)
592 fprintf(debug,
593 "Pull ax^2+bx+c=0: a=%e b=%e c=%e lambda=%e\n",
594 a,b,c,lambda);
597 /* The position corrections dr due to the constraints */
598 dsvmul(-lambda*rm*pgrp->invtm, r_ij[g], dr[g]);
599 dsvmul( lambda*rm*pref->invtm, r_ij[g], ref_dr);
600 break;
601 case epullgDIR:
602 case epullgDIRPBC:
603 case epullgCYL:
604 /* A 1-dimensional constraint along a vector */
605 a = 0;
606 for(m=0; m<DIM; m++)
608 vec[m] = pgrp->vec[m];
609 a += unc_ij[m]*vec[m];
611 /* Select only the component along the vector */
612 dsvmul(a,vec,unc_ij);
613 lambda = a - ref[0];
614 if (debug)
616 fprintf(debug,"Pull inpr %e lambda: %e\n",a,lambda);
619 /* The position corrections dr due to the constraints */
620 dsvmul(-lambda*rm*pull->grp[g].invtm, vec, dr[g]);
621 dsvmul( lambda*rm* pref->invtm, vec,ref_dr);
622 break;
623 case epullgPOS:
624 for(m=0; m<DIM; m++)
626 if (pull->dim[m])
628 lambda = r_ij[g][m] - ref[m];
629 /* The position corrections dr due to the constraints */
630 dr[g][m] = -lambda*rm*pull->grp[g].invtm;
631 ref_dr[m] = lambda*rm*pref->invtm;
633 else
635 dr[g][m] = 0;
636 ref_dr[m] = 0;
639 break;
642 /* DEBUG */
643 if (debug)
645 j = (PULL_CYL(pull) ? g : 0);
646 get_pullgrps_dr(pull,pbc,g,t,rinew[g],rjnew[j],-1,tmp);
647 get_pullgrps_dr(pull,pbc,g,t,dr[g] ,ref_dr ,-1,tmp3);
648 fprintf(debug,
649 "Pull cur %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
650 rinew[g][0],rinew[g][1],rinew[g][2],
651 rjnew[j][0],rjnew[j][1],rjnew[j][2], dnorm(tmp));
652 if (pull->eGeom == epullgPOS)
654 fprintf(debug,
655 "Pull ref %8.5f %8.5f %8.5f\n",
656 pgrp->vec[0],pgrp->vec[1],pgrp->vec[2]);
658 else
660 fprintf(debug,
661 "Pull ref %8s %8s %8s %8s %8s %8s d: %8.5f %8.5f %8.5f\n",
662 "","","","","","",ref[0],ref[1],ref[2]);
664 fprintf(debug,
665 "Pull cor %8.5f %8.5f %8.5f j:%8.5f %8.5f %8.5f d: %8.5f\n",
666 dr[g][0],dr[g][1],dr[g][2],
667 ref_dr[0],ref_dr[1],ref_dr[2],
668 dnorm(tmp3));
669 fprintf(debug,
670 "Pull cor %10.7f %10.7f %10.7f\n",
671 dr[g][0],dr[g][1],dr[g][2]);
672 } /* END DEBUG */
674 /* Update the COMs with dr */
675 dvec_inc(rinew[g], dr[g]);
676 dvec_inc(rjnew[PULL_CYL(pull) ? g : 0],ref_dr);
679 /* Check if all constraints are fullfilled now */
680 for(g=1; g<1+pull->ngrp; g++)
682 pgrp = &pull->grp[g];
684 get_pullgrps_dr(pull,pbc,g,t,rinew[g],rjnew[PULL_CYL(pull) ? g : 0],
685 -1,unc_ij);
687 switch (pull->eGeom)
689 case epullgDIST:
690 bConverged = fabs(dnorm(unc_ij) - ref[0]) < pull->constr_tol;
691 break;
692 case epullgDIR:
693 case epullgDIRPBC:
694 case epullgCYL:
695 for(m=0; m<DIM; m++)
697 vec[m] = pgrp->vec[m];
699 inpr = diprod(unc_ij,vec);
700 dsvmul(inpr,vec,unc_ij);
701 bConverged =
702 fabs(diprod(unc_ij,vec) - ref[0]) < pull->constr_tol;
703 break;
704 case epullgPOS:
705 bConverged = TRUE;
706 for(m=0; m<DIM; m++)
708 if (pull->dim[m] &&
709 fabs(unc_ij[m] - ref[m]) >= pull->constr_tol)
711 bConverged = FALSE;
714 break;
717 if (!bConverged)
719 if (debug)
721 fprintf(debug,"NOT CONVERGED YET: Group %d:"
722 "d_ref = %f %f %f, current d = %f\n",
723 g,ref[0],ref[1],ref[2],dnorm(unc_ij));
726 bConverged_all = FALSE;
730 niter++;
731 /* if after all constraints are dealt with and bConverged is still TRUE
732 we're finished, if not we do another iteration */
734 if (niter > max_iter)
736 gmx_fatal(FARGS,"Too many iterations for constraint run: %d",niter);
739 /* DONE ITERATING, NOW UPDATE COORDINATES AND CALC. CONSTRAINT FORCES */
741 if (v)
743 invdt = 1/dt;
746 /* update the normal groups */
747 for(g=1; g<1+pull->ngrp; g++)
749 pgrp = &pull->grp[g];
750 /* get the final dr and constraint force for group i */
751 dvec_sub(rinew[g],pgrp->xp,dr[g]);
752 /* select components of dr */
753 for(m=0; m<DIM; m++)
755 dr[g][m] *= pull->dim[m];
757 dsvmul(1.0/(pgrp->invtm*dt*dt),dr[g],f);
758 dvec_inc(pgrp->f,f);
759 switch (pull->eGeom)
761 case epullgDIST:
762 for(m=0; m<DIM; m++)
764 pgrp->f_scal += r_ij[g][m]*f[m]/dnorm(r_ij[g]);
766 break;
767 case epullgDIR:
768 case epullgDIRPBC:
769 case epullgCYL:
770 for(m=0; m<DIM; m++)
772 pgrp->f_scal += pgrp->vec[m]*f[m];
774 break;
775 case epullgPOS:
776 break;
779 if (vir && bMaster) {
780 /* Add the pull contribution to the virial */
781 for(j=0; j<DIM; j++)
783 for(m=0; m<DIM; m++)
785 vir[j][m] -= 0.5*f[j]*r_ij[g][m];
790 /* update the atom positions */
791 copy_dvec(dr[g],tmp);
792 for(j=0;j<pgrp->nat_loc;j++)
794 ii = pgrp->ind_loc[j];
795 if (pgrp->weight_loc)
797 dsvmul(pgrp->wscale*pgrp->weight_loc[j],dr[g],tmp);
799 for(m=0; m<DIM; m++)
801 x[ii][m] += tmp[m];
803 if (v)
805 for(m=0; m<DIM; m++)
807 v[ii][m] += invdt*tmp[m];
813 /* update the reference groups */
814 if (PULL_CYL(pull))
816 /* update the dynamic reference groups */
817 for(g=1; g<1+pull->ngrp; g++)
819 pdyna = &pull->dyna[g];
820 dvec_sub(rjnew[g],pdyna->xp,ref_dr);
821 /* select components of ref_dr */
822 for(m=0; m<DIM; m++)
824 ref_dr[m] *= pull->dim[m];
827 for(j=0;j<pdyna->nat_loc;j++)
829 /* reset the atoms with dr, weighted by w_i */
830 dsvmul(pdyna->wscale*pdyna->weight_loc[j],ref_dr,tmp);
831 ii = pdyna->ind_loc[j];
832 for(m=0; m<DIM; m++)
834 x[ii][m] += tmp[m];
836 if (v)
838 for(m=0; m<DIM; m++)
840 v[ii][m] += invdt*tmp[m];
846 else
848 pgrp = &pull->grp[0];
849 /* update the reference group */
850 dvec_sub(rjnew[0],pgrp->xp, ref_dr);
851 /* select components of ref_dr */
852 for(m=0;m<DIM;m++)
854 ref_dr[m] *= pull->dim[m];
857 copy_dvec(ref_dr,tmp);
858 for(j=0; j<pgrp->nat_loc;j++)
860 ii = pgrp->ind_loc[j];
861 if (pgrp->weight_loc)
863 dsvmul(pgrp->wscale*pgrp->weight_loc[j],ref_dr,tmp);
865 for(m=0; m<DIM; m++)
867 x[ii][m] += tmp[m];
869 if (v)
871 for(m=0; m<DIM; m++)
873 v[ii][m] += invdt*tmp[m];
879 /* finished! I hope. Give back some memory */
880 sfree(r_ij);
881 sfree(rinew);
882 sfree(rjnew);
883 sfree(dr);
886 /* Pulling with a harmonic umbrella potential or constant force */
887 static void do_pull_pot(int ePull,
888 t_pull *pull, t_pbc *pbc, double t, real lambda,
889 real *V, tensor vir, real *dVdl)
891 int g,j,m;
892 dvec dev;
893 double ndr,invdr;
894 real k,dkdl;
895 t_pullgrp *pgrp;
897 /* loop over the groups that are being pulled */
898 *V = 0;
899 *dVdl = 0;
900 for(g=1; g<1+pull->ngrp; g++)
902 pgrp = &pull->grp[g];
903 get_pullgrp_distance(pull,pbc,g,t,pgrp->dr,dev);
905 k = (1.0 - lambda)*pgrp->k + lambda*pgrp->kB;
906 dkdl = pgrp->kB - pgrp->k;
908 switch (pull->eGeom)
910 case epullgDIST:
911 ndr = dnorm(pgrp->dr);
912 invdr = 1/ndr;
913 if (ePull == epullUMBRELLA)
915 pgrp->f_scal = -k*dev[0];
916 *V += 0.5* k*dsqr(dev[0]);
917 *dVdl += 0.5*dkdl*dsqr(dev[0]);
919 else
921 pgrp->f_scal = -k;
922 *V += k*ndr;
923 *dVdl += dkdl*ndr;
925 for(m=0; m<DIM; m++)
927 pgrp->f[m] = pgrp->f_scal*pgrp->dr[m]*invdr;
929 break;
930 case epullgDIR:
931 case epullgDIRPBC:
932 case epullgCYL:
933 if (ePull == epullUMBRELLA)
935 pgrp->f_scal = -k*dev[0];
936 *V += 0.5* k*dsqr(dev[0]);
937 *dVdl += 0.5*dkdl*dsqr(dev[0]);
939 else
941 ndr = 0;
942 for(m=0; m<DIM; m++)
944 ndr += pgrp->vec[m]*pgrp->dr[m];
946 pgrp->f_scal = -k;
947 *V += k*ndr;
948 *dVdl += dkdl*ndr;
950 for(m=0; m<DIM; m++)
952 pgrp->f[m] = pgrp->f_scal*pgrp->vec[m];
954 break;
955 case epullgPOS:
956 for(m=0; m<DIM; m++)
958 if (ePull == epullUMBRELLA)
960 pgrp->f[m] = -k*dev[m];
961 *V += 0.5* k*dsqr(dev[m]);
962 *dVdl += 0.5*dkdl*dsqr(dev[m]);
964 else
966 pgrp->f[m] = -k*pull->dim[m];
967 *V += k*pgrp->dr[m]*pull->dim[m];
968 *dVdl += dkdl*pgrp->dr[m]*pull->dim[m];
971 break;
974 if (vir)
976 /* Add the pull contribution to the virial */
977 for(j=0; j<DIM; j++)
979 for(m=0;m<DIM;m++)
981 vir[j][m] -= 0.5*pgrp->f[j]*pgrp->dr[m];
988 real pull_potential(int ePull,t_pull *pull, t_mdatoms *md, t_pbc *pbc,
989 t_commrec *cr, double t, real lambda,
990 rvec *x, rvec *f, tensor vir, real *dvdlambda)
992 real V,dVdl;
994 pull_calc_coms(cr,pull,md,pbc,t,x,NULL);
996 do_pull_pot(ePull,pull,pbc,t,lambda,
997 &V,pull->bVirial && MASTER(cr) ? vir : NULL,&dVdl);
999 /* Distribute forces over pulled groups */
1000 apply_forces(pull, md, DOMAINDECOMP(cr) ? cr->dd->ga2la : NULL, f);
1002 if (MASTER(cr)) {
1003 *dvdlambda += dVdl;
1006 return (MASTER(cr) ? V : 0.0);
1009 void pull_constraint(t_pull *pull, t_mdatoms *md, t_pbc *pbc,
1010 t_commrec *cr, double dt, double t,
1011 rvec *x, rvec *xp, rvec *v, tensor vir)
1013 pull_calc_coms(cr,pull,md,pbc,t,x,xp);
1015 do_constraint(pull,md,pbc,xp,v,pull->bVirial && MASTER(cr),vir,dt,t);
1018 static void make_local_pull_group(gmx_ga2la_t ga2la,
1019 t_pullgrp *pg,int start,int end)
1021 int i,ii;
1023 pg->nat_loc = 0;
1024 for(i=0; i<pg->nat; i++) {
1025 ii = pg->ind[i];
1026 if (ga2la) {
1027 if (!ga2la_get_home(ga2la,ii,&ii)) {
1028 ii = -1;
1031 if (ii >= start && ii < end) {
1032 /* This is a home atom, add it to the local pull group */
1033 if (pg->nat_loc >= pg->nalloc_loc) {
1034 pg->nalloc_loc = over_alloc_dd(pg->nat_loc+1);
1035 srenew(pg->ind_loc,pg->nalloc_loc);
1036 if (pg->epgrppbc == epgrppbcCOS || pg->weight) {
1037 srenew(pg->weight_loc,pg->nalloc_loc);
1040 pg->ind_loc[pg->nat_loc] = ii;
1041 if (pg->weight) {
1042 pg->weight_loc[pg->nat_loc] = pg->weight[i];
1044 pg->nat_loc++;
1049 void dd_make_local_pull_groups(gmx_domdec_t *dd,t_pull *pull,t_mdatoms *md)
1051 gmx_ga2la_t ga2la;
1052 int g;
1054 if (dd) {
1055 ga2la = dd->ga2la;
1056 } else {
1057 ga2la = NULL;
1060 if (pull->grp[0].nat > 0)
1061 make_local_pull_group(ga2la,&pull->grp[0],md->start,md->start+md->homenr);
1062 for(g=1; g<1+pull->ngrp; g++)
1063 make_local_pull_group(ga2la,&pull->grp[g],md->start,md->start+md->homenr);
1066 static void init_pull_group_index(FILE *fplog,t_commrec *cr,
1067 int start,int end,
1068 int g,t_pullgrp *pg,ivec pulldims,
1069 gmx_mtop_t *mtop,t_inputrec *ir, real lambda)
1071 int i,ii,d,nfrozen,ndim;
1072 real m,w,mbd;
1073 double tmass,wmass,wwmass;
1074 gmx_bool bDomDec;
1075 gmx_ga2la_t ga2la=NULL;
1076 gmx_groups_t *groups;
1077 t_atom *atom;
1079 bDomDec = (cr && DOMAINDECOMP(cr));
1080 if (bDomDec) {
1081 ga2la = cr->dd->ga2la;
1084 if (EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD) {
1085 /* There are no masses in the integrator.
1086 * But we still want to have the correct mass-weighted COMs.
1087 * So we store the real masses in the weights.
1088 * We do not set nweight, so these weights do not end up in the tpx file.
1090 if (pg->nweight == 0) {
1091 snew(pg->weight,pg->nat);
1095 if (cr && PAR(cr)) {
1096 pg->nat_loc = 0;
1097 pg->nalloc_loc = 0;
1098 pg->ind_loc = NULL;
1099 pg->weight_loc = NULL;
1100 } else {
1101 pg->nat_loc = pg->nat;
1102 pg->ind_loc = pg->ind;
1103 if (pg->epgrppbc == epgrppbcCOS) {
1104 snew(pg->weight_loc,pg->nat);
1105 } else {
1106 pg->weight_loc = pg->weight;
1110 groups = &mtop->groups;
1112 nfrozen = 0;
1113 tmass = 0;
1114 wmass = 0;
1115 wwmass = 0;
1116 for(i=0; i<pg->nat; i++) {
1117 ii = pg->ind[i];
1118 gmx_mtop_atomnr_to_atom(mtop,ii,&atom);
1119 if (cr && PAR(cr) && !bDomDec && ii >= start && ii < end)
1120 pg->ind_loc[pg->nat_loc++] = ii;
1121 if (ir->opts.nFreeze) {
1122 for(d=0; d<DIM; d++)
1123 if (pulldims[d] && ir->opts.nFreeze[ggrpnr(groups,egcFREEZE,ii)][d])
1124 nfrozen++;
1126 if (ir->efep == efepNO) {
1127 m = atom->m;
1128 } else {
1129 m = (1 - lambda)*atom->m + lambda*atom->mB;
1131 if (pg->nweight > 0) {
1132 w = pg->weight[i];
1133 } else {
1134 w = 1;
1136 if (EI_ENERGY_MINIMIZATION(ir->eI)) {
1137 /* Move the mass to the weight */
1138 w *= m;
1139 m = 1;
1140 pg->weight[i] = w;
1141 } else if (ir->eI == eiBD) {
1142 if (ir->bd_fric) {
1143 mbd = ir->bd_fric*ir->delta_t;
1144 } else {
1145 if (groups->grpnr[egcTC] == NULL) {
1146 mbd = ir->delta_t/ir->opts.tau_t[0];
1147 } else {
1148 mbd = ir->delta_t/ir->opts.tau_t[groups->grpnr[egcTC][ii]];
1151 w *= m/mbd;
1152 m = mbd;
1153 pg->weight[i] = w;
1155 tmass += m;
1156 wmass += m*w;
1157 wwmass += m*w*w;
1160 if (wmass == 0) {
1161 gmx_fatal(FARGS,"The total%s mass of pull group %d is zero",
1162 pg->weight ? " weighted" : "",g);
1164 if (fplog) {
1165 fprintf(fplog,
1166 "Pull group %d: %5d atoms, mass %9.3f",g,pg->nat,tmass);
1167 if (pg->weight || EI_ENERGY_MINIMIZATION(ir->eI) || ir->eI == eiBD) {
1168 fprintf(fplog,", weighted mass %9.3f",wmass*wmass/wwmass);
1170 if (pg->epgrppbc == epgrppbcCOS) {
1171 fprintf(fplog,", cosine weighting will be used");
1173 fprintf(fplog,"\n");
1176 if (nfrozen == 0) {
1177 /* A value > 0 signals not frozen, it is updated later */
1178 pg->invtm = 1.0;
1179 } else {
1180 ndim = 0;
1181 for(d=0; d<DIM; d++)
1182 ndim += pulldims[d]*pg->nat;
1183 if (fplog && nfrozen > 0 && nfrozen < ndim) {
1184 fprintf(fplog,
1185 "\nWARNING: In pull group %d some, but not all of the degrees of freedom\n"
1186 " that are subject to pulling are frozen.\n"
1187 " For pulling the whole group will be frozen.\n\n",
1190 pg->invtm = 0.0;
1191 pg->wscale = 1.0;
1195 void init_pull(FILE *fplog,t_inputrec *ir,int nfile,const t_filenm fnm[],
1196 gmx_mtop_t *mtop,t_commrec *cr,const output_env_t oenv, real lambda,
1197 gmx_bool bOutFile, unsigned long Flags)
1199 t_pull *pull;
1200 t_pullgrp *pgrp;
1201 int g,start=0,end=0,m;
1202 gmx_bool bCite;
1204 pull = ir->pull;
1206 pull->ePBC = ir->ePBC;
1207 switch (pull->ePBC)
1209 case epbcNONE: pull->npbcdim = 0; break;
1210 case epbcXY: pull->npbcdim = 2; break;
1211 default: pull->npbcdim = 3; break;
1214 if (fplog)
1216 fprintf(fplog,"\nWill apply %s COM pulling in geometry '%s'\n",
1217 EPULLTYPE(ir->ePull),EPULLGEOM(pull->eGeom));
1218 if (pull->grp[0].nat > 0)
1220 fprintf(fplog,"between a reference group and %d group%s\n",
1221 pull->ngrp,pull->ngrp==1 ? "" : "s");
1223 else
1225 fprintf(fplog,"with an absolute reference on %d group%s\n",
1226 pull->ngrp,pull->ngrp==1 ? "" : "s");
1228 bCite = FALSE;
1229 for(g=0; g<pull->ngrp+1; g++)
1231 if (pull->grp[g].nat > 1 &&
1232 pull->grp[g].pbcatom < 0)
1234 /* We are using cosine weighting */
1235 fprintf(fplog,"Cosine weighting is used for group %d\n",g);
1236 bCite = TRUE;
1239 if (bCite)
1241 please_cite(fplog,"Engin2010");
1245 /* We always add the virial contribution,
1246 * except for geometry = direction_periodic where this is impossible.
1248 pull->bVirial = (pull->eGeom != epullgDIRPBC);
1249 if (getenv("GMX_NO_PULLVIR") != NULL)
1251 if (fplog)
1253 fprintf(fplog,"Found env. var., will not add the virial contribution of the COM pull forces\n");
1255 pull->bVirial = FALSE;
1258 if (cr && PARTDECOMP(cr))
1260 pd_at_range(cr,&start,&end);
1262 pull->rbuf=NULL;
1263 pull->dbuf=NULL;
1264 pull->dbuf_cyl=NULL;
1265 pull->bRefAt = FALSE;
1266 pull->cosdim = -1;
1267 for(g=0; g<pull->ngrp+1; g++)
1269 pgrp = &pull->grp[g];
1270 pgrp->epgrppbc = epgrppbcNONE;
1271 if (pgrp->nat > 0)
1273 /* Determine if we need to take PBC into account for calculating
1274 * the COM's of the pull groups.
1276 for(m=0; m<pull->npbcdim; m++)
1278 if (pull->dim[m] && pgrp->nat > 1)
1280 if (pgrp->pbcatom >= 0)
1282 pgrp->epgrppbc = epgrppbcREFAT;
1283 pull->bRefAt = TRUE;
1285 else
1287 if (pgrp->weight)
1289 gmx_fatal(FARGS,"Pull groups can not have relative weights and cosine weighting at same time");
1291 pgrp->epgrppbc = epgrppbcCOS;
1292 if (pull->cosdim >= 0 && pull->cosdim != m)
1294 gmx_fatal(FARGS,"Can only use cosine weighting with pulling in one dimension (use mdp option pull_dim)");
1296 pull->cosdim = m;
1300 /* Set the indices */
1301 init_pull_group_index(fplog,cr,start,end,g,pgrp,pull->dim,mtop,ir,lambda);
1302 if (PULL_CYL(pull) && pgrp->invtm == 0)
1304 gmx_fatal(FARGS,"Can not have frozen atoms in a cylinder pull group");
1307 else
1309 /* Absolute reference, set the inverse mass to zero */
1310 pgrp->invtm = 0;
1311 pgrp->wscale = 1;
1315 /* if we use dynamic reference groups, do some initialising for them */
1316 if (PULL_CYL(pull))
1318 if (pull->grp[0].nat == 0)
1320 gmx_fatal(FARGS, "Dynamic reference groups are not supported when using absolute reference!\n");
1322 snew(pull->dyna,pull->ngrp+1);
1325 /* Only do I/O when we are doing dynamics and if we are the MASTER */
1326 pull->out_x = NULL;
1327 pull->out_f = NULL;
1328 if (bOutFile)
1330 if (pull->nstxout > 0)
1332 pull->out_x = open_pull_out(opt2fn("-px",nfile,fnm),pull,oenv,TRUE,Flags);
1334 if (pull->nstfout > 0)
1336 pull->out_f = open_pull_out(opt2fn("-pf",nfile,fnm),pull,oenv,
1337 FALSE,Flags);
1342 void finish_pull(FILE *fplog,t_pull *pull)
1344 if (pull->out_x)
1346 gmx_fio_fclose(pull->out_x);
1348 if (pull->out_f)
1350 gmx_fio_fclose(pull->out_f);