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