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