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,2015,2016,2017,2018,2019, 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.
39 * \brief This file defines the integrator for test particle insertion
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_mdrun
55 #include "gromacs/commandline/filenm.h"
56 #include "gromacs/domdec/dlbtiming.h"
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/ewald/pme.h"
59 #include "gromacs/fileio/confio.h"
60 #include "gromacs/fileio/trxio.h"
61 #include "gromacs/fileio/xvgr.h"
62 #include "gromacs/gmxlib/chargegroup.h"
63 #include "gromacs/gmxlib/conformation_utilities.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/gmxlib/nrnb.h"
66 #include "gromacs/math/units.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/mdlib/constr.h"
69 #include "gromacs/mdlib/dispersioncorrection.h"
70 #include "gromacs/mdlib/energyoutput.h"
71 #include "gromacs/mdlib/force.h"
72 #include "gromacs/mdlib/force_flags.h"
73 #include "gromacs/mdlib/mdatoms.h"
74 #include "gromacs/mdlib/ns.h"
75 #include "gromacs/mdlib/tgroup.h"
76 #include "gromacs/mdlib/update.h"
77 #include "gromacs/mdlib/vsite.h"
78 #include "gromacs/mdrunutility/printtime.h"
79 #include "gromacs/mdtypes/commrec.h"
80 #include "gromacs/mdtypes/forcerec.h"
81 #include "gromacs/mdtypes/group.h"
82 #include "gromacs/mdtypes/inputrec.h"
83 #include "gromacs/mdtypes/md_enums.h"
84 #include "gromacs/mdtypes/mdrunoptions.h"
85 #include "gromacs/mdtypes/state.h"
86 #include "gromacs/pbcutil/pbc.h"
87 #include "gromacs/random/threefry.h"
88 #include "gromacs/random/uniformrealdistribution.h"
89 #include "gromacs/timing/wallcycle.h"
90 #include "gromacs/timing/walltime_accounting.h"
91 #include "gromacs/topology/mtop_util.h"
92 #include "gromacs/trajectory/trajectoryframe.h"
93 #include "gromacs/utility/cstringutil.h"
94 #include "gromacs/utility/fatalerror.h"
95 #include "gromacs/utility/gmxassert.h"
96 #include "gromacs/utility/logger.h"
97 #include "gromacs/utility/smalloc.h"
99 #include "integrator.h"
101 //! Global max algorithm
102 static void global_max(t_commrec *cr, int *n)
106 snew(sum, cr->nnodes);
107 sum[cr->nodeid] = *n;
108 gmx_sumi(cr->nnodes, sum, cr);
109 for (i = 0; i < cr->nnodes; i++)
111 *n = std::max(*n, sum[i]);
117 //! Reallocate arrays.
118 static void realloc_bins(double **bin, int *nbin, int nbin_new)
122 if (nbin_new != *nbin)
124 srenew(*bin, nbin_new);
125 for (i = *nbin; i < nbin_new; i++)
140 gmx_groups_t *groups;
141 gmx_enerdata_t *enerd;
142 PaddedVector<gmx::RVec> f {};
143 real lambda, t, temp, beta, drmax, epot;
144 double embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
147 gmx_bool bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
148 tensor force_vir, shake_vir, vir, pres;
149 int cg_tp, a_tp0, a_tp1, ngid, gid_tp, nener, e;
151 rvec mu_tot, x_init, dx, x_tp;
153 int64_t frame_step_prev, frame_step;
154 int64_t nsteps, stepblocksize = 0, step;
157 FILE *fp_tpi = nullptr;
158 char *ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
159 double dbl, dump_ener;
161 int nat_cavity = 0, d;
162 real *mass_cavity = nullptr, mass_tot;
164 double invbinw, *bin, refvolshift, logV, bUlogV;
165 real prescorr, enercorr, dvdlcorr;
166 gmx_bool bEnergyOutOfBounds;
167 const char *tpid_leg[2] = {"direct", "reweighted"};
168 auto mdatoms = mdAtoms->mdatoms();
170 GMX_UNUSED_VALUE(outputProvider);
172 GMX_LOG(mdlog.info).asParagraph().
173 appendText("Note that it is planned to change the command gmx mdrun -tpi "
174 "(and -tpic) to make the functionality available in a different "
175 "form in a future version of GROMACS, e.g. gmx test-particle-insertion.");
177 /* Since there is no upper limit to the insertion energies,
178 * we need to set an upper limit for the distribution output.
180 real bU_bin_limit = 50;
181 real bU_logV_bin_limit = bU_bin_limit + 10;
183 if (inputrec->cutoff_scheme == ecutsVERLET)
185 gmx_fatal(FARGS, "TPI does not work (yet) with the Verlet cut-off scheme");
190 gmx_mtop_generate_local_top(*top_global, &top, inputrec->efep != efepNO);
192 groups = &top_global->groups;
194 bCavity = (inputrec->eI == eiTPIC);
197 ptr = getenv("GMX_TPIC_MASSES");
204 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
205 * The center of mass of the last atoms is then used for TPIC.
208 while (sscanf(ptr, "%20lf%n", &dbl, &i) > 0)
210 srenew(mass_cavity, nat_cavity+1);
211 mass_cavity[nat_cavity] = dbl;
212 fprintf(fplog, "mass[%d] = %f\n",
213 nat_cavity+1, mass_cavity[nat_cavity]);
219 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
225 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
226 state_global->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
227 /* We never need full pbc for TPI */
229 /* Determine the temperature for the Boltzmann weighting */
230 temp = inputrec->opts.ref_t[0];
233 for (i = 1; (i < inputrec->opts.ngtc); i++)
235 if (inputrec->opts.ref_t[i] != temp)
237 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
238 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
242 "\n The temperature for test particle insertion is %.3f K\n\n",
245 beta = 1.0/(BOLTZ*temp);
247 /* Number of insertions per frame */
248 nsteps = inputrec->nsteps;
250 /* Use the same neighborlist with more insertions points
251 * in a sphere of radius drmax around the initial point
253 /* This should be a proper mdp parameter */
254 drmax = inputrec->rtpi;
256 /* An environment variable can be set to dump all configurations
257 * to pdb with an insertion energy <= this value.
259 dump_pdb = getenv("GMX_TPI_DUMP");
263 sscanf(dump_pdb, "%20lf", &dump_ener);
266 atoms2md(top_global, inputrec, -1, nullptr, top_global->natoms, mdAtoms);
267 update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
270 init_enerdata(groups->grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
271 f.resizeWithPadding(top_global->natoms);
273 /* Print to log file */
274 walltime_accounting_start_time(walltime_accounting);
275 wallcycle_start(wcycle, ewcRUN);
276 print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
278 /* The last charge group is the group to be inserted */
279 cg_tp = top.cgs.nr - 1;
280 a_tp0 = top.cgs.index[cg_tp];
281 a_tp1 = top.cgs.index[cg_tp+1];
284 fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
287 GMX_RELEASE_ASSERT(inputrec->rcoulomb <= inputrec->rlist && inputrec->rvdw <= inputrec->rlist, "Twin-range interactions are not supported with TPI");
289 snew(x_mol, a_tp1-a_tp0);
291 bDispCorr = (inputrec->eDispCorr != edispcNO);
293 auto x = makeArrayRef(state_global->x);
294 for (i = a_tp0; i < a_tp1; i++)
296 /* Copy the coordinates of the molecule to be insterted */
297 copy_rvec(x[i], x_mol[i-a_tp0]);
298 /* Check if we need to print electrostatic energies */
299 bCharge |= (mdatoms->chargeA[i] != 0 ||
300 ((mdatoms->chargeB != nullptr) && mdatoms->chargeB[i] != 0));
302 bRFExcl = (bCharge && EEL_RF(fr->ic->eeltype));
304 calc_cgcm(fplog, cg_tp, cg_tp+1, &(top.cgs), state_global->x.rvec_array(), fr->cg_cm);
307 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
309 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
310 fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
315 /* Center the molecule to be inserted at zero */
316 for (i = 0; i < a_tp1-a_tp0; i++)
318 rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
324 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
325 a_tp1-a_tp0, bCharge ? "with" : "without");
327 fprintf(fplog, "\nWill insert %" PRId64 " times in each frame of %s\n",
328 nsteps, opt2fn("-rerun", nfile, fnm));
333 if (inputrec->nstlist > 1)
335 if (drmax == 0 && a_tp1-a_tp0 == 1)
337 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);
341 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
349 fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
353 ngid = groups->grps[egcENER].nr;
354 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
367 if (EEL_FULL(fr->ic->eeltype))
372 snew(sum_UgembU, nener);
374 /* Copy the random seed set by the user */
375 seed = inputrec->ld_seed;
377 gmx::ThreeFry2x64<16> rng(seed, gmx::RandomDomain::TestParticleInsertion); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
378 gmx::UniformRealDistribution<real> dist;
382 fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
383 "TPI energies", "Time (ps)",
384 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
385 xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
388 sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
389 leg[e++] = gmx_strdup(str);
390 sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
391 leg[e++] = gmx_strdup(str);
392 sprintf(str, "f. <e\\S-\\betaU\\N>");
393 leg[e++] = gmx_strdup(str);
394 sprintf(str, "f. V");
395 leg[e++] = gmx_strdup(str);
396 sprintf(str, "f. <Ue\\S-\\betaU\\N>");
397 leg[e++] = gmx_strdup(str);
398 for (i = 0; i < ngid; i++)
400 sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
401 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
402 leg[e++] = gmx_strdup(str);
406 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
407 leg[e++] = gmx_strdup(str);
411 for (i = 0; i < ngid; i++)
413 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
414 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
415 leg[e++] = gmx_strdup(str);
419 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
420 leg[e++] = gmx_strdup(str);
422 if (EEL_FULL(fr->ic->eeltype))
424 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
425 leg[e++] = gmx_strdup(str);
428 xvgr_legend(fp_tpi, 4+nener, leg, oenv);
429 for (i = 0; i < 4+nener; i++)
443 /* Avoid frame step numbers <= -1 */
444 frame_step_prev = -1;
446 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
447 &rerun_fr, TRX_NEED_X);
450 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
451 mdatoms->nr - (a_tp1 - a_tp0))
453 gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
454 "is not equal the number in the run input file (%d) "
455 "minus the number of atoms to insert (%d)\n",
456 rerun_fr.natoms, bCavity ? " minus one" : "",
457 mdatoms->nr, a_tp1-a_tp0);
460 refvolshift = log(det(rerun_fr.box));
462 switch (inputrec->eI)
465 stepblocksize = inputrec->nstlist;
471 gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
474 while (bNotLastFrame)
476 frame_step = rerun_fr.step;
477 if (frame_step <= frame_step_prev)
479 /* We don't have step number in the trajectory file,
480 * or we have constant or decreasing step numbers.
481 * Ensure we have increasing step numbers, since we use
482 * the step numbers as a counter for random numbers.
484 frame_step = frame_step_prev + 1;
486 frame_step_prev = frame_step;
488 lambda = rerun_fr.lambda;
492 for (e = 0; e < nener; e++)
497 /* Copy the coordinates from the input trajectory */
498 auto x = makeArrayRef(state_global->x);
499 for (i = 0; i < rerun_fr.natoms; i++)
501 copy_rvec(rerun_fr.x[i], x[i]);
503 copy_mat(rerun_fr.box, state_global->box);
505 V = det(state_global->box);
508 bStateChanged = TRUE;
511 step = cr->nodeid*stepblocksize;
512 while (step < nsteps)
514 /* Restart random engine using the frame and insertion step
516 * Note that we need to draw several random values per iteration,
517 * but by using the internal subcounter functionality of ThreeFry2x64
518 * we can draw 131072 unique 64-bit values before exhausting
519 * the stream. This is a huge margin, and if something still goes
520 * wrong you will get an exception when the stream is exhausted.
522 rng.restart(frame_step, step);
523 dist.reset(); // erase any memory in the distribution
527 /* Random insertion in the whole volume */
528 bNS = (step % inputrec->nstlist == 0);
531 /* Generate a random position in the box */
532 for (d = 0; d < DIM; d++)
534 x_init[d] = dist(rng)*state_global->box[d][d];
538 if (inputrec->nstlist == 1)
540 copy_rvec(x_init, x_tp);
544 /* Generate coordinates within |dx|=drmax of x_init */
547 for (d = 0; d < DIM; d++)
549 dx[d] = (2*dist(rng) - 1)*drmax;
552 while (norm2(dx) > drmax*drmax);
553 rvec_add(x_init, dx, x_tp);
558 /* Random insertion around a cavity location
559 * given by the last coordinate of the trajectory.
565 /* Copy the location of the cavity */
566 copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
570 /* Determine the center of mass of the last molecule */
573 for (i = 0; i < nat_cavity; i++)
575 for (d = 0; d < DIM; d++)
578 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
580 mass_tot += mass_cavity[i];
582 for (d = 0; d < DIM; d++)
584 x_init[d] /= mass_tot;
588 /* Generate coordinates within |dx|=drmax of x_init */
591 for (d = 0; d < DIM; d++)
593 dx[d] = (2*dist(rng) - 1)*drmax;
596 while (norm2(dx) > drmax*drmax);
597 rvec_add(x_init, dx, x_tp);
600 if (a_tp1 - a_tp0 == 1)
602 /* Insert a single atom, just copy the insertion location */
603 copy_rvec(x_tp, x[a_tp0]);
607 /* Copy the coordinates from the top file */
608 for (i = a_tp0; i < a_tp1; i++)
610 copy_rvec(x_mol[i-a_tp0], x[i]);
612 /* Rotate the molecule randomly */
613 real angleX = 2*M_PI*dist(rng);
614 real angleY = 2*M_PI*dist(rng);
615 real angleZ = 2*M_PI*dist(rng);
616 rotate_conf(a_tp1-a_tp0, state_global->x.rvec_array()+a_tp0, nullptr,
617 angleX, angleY, angleZ);
618 /* Shift to the insertion location */
619 for (i = a_tp0; i < a_tp1; i++)
621 rvec_inc(x[i], x_tp);
625 /* Clear some matrix variables */
626 clear_mat(force_vir);
627 clear_mat(shake_vir);
631 /* Set the charge group center of mass of the test particle */
632 copy_rvec(x_init, fr->cg_cm[top.cgs.nr-1]);
634 /* Calc energy (no forces) on new positions.
635 * Since we only need the intermolecular energy
636 * and the RF exclusion terms of the inserted molecule occur
637 * within a single charge group we can pass NULL for the graph.
638 * This also avoids shifts that would move charge groups
640 /* Make do_force do a single node force calculation */
643 // TPI might place a particle so close that the potential
644 // is infinite. Since this is intended to happen, we
645 // temporarily suppress any exceptions that the processor
646 // might raise, then restore the old behaviour.
647 std::fenv_t floatingPointEnvironment;
648 std::feholdexcept(&floatingPointEnvironment);
649 do_force(fplog, cr, ms, inputrec, nullptr, nullptr,
650 step, nrnb, wcycle, &top, &top_global->groups,
651 state_global->box, state_global->x.arrayRefWithPadding(), &state_global->hist,
652 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
653 state_global->lambda,
654 nullptr, fr, ppForceWorkload, nullptr, mu_tot, t, nullptr,
655 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
656 (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS : 0) |
657 (bStateChanged ? GMX_FORCE_STATECHANGED : 0),
658 DDBalanceRegionHandler(nullptr));
659 std::feclearexcept(FE_DIVBYZERO | FE_INVALID | FE_OVERFLOW);
660 std::feupdateenv(&floatingPointEnvironment);
663 bStateChanged = FALSE;
666 /* Calculate long range corrections to pressure and energy */
667 calc_dispcorr(inputrec, fr, state_global->box,
668 lambda, pres, vir, &prescorr, &enercorr, &dvdlcorr);
669 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
670 enerd->term[F_DISPCORR] = enercorr;
671 enerd->term[F_EPOT] += enercorr;
672 enerd->term[F_PRES] += prescorr;
673 enerd->term[F_DVDL_VDW] += dvdlcorr;
675 epot = enerd->term[F_EPOT];
676 bEnergyOutOfBounds = FALSE;
678 /* If the compiler doesn't optimize this check away
679 * we catch the NAN energies.
680 * The epot>GMX_REAL_MAX check catches inf values,
681 * which should nicely result in embU=0 through the exp below,
682 * but it does not hurt to check anyhow.
684 /* Non-bonded Interaction usually diverge at r=0.
685 * With tabulated interaction functions the first few entries
686 * should be capped in a consistent fashion between
687 * repulsion, dispersion and Coulomb to avoid accidental
688 * negative values in the total energy.
689 * The table generation code in tables.c does this.
690 * With user tbales the user should take care of this.
692 if (epot != epot || epot > GMX_REAL_MAX)
694 bEnergyOutOfBounds = TRUE;
696 if (bEnergyOutOfBounds)
700 fprintf(debug, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, static_cast<int>(step), epot);
706 // Exponent argument is fine in SP range, but output can be in DP range
707 embU = exp(static_cast<double>(-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)]*embU;
722 for (i = 0; i < ngid; i++)
725 enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)]*embU;
730 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
734 for (i = 0; i < ngid; i++)
736 sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
740 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
742 if (EEL_FULL(fr->ic->eeltype))
744 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
749 if (embU == 0 || beta*epot > bU_bin_limit)
755 i = gmx::roundToInt((bU_logV_bin_limit
756 - (beta*epot - logV + refvolshift))*invbinw);
763 realloc_bins(&bin, &nbin, i+10);
770 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
771 static_cast<int>(step), epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
774 if (dump_pdb && epot <= dump_ener)
776 sprintf(str, "t%g_step%d.pdb", t, static_cast<int>(step));
777 sprintf(str2, "t: %f step %d ener: %f", t, static_cast<int>(step), epot);
778 write_sto_conf_mtop(str, str2, top_global, state_global->x.rvec_array(), state_global->v.rvec_array(),
779 inputrec->ePBC, state_global->box);
783 if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
785 /* Skip all steps assigned to the other MPI ranks */
786 step += (cr->nnodes - 1)*stepblocksize;
792 /* When running in parallel sum the energies over the processes */
793 gmx_sumd(1, &sum_embU, cr);
794 gmx_sumd(nener, sum_UgembU, cr);
799 VembU_all += V*sum_embU/nsteps;
803 if (mdrunOptions.verbose || frame%10 == 0 || frame < 10)
805 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
806 -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
809 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
811 VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
812 sum_embU == 0 ? 20/beta : -log(sum_embU/nsteps)/beta,
814 for (e = 0; e < nener; e++)
816 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
818 fprintf(fp_tpi, "\n");
822 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
823 } /* End of the loop */
824 walltime_accounting_end_time(walltime_accounting);
828 if (fp_tpi != nullptr)
833 if (fplog != nullptr)
835 fprintf(fplog, "\n");
836 fprintf(fplog, " <V> = %12.5e nm^3\n", V_all/frame);
837 const double mu = -log(VembU_all/V_all)/beta;
838 fprintf(fplog, " <mu> = %12.5e kJ/mol\n", mu);
840 if (!std::isfinite(mu))
842 fprintf(fplog, "\nThe computed chemical potential is not finite - consider increasing the number of steps and/or the number of frames to insert into.\n");
846 /* Write the Boltzmann factor histogram */
849 /* When running in parallel sum the bins over the processes */
852 realloc_bins(&bin, &nbin, i);
853 gmx_sumd(nbin, bin, cr);
857 fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
858 "TPI energy distribution",
859 "\\betaU - log(V/<V>)", "count", oenv);
860 sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
861 xvgr_subtitle(fp_tpi, str, oenv);
862 xvgr_legend(fp_tpi, 2, tpid_leg, oenv);
863 for (i = nbin-1; i > 0; i--)
865 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
866 fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
869 bin[i]*exp(-bUlogV)*V_all/VembU_all);
877 walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);