2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
47 #include "gromacs/utility/smalloc.h"
49 #include "chargegroup.h"
53 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/math/vec.h"
65 #include "gromacs/random/random.h"
67 #include "gromacs/fileio/xvgr.h"
70 #include "mtop_util.h"
72 #include "gromacs/gmxlib/conformation-utilities.h"
74 #include "gromacs/fileio/confio.h"
75 #include "gromacs/fileio/gmxfio.h"
76 #include "gromacs/fileio/trxio.h"
77 #include "gromacs/timing/wallcycle.h"
78 #include "gromacs/timing/walltime_accounting.h"
80 static void global_max(t_commrec *cr, int *n)
84 snew(sum, cr->nnodes);
86 gmx_sumi(cr->nnodes, sum, cr);
87 for (i = 0; i < cr->nnodes; i++)
95 static void realloc_bins(double **bin, int *nbin, int nbin_new)
99 if (nbin_new != *nbin)
101 srenew(*bin, nbin_new);
102 for (i = *nbin; i < nbin_new; i++)
110 double do_tpi(FILE *fplog, t_commrec *cr,
111 int nfile, const t_filenm fnm[],
112 const output_env_t oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
113 int gmx_unused nstglobalcomm,
114 gmx_vsite_t gmx_unused *vsite, gmx_constr_t gmx_unused constr,
115 int gmx_unused stepout,
116 t_inputrec *inputrec,
117 gmx_mtop_t *top_global, t_fcdata *fcd,
120 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
121 gmx_edsam_t gmx_unused ed,
123 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
124 gmx_membed_t gmx_unused membed,
125 real gmx_unused cpt_period, real gmx_unused max_hours,
126 const char gmx_unused *deviceOptions,
127 int gmx_unused imdport,
128 unsigned long gmx_unused Flags,
129 gmx_walltime_accounting_t walltime_accounting)
131 const char *TPI = "Test Particle Insertion";
133 gmx_groups_t *groups;
134 gmx_enerdata_t *enerd;
136 real lambda, t, temp, beta, drmax, epot;
137 double embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
140 gmx_bool bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
141 tensor force_vir, shake_vir, vir, pres;
142 int cg_tp, a_tp0, a_tp1, ngid, gid_tp, nener, e;
144 rvec mu_tot, x_init, dx, x_tp;
146 gmx_int64_t frame_step_prev, frame_step;
147 gmx_int64_t nsteps, stepblocksize = 0, step;
148 gmx_int64_t rnd_count_stride, rnd_count;
153 char *ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
154 double dbl, dump_ener;
156 int nat_cavity = 0, d;
157 real *mass_cavity = NULL, mass_tot;
159 double invbinw, *bin, refvolshift, logV, bUlogV;
160 real dvdl, prescorr, enercorr, dvdlcorr;
161 gmx_bool bEnergyOutOfBounds;
162 const char *tpid_leg[2] = {"direct", "reweighted"};
164 /* Since there is no upper limit to the insertion energies,
165 * we need to set an upper limit for the distribution output.
167 real bU_bin_limit = 50;
168 real bU_logV_bin_limit = bU_bin_limit + 10;
172 top = gmx_mtop_generate_local_top(top_global, inputrec);
174 groups = &top_global->groups;
176 bCavity = (inputrec->eI == eiTPIC);
179 ptr = getenv("GMX_TPIC_MASSES");
186 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
187 * The center of mass of the last atoms is then used for TPIC.
190 while (sscanf(ptr, "%lf%n", &dbl, &i) > 0)
192 srenew(mass_cavity, nat_cavity+1);
193 mass_cavity[nat_cavity] = dbl;
194 fprintf(fplog, "mass[%d] = %f\n",
195 nat_cavity+1, mass_cavity[nat_cavity]);
201 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
207 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
208 state->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
209 /* We never need full pbc for TPI */
211 /* Determine the temperature for the Boltzmann weighting */
212 temp = inputrec->opts.ref_t[0];
215 for (i = 1; (i < inputrec->opts.ngtc); i++)
217 if (inputrec->opts.ref_t[i] != temp)
219 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
220 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
224 "\n The temperature for test particle insertion is %.3f K\n\n",
227 beta = 1.0/(BOLTZ*temp);
229 /* Number of insertions per frame */
230 nsteps = inputrec->nsteps;
232 /* Use the same neighborlist with more insertions points
233 * in a sphere of radius drmax around the initial point
235 /* This should be a proper mdp parameter */
236 drmax = inputrec->rtpi;
238 /* An environment variable can be set to dump all configurations
239 * to pdb with an insertion energy <= this value.
241 dump_pdb = getenv("GMX_TPI_DUMP");
245 sscanf(dump_pdb, "%lf", &dump_ener);
248 atoms2md(top_global, inputrec, 0, NULL, top_global->natoms, mdatoms);
249 update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
252 init_enerdata(groups->grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
253 snew(f, top_global->natoms);
255 /* Print to log file */
256 walltime_accounting_start(walltime_accounting);
257 wallcycle_start(wcycle, ewcRUN);
258 print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
260 /* The last charge group is the group to be inserted */
261 cg_tp = top->cgs.nr - 1;
262 a_tp0 = top->cgs.index[cg_tp];
263 a_tp1 = top->cgs.index[cg_tp+1];
266 fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
268 if (a_tp1 - a_tp0 > 1 &&
269 (inputrec->rlist < inputrec->rcoulomb ||
270 inputrec->rlist < inputrec->rvdw))
272 gmx_fatal(FARGS, "Can not do TPI for multi-atom molecule with a twin-range cut-off");
274 snew(x_mol, a_tp1-a_tp0);
276 bDispCorr = (inputrec->eDispCorr != edispcNO);
278 for (i = a_tp0; i < a_tp1; i++)
280 /* Copy the coordinates of the molecule to be insterted */
281 copy_rvec(state->x[i], x_mol[i-a_tp0]);
282 /* Check if we need to print electrostatic energies */
283 bCharge |= (mdatoms->chargeA[i] != 0 ||
284 (mdatoms->chargeB && mdatoms->chargeB[i] != 0));
286 bRFExcl = (bCharge && EEL_RF(fr->eeltype) && fr->eeltype != eelRF_NEC);
288 calc_cgcm(fplog, cg_tp, cg_tp+1, &(top->cgs), state->x, fr->cg_cm);
291 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
293 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
294 fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
299 /* Center the molecule to be inserted at zero */
300 for (i = 0; i < a_tp1-a_tp0; i++)
302 rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
308 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
309 a_tp1-a_tp0, bCharge ? "with" : "without");
311 fprintf(fplog, "\nWill insert %d times in each frame of %s\n",
312 (int)nsteps, opt2fn("-rerun", nfile, fnm));
317 if (inputrec->nstlist > 1)
319 if (drmax == 0 && a_tp1-a_tp0 == 1)
321 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);
325 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
333 fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
337 ngid = groups->grps[egcENER].nr;
338 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
351 if (EEL_FULL(fr->eeltype))
356 snew(sum_UgembU, nener);
358 /* Copy the random seed set by the user */
359 seed = inputrec->ld_seed;
360 /* We use the frame step number as one random counter.
361 * The second counter use the insertion (step) count. But we
362 * need multiple random numbers per insertion. This number is
363 * not fixed, since we generate random locations in a sphere
364 * by putting locations in a cube and some of these fail.
365 * A count of 20 is already extremely unlikely, so 10000 is
366 * a safe margin for random numbers per insertion.
368 rnd_count_stride = 10000;
372 fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
373 "TPI energies", "Time (ps)",
374 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
375 xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
378 sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
379 leg[e++] = strdup(str);
380 sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
381 leg[e++] = strdup(str);
382 sprintf(str, "f. <e\\S-\\betaU\\N>");
383 leg[e++] = strdup(str);
384 sprintf(str, "f. V");
385 leg[e++] = strdup(str);
386 sprintf(str, "f. <Ue\\S-\\betaU\\N>");
387 leg[e++] = strdup(str);
388 for (i = 0; i < ngid; i++)
390 sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
391 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
392 leg[e++] = strdup(str);
396 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
397 leg[e++] = strdup(str);
401 for (i = 0; i < ngid; i++)
403 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
404 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
405 leg[e++] = strdup(str);
409 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
410 leg[e++] = strdup(str);
412 if (EEL_FULL(fr->eeltype))
414 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
415 leg[e++] = strdup(str);
418 xvgr_legend(fp_tpi, 4+nener, (const char**)leg, oenv);
419 for (i = 0; i < 4+nener; i++)
433 /* Avoid frame step numbers <= -1 */
434 frame_step_prev = -1;
436 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
437 &rerun_fr, TRX_NEED_X);
440 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
441 mdatoms->nr - (a_tp1 - a_tp0))
443 gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
444 "is not equal the number in the run input file (%d) "
445 "minus the number of atoms to insert (%d)\n",
446 rerun_fr.natoms, bCavity ? " minus one" : "",
447 mdatoms->nr, a_tp1-a_tp0);
450 refvolshift = log(det(rerun_fr.box));
452 switch (inputrec->eI)
455 stepblocksize = inputrec->nstlist;
461 gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
465 /* Make sure we don't detect SIMD overflow generated before this point */
466 gmx_simd_check_and_reset_overflow();
469 while (bNotLastFrame)
471 frame_step = rerun_fr.step;
472 if (frame_step <= frame_step_prev)
474 /* We don't have step number in the trajectory file,
475 * or we have constant or decreasing step numbers.
476 * Ensure we have increasing step numbers, since we use
477 * the step numbers as a counter for random numbers.
479 frame_step = frame_step_prev + 1;
481 frame_step_prev = frame_step;
483 lambda = rerun_fr.lambda;
487 for (e = 0; e < nener; e++)
492 /* Copy the coordinates from the input trajectory */
493 for (i = 0; i < rerun_fr.natoms; i++)
495 copy_rvec(rerun_fr.x[i], state->x[i]);
497 copy_mat(rerun_fr.box, state->box);
502 bStateChanged = TRUE;
505 step = cr->nodeid*stepblocksize;
506 while (step < nsteps)
508 /* Initialize the second counter for random numbers using
509 * the insertion step index. This ensures that we get
510 * the same random numbers independently of how many
511 * MPI ranks we use. Also for the same seed, we get
512 * the same initial random sequence for different nsteps.
514 rnd_count = step*rnd_count_stride;
518 /* Random insertion in the whole volume */
519 bNS = (step % inputrec->nstlist == 0);
522 /* Generate a random position in the box */
523 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd);
524 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd+2);
525 for (d = 0; d < DIM; d++)
527 x_init[d] = rnd[d]*state->box[d][d];
530 if (inputrec->nstlist == 1)
532 copy_rvec(x_init, x_tp);
536 /* Generate coordinates within |dx|=drmax of x_init */
539 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd);
540 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd+2);
541 for (d = 0; d < DIM; d++)
543 dx[d] = (2*rnd[d] - 1)*drmax;
546 while (norm2(dx) > drmax*drmax);
547 rvec_add(x_init, dx, x_tp);
552 /* Random insertion around a cavity location
553 * given by the last coordinate of the trajectory.
559 /* Copy the location of the cavity */
560 copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
564 /* Determine the center of mass of the last molecule */
567 for (i = 0; i < nat_cavity; i++)
569 for (d = 0; d < DIM; d++)
572 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
574 mass_tot += mass_cavity[i];
576 for (d = 0; d < DIM; d++)
578 x_init[d] /= mass_tot;
582 /* Generate coordinates within |dx|=drmax of x_init */
585 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd);
586 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd+2);
587 for (d = 0; d < DIM; d++)
589 dx[d] = (2*rnd[d] - 1)*drmax;
592 while (norm2(dx) > drmax*drmax);
593 rvec_add(x_init, dx, x_tp);
596 if (a_tp1 - a_tp0 == 1)
598 /* Insert a single atom, just copy the insertion location */
599 copy_rvec(x_tp, state->x[a_tp0]);
603 /* Copy the coordinates from the top file */
604 for (i = a_tp0; i < a_tp1; i++)
606 copy_rvec(x_mol[i-a_tp0], state->x[i]);
608 /* Rotate the molecule randomly */
609 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd);
610 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd+2);
611 rotate_conf(a_tp1-a_tp0, state->x+a_tp0, NULL,
615 /* Shift to the insertion location */
616 for (i = a_tp0; i < a_tp1; i++)
618 rvec_inc(state->x[i], x_tp);
622 /* Clear some matrix variables */
623 clear_mat(force_vir);
624 clear_mat(shake_vir);
628 /* Set the charge group center of mass of the test particle */
629 copy_rvec(x_init, fr->cg_cm[top->cgs.nr-1]);
631 /* Calc energy (no forces) on new positions.
632 * Since we only need the intermolecular energy
633 * and the RF exclusion terms of the inserted molecule occur
634 * within a single charge group we can pass NULL for the graph.
635 * This also avoids shifts that would move charge groups
638 * Some checks above ensure than we can not have
639 * twin-range interactions together with nstlist > 1,
640 * therefore we do not need to remember the LR energies.
642 /* Make do_force do a single node force calculation */
644 do_force(fplog, cr, inputrec,
645 step, nrnb, wcycle, top, &top_global->groups,
646 state->box, state->x, &state->hist,
647 f, force_vir, mdatoms, enerd, fcd,
649 NULL, fr, NULL, mu_tot, t, NULL, NULL, FALSE,
650 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
651 (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS | GMX_FORCE_DO_LR : 0) |
652 (bStateChanged ? GMX_FORCE_STATECHANGED : 0));
654 bStateChanged = FALSE;
657 /* Calculate long range corrections to pressure and energy */
658 calc_dispcorr(fplog, inputrec, fr, step, top_global->natoms, state->box,
659 lambda, pres, vir, &prescorr, &enercorr, &dvdlcorr);
660 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
661 enerd->term[F_DISPCORR] = enercorr;
662 enerd->term[F_EPOT] += enercorr;
663 enerd->term[F_PRES] += prescorr;
664 enerd->term[F_DVDL_VDW] += dvdlcorr;
666 epot = enerd->term[F_EPOT];
667 bEnergyOutOfBounds = FALSE;
668 #ifdef GMX_SIMD_X86_SSE2_OR_HIGHER
669 /* With SSE the energy can overflow, check for this */
670 if (gmx_mm_check_and_reset_overflow())
674 fprintf(debug, "Found an SSE overflow, assuming the energy is out of bounds\n");
676 bEnergyOutOfBounds = TRUE;
679 /* If the compiler doesn't optimize this check away
680 * we catch the NAN energies.
681 * The epot>GMX_REAL_MAX check catches inf values,
682 * which should nicely result in embU=0 through the exp below,
683 * but it does not hurt to check anyhow.
685 /* Non-bonded Interaction usually diverge at r=0.
686 * With tabulated interaction functions the first few entries
687 * should be capped in a consistent fashion between
688 * repulsion, dispersion and Coulomb to avoid accidental
689 * negative values in the total energy.
690 * The table generation code in tables.c does this.
691 * With user tbales the user should take care of this.
693 if (epot != epot || epot > GMX_REAL_MAX)
695 bEnergyOutOfBounds = TRUE;
697 if (bEnergyOutOfBounds)
701 fprintf(debug, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, (int)step, epot);
707 embU = exp(-beta*epot);
709 /* Determine the weighted energy contributions of each energy group */
711 sum_UgembU[e++] += epot*embU;
714 for (i = 0; i < ngid; i++)
717 (enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)] +
718 enerd->grpp.ener[egBHAMLR][GID(i, gid_tp, ngid)])*embU;
723 for (i = 0; i < ngid; i++)
726 (enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)] +
727 enerd->grpp.ener[egLJLR][GID(i, gid_tp, ngid)])*embU;
732 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
736 for (i = 0; i < ngid; i++)
739 (enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] +
740 enerd->grpp.ener[egCOULLR][GID(i, gid_tp, ngid)])*embU;
744 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
746 if (EEL_FULL(fr->eeltype))
748 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
753 if (embU == 0 || beta*epot > bU_bin_limit)
759 i = (int)((bU_logV_bin_limit
760 - (beta*epot - logV + refvolshift))*invbinw
768 realloc_bins(&bin, &nbin, i+10);
775 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
776 (int)step, epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
779 if (dump_pdb && epot <= dump_ener)
781 sprintf(str, "t%g_step%d.pdb", t, (int)step);
782 sprintf(str2, "t: %f step %d ener: %f", t, (int)step, epot);
783 write_sto_conf_mtop(str, str2, top_global, state->x, state->v,
784 inputrec->ePBC, state->box);
788 if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
790 /* Skip all steps assigned to the other MPI ranks */
791 step += (cr->nnodes - 1)*stepblocksize;
797 /* When running in parallel sum the energies over the processes */
798 gmx_sumd(1, &sum_embU, cr);
799 gmx_sumd(nener, sum_UgembU, cr);
804 VembU_all += V*sum_embU/nsteps;
808 if (bVerbose || frame%10 == 0 || frame < 10)
810 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
811 -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
814 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
816 VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
817 sum_embU == 0 ? 20/beta : -log(sum_embU/nsteps)/beta,
819 for (e = 0; e < nener; e++)
821 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
823 fprintf(fp_tpi, "\n");
827 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
828 } /* End of the loop */
829 walltime_accounting_end(walltime_accounting);
835 gmx_fio_fclose(fp_tpi);
840 fprintf(fplog, "\n");
841 fprintf(fplog, " <V> = %12.5e nm^3\n", V_all/frame);
842 fprintf(fplog, " <mu> = %12.5e kJ/mol\n", -log(VembU_all/V_all)/beta);
845 /* Write the Boltzmann factor histogram */
848 /* When running in parallel sum the bins over the processes */
851 realloc_bins(&bin, &nbin, i);
852 gmx_sumd(nbin, bin, cr);
856 fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
857 "TPI energy distribution",
858 "\\betaU - log(V/<V>)", "count", oenv);
859 sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
860 xvgr_subtitle(fp_tpi, str, oenv);
861 xvgr_legend(fp_tpi, 2, (const char **)tpid_leg, oenv);
862 for (i = nbin-1; i > 0; i--)
864 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
865 fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
868 bin[i]*exp(-bUlogV)*V_all/VembU_all);
870 gmx_fio_fclose(fp_tpi);
876 walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);