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/mdrun.h"
75 #include "gromacs/mdlib/ns.h"
76 #include "gromacs/mdlib/sim_util.h"
77 #include "gromacs/mdlib/tgroup.h"
78 #include "gromacs/mdlib/update.h"
79 #include "gromacs/mdlib/vsite.h"
80 #include "gromacs/mdrunutility/printtime.h"
81 #include "gromacs/mdtypes/commrec.h"
82 #include "gromacs/mdtypes/forcerec.h"
83 #include "gromacs/mdtypes/group.h"
84 #include "gromacs/mdtypes/inputrec.h"
85 #include "gromacs/mdtypes/md_enums.h"
86 #include "gromacs/mdtypes/state.h"
87 #include "gromacs/pbcutil/pbc.h"
88 #include "gromacs/random/threefry.h"
89 #include "gromacs/random/uniformrealdistribution.h"
90 #include "gromacs/timing/wallcycle.h"
91 #include "gromacs/timing/walltime_accounting.h"
92 #include "gromacs/topology/mtop_util.h"
93 #include "gromacs/trajectory/trajectoryframe.h"
94 #include "gromacs/utility/cstringutil.h"
95 #include "gromacs/utility/fatalerror.h"
96 #include "gromacs/utility/gmxassert.h"
97 #include "gromacs/utility/logger.h"
98 #include "gromacs/utility/smalloc.h"
100 #include "integrator.h"
102 //! Global max algorithm
103 static void global_max(t_commrec *cr, int *n)
107 snew(sum, cr->nnodes);
108 sum[cr->nodeid] = *n;
109 gmx_sumi(cr->nnodes, sum, cr);
110 for (i = 0; i < cr->nnodes; i++)
112 *n = std::max(*n, sum[i]);
118 //! Reallocate arrays.
119 static void realloc_bins(double **bin, int *nbin, int nbin_new)
123 if (nbin_new != *nbin)
125 srenew(*bin, nbin_new);
126 for (i = *nbin; i < nbin_new; i++)
141 gmx_groups_t *groups;
142 gmx_enerdata_t *enerd;
143 PaddedVector<gmx::RVec> f {};
144 real lambda, t, temp, beta, drmax, epot;
145 double embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
148 gmx_bool bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
149 tensor force_vir, shake_vir, vir, pres;
150 int cg_tp, a_tp0, a_tp1, ngid, gid_tp, nener, e;
152 rvec mu_tot, x_init, dx, x_tp;
154 int64_t frame_step_prev, frame_step;
155 int64_t nsteps, stepblocksize = 0, step;
158 FILE *fp_tpi = nullptr;
159 char *ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
160 double dbl, dump_ener;
162 int nat_cavity = 0, d;
163 real *mass_cavity = nullptr, mass_tot;
165 double invbinw, *bin, refvolshift, logV, bUlogV;
166 real prescorr, enercorr, dvdlcorr;
167 gmx_bool bEnergyOutOfBounds;
168 const char *tpid_leg[2] = {"direct", "reweighted"};
169 auto mdatoms = mdAtoms->mdatoms();
171 GMX_UNUSED_VALUE(outputProvider);
173 GMX_LOG(mdlog.info).asParagraph().
174 appendText("Note that it is planned to change the command gmx mdrun -tpi "
175 "(and -tpic) to make the functionality available in a different "
176 "form in a future version of GROMACS, e.g. gmx test-particle-insertion.");
178 /* Since there is no upper limit to the insertion energies,
179 * we need to set an upper limit for the distribution output.
181 real bU_bin_limit = 50;
182 real bU_logV_bin_limit = bU_bin_limit + 10;
184 if (inputrec->cutoff_scheme == ecutsVERLET)
186 gmx_fatal(FARGS, "TPI does not work (yet) with the Verlet cut-off scheme");
191 gmx_mtop_generate_local_top(*top_global, &top, inputrec->efep != efepNO);
193 groups = &top_global->groups;
195 bCavity = (inputrec->eI == eiTPIC);
198 ptr = getenv("GMX_TPIC_MASSES");
205 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
206 * The center of mass of the last atoms is then used for TPIC.
209 while (sscanf(ptr, "%20lf%n", &dbl, &i) > 0)
211 srenew(mass_cavity, nat_cavity+1);
212 mass_cavity[nat_cavity] = dbl;
213 fprintf(fplog, "mass[%d] = %f\n",
214 nat_cavity+1, mass_cavity[nat_cavity]);
220 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
226 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
227 state_global->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
228 /* We never need full pbc for TPI */
230 /* Determine the temperature for the Boltzmann weighting */
231 temp = inputrec->opts.ref_t[0];
234 for (i = 1; (i < inputrec->opts.ngtc); i++)
236 if (inputrec->opts.ref_t[i] != temp)
238 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
239 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
243 "\n The temperature for test particle insertion is %.3f K\n\n",
246 beta = 1.0/(BOLTZ*temp);
248 /* Number of insertions per frame */
249 nsteps = inputrec->nsteps;
251 /* Use the same neighborlist with more insertions points
252 * in a sphere of radius drmax around the initial point
254 /* This should be a proper mdp parameter */
255 drmax = inputrec->rtpi;
257 /* An environment variable can be set to dump all configurations
258 * to pdb with an insertion energy <= this value.
260 dump_pdb = getenv("GMX_TPI_DUMP");
264 sscanf(dump_pdb, "%20lf", &dump_ener);
267 atoms2md(top_global, inputrec, -1, nullptr, top_global->natoms, mdAtoms);
268 update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
271 init_enerdata(groups->grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
272 f.resizeWithPadding(top_global->natoms);
274 /* Print to log file */
275 walltime_accounting_start_time(walltime_accounting);
276 wallcycle_start(wcycle, ewcRUN);
277 print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
279 /* The last charge group is the group to be inserted */
280 cg_tp = top.cgs.nr - 1;
281 a_tp0 = top.cgs.index[cg_tp];
282 a_tp1 = top.cgs.index[cg_tp+1];
285 fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
288 GMX_RELEASE_ASSERT(inputrec->rcoulomb <= inputrec->rlist && inputrec->rvdw <= inputrec->rlist, "Twin-range interactions are not supported with TPI");
290 snew(x_mol, a_tp1-a_tp0);
292 bDispCorr = (inputrec->eDispCorr != edispcNO);
294 auto x = makeArrayRef(state_global->x);
295 for (i = a_tp0; i < a_tp1; i++)
297 /* Copy the coordinates of the molecule to be insterted */
298 copy_rvec(x[i], x_mol[i-a_tp0]);
299 /* Check if we need to print electrostatic energies */
300 bCharge |= (mdatoms->chargeA[i] != 0 ||
301 ((mdatoms->chargeB != nullptr) && mdatoms->chargeB[i] != 0));
303 bRFExcl = (bCharge && EEL_RF(fr->ic->eeltype));
305 calc_cgcm(fplog, cg_tp, cg_tp+1, &(top.cgs), state_global->x.rvec_array(), fr->cg_cm);
308 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
310 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
311 fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
316 /* Center the molecule to be inserted at zero */
317 for (i = 0; i < a_tp1-a_tp0; i++)
319 rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
325 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
326 a_tp1-a_tp0, bCharge ? "with" : "without");
328 fprintf(fplog, "\nWill insert %" PRId64 " times in each frame of %s\n",
329 nsteps, opt2fn("-rerun", nfile, fnm));
334 if (inputrec->nstlist > 1)
336 if (drmax == 0 && a_tp1-a_tp0 == 1)
338 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);
342 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
350 fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
354 ngid = groups->grps[egcENER].nr;
355 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
368 if (EEL_FULL(fr->ic->eeltype))
373 snew(sum_UgembU, nener);
375 /* Copy the random seed set by the user */
376 seed = inputrec->ld_seed;
378 gmx::ThreeFry2x64<16> rng(seed, gmx::RandomDomain::TestParticleInsertion); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
379 gmx::UniformRealDistribution<real> dist;
383 fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
384 "TPI energies", "Time (ps)",
385 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
386 xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
389 sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
390 leg[e++] = gmx_strdup(str);
391 sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
392 leg[e++] = gmx_strdup(str);
393 sprintf(str, "f. <e\\S-\\betaU\\N>");
394 leg[e++] = gmx_strdup(str);
395 sprintf(str, "f. V");
396 leg[e++] = gmx_strdup(str);
397 sprintf(str, "f. <Ue\\S-\\betaU\\N>");
398 leg[e++] = gmx_strdup(str);
399 for (i = 0; i < ngid; i++)
401 sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
402 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
403 leg[e++] = gmx_strdup(str);
407 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
408 leg[e++] = gmx_strdup(str);
412 for (i = 0; i < ngid; i++)
414 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
415 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
416 leg[e++] = gmx_strdup(str);
420 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
421 leg[e++] = gmx_strdup(str);
423 if (EEL_FULL(fr->ic->eeltype))
425 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
426 leg[e++] = gmx_strdup(str);
429 xvgr_legend(fp_tpi, 4+nener, leg, oenv);
430 for (i = 0; i < 4+nener; i++)
444 /* Avoid frame step numbers <= -1 */
445 frame_step_prev = -1;
447 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
448 &rerun_fr, TRX_NEED_X);
451 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
452 mdatoms->nr - (a_tp1 - a_tp0))
454 gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
455 "is not equal the number in the run input file (%d) "
456 "minus the number of atoms to insert (%d)\n",
457 rerun_fr.natoms, bCavity ? " minus one" : "",
458 mdatoms->nr, a_tp1-a_tp0);
461 refvolshift = log(det(rerun_fr.box));
463 switch (inputrec->eI)
466 stepblocksize = inputrec->nstlist;
472 gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
475 while (bNotLastFrame)
477 frame_step = rerun_fr.step;
478 if (frame_step <= frame_step_prev)
480 /* We don't have step number in the trajectory file,
481 * or we have constant or decreasing step numbers.
482 * Ensure we have increasing step numbers, since we use
483 * the step numbers as a counter for random numbers.
485 frame_step = frame_step_prev + 1;
487 frame_step_prev = frame_step;
489 lambda = rerun_fr.lambda;
493 for (e = 0; e < nener; e++)
498 /* Copy the coordinates from the input trajectory */
499 auto x = makeArrayRef(state_global->x);
500 for (i = 0; i < rerun_fr.natoms; i++)
502 copy_rvec(rerun_fr.x[i], x[i]);
504 copy_mat(rerun_fr.box, state_global->box);
506 V = det(state_global->box);
509 bStateChanged = TRUE;
512 step = cr->nodeid*stepblocksize;
513 while (step < nsteps)
515 /* Restart random engine using the frame and insertion step
517 * Note that we need to draw several random values per iteration,
518 * but by using the internal subcounter functionality of ThreeFry2x64
519 * we can draw 131072 unique 64-bit values before exhausting
520 * the stream. This is a huge margin, and if something still goes
521 * wrong you will get an exception when the stream is exhausted.
523 rng.restart(frame_step, step);
524 dist.reset(); // erase any memory in the distribution
528 /* Random insertion in the whole volume */
529 bNS = (step % inputrec->nstlist == 0);
532 /* Generate a random position in the box */
533 for (d = 0; d < DIM; d++)
535 x_init[d] = dist(rng)*state_global->box[d][d];
539 if (inputrec->nstlist == 1)
541 copy_rvec(x_init, x_tp);
545 /* Generate coordinates within |dx|=drmax of x_init */
548 for (d = 0; d < DIM; d++)
550 dx[d] = (2*dist(rng) - 1)*drmax;
553 while (norm2(dx) > drmax*drmax);
554 rvec_add(x_init, dx, x_tp);
559 /* Random insertion around a cavity location
560 * given by the last coordinate of the trajectory.
566 /* Copy the location of the cavity */
567 copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
571 /* Determine the center of mass of the last molecule */
574 for (i = 0; i < nat_cavity; i++)
576 for (d = 0; d < DIM; d++)
579 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
581 mass_tot += mass_cavity[i];
583 for (d = 0; d < DIM; d++)
585 x_init[d] /= mass_tot;
589 /* Generate coordinates within |dx|=drmax of x_init */
592 for (d = 0; d < DIM; d++)
594 dx[d] = (2*dist(rng) - 1)*drmax;
597 while (norm2(dx) > drmax*drmax);
598 rvec_add(x_init, dx, x_tp);
601 if (a_tp1 - a_tp0 == 1)
603 /* Insert a single atom, just copy the insertion location */
604 copy_rvec(x_tp, x[a_tp0]);
608 /* Copy the coordinates from the top file */
609 for (i = a_tp0; i < a_tp1; i++)
611 copy_rvec(x_mol[i-a_tp0], x[i]);
613 /* Rotate the molecule randomly */
614 real angleX = 2*M_PI*dist(rng);
615 real angleY = 2*M_PI*dist(rng);
616 real angleZ = 2*M_PI*dist(rng);
617 rotate_conf(a_tp1-a_tp0, state_global->x.rvec_array()+a_tp0, nullptr,
618 angleX, angleY, angleZ);
619 /* Shift to the insertion location */
620 for (i = a_tp0; i < a_tp1; i++)
622 rvec_inc(x[i], x_tp);
626 /* Clear some matrix variables */
627 clear_mat(force_vir);
628 clear_mat(shake_vir);
632 /* Set the charge group center of mass of the test particle */
633 copy_rvec(x_init, fr->cg_cm[top.cgs.nr-1]);
635 /* Calc energy (no forces) on new positions.
636 * Since we only need the intermolecular energy
637 * and the RF exclusion terms of the inserted molecule occur
638 * within a single charge group we can pass NULL for the graph.
639 * This also avoids shifts that would move charge groups
641 /* Make do_force do a single node force calculation */
644 // TPI might place a particle so close that the potential
645 // is infinite. Since this is intended to happen, we
646 // temporarily suppress any exceptions that the processor
647 // might raise, then restore the old behaviour.
648 std::fenv_t floatingPointEnvironment;
649 std::feholdexcept(&floatingPointEnvironment);
650 do_force(fplog, cr, ms, inputrec, nullptr, nullptr,
651 step, nrnb, wcycle, &top, &top_global->groups,
652 state_global->box, state_global->x.arrayRefWithPadding(), &state_global->hist,
653 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
654 state_global->lambda,
655 nullptr, fr, ppForceWorkload, nullptr, mu_tot, t, nullptr,
656 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
657 (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS : 0) |
658 (bStateChanged ? GMX_FORCE_STATECHANGED : 0),
659 DDBalanceRegionHandler(nullptr));
660 std::feclearexcept(FE_DIVBYZERO | FE_INVALID | FE_OVERFLOW);
661 std::feupdateenv(&floatingPointEnvironment);
664 bStateChanged = FALSE;
667 /* Calculate long range corrections to pressure and energy */
668 calc_dispcorr(inputrec, fr, state_global->box,
669 lambda, pres, vir, &prescorr, &enercorr, &dvdlcorr);
670 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
671 enerd->term[F_DISPCORR] = enercorr;
672 enerd->term[F_EPOT] += enercorr;
673 enerd->term[F_PRES] += prescorr;
674 enerd->term[F_DVDL_VDW] += dvdlcorr;
676 epot = enerd->term[F_EPOT];
677 bEnergyOutOfBounds = FALSE;
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, static_cast<int>(step), epot);
707 // Exponent argument is fine in SP range, but output can be in DP range
708 embU = exp(static_cast<double>(-beta*epot));
710 /* Determine the weighted energy contributions of each energy group */
712 sum_UgembU[e++] += epot*embU;
715 for (i = 0; i < ngid; i++)
718 enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)]*embU;
723 for (i = 0; i < ngid; i++)
726 enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)]*embU;
731 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
735 for (i = 0; i < ngid; i++)
737 sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
741 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
743 if (EEL_FULL(fr->ic->eeltype))
745 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
750 if (embU == 0 || beta*epot > bU_bin_limit)
756 i = gmx::roundToInt((bU_logV_bin_limit
757 - (beta*epot - logV + refvolshift))*invbinw);
764 realloc_bins(&bin, &nbin, i+10);
771 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
772 static_cast<int>(step), epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
775 if (dump_pdb && epot <= dump_ener)
777 sprintf(str, "t%g_step%d.pdb", t, static_cast<int>(step));
778 sprintf(str2, "t: %f step %d ener: %f", t, static_cast<int>(step), epot);
779 write_sto_conf_mtop(str, str2, top_global, state_global->x.rvec_array(), state_global->v.rvec_array(),
780 inputrec->ePBC, state_global->box);
784 if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
786 /* Skip all steps assigned to the other MPI ranks */
787 step += (cr->nnodes - 1)*stepblocksize;
793 /* When running in parallel sum the energies over the processes */
794 gmx_sumd(1, &sum_embU, cr);
795 gmx_sumd(nener, sum_UgembU, cr);
800 VembU_all += V*sum_embU/nsteps;
804 if (mdrunOptions.verbose || frame%10 == 0 || frame < 10)
806 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
807 -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
810 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
812 VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
813 sum_embU == 0 ? 20/beta : -log(sum_embU/nsteps)/beta,
815 for (e = 0; e < nener; e++)
817 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
819 fprintf(fp_tpi, "\n");
823 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
824 } /* End of the loop */
825 walltime_accounting_end_time(walltime_accounting);
829 if (fp_tpi != nullptr)
834 if (fplog != nullptr)
836 fprintf(fplog, "\n");
837 fprintf(fplog, " <V> = %12.5e nm^3\n", V_all/frame);
838 const double mu = -log(VembU_all/V_all)/beta;
839 fprintf(fplog, " <mu> = %12.5e kJ/mol\n", mu);
841 if (!std::isfinite(mu))
843 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");
847 /* Write the Boltzmann factor histogram */
850 /* When running in parallel sum the bins over the processes */
853 realloc_bins(&bin, &nbin, i);
854 gmx_sumd(nbin, bin, cr);
858 fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
859 "TPI energy distribution",
860 "\\betaU - log(V/<V>)", "count", oenv);
861 sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
862 xvgr_subtitle(fp_tpi, str, oenv);
863 xvgr_legend(fp_tpi, 2, tpid_leg, oenv);
864 for (i = nbin-1; i > 0; i--)
866 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
867 fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
870 bin[i]*exp(-bUlogV)*V_all/VembU_all);
878 walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);