6be238a14baf76144e9862c5c30f49cef8355b64
[alexxy/gromacs.git] / src / mdlib / tpi.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
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.
15
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.
20  * 
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.
27  * 
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.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * GROwing Monsters And Cloning Shrimps
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <string.h>
41 #include <time.h>
42 #include <math.h>
43 #include "sysstuff.h"
44 #include "string2.h"
45 #include "network.h"
46 #include "confio.h"
47 #include "copyrite.h"
48 #include "smalloc.h"
49 #include "nrnb.h"
50 #include "main.h"
51 #include "chargegroup.h"
52 #include "force.h"
53 #include "macros.h"
54 #include "random.h"
55 #include "names.h"
56 #include "gmx_fatal.h"
57 #include "txtdump.h"
58 #include "typedefs.h"
59 #include "update.h"
60 #include "random.h"
61 #include "constr.h"
62 #include "vec.h"
63 #include "statutil.h"
64 #include "tgroup.h"
65 #include "mdebin.h"
66 #include "vsite.h"
67 #include "force.h"
68 #include "mdrun.h"
69 #include "domdec.h"
70 #include "partdec.h"
71 #include "gmx_random.h"
72 #include "physics.h"
73 #include "xvgr.h"
74 #include "mdatoms.h"
75 #include "ns.h"
76 #include "gmx_wallcycle.h"
77 #include "mtop_util.h"
78 #include "gmxfio.h"
79 #include "pme.h"
80 #include "gbutil.h"
81
82
83 static void global_max(t_commrec *cr,int *n)
84 {
85   int *sum,i;
86
87   snew(sum,cr->nnodes);
88   sum[cr->nodeid] = *n;
89   gmx_sumi(cr->nnodes,sum,cr);
90   for(i=0; i<cr->nnodes; i++)
91     *n = max(*n,sum[i]);
92   
93   sfree(sum);
94 }
95
96 static void realloc_bins(double **bin,int *nbin,int nbin_new)
97 {
98   int i;
99
100   if (nbin_new != *nbin) {
101     srenew(*bin,nbin_new);
102     for(i=*nbin; i<nbin_new; i++)
103       (*bin)[i] = 0;
104     *nbin = nbin_new;
105   }
106 }
107
108 double do_tpi(FILE *fplog,t_commrec *cr,
109               int nfile, const t_filenm fnm[],
110               const output_env_t oenv, bool bVerbose,bool bCompact,
111               int nstglobalcomm,
112               gmx_vsite_t *vsite,gmx_constr_t constr,
113               int stepout,
114               t_inputrec *inputrec,
115               gmx_mtop_t *top_global,t_fcdata *fcd,
116               t_state *state,
117               t_mdatoms *mdatoms,
118               t_nrnb *nrnb,gmx_wallcycle_t wcycle,
119               gmx_edsam_t ed,
120               t_forcerec *fr,
121               int repl_ex_nst,int repl_ex_seed,
122               real cpt_period,real max_hours,
123               const char *deviceOptions,
124               unsigned long Flags,
125               gmx_runtime_t *runtime)
126 {
127   const char *TPI="Test Particle Insertion"; 
128   gmx_localtop_t *top;
129   gmx_groups_t *groups;
130   gmx_enerdata_t *enerd;
131   rvec   *f;
132   real   lambda,t,temp,beta,drmax,epot;
133   double embU,sum_embU,*sum_UgembU,V,V_all,VembU_all;
134   t_trxstatus   *status;
135   t_trxframe rerun_fr;
136   bool   bDispCorr,bCharge,bRFExcl,bNotLastFrame,bStateChanged,bNS,bOurStep;
137   tensor force_vir,shake_vir,vir,pres;
138   int    cg_tp,a_tp0,a_tp1,ngid,gid_tp,nener,e;
139   rvec   *x_mol;
140   rvec   mu_tot,x_init,dx,x_tp;
141   int    nnodes,frame,nsteps,step;
142   int    i,start,end;
143   gmx_rng_t tpi_rand;
144   FILE   *fp_tpi=NULL;
145   char   *ptr,*dump_pdb,**leg,str[STRLEN],str2[STRLEN];
146   double dbl,dump_ener;
147   bool   bCavity;
148   int    nat_cavity=0,d;
149   real   *mass_cavity=NULL,mass_tot;
150   int    nbin;
151   double invbinw,*bin,refvolshift,logV,bUlogV;
152   real dvdl,prescorr,enercorr,dvdlcorr;
153   const char *tpid_leg[2]={"direct","reweighted"};
154
155   /* Since numerical problems can lead to extreme negative energies
156    * when atoms overlap, we need to set a lower limit for beta*U.
157    */
158   real bU_neg_limit = -50;
159
160   /* Since there is no upper limit to the insertion energies,
161    * we need to set an upper limit for the distribution output.
162    */
163   real bU_bin_limit      = 50;
164   real bU_logV_bin_limit = bU_bin_limit + 10;
165
166   nnodes = cr->nnodes;
167
168   top = gmx_mtop_generate_local_top(top_global,inputrec);
169
170   groups = &top_global->groups;
171
172   bCavity = (inputrec->eI == eiTPIC);
173   if (bCavity) {
174     ptr = getenv("GMX_TPIC_MASSES");
175     if (ptr == NULL) {
176       nat_cavity = 1;
177     } else {
178       /* Read (multiple) masses from env var GMX_TPIC_MASSES,
179        * The center of mass of the last atoms is then used for TPIC.
180        */
181       nat_cavity = 0;
182       while (sscanf(ptr,"%lf%n",&dbl,&i) > 0) {
183         srenew(mass_cavity,nat_cavity+1);
184         mass_cavity[nat_cavity] = dbl;
185         fprintf(fplog,"mass[%d] = %f\n",
186                 nat_cavity+1,mass_cavity[nat_cavity]);
187         nat_cavity++;
188         ptr += i;
189       }
190       if (nat_cavity == 0)
191         gmx_fatal(FARGS,"Found %d masses in GMX_TPIC_MASSES",nat_cavity);
192     }
193   }
194
195   /*
196   init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
197   state->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
198   /* We never need full pbc for TPI */
199   fr->ePBC = epbcXYZ;
200   /* Determine the temperature for the Boltzmann weighting */
201   temp = inputrec->opts.ref_t[0];
202   if (fplog) {
203     for(i=1; (i<inputrec->opts.ngtc); i++) {
204       if (inputrec->opts.ref_t[i] != temp) {
205         fprintf(fplog,"\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
206         fprintf(stderr,"\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
207       }
208     }
209     fprintf(fplog,
210             "\n  The temperature for test particle insertion is %.3f K\n\n",
211             temp);
212   }
213   beta = 1.0/(BOLTZ*temp);
214
215   /* Number of insertions per frame */
216   nsteps = inputrec->nsteps; 
217
218   /* Use the same neighborlist with more insertions points
219    * in a sphere of radius drmax around the initial point
220    */
221   /* This should be a proper mdp parameter */
222   drmax = inputrec->rtpi;
223
224   /* An environment variable can be set to dump all configurations
225    * to pdb with an insertion energy <= this value.
226    */
227   dump_pdb = getenv("GMX_TPI_DUMP");
228   dump_ener = 0;
229   if (dump_pdb)
230     sscanf(dump_pdb,"%lf",&dump_ener);
231
232   atoms2md(top_global,inputrec,0,NULL,0,top_global->natoms,mdatoms);
233   update_mdatoms(mdatoms,inputrec->init_lambda);
234
235   snew(enerd,1);
236   init_enerdata(groups->grps[egcENER].nr,inputrec->n_flambda,enerd);
237   snew(f,top_global->natoms);
238
239   /* Print to log file  */
240   runtime_start(runtime);
241   print_date_and_time(fplog,cr->nodeid,
242                       "Started Test Particle Insertion",runtime); 
243   wallcycle_start(wcycle,ewcRUN);
244
245   /* The last charge group is the group to be inserted */
246   cg_tp = top->cgs.nr - 1;
247   a_tp0 = top->cgs.index[cg_tp];
248   a_tp1 = top->cgs.index[cg_tp+1];
249   if (debug)
250     fprintf(debug,"TPI cg %d, atoms %d-%d\n",cg_tp,a_tp0,a_tp1);
251   if (a_tp1 - a_tp0 > 1 &&
252       (inputrec->rlist < inputrec->rcoulomb ||
253        inputrec->rlist < inputrec->rvdw))
254     gmx_fatal(FARGS,"Can not do TPI for multi-atom molecule with a twin-range cut-off");
255   snew(x_mol,a_tp1-a_tp0);
256
257   bDispCorr = (inputrec->eDispCorr != edispcNO);
258   bCharge = FALSE;
259   for(i=a_tp0; i<a_tp1; i++) {
260     /* Copy the coordinates of the molecule to be insterted */
261     copy_rvec(state->x[i],x_mol[i-a_tp0]);
262     /* Check if we need to print electrostatic energies */
263     bCharge |= (mdatoms->chargeA[i] != 0 ||
264                 (mdatoms->chargeB && mdatoms->chargeB[i] != 0));
265   }
266   bRFExcl = (bCharge && EEL_RF(fr->eeltype) && fr->eeltype!=eelRF_NEC);
267
268   calc_cgcm(fplog,cg_tp,cg_tp+1,&(top->cgs),state->x,fr->cg_cm);
269   if (bCavity) {
270     if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog) {
271       fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
272       fprintf(stderr,"WARNING: Your TPI molecule is not centered at 0,0,0\n");
273     }
274   } else {
275     /* Center the molecule to be inserted at zero */
276      for(i=0; i<a_tp1-a_tp0; i++)
277       rvec_dec(x_mol[i],fr->cg_cm[cg_tp]);
278   }
279
280   if (fplog) {
281     fprintf(fplog,"\nWill insert %d atoms %s partial charges\n",
282             a_tp1-a_tp0,bCharge ? "with" : "without");
283     
284     fprintf(fplog,"\nWill insert %d times in each frame of %s\n",
285             nsteps,opt2fn("-rerun",nfile,fnm));
286   }
287   
288     if (!bCavity)
289     {
290         if (inputrec->nstlist > 1)
291         {
292             if (drmax==0 && a_tp1-a_tp0==1)
293             {
294                 gmx_fatal(FARGS,"Re-using the neighborlist %d times for insertions of a single atom in a sphere of radius %f does not make sense",inputrec->nstlist,drmax);
295             }
296             if (fplog)
297             {
298                 fprintf(fplog,"Will use the same neighborlist for %d insertions in a sphere of radius %f\n",inputrec->nstlist,drmax);
299             }
300         }
301     }
302     else
303     {
304         if (fplog)
305         {
306             fprintf(fplog,"Will insert randomly in a sphere of radius %f around the center of the cavity\n",drmax);
307         }
308     }
309
310   ngid = groups->grps[egcENER].nr;
311   gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
312   nener = 1 + ngid;
313   if (bDispCorr)
314     nener += 1;
315   if (bCharge) {
316     nener += ngid;
317     if (bRFExcl)
318       nener += 1;
319     if (EEL_FULL(fr->eeltype))
320       nener += 1;
321   }
322   snew(sum_UgembU,nener);
323
324   /* Initialize random generator */
325   tpi_rand = gmx_rng_init(inputrec->ld_seed);
326
327   if (MASTER(cr)) {
328     fp_tpi = xvgropen(opt2fn("-tpi",nfile,fnm),
329                       "TPI energies","Time (ps)",
330                       "(kJ mol\\S-1\\N) / (nm\\S3\\N)",oenv);
331     xvgr_subtitle(fp_tpi,"f. are averages over one frame",oenv);
332     snew(leg,4+nener);
333     e = 0;
334     sprintf(str,"-kT log(<Ve\\S-\\betaU\\N>/<V>)");
335     leg[e++] = strdup(str);
336     sprintf(str,"f. -kT log<e\\S-\\betaU\\N>");
337     leg[e++] = strdup(str);
338     sprintf(str,"f. <e\\S-\\betaU\\N>");
339     leg[e++] = strdup(str);
340     sprintf(str,"f. V");
341     leg[e++] = strdup(str);
342     sprintf(str,"f. <Ue\\S-\\betaU\\N>");
343     leg[e++] = strdup(str);
344     for(i=0; i<ngid; i++) {
345       sprintf(str,"f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
346               *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
347       leg[e++] = strdup(str);
348     }
349     if (bDispCorr) {
350       sprintf(str,"f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
351       leg[e++] = strdup(str);
352     }
353     if (bCharge) {
354       for(i=0; i<ngid; i++) {
355         sprintf(str,"f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
356                 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
357         leg[e++] = strdup(str);
358       }
359       if (bRFExcl) {
360         sprintf(str,"f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
361         leg[e++] = strdup(str);
362       }
363       if (EEL_FULL(fr->eeltype)) {
364         sprintf(str,"f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
365         leg[e++] = strdup(str);
366       }
367     }
368     xvgr_legend(fp_tpi,4+nener,(const char**)leg,oenv);
369     for(i=0; i<4+nener; i++)
370       sfree(leg[i]);
371     sfree(leg);
372   }
373   clear_rvec(x_init);
374   V_all = 0;
375   VembU_all = 0;
376
377   invbinw = 10;
378   nbin = 10;
379   snew(bin,nbin);
380
381   bNotLastFrame = read_first_frame(oenv,&status,opt2fn("-rerun",nfile,fnm),
382                                    &rerun_fr,TRX_NEED_X);
383   frame = 0;
384
385   if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
386       mdatoms->nr - (a_tp1 - a_tp0))
387     gmx_fatal(FARGS,"Number of atoms in trajectory (%d)%s "
388               "is not equal the number in the run input file (%d) "
389               "minus the number of atoms to insert (%d)\n",
390               rerun_fr.natoms,bCavity ? " minus one" : "",
391               mdatoms->nr,a_tp1-a_tp0);
392
393   refvolshift = log(det(rerun_fr.box));
394   
395     while (bNotLastFrame)
396     {
397         lambda = rerun_fr.lambda;
398         t = rerun_fr.time;
399         
400         sum_embU = 0;
401         for(e=0; e<nener; e++)
402         {
403             sum_UgembU[e] = 0;
404         }
405         
406         /* Copy the coordinates from the input trajectory */
407         for(i=0; i<rerun_fr.natoms; i++)
408         {
409             copy_rvec(rerun_fr.x[i],state->x[i]);
410         }
411         
412         V = det(rerun_fr.box);
413         logV = log(V);
414         
415         bStateChanged = TRUE;
416         bNS = TRUE;
417         for(step=0; step<nsteps; step++)
418         {
419             /* In parallel all nodes generate all random configurations.
420              * In that way the result is identical to a single cpu tpi run.
421              */
422             if (!bCavity)
423             {
424                 /* Random insertion in the whole volume */
425                 bNS = (step % inputrec->nstlist == 0);
426                 if (bNS)
427                 {
428                     /* Generate a random position in the box */
429                     x_init[XX] = gmx_rng_uniform_real(tpi_rand)*state->box[XX][XX];
430                     x_init[YY] = gmx_rng_uniform_real(tpi_rand)*state->box[YY][YY];
431                     x_init[ZZ] = gmx_rng_uniform_real(tpi_rand)*state->box[ZZ][ZZ];
432                 }
433                 if (inputrec->nstlist == 1)
434                 {
435                     copy_rvec(x_init,x_tp);
436                 }
437                 else
438                 {
439                     /* Generate coordinates within |dx|=drmax of x_init */
440                     do
441                     {
442                         dx[XX] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
443                         dx[YY] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
444                         dx[ZZ] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
445                     }
446                     while (norm2(dx) > drmax*drmax);
447                     rvec_add(x_init,dx,x_tp);
448                 }
449             }
450             else
451             {
452                 /* Random insertion around a cavity location
453                  * given by the last coordinate of the trajectory.
454                  */
455                 if (step == 0)
456                 {
457                     if (nat_cavity == 1)
458                     {
459                         /* Copy the location of the cavity */
460                         copy_rvec(rerun_fr.x[rerun_fr.natoms-1],x_init);
461                     }
462                     else
463                     {
464                         /* Determine the center of mass of the last molecule */
465                         clear_rvec(x_init);
466                         mass_tot = 0;
467                         for(i=0; i<nat_cavity; i++)
468                         {
469                             for(d=0; d<DIM; d++)
470                             {
471                                 x_init[d] +=
472                                     mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
473                             }
474                             mass_tot += mass_cavity[i];
475                         }
476                         for(d=0; d<DIM; d++)
477                         {
478                             x_init[d] /= mass_tot;
479                         }
480                     }
481                 }
482                 /* Generate coordinates within |dx|=drmax of x_init */
483                 do
484                 {
485                     dx[XX] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
486                     dx[YY] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
487                     dx[ZZ] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
488                 }
489                 while (norm2(dx) > drmax*drmax);
490                 rvec_add(x_init,dx,x_tp);
491             }
492             
493             if (a_tp1 - a_tp0 == 1)
494             {
495                 /* Insert a single atom, just copy the insertion location */
496         copy_rvec(x_tp,state->x[a_tp0]);
497             }
498             else
499             {
500                 /* Copy the coordinates from the top file */
501                 for(i=a_tp0; i<a_tp1; i++)
502                 {
503                     copy_rvec(x_mol[i-a_tp0],state->x[i]);
504                 }
505                 /* Rotate the molecule randomly */
506                 rotate_conf(a_tp1-a_tp0,state->x+a_tp0,NULL,
507                             2*M_PI*gmx_rng_uniform_real(tpi_rand),
508                             2*M_PI*gmx_rng_uniform_real(tpi_rand),
509                             2*M_PI*gmx_rng_uniform_real(tpi_rand));
510                 /* Shift to the insertion location */
511                 for(i=a_tp0; i<a_tp1; i++)
512                 {
513                     rvec_inc(state->x[i],x_tp);
514                 }
515             }
516             
517             /* Check if this insertion belongs to this node */
518             bOurStep = TRUE;
519             if (PAR(cr))
520             {
521                 switch (inputrec->eI)
522                 {
523                 case eiTPI:
524                     bOurStep = ((step / inputrec->nstlist) % nnodes == cr->nodeid);
525                     break;
526                 case eiTPIC:
527                     bOurStep = (step % nnodes == cr->nodeid);
528                     break;
529                 default:
530                     gmx_fatal(FARGS,"Unknown integrator %s",ei_names[inputrec->eI]);
531                 }
532             }
533             if (bOurStep)
534             {
535                 /* Clear some matrix variables  */
536                 clear_mat(force_vir); 
537                 clear_mat(shake_vir);
538                 clear_mat(vir);
539                 clear_mat(pres);
540                 
541                 /* Set the charge group center of mass of the test particle */
542                 copy_rvec(x_init,fr->cg_cm[top->cgs.nr-1]);
543                 
544                 /* Calc energy (no forces) on new positions.
545                  * Since we only need the intermolecular energy
546                  * and the RF exclusion terms of the inserted molecule occur
547                  * within a single charge group we can pass NULL for the graph.
548                  * This also avoids shifts that would move charge groups
549                  * out of the box.
550                  *
551                  * Some checks above ensure than we can not have
552                  * twin-range interactions together with nstlist > 1,
553                  * therefore we do not need to remember the LR energies.
554                  */
555                 /* Make do_force do a single node force calculation */
556                 cr->nnodes = 1;
557                 do_force(fplog,cr,inputrec,
558                          step,nrnb,wcycle,top,top_global,&top_global->groups,
559                          rerun_fr.box,state->x,&state->hist,
560                          f,force_vir,mdatoms,enerd,fcd,
561                          lambda,NULL,fr,NULL,mu_tot,t,NULL,NULL,FALSE,
562                          GMX_FORCE_NONBONDED |
563                          (bNS ? GMX_FORCE_NS | GMX_FORCE_DOLR : 0) |
564                          (bStateChanged ? GMX_FORCE_STATECHANGED : 0)); 
565                 cr->nnodes = nnodes;
566                 bStateChanged = FALSE;
567                 bNS = FALSE;
568                 
569                 /* Calculate long range corrections to pressure and energy */
570                 calc_dispcorr(fplog,inputrec,fr,step,top_global->natoms,rerun_fr.box,
571                               lambda,pres,vir,&prescorr,&enercorr,&dvdlcorr);
572                 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
573                 enerd->term[F_DISPCORR] = enercorr;
574                 enerd->term[F_EPOT] += enercorr;
575                 enerd->term[F_PRES] += prescorr;
576                 enerd->term[F_DVDL] += dvdlcorr;        
577                 
578                 /* If the compiler doesn't optimize this check away
579                  * we catch the NAN energies. With tables extreme negative
580                  * energies might occur close to r=0.
581                  */
582                 epot = enerd->term[F_EPOT];
583                 if (epot != epot || epot*beta < bU_neg_limit)
584                 {
585                     if (debug)
586                     {
587                         fprintf(debug,"\n  time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n",t,step,epot);
588                     }
589                     embU = 0;
590                 }
591                 else
592                 {
593                     embU = exp(-beta*epot);
594                     sum_embU += embU;
595                     /* Determine the weighted energy contributions of each energy group */
596                     e = 0;
597                     sum_UgembU[e++] += epot*embU;
598                     if (fr->bBHAM)
599                     {
600                         for(i=0; i<ngid; i++)
601                         {
602                             sum_UgembU[e++] +=
603                                 (enerd->grpp.ener[egBHAMSR][GID(i,gid_tp,ngid)] +
604                                  enerd->grpp.ener[egBHAMLR][GID(i,gid_tp,ngid)])*embU;
605                         }
606                     }
607                     else
608                     {
609                         for(i=0; i<ngid; i++)
610                         {
611                             sum_UgembU[e++] +=
612                                 (enerd->grpp.ener[egLJSR][GID(i,gid_tp,ngid)] +
613                                  enerd->grpp.ener[egLJLR][GID(i,gid_tp,ngid)])*embU;
614                         }
615                     }
616                     if (bDispCorr)
617                     {
618                         sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
619                     }
620                     if (bCharge)
621                     {
622                         for(i=0; i<ngid; i++)
623                         {
624                             sum_UgembU[e++] +=
625                                 (enerd->grpp.ener[egCOULSR][GID(i,gid_tp,ngid)] +
626                                  enerd->grpp.ener[egCOULLR][GID(i,gid_tp,ngid)])*embU;
627                         }
628                         if (bRFExcl)
629                         {
630                             sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
631                         }
632                         if (EEL_FULL(fr->eeltype))
633                         {
634                             sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
635                         }
636                     }
637                 }
638                 
639                 if (embU == 0 || beta*epot > bU_bin_limit)
640                 {
641                     bin[0]++;
642                 }
643                 else
644                 {
645                     i = (int)((bU_logV_bin_limit
646                                - (beta*epot - logV + refvolshift))*invbinw
647                               + 0.5);
648                     if (i < 0)
649                     {
650                         i = 0;
651                     }
652                     if (i >= nbin)
653                     {
654                         realloc_bins(&bin,&nbin,i+10);
655                     }
656                     bin[i]++;
657                 }
658                 
659                 if (debug)
660                 {
661                     fprintf(debug,"TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
662                             step,epot,x_tp[XX],x_tp[YY],x_tp[ZZ]);
663                 }
664
665                 if (dump_pdb && epot <= dump_ener)
666                 {
667                     sprintf(str,"t%g_step%d.pdb",t,step);
668                     sprintf(str2,"t: %f step %d ener: %f",t,step,epot);
669                     write_sto_conf_mtop(str,str2,top_global,state->x,state->v,
670                                         inputrec->ePBC,state->box);
671                 }
672             }
673         }
674         
675         if (PAR(cr))
676         {
677             /* When running in parallel sum the energies over the processes */
678             gmx_sumd(1,    &sum_embU, cr);
679             gmx_sumd(nener,sum_UgembU,cr);
680         }
681
682         frame++;
683         V_all += V;
684         VembU_all += V*sum_embU/nsteps;
685         
686         if (fp_tpi)
687         {
688             if (bVerbose || frame%10==0 || frame<10)
689             {
690                 fprintf(stderr,"mu %10.3e <mu> %10.3e\n",
691                         -log(sum_embU/nsteps)/beta,-log(VembU_all/V_all)/beta);
692             }
693             
694             fprintf(fp_tpi,"%10.3f %12.5e %12.5e %12.5e %12.5e",
695                     t,
696                     VembU_all==0 ? 20/beta : -log(VembU_all/V_all)/beta,
697                     sum_embU==0  ? 20/beta : -log(sum_embU/nsteps)/beta,
698                     sum_embU/nsteps,V);
699             for(e=0; e<nener; e++)
700             {
701                 fprintf(fp_tpi," %12.5e",sum_UgembU[e]/nsteps);
702             }
703             fprintf(fp_tpi,"\n");
704             fflush(fp_tpi);
705         }
706         
707         bNotLastFrame = read_next_frame(oenv, status,&rerun_fr);
708     } /* End of the loop  */
709     runtime_end(runtime);
710
711     close_trj(status);
712
713     if (fp_tpi != NULL)
714     {
715         gmx_fio_fclose(fp_tpi);
716     }
717
718     if (fplog != NULL)
719     {
720         fprintf(fplog,"\n");
721         fprintf(fplog,"  <V>  = %12.5e nm^3\n",V_all/frame);
722         fprintf(fplog,"  <mu> = %12.5e kJ/mol\n",-log(VembU_all/V_all)/beta);
723     }
724   
725     /* Write the Boltzmann factor histogram */
726     if (PAR(cr))
727     {
728         /* When running in parallel sum the bins over the processes */
729         i = nbin;
730         global_max(cr,&i);
731         realloc_bins(&bin,&nbin,i);
732         gmx_sumd(nbin,bin,cr);
733     }
734     if (MASTER(cr))
735     {
736         fp_tpi = xvgropen(opt2fn("-tpid",nfile,fnm),
737                           "TPI energy distribution",
738                           "\\betaU - log(V/<V>)","count",oenv);
739         sprintf(str,"number \\betaU > %g: %9.3e",bU_bin_limit,bin[0]);
740         xvgr_subtitle(fp_tpi,str,oenv);
741         xvgr_legend(fp_tpi,2,(const char **)tpid_leg,oenv);
742         for(i=nbin-1; i>0; i--)
743         {
744             bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
745             fprintf(fp_tpi,"%6.2f %10d %12.5e\n",
746                     bUlogV,
747                     (int)(bin[i]+0.5),
748                     bin[i]*exp(-bUlogV)*V_all/VembU_all);
749         }
750         gmx_fio_fclose(fp_tpi);
751     }
752     sfree(bin);
753
754     sfree(sum_UgembU);
755
756     runtime->nsteps_done = frame*inputrec->nsteps;
757
758     return 0;
759 }