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