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