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