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/conformation_utilities.h"
63 #include "gromacs/gmxlib/network.h"
64 #include "gromacs/gmxlib/nrnb.h"
65 #include "gromacs/math/units.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/mdlib/constr.h"
68 #include "gromacs/mdlib/dispersioncorrection.h"
69 #include "gromacs/mdlib/energyoutput.h"
70 #include "gromacs/mdlib/force.h"
71 #include "gromacs/mdlib/force_flags.h"
72 #include "gromacs/mdlib/gmx_omp_nthreads.h"
73 #include "gromacs/mdlib/mdatoms.h"
74 #include "gromacs/mdlib/tgroup.h"
75 #include "gromacs/mdlib/update.h"
76 #include "gromacs/mdlib/vsite.h"
77 #include "gromacs/mdrunutility/printtime.h"
78 #include "gromacs/mdtypes/commrec.h"
79 #include "gromacs/mdtypes/forcerec.h"
80 #include "gromacs/mdtypes/group.h"
81 #include "gromacs/mdtypes/inputrec.h"
82 #include "gromacs/mdtypes/md_enums.h"
83 #include "gromacs/mdtypes/mdrunoptions.h"
84 #include "gromacs/mdtypes/state.h"
85 #include "gromacs/nbnxm/nbnxm.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 "legacysimulator.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++)
133 //! Computes and returns the RF exclusion energy for the last molecule starting at \p beginAtom
134 static real reactionFieldExclusionCorrection(gmx::ArrayRef<const gmx::RVec> x,
135 const t_mdatoms& mdatoms,
136 const interaction_const_t& ic,
141 for (int i = beginAtom; i < mdatoms.homenr; i++)
143 const real qi = mdatoms.chargeA[i];
144 energy -= 0.5 * qi * qi * ic.c_rf;
146 for (int j = i + 1; j < mdatoms.homenr; j++)
148 const real qj = mdatoms.chargeA[j];
149 const real rsq = distance2(x[i], x[j]);
150 energy += qi * qj * (ic.k_rf * rsq - ic.c_rf);
154 return ic.epsfac * energy;
160 // TODO: Convert to use the nbnxm kernels by putting the system and the teset molecule on two separate search grids
161 void LegacySimulator::do_tpi()
163 GMX_RELEASE_ASSERT(gmx_omp_nthreads_get(emntDefault) == 1, "TPI does not support OpenMP");
166 PaddedHostVector<gmx::RVec> f{};
167 real lambda, t, temp, beta, drmax, epot;
168 double embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
171 gmx_bool bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
172 tensor force_vir, shake_vir, vir, pres;
173 int a_tp0, a_tp1, ngid, gid_tp, nener, e;
175 rvec mu_tot, x_init, dx;
177 int64_t frame_step_prev, frame_step;
178 int64_t nsteps, stepblocksize = 0, step;
181 FILE* fp_tpi = nullptr;
182 char * ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
183 double dbl, dump_ener;
185 int nat_cavity = 0, d;
186 real * mass_cavity = nullptr, mass_tot;
188 double invbinw, *bin, refvolshift, logV, bUlogV;
189 gmx_bool bEnergyOutOfBounds;
190 const char* tpid_leg[2] = { "direct", "reweighted" };
191 auto mdatoms = mdAtoms->mdatoms();
193 GMX_UNUSED_VALUE(outputProvider);
198 "Note that it is planned to change the command gmx mdrun -tpi "
199 "(and -tpic) to make the functionality available in a different "
200 "form in a future version of GROMACS, e.g. gmx test-particle-insertion.");
202 /* Since there is no upper limit to the insertion energies,
203 * we need to set an upper limit for the distribution output.
205 real bU_bin_limit = 50;
206 real bU_logV_bin_limit = bU_bin_limit + 10;
210 gmx_mtop_generate_local_top(*top_global, &top, inputrec->efep != efepNO);
212 SimulationGroups* groups = &top_global->groups;
214 bCavity = (inputrec->eI == eiTPIC);
217 ptr = getenv("GMX_TPIC_MASSES");
224 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
225 * The center of mass of the last atoms is then used for TPIC.
228 while (sscanf(ptr, "%20lf%n", &dbl, &i) > 0)
230 srenew(mass_cavity, nat_cavity + 1);
231 mass_cavity[nat_cavity] = dbl;
232 fprintf(fplog, "mass[%d] = %f\n", nat_cavity + 1, mass_cavity[nat_cavity]);
238 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
244 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
245 state_global->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
246 /* We never need full pbc for TPI */
248 /* Determine the temperature for the Boltzmann weighting */
249 temp = inputrec->opts.ref_t[0];
252 for (i = 1; (i < inputrec->opts.ngtc); i++)
254 if (inputrec->opts.ref_t[i] != temp)
257 "\nWARNING: The temperatures of the different temperature coupling groups "
258 "are not identical\n\n");
260 "\nWARNING: The temperatures of the different temperature coupling groups "
261 "are not identical\n\n");
264 fprintf(fplog, "\n The temperature for test particle insertion is %.3f K\n\n", temp);
266 beta = 1.0 / (BOLTZ * temp);
268 /* Number of insertions per frame */
269 nsteps = inputrec->nsteps;
271 /* Use the same neighborlist with more insertions points
272 * in a sphere of radius drmax around the initial point
274 /* This should be a proper mdp parameter */
275 drmax = inputrec->rtpi;
277 /* An environment variable can be set to dump all configurations
278 * to pdb with an insertion energy <= this value.
280 dump_pdb = getenv("GMX_TPI_DUMP");
284 sscanf(dump_pdb, "%20lf", &dump_ener);
287 atoms2md(top_global, inputrec, -1, nullptr, top_global->natoms, mdAtoms);
288 update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
290 f.resizeWithPadding(top_global->natoms);
292 /* Print to log file */
293 walltime_accounting_start_time(walltime_accounting);
294 wallcycle_start(wcycle, ewcRUN);
295 print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
297 /* The last charge group is the group to be inserted */
298 const t_atoms& atomsToInsert = top_global->moltype[top_global->molblock.back().type].atoms;
299 a_tp0 = top_global->natoms - atomsToInsert.nr;
300 a_tp1 = top_global->natoms;
303 fprintf(debug, "TPI atoms %d-%d\n", a_tp0, a_tp1);
306 auto x = makeArrayRef(state_global->x);
308 if (EEL_PME(fr->ic->eeltype))
310 gmx_pme_reinit_atoms(fr->pmedata, a_tp0, nullptr);
313 /* With reacion-field we have distance dependent potentials
314 * between excluded atoms, we need to add these separately
315 * for the inserted molecule.
317 real rfExclusionEnergy = 0;
318 if (EEL_RF(fr->ic->eeltype))
320 rfExclusionEnergy = reactionFieldExclusionCorrection(x, *mdatoms, *fr->ic, a_tp0);
323 fprintf(debug, "RF exclusion correction for inserted molecule: %f kJ/mol\n", rfExclusionEnergy);
327 snew(x_mol, a_tp1 - a_tp0);
329 bDispCorr = (inputrec->eDispCorr != edispcNO);
331 for (i = a_tp0; i < a_tp1; i++)
333 /* Copy the coordinates of the molecule to be insterted */
334 copy_rvec(x[i], x_mol[i - a_tp0]);
335 /* Check if we need to print electrostatic energies */
336 bCharge |= (mdatoms->chargeA[i] != 0
337 || ((mdatoms->chargeB != nullptr) && mdatoms->chargeB[i] != 0));
339 bRFExcl = (bCharge && EEL_RF(fr->ic->eeltype));
341 // Calculate the center of geometry of the molecule to insert
342 rvec cog = { 0, 0, 0 };
343 for (int a = a_tp0; a < a_tp1; a++)
347 svmul(1.0_real / (a_tp1 - a_tp0), cog, cog);
349 for (int a = a_tp0; a < a_tp1; a++)
351 molRadius = std::max(molRadius, distance2(x[a], cog));
353 molRadius = std::sqrt(molRadius);
355 const real maxCutoff = std::max(inputrec->rvdw, inputrec->rcoulomb);
358 if (norm(cog) > 0.5 * maxCutoff && fplog)
360 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
361 fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
366 /* Center the molecule to be inserted at zero */
367 for (i = 0; i < a_tp1 - a_tp0; i++)
369 rvec_dec(x_mol[i], cog);
375 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n", a_tp1 - a_tp0,
376 bCharge ? "with" : "without");
378 fprintf(fplog, "\nWill insert %" PRId64 " times in each frame of %s\n", nsteps,
379 opt2fn("-rerun", nfile, fnm));
384 if (inputrec->nstlist > 1)
387 /* With the same pair list we insert in a sphere of radius rtpi in different orientations */
388 if (drmax == 0 && a_tp1 - a_tp0 == 1)
391 "Re-using the neighborlist %d times for insertions of a single atom in a "
392 "sphere of radius %f does not make sense",
393 inputrec->nstlist, drmax);
398 "Will use the same neighborlist for %d insertions in a sphere of radius "
400 inputrec->nstlist, drmax);
409 "Will insert randomly in a sphere of radius %f around the center of the "
415 /* With the same pair list we insert in a sphere of radius rtpi
416 * in different orientations. We generate the pairlist with all
417 * inserted atoms located in the center of the sphere, so we need
418 * a buffer of size of the sphere and molecule radius.
420 inputrec->rlist = maxCutoff + 2 * inputrec->rtpi + 2 * molRadius;
421 fr->rlist = inputrec->rlist;
422 fr->nbv->changePairlistRadii(inputrec->rlist, inputrec->rlist);
424 ngid = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
425 gid_tp = GET_CGINFO_GID(fr->cginfo[a_tp0]);
426 for (int a = a_tp0 + 1; a < a_tp1; a++)
428 if (GET_CGINFO_GID(fr->cginfo[a]) != gid_tp)
431 "NOTE: Atoms in the molecule to insert belong to different energy groups.\n"
432 " Only contributions to the group of the first atom will be reported.\n");
448 if (EEL_FULL(fr->ic->eeltype))
453 snew(sum_UgembU, nener);
455 /* Copy the random seed set by the user */
456 seed = inputrec->ld_seed;
458 gmx::ThreeFry2x64<16> rng(
459 seed, gmx::RandomDomain::TestParticleInsertion); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
460 gmx::UniformRealDistribution<real> dist;
464 fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm), "TPI energies", "Time (ps)",
465 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
466 xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
467 snew(leg, 4 + nener);
469 sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
470 leg[e++] = gmx_strdup(str);
471 sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
472 leg[e++] = gmx_strdup(str);
473 sprintf(str, "f. <e\\S-\\betaU\\N>");
474 leg[e++] = gmx_strdup(str);
475 sprintf(str, "f. V");
476 leg[e++] = gmx_strdup(str);
477 sprintf(str, "f. <Ue\\S-\\betaU\\N>");
478 leg[e++] = gmx_strdup(str);
479 for (i = 0; i < ngid; i++)
481 sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
482 *(groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput][i]]));
483 leg[e++] = gmx_strdup(str);
487 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
488 leg[e++] = gmx_strdup(str);
492 for (i = 0; i < ngid; i++)
494 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
495 *(groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput][i]]));
496 leg[e++] = gmx_strdup(str);
500 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
501 leg[e++] = gmx_strdup(str);
503 if (EEL_FULL(fr->ic->eeltype))
505 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
506 leg[e++] = gmx_strdup(str);
509 xvgr_legend(fp_tpi, 4 + nener, leg, oenv);
510 for (i = 0; i < 4 + nener; i++)
524 /* Avoid frame step numbers <= -1 */
525 frame_step_prev = -1;
527 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
530 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) != mdatoms->nr - (a_tp1 - a_tp0))
533 "Number of atoms in trajectory (%d)%s "
534 "is not equal the number in the run input file (%d) "
535 "minus the number of atoms to insert (%d)\n",
536 rerun_fr.natoms, bCavity ? " minus one" : "", mdatoms->nr, a_tp1 - a_tp0);
539 refvolshift = log(det(rerun_fr.box));
541 switch (inputrec->eI)
543 case eiTPI: stepblocksize = inputrec->nstlist; break;
544 case eiTPIC: stepblocksize = 1; break;
545 default: gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
548 while (bNotLastFrame)
550 frame_step = rerun_fr.step;
551 if (frame_step <= frame_step_prev)
553 /* We don't have step number in the trajectory file,
554 * or we have constant or decreasing step numbers.
555 * Ensure we have increasing step numbers, since we use
556 * the step numbers as a counter for random numbers.
558 frame_step = frame_step_prev + 1;
560 frame_step_prev = frame_step;
562 lambda = rerun_fr.lambda;
566 for (e = 0; e < nener; e++)
571 /* Copy the coordinates from the input trajectory */
572 auto x = makeArrayRef(state_global->x);
573 for (i = 0; i < rerun_fr.natoms; i++)
575 copy_rvec(rerun_fr.x[i], x[i]);
577 copy_mat(rerun_fr.box, state_global->box);
578 const matrix& box = state_global->box;
583 bStateChanged = TRUE;
586 put_atoms_in_box(fr->ePBC, box, x);
588 /* Put all atoms except for the inserted ones on the grid */
589 rvec vzero = { 0, 0, 0 };
590 rvec boxDiagonal = { box[XX][XX], box[YY][YY], box[ZZ][ZZ] };
591 nbnxn_put_on_grid(fr->nbv.get(), box, 0, vzero, boxDiagonal, nullptr, { 0, a_tp0 }, -1,
592 fr->cginfo, x, 0, nullptr);
594 step = cr->nodeid * stepblocksize;
595 while (step < nsteps)
597 /* Restart random engine using the frame and insertion step
599 * Note that we need to draw several random values per iteration,
600 * but by using the internal subcounter functionality of ThreeFry2x64
601 * we can draw 131072 unique 64-bit values before exhausting
602 * the stream. This is a huge margin, and if something still goes
603 * wrong you will get an exception when the stream is exhausted.
605 rng.restart(frame_step, step);
606 dist.reset(); // erase any memory in the distribution
610 /* Random insertion in the whole volume */
611 bNS = (step % inputrec->nstlist == 0);
614 /* Generate a random position in the box */
615 for (d = 0; d < DIM; d++)
617 x_init[d] = dist(rng) * state_global->box[d][d];
623 /* Random insertion around a cavity location
624 * given by the last coordinate of the trajectory.
630 /* Copy the location of the cavity */
631 copy_rvec(rerun_fr.x[rerun_fr.natoms - 1], x_init);
635 /* Determine the center of mass of the last molecule */
638 for (i = 0; i < nat_cavity; i++)
640 for (d = 0; d < DIM; d++)
642 x_init[d] += mass_cavity[i]
643 * rerun_fr.x[rerun_fr.natoms - nat_cavity + i][d];
645 mass_tot += mass_cavity[i];
647 for (d = 0; d < DIM; d++)
649 x_init[d] /= mass_tot;
657 for (int a = a_tp0; a < a_tp1; a++)
662 /* Put the inserted molecule on it's own search grid */
663 nbnxn_put_on_grid(fr->nbv.get(), box, 1, x_init, x_init, nullptr, { a_tp0, a_tp1 },
664 -1, fr->cginfo, x, 0, nullptr);
666 /* TODO: Avoid updating all atoms at every bNS step */
667 fr->nbv->setAtomProperties(*mdatoms, fr->cginfo);
669 fr->nbv->constructPairlist(InteractionLocality::Local, &top.excls, step, nrnb);
674 /* Add random displacement uniformly distributed in a sphere
675 * of radius rtpi. We don't need to do this is we generate
676 * a new center location every step.
679 if (bCavity || inputrec->nstlist > 1)
681 /* Generate coordinates within |dx|=drmax of x_init */
684 for (d = 0; d < DIM; d++)
686 dx[d] = (2 * dist(rng) - 1) * drmax;
688 } while (norm2(dx) > drmax * drmax);
689 rvec_add(x_init, dx, x_tp);
693 copy_rvec(x_init, x_tp);
696 if (a_tp1 - a_tp0 == 1)
698 /* Insert a single atom, just copy the insertion location */
699 copy_rvec(x_tp, x[a_tp0]);
703 /* Copy the coordinates from the top file */
704 for (i = a_tp0; i < a_tp1; i++)
706 copy_rvec(x_mol[i - a_tp0], x[i]);
708 /* Rotate the molecule randomly */
709 real angleX = 2 * M_PI * dist(rng);
710 real angleY = 2 * M_PI * dist(rng);
711 real angleZ = 2 * M_PI * dist(rng);
712 rotate_conf(a_tp1 - a_tp0, state_global->x.rvec_array() + a_tp0, nullptr, angleX,
714 /* Shift to the insertion location */
715 for (i = a_tp0; i < a_tp1; i++)
717 rvec_inc(x[i], x_tp);
721 /* Note: NonLocal refers to the inserted molecule */
722 fr->nbv->convertCoordinates(AtomLocality::NonLocal, false, x);
724 /* Clear some matrix variables */
725 clear_mat(force_vir);
726 clear_mat(shake_vir);
730 /* Calc energy (no forces) on new positions.
731 * Since we only need the intermolecular energy
732 * and the RF exclusion terms of the inserted molecule occur
733 * within a single charge group we can pass NULL for the graph.
734 * This also avoids shifts that would move charge groups
736 /* Make do_force do a single node force calculation */
739 // TPI might place a particle so close that the potential
740 // is infinite. Since this is intended to happen, we
741 // temporarily suppress any exceptions that the processor
742 // might raise, then restore the old behaviour.
743 std::fenv_t floatingPointEnvironment;
744 std::feholdexcept(&floatingPointEnvironment);
745 do_force(fplog, cr, ms, inputrec, nullptr, nullptr, imdSession, pull_work, step, nrnb,
746 wcycle, &top, state_global->box, state_global->x.arrayRefWithPadding(),
747 &state_global->hist, f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
748 state_global->lambda, nullptr, fr, runScheduleWork, nullptr, mu_tot, t, nullptr,
749 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY | (bStateChanged ? GMX_FORCE_STATECHANGED : 0),
750 DDBalanceRegionHandler(nullptr));
751 std::feclearexcept(FE_DIVBYZERO | FE_INVALID | FE_OVERFLOW);
752 std::feupdateenv(&floatingPointEnvironment);
755 bStateChanged = FALSE;
757 if (fr->dispersionCorrection)
759 /* Calculate long range corrections to pressure and energy */
760 const DispersionCorrection::Correction correction =
761 fr->dispersionCorrection->calculate(state_global->box, lambda);
762 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
763 enerd->term[F_DISPCORR] = correction.energy;
764 enerd->term[F_EPOT] += correction.energy;
765 enerd->term[F_PRES] += correction.pressure;
766 enerd->term[F_DVDL] += correction.dvdl;
770 enerd->term[F_DISPCORR] = 0;
772 if (EEL_RF(fr->ic->eeltype))
774 enerd->term[F_EPOT] += rfExclusionEnergy;
777 epot = enerd->term[F_EPOT];
778 bEnergyOutOfBounds = FALSE;
780 /* If the compiler doesn't optimize this check away
781 * we catch the NAN energies.
782 * The epot>GMX_REAL_MAX check catches inf values,
783 * which should nicely result in embU=0 through the exp below,
784 * but it does not hurt to check anyhow.
786 /* Non-bonded Interaction usually diverge at r=0.
787 * With tabulated interaction functions the first few entries
788 * should be capped in a consistent fashion between
789 * repulsion, dispersion and Coulomb to avoid accidental
790 * negative values in the total energy.
791 * The table generation code in tables.c does this.
792 * With user tbales the user should take care of this.
794 if (epot != epot || epot > GMX_REAL_MAX)
796 bEnergyOutOfBounds = TRUE;
798 if (bEnergyOutOfBounds)
803 "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t,
804 static_cast<int>(step), epot);
810 // Exponent argument is fine in SP range, but output can be in DP range
811 embU = exp(static_cast<double>(-beta * epot));
813 /* Determine the weighted energy contributions of each energy group */
815 sum_UgembU[e++] += epot * embU;
818 for (i = 0; i < ngid; i++)
820 sum_UgembU[e++] += enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)] * embU;
825 for (i = 0; i < ngid; i++)
827 sum_UgembU[e++] += enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)] * embU;
832 sum_UgembU[e++] += enerd->term[F_DISPCORR] * embU;
836 for (i = 0; i < ngid; i++)
838 sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
842 sum_UgembU[e++] += rfExclusionEnergy * embU;
844 if (EEL_FULL(fr->ic->eeltype))
846 sum_UgembU[e++] += enerd->term[F_COUL_RECIP] * embU;
851 if (embU == 0 || beta * epot > bU_bin_limit)
857 i = gmx::roundToInt((bU_logV_bin_limit - (beta * epot - logV + refvolshift)) * invbinw);
864 realloc_bins(&bin, &nbin, i + 10);
871 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n", static_cast<int>(step),
872 epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
875 if (dump_pdb && epot <= dump_ener)
877 sprintf(str, "t%g_step%d.pdb", t, static_cast<int>(step));
878 sprintf(str2, "t: %f step %d ener: %f", t, static_cast<int>(step), epot);
879 write_sto_conf_mtop(str, str2, top_global, state_global->x.rvec_array(),
880 state_global->v.rvec_array(), inputrec->ePBC, state_global->box);
884 if ((step / stepblocksize) % cr->nnodes != cr->nodeid)
886 /* Skip all steps assigned to the other MPI ranks */
887 step += (cr->nnodes - 1) * stepblocksize;
893 /* When running in parallel sum the energies over the processes */
894 gmx_sumd(1, &sum_embU, cr);
895 gmx_sumd(nener, sum_UgembU, cr);
900 VembU_all += V * sum_embU / nsteps;
904 if (mdrunOptions.verbose || frame % 10 == 0 || frame < 10)
906 fprintf(stderr, "mu %10.3e <mu> %10.3e\n", -log(sum_embU / nsteps) / beta,
907 -log(VembU_all / V_all) / beta);
910 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e", t,
911 VembU_all == 0 ? 20 / beta : -log(VembU_all / V_all) / beta,
912 sum_embU == 0 ? 20 / beta : -log(sum_embU / nsteps) / beta, sum_embU / nsteps, V);
913 for (e = 0; e < nener; e++)
915 fprintf(fp_tpi, " %12.5e", sum_UgembU[e] / nsteps);
917 fprintf(fp_tpi, "\n");
921 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
922 } /* End of the loop */
923 walltime_accounting_end_time(walltime_accounting);
927 if (fp_tpi != nullptr)
932 if (fplog != nullptr)
934 fprintf(fplog, "\n");
935 fprintf(fplog, " <V> = %12.5e nm^3\n", V_all / frame);
936 const double mu = -log(VembU_all / V_all) / beta;
937 fprintf(fplog, " <mu> = %12.5e kJ/mol\n", mu);
939 if (!std::isfinite(mu))
942 "\nThe computed chemical potential is not finite - consider increasing the "
943 "number of steps and/or the number of frames to insert into.\n");
947 /* Write the Boltzmann factor histogram */
950 /* When running in parallel sum the bins over the processes */
953 realloc_bins(&bin, &nbin, i);
954 gmx_sumd(nbin, bin, cr);
958 fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm), "TPI energy distribution",
959 "\\betaU - log(V/<V>)", "count", oenv);
960 sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
961 xvgr_subtitle(fp_tpi, str, oenv);
962 xvgr_legend(fp_tpi, 2, tpid_leg, oenv);
963 for (i = nbin - 1; i > 0; i--)
965 bUlogV = -i / invbinw + bU_logV_bin_limit - refvolshift + log(V_all / frame);
966 fprintf(fp_tpi, "%6.2f %10d %12.5e\n", bUlogV, roundToInt(bin[i]),
967 bin[i] * exp(-bUlogV) * V_all / VembU_all);
975 walltime_accounting_set_nsteps_done(walltime_accounting, frame * inputrec->nsteps);