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 PaddedVector<gmx::RVec> f {};
141 real lambda, t, temp, beta, drmax, epot;
142 double embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
145 gmx_bool bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
146 tensor force_vir, shake_vir, vir, pres;
147 int cg_tp, a_tp0, a_tp1, ngid, gid_tp, nener, e;
149 rvec mu_tot, x_init, dx, x_tp;
151 int64_t frame_step_prev, frame_step;
152 int64_t nsteps, stepblocksize = 0, step;
155 FILE *fp_tpi = nullptr;
156 char *ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
157 double dbl, dump_ener;
159 int nat_cavity = 0, d;
160 real *mass_cavity = nullptr, mass_tot;
162 double invbinw, *bin, refvolshift, logV, bUlogV;
163 real prescorr, enercorr, dvdlcorr;
164 gmx_bool bEnergyOutOfBounds;
165 const char *tpid_leg[2] = {"direct", "reweighted"};
166 auto mdatoms = mdAtoms->mdatoms();
168 GMX_UNUSED_VALUE(outputProvider);
170 GMX_LOG(mdlog.info).asParagraph().
171 appendText("Note that it is planned to change the command gmx mdrun -tpi "
172 "(and -tpic) to make the functionality available in a different "
173 "form in a future version of GROMACS, e.g. gmx test-particle-insertion.");
175 /* Since there is no upper limit to the insertion energies,
176 * we need to set an upper limit for the distribution output.
178 real bU_bin_limit = 50;
179 real bU_logV_bin_limit = bU_bin_limit + 10;
181 if (inputrec->cutoff_scheme == ecutsVERLET)
183 gmx_fatal(FARGS, "TPI does not work (yet) with the Verlet cut-off scheme");
188 gmx_mtop_generate_local_top(*top_global, &top, inputrec->efep != efepNO);
190 SimulationGroups *groups = &top_global->groups;
192 bCavity = (inputrec->eI == eiTPIC);
195 ptr = getenv("GMX_TPIC_MASSES");
202 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
203 * The center of mass of the last atoms is then used for TPIC.
206 while (sscanf(ptr, "%20lf%n", &dbl, &i) > 0)
208 srenew(mass_cavity, nat_cavity+1);
209 mass_cavity[nat_cavity] = dbl;
210 fprintf(fplog, "mass[%d] = %f\n",
211 nat_cavity+1, mass_cavity[nat_cavity]);
217 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
223 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
224 state_global->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
225 /* We never need full pbc for TPI */
227 /* Determine the temperature for the Boltzmann weighting */
228 temp = inputrec->opts.ref_t[0];
231 for (i = 1; (i < inputrec->opts.ngtc); i++)
233 if (inputrec->opts.ref_t[i] != temp)
235 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
236 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
240 "\n The temperature for test particle insertion is %.3f K\n\n",
243 beta = 1.0/(BOLTZ*temp);
245 /* Number of insertions per frame */
246 nsteps = inputrec->nsteps;
248 /* Use the same neighborlist with more insertions points
249 * in a sphere of radius drmax around the initial point
251 /* This should be a proper mdp parameter */
252 drmax = inputrec->rtpi;
254 /* An environment variable can be set to dump all configurations
255 * to pdb with an insertion energy <= this value.
257 dump_pdb = getenv("GMX_TPI_DUMP");
261 sscanf(dump_pdb, "%20lf", &dump_ener);
264 atoms2md(top_global, inputrec, -1, nullptr, top_global->natoms, mdAtoms);
265 update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
267 f.resizeWithPadding(top_global->natoms);
269 /* Print to log file */
270 walltime_accounting_start_time(walltime_accounting);
271 wallcycle_start(wcycle, ewcRUN);
272 print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
274 /* The last charge group is the group to be inserted */
275 cg_tp = top.cgs.nr - 1;
276 a_tp0 = top.cgs.index[cg_tp];
277 a_tp1 = top.cgs.index[cg_tp+1];
280 fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
283 GMX_RELEASE_ASSERT(inputrec->rcoulomb <= inputrec->rlist && inputrec->rvdw <= inputrec->rlist, "Twin-range interactions are not supported with TPI");
285 snew(x_mol, a_tp1-a_tp0);
287 bDispCorr = (inputrec->eDispCorr != edispcNO);
289 auto x = makeArrayRef(state_global->x);
290 for (i = a_tp0; i < a_tp1; i++)
292 /* Copy the coordinates of the molecule to be insterted */
293 copy_rvec(x[i], x_mol[i-a_tp0]);
294 /* Check if we need to print electrostatic energies */
295 bCharge |= (mdatoms->chargeA[i] != 0 ||
296 ((mdatoms->chargeB != nullptr) && mdatoms->chargeB[i] != 0));
298 bRFExcl = (bCharge && EEL_RF(fr->ic->eeltype));
300 calc_cgcm(fplog, cg_tp, cg_tp+1, &(top.cgs), state_global->x.rvec_array(), fr->cg_cm);
303 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
305 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
306 fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
311 /* Center the molecule to be inserted at zero */
312 for (i = 0; i < a_tp1-a_tp0; i++)
314 rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
320 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
321 a_tp1-a_tp0, bCharge ? "with" : "without");
323 fprintf(fplog, "\nWill insert %" PRId64 " times in each frame of %s\n",
324 nsteps, opt2fn("-rerun", nfile, fnm));
329 if (inputrec->nstlist > 1)
331 if (drmax == 0 && a_tp1-a_tp0 == 1)
333 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);
337 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
345 fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
349 ngid = groups->groups[SimulationAtomGroupType::EnergyOutput].nr;
350 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
363 if (EEL_FULL(fr->ic->eeltype))
368 snew(sum_UgembU, nener);
370 /* Copy the random seed set by the user */
371 seed = inputrec->ld_seed;
373 gmx::ThreeFry2x64<16> rng(seed, gmx::RandomDomain::TestParticleInsertion); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
374 gmx::UniformRealDistribution<real> dist;
378 fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
379 "TPI energies", "Time (ps)",
380 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
381 xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
384 sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
385 leg[e++] = gmx_strdup(str);
386 sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
387 leg[e++] = gmx_strdup(str);
388 sprintf(str, "f. <e\\S-\\betaU\\N>");
389 leg[e++] = gmx_strdup(str);
390 sprintf(str, "f. V");
391 leg[e++] = gmx_strdup(str);
392 sprintf(str, "f. <Ue\\S-\\betaU\\N>");
393 leg[e++] = gmx_strdup(str);
394 for (i = 0; i < ngid; i++)
396 sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
397 *(groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[i]]));
398 leg[e++] = gmx_strdup(str);
402 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
403 leg[e++] = gmx_strdup(str);
407 for (i = 0; i < ngid; i++)
409 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
410 *(groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind[i]]));
411 leg[e++] = gmx_strdup(str);
415 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
416 leg[e++] = gmx_strdup(str);
418 if (EEL_FULL(fr->ic->eeltype))
420 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
421 leg[e++] = gmx_strdup(str);
424 xvgr_legend(fp_tpi, 4+nener, leg, oenv);
425 for (i = 0; i < 4+nener; i++)
439 /* Avoid frame step numbers <= -1 */
440 frame_step_prev = -1;
442 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
443 &rerun_fr, TRX_NEED_X);
446 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
447 mdatoms->nr - (a_tp1 - a_tp0))
449 gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
450 "is not equal the number in the run input file (%d) "
451 "minus the number of atoms to insert (%d)\n",
452 rerun_fr.natoms, bCavity ? " minus one" : "",
453 mdatoms->nr, a_tp1-a_tp0);
456 refvolshift = log(det(rerun_fr.box));
458 switch (inputrec->eI)
461 stepblocksize = inputrec->nstlist;
467 gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
470 while (bNotLastFrame)
472 frame_step = rerun_fr.step;
473 if (frame_step <= frame_step_prev)
475 /* We don't have step number in the trajectory file,
476 * or we have constant or decreasing step numbers.
477 * Ensure we have increasing step numbers, since we use
478 * the step numbers as a counter for random numbers.
480 frame_step = frame_step_prev + 1;
482 frame_step_prev = frame_step;
484 lambda = rerun_fr.lambda;
488 for (e = 0; e < nener; e++)
493 /* Copy the coordinates from the input trajectory */
494 auto x = makeArrayRef(state_global->x);
495 for (i = 0; i < rerun_fr.natoms; i++)
497 copy_rvec(rerun_fr.x[i], x[i]);
499 copy_mat(rerun_fr.box, state_global->box);
501 V = det(state_global->box);
504 bStateChanged = TRUE;
507 step = cr->nodeid*stepblocksize;
508 while (step < nsteps)
510 /* Restart random engine using the frame and insertion step
512 * Note that we need to draw several random values per iteration,
513 * but by using the internal subcounter functionality of ThreeFry2x64
514 * we can draw 131072 unique 64-bit values before exhausting
515 * the stream. This is a huge margin, and if something still goes
516 * wrong you will get an exception when the stream is exhausted.
518 rng.restart(frame_step, step);
519 dist.reset(); // erase any memory in the distribution
523 /* Random insertion in the whole volume */
524 bNS = (step % inputrec->nstlist == 0);
527 /* Generate a random position in the box */
528 for (d = 0; d < DIM; d++)
530 x_init[d] = dist(rng)*state_global->box[d][d];
534 if (inputrec->nstlist == 1)
536 copy_rvec(x_init, x_tp);
540 /* Generate coordinates within |dx|=drmax of x_init */
543 for (d = 0; d < DIM; d++)
545 dx[d] = (2*dist(rng) - 1)*drmax;
548 while (norm2(dx) > drmax*drmax);
549 rvec_add(x_init, dx, x_tp);
554 /* Random insertion around a cavity location
555 * given by the last coordinate of the trajectory.
561 /* Copy the location of the cavity */
562 copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
566 /* Determine the center of mass of the last molecule */
569 for (i = 0; i < nat_cavity; i++)
571 for (d = 0; d < DIM; d++)
574 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
576 mass_tot += mass_cavity[i];
578 for (d = 0; d < DIM; d++)
580 x_init[d] /= mass_tot;
584 /* Generate coordinates within |dx|=drmax of x_init */
587 for (d = 0; d < DIM; d++)
589 dx[d] = (2*dist(rng) - 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, 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], x[i]);
608 /* Rotate the molecule randomly */
609 real angleX = 2*M_PI*dist(rng);
610 real angleY = 2*M_PI*dist(rng);
611 real angleZ = 2*M_PI*dist(rng);
612 rotate_conf(a_tp1-a_tp0, state_global->x.rvec_array()+a_tp0, nullptr,
613 angleX, angleY, angleZ);
614 /* Shift to the insertion location */
615 for (i = a_tp0; i < a_tp1; i++)
617 rvec_inc(x[i], x_tp);
621 /* Clear some matrix variables */
622 clear_mat(force_vir);
623 clear_mat(shake_vir);
627 /* Set the charge group center of mass of the test particle */
628 copy_rvec(x_init, fr->cg_cm[top.cgs.nr-1]);
630 /* Calc energy (no forces) on new positions.
631 * Since we only need the intermolecular energy
632 * and the RF exclusion terms of the inserted molecule occur
633 * within a single charge group we can pass NULL for the graph.
634 * This also avoids shifts that would move charge groups
636 /* Make do_force do a single node force calculation */
639 // TPI might place a particle so close that the potential
640 // is infinite. Since this is intended to happen, we
641 // temporarily suppress any exceptions that the processor
642 // might raise, then restore the old behaviour.
643 std::fenv_t floatingPointEnvironment;
644 std::feholdexcept(&floatingPointEnvironment);
645 do_force(fplog, cr, ms, inputrec, nullptr, nullptr, imdSession,
646 step, nrnb, wcycle, &top,
647 state_global->box, state_global->x.arrayRefWithPadding(), &state_global->hist,
648 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
649 state_global->lambda,
650 nullptr, fr, ppForceWorkload, nullptr, mu_tot, t, nullptr,
651 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
652 (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS : 0) |
653 (bStateChanged ? GMX_FORCE_STATECHANGED : 0),
654 DDBalanceRegionHandler(nullptr));
655 std::feclearexcept(FE_DIVBYZERO | FE_INVALID | FE_OVERFLOW);
656 std::feupdateenv(&floatingPointEnvironment);
659 bStateChanged = FALSE;
662 /* Calculate long range corrections to pressure and energy */
663 calc_dispcorr(inputrec, fr, state_global->box,
664 lambda, pres, vir, &prescorr, &enercorr, &dvdlcorr);
665 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
666 enerd->term[F_DISPCORR] = enercorr;
667 enerd->term[F_EPOT] += enercorr;
668 enerd->term[F_PRES] += prescorr;
669 enerd->term[F_DVDL_VDW] += dvdlcorr;
671 epot = enerd->term[F_EPOT];
672 bEnergyOutOfBounds = FALSE;
674 /* If the compiler doesn't optimize this check away
675 * we catch the NAN energies.
676 * The epot>GMX_REAL_MAX check catches inf values,
677 * which should nicely result in embU=0 through the exp below,
678 * but it does not hurt to check anyhow.
680 /* Non-bonded Interaction usually diverge at r=0.
681 * With tabulated interaction functions the first few entries
682 * should be capped in a consistent fashion between
683 * repulsion, dispersion and Coulomb to avoid accidental
684 * negative values in the total energy.
685 * The table generation code in tables.c does this.
686 * With user tbales the user should take care of this.
688 if (epot != epot || epot > GMX_REAL_MAX)
690 bEnergyOutOfBounds = TRUE;
692 if (bEnergyOutOfBounds)
696 fprintf(debug, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, static_cast<int>(step), epot);
702 // Exponent argument is fine in SP range, but output can be in DP range
703 embU = exp(static_cast<double>(-beta*epot));
705 /* Determine the weighted energy contributions of each energy group */
707 sum_UgembU[e++] += epot*embU;
710 for (i = 0; i < ngid; i++)
713 enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)]*embU;
718 for (i = 0; i < ngid; i++)
721 enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)]*embU;
726 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
730 for (i = 0; i < ngid; i++)
732 sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
736 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
738 if (EEL_FULL(fr->ic->eeltype))
740 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
745 if (embU == 0 || beta*epot > bU_bin_limit)
751 i = gmx::roundToInt((bU_logV_bin_limit
752 - (beta*epot - logV + refvolshift))*invbinw);
759 realloc_bins(&bin, &nbin, i+10);
766 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
767 static_cast<int>(step), epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
770 if (dump_pdb && epot <= dump_ener)
772 sprintf(str, "t%g_step%d.pdb", t, static_cast<int>(step));
773 sprintf(str2, "t: %f step %d ener: %f", t, static_cast<int>(step), epot);
774 write_sto_conf_mtop(str, str2, top_global, state_global->x.rvec_array(), state_global->v.rvec_array(),
775 inputrec->ePBC, state_global->box);
779 if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
781 /* Skip all steps assigned to the other MPI ranks */
782 step += (cr->nnodes - 1)*stepblocksize;
788 /* When running in parallel sum the energies over the processes */
789 gmx_sumd(1, &sum_embU, cr);
790 gmx_sumd(nener, sum_UgembU, cr);
795 VembU_all += V*sum_embU/nsteps;
799 if (mdrunOptions.verbose || frame%10 == 0 || frame < 10)
801 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
802 -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
805 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
807 VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
808 sum_embU == 0 ? 20/beta : -log(sum_embU/nsteps)/beta,
810 for (e = 0; e < nener; e++)
812 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
814 fprintf(fp_tpi, "\n");
818 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
819 } /* End of the loop */
820 walltime_accounting_end_time(walltime_accounting);
824 if (fp_tpi != nullptr)
829 if (fplog != nullptr)
831 fprintf(fplog, "\n");
832 fprintf(fplog, " <V> = %12.5e nm^3\n", V_all/frame);
833 const double mu = -log(VembU_all/V_all)/beta;
834 fprintf(fplog, " <mu> = %12.5e kJ/mol\n", mu);
836 if (!std::isfinite(mu))
838 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");
842 /* Write the Boltzmann factor histogram */
845 /* When running in parallel sum the bins over the processes */
848 realloc_bins(&bin, &nbin, i);
849 gmx_sumd(nbin, bin, cr);
853 fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
854 "TPI energy distribution",
855 "\\betaU - log(V/<V>)", "count", oenv);
856 sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
857 xvgr_subtitle(fp_tpi, str, oenv);
858 xvgr_legend(fp_tpi, 2, tpid_leg, oenv);
859 for (i = nbin-1; i > 0; i--)
861 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
862 fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
865 bin[i]*exp(-bUlogV)*V_all/VembU_all);
873 walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);