Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / mdrun / tpi.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /*! \internal \file
38  *
39  * \brief This file defines the integrator for test particle insertion
40  *
41  * \author Berk Hess <hess@kth.se>
42  * \ingroup module_mdrun
43  */
44 #include "gmxpre.h"
45
46 #include <cmath>
47 #include <cstdlib>
48 #include <cstring>
49 #include <ctime>
50
51 #include <algorithm>
52
53 #include <cfenv>
54
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"
98
99 #include "legacysimulator.h"
100
101 //! Global max algorithm
102 static void global_max(t_commrec* cr, int* n)
103 {
104     int *sum, i;
105
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++)
110     {
111         *n = std::max(*n, sum[i]);
112     }
113
114     sfree(sum);
115 }
116
117 //! Reallocate arrays.
118 static void realloc_bins(double** bin, int* nbin, int nbin_new)
119 {
120     int i;
121
122     if (nbin_new != *nbin)
123     {
124         srenew(*bin, nbin_new);
125         for (i = *nbin; i < nbin_new; i++)
126         {
127             (*bin)[i] = 0;
128         }
129         *nbin = nbin_new;
130     }
131 }
132
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,
137                                              const int                      beginAtom)
138 {
139     real energy = 0;
140
141     for (int i = beginAtom; i < mdatoms.homenr; i++)
142     {
143         const real qi = mdatoms.chargeA[i];
144         energy -= 0.5 * qi * qi * ic.c_rf;
145
146         for (int j = i + 1; j < mdatoms.homenr; j++)
147         {
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);
151         }
152     }
153
154     return ic.epsfac * energy;
155 }
156
157 namespace gmx
158 {
159
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()
162 {
163     GMX_RELEASE_ASSERT(gmx_omp_nthreads_get(emntDefault) == 1, "TPI does not support OpenMP");
164
165     gmx_localtop_t              top;
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;
169     t_trxstatus*                status;
170     t_trxframe                  rerun_fr;
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;
174     rvec*                       x_mol;
175     rvec                        mu_tot, x_init, dx;
176     int                         nnodes, frame;
177     int64_t                     frame_step_prev, frame_step;
178     int64_t                     nsteps, stepblocksize = 0, step;
179     int64_t                     seed;
180     int                         i;
181     FILE*                       fp_tpi = nullptr;
182     char *                      ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
183     double                      dbl, dump_ener;
184     gmx_bool                    bCavity;
185     int                         nat_cavity  = 0, d;
186     real *                      mass_cavity = nullptr, mass_tot;
187     int                         nbin;
188     double                      invbinw, *bin, refvolshift, logV, bUlogV;
189     gmx_bool                    bEnergyOutOfBounds;
190     const char*                 tpid_leg[2] = { "direct", "reweighted" };
191     auto                        mdatoms     = mdAtoms->mdatoms();
192
193     GMX_UNUSED_VALUE(outputProvider);
194
195     GMX_LOG(mdlog.info)
196             .asParagraph()
197             .appendText(
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.");
201
202     /* Since there is no upper limit to the insertion energies,
203      * we need to set an upper limit for the distribution output.
204      */
205     real bU_bin_limit      = 50;
206     real bU_logV_bin_limit = bU_bin_limit + 10;
207
208     nnodes = cr->nnodes;
209
210     gmx_mtop_generate_local_top(*top_global, &top, inputrec->efep != efepNO);
211
212     SimulationGroups* groups = &top_global->groups;
213
214     bCavity = (inputrec->eI == eiTPIC);
215     if (bCavity)
216     {
217         ptr = getenv("GMX_TPIC_MASSES");
218         if (ptr == nullptr)
219         {
220             nat_cavity = 1;
221         }
222         else
223         {
224             /* Read (multiple) masses from env var GMX_TPIC_MASSES,
225              * The center of mass of the last atoms is then used for TPIC.
226              */
227             nat_cavity = 0;
228             while (sscanf(ptr, "%20lf%n", &dbl, &i) > 0)
229             {
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]);
233                 nat_cavity++;
234                 ptr += i;
235             }
236             if (nat_cavity == 0)
237             {
238                 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
239             }
240         }
241     }
242
243     /*
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 */
247     fr->ePBC = epbcXYZ;
248     /* Determine the temperature for the Boltzmann weighting */
249     temp = inputrec->opts.ref_t[0];
250     if (fplog)
251     {
252         for (i = 1; (i < inputrec->opts.ngtc); i++)
253         {
254             if (inputrec->opts.ref_t[i] != temp)
255             {
256                 fprintf(fplog,
257                         "\nWARNING: The temperatures of the different temperature coupling groups "
258                         "are not identical\n\n");
259                 fprintf(stderr,
260                         "\nWARNING: The temperatures of the different temperature coupling groups "
261                         "are not identical\n\n");
262             }
263         }
264         fprintf(fplog, "\n  The temperature for test particle insertion is %.3f K\n\n", temp);
265     }
266     beta = 1.0 / (BOLTZ * temp);
267
268     /* Number of insertions per frame */
269     nsteps = inputrec->nsteps;
270
271     /* Use the same neighborlist with more insertions points
272      * in a sphere of radius drmax around the initial point
273      */
274     /* This should be a proper mdp parameter */
275     drmax = inputrec->rtpi;
276
277     /* An environment variable can be set to dump all configurations
278      * to pdb with an insertion energy <= this value.
279      */
280     dump_pdb  = getenv("GMX_TPI_DUMP");
281     dump_ener = 0;
282     if (dump_pdb)
283     {
284         sscanf(dump_pdb, "%20lf", &dump_ener);
285     }
286
287     atoms2md(top_global, inputrec, -1, nullptr, top_global->natoms, mdAtoms);
288     update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
289
290     f.resizeWithPadding(top_global->natoms);
291
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");
296
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;
301     if (debug)
302     {
303         fprintf(debug, "TPI atoms %d-%d\n", a_tp0, a_tp1);
304     }
305
306     auto x = makeArrayRef(state_global->x);
307
308     if (EEL_PME(fr->ic->eeltype))
309     {
310         gmx_pme_reinit_atoms(fr->pmedata, a_tp0, nullptr);
311     }
312
313     /* With reacion-field we have distance dependent potentials
314      * between excluded atoms, we need to add these separately
315      * for the inserted molecule.
316      */
317     real rfExclusionEnergy = 0;
318     if (EEL_RF(fr->ic->eeltype))
319     {
320         rfExclusionEnergy = reactionFieldExclusionCorrection(x, *mdatoms, *fr->ic, a_tp0);
321         if (debug)
322         {
323             fprintf(debug, "RF exclusion correction for inserted molecule: %f kJ/mol\n", rfExclusionEnergy);
324         }
325     }
326
327     snew(x_mol, a_tp1 - a_tp0);
328
329     bDispCorr = (inputrec->eDispCorr != edispcNO);
330     bCharge   = FALSE;
331     for (i = a_tp0; i < a_tp1; i++)
332     {
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));
338     }
339     bRFExcl = (bCharge && EEL_RF(fr->ic->eeltype));
340
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++)
344     {
345         rvec_inc(cog, x[a]);
346     }
347     svmul(1.0_real / (a_tp1 - a_tp0), cog, cog);
348     real molRadius = 0;
349     for (int a = a_tp0; a < a_tp1; a++)
350     {
351         molRadius = std::max(molRadius, distance2(x[a], cog));
352     }
353     molRadius = std::sqrt(molRadius);
354
355     const real maxCutoff = std::max(inputrec->rvdw, inputrec->rcoulomb);
356     if (bCavity)
357     {
358         if (norm(cog) > 0.5 * maxCutoff && fplog)
359         {
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");
362         }
363     }
364     else
365     {
366         /* Center the molecule to be inserted at zero */
367         for (i = 0; i < a_tp1 - a_tp0; i++)
368         {
369             rvec_dec(x_mol[i], cog);
370         }
371     }
372
373     if (fplog)
374     {
375         fprintf(fplog, "\nWill insert %d atoms %s partial charges\n", a_tp1 - a_tp0,
376                 bCharge ? "with" : "without");
377
378         fprintf(fplog, "\nWill insert %" PRId64 " times in each frame of %s\n", nsteps,
379                 opt2fn("-rerun", nfile, fnm));
380     }
381
382     if (!bCavity)
383     {
384         if (inputrec->nstlist > 1)
385         {
386
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)
389             {
390                 gmx_fatal(FARGS,
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);
394             }
395             if (fplog)
396             {
397                 fprintf(fplog,
398                         "Will use the same neighborlist for %d insertions in a sphere of radius "
399                         "%f\n",
400                         inputrec->nstlist, drmax);
401             }
402         }
403     }
404     else
405     {
406         if (fplog)
407         {
408             fprintf(fplog,
409                     "Will insert randomly in a sphere of radius %f around the center of the "
410                     "cavity\n",
411                     drmax);
412         }
413     }
414
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.
419      */
420     inputrec->rlist = maxCutoff + 2 * inputrec->rtpi + 2 * molRadius;
421     fr->rlist       = inputrec->rlist;
422     fr->nbv->changePairlistRadii(inputrec->rlist, inputrec->rlist);
423
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++)
427     {
428         if (GET_CGINFO_GID(fr->cginfo[a]) != gid_tp)
429         {
430             fprintf(fplog,
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");
433             break;
434         }
435     }
436     nener = 1 + ngid;
437     if (bDispCorr)
438     {
439         nener += 1;
440     }
441     if (bCharge)
442     {
443         nener += ngid;
444         if (bRFExcl)
445         {
446             nener += 1;
447         }
448         if (EEL_FULL(fr->ic->eeltype))
449         {
450             nener += 1;
451         }
452     }
453     snew(sum_UgembU, nener);
454
455     /* Copy the random seed set by the user */
456     seed = inputrec->ld_seed;
457
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;
461
462     if (MASTER(cr))
463     {
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);
468         e = 0;
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++)
480         {
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);
484         }
485         if (bDispCorr)
486         {
487             sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
488             leg[e++] = gmx_strdup(str);
489         }
490         if (bCharge)
491         {
492             for (i = 0; i < ngid; i++)
493             {
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);
497             }
498             if (bRFExcl)
499             {
500                 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
501                 leg[e++] = gmx_strdup(str);
502             }
503             if (EEL_FULL(fr->ic->eeltype))
504             {
505                 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
506                 leg[e++] = gmx_strdup(str);
507             }
508         }
509         xvgr_legend(fp_tpi, 4 + nener, leg, oenv);
510         for (i = 0; i < 4 + nener; i++)
511         {
512             sfree(leg[i]);
513         }
514         sfree(leg);
515     }
516     clear_rvec(x_init);
517     V_all     = 0;
518     VembU_all = 0;
519
520     invbinw = 10;
521     nbin    = 10;
522     snew(bin, nbin);
523
524     /* Avoid frame step numbers <= -1 */
525     frame_step_prev = -1;
526
527     bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
528     frame         = 0;
529
530     if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) != mdatoms->nr - (a_tp1 - a_tp0))
531     {
532         gmx_fatal(FARGS,
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);
537     }
538
539     refvolshift = log(det(rerun_fr.box));
540
541     switch (inputrec->eI)
542     {
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]);
546     }
547
548     while (bNotLastFrame)
549     {
550         frame_step = rerun_fr.step;
551         if (frame_step <= frame_step_prev)
552         {
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.
557              */
558             frame_step = frame_step_prev + 1;
559         }
560         frame_step_prev = frame_step;
561
562         lambda = rerun_fr.lambda;
563         t      = rerun_fr.time;
564
565         sum_embU = 0;
566         for (e = 0; e < nener; e++)
567         {
568             sum_UgembU[e] = 0;
569         }
570
571         /* Copy the coordinates from the input trajectory */
572         auto x = makeArrayRef(state_global->x);
573         for (i = 0; i < rerun_fr.natoms; i++)
574         {
575             copy_rvec(rerun_fr.x[i], x[i]);
576         }
577         copy_mat(rerun_fr.box, state_global->box);
578         const matrix& box = state_global->box;
579
580         V    = det(box);
581         logV = log(V);
582
583         bStateChanged = TRUE;
584         bNS           = TRUE;
585
586         put_atoms_in_box(fr->ePBC, box, x);
587
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);
593
594         step = cr->nodeid * stepblocksize;
595         while (step < nsteps)
596         {
597             /* Restart random engine using the frame and insertion step
598              * as counters.
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.
604              */
605             rng.restart(frame_step, step);
606             dist.reset(); // erase any memory in the distribution
607
608             if (!bCavity)
609             {
610                 /* Random insertion in the whole volume */
611                 bNS = (step % inputrec->nstlist == 0);
612                 if (bNS)
613                 {
614                     /* Generate a random position in the box */
615                     for (d = 0; d < DIM; d++)
616                     {
617                         x_init[d] = dist(rng) * state_global->box[d][d];
618                     }
619                 }
620             }
621             else
622             {
623                 /* Random insertion around a cavity location
624                  * given by the last coordinate of the trajectory.
625                  */
626                 if (step == 0)
627                 {
628                     if (nat_cavity == 1)
629                     {
630                         /* Copy the location of the cavity */
631                         copy_rvec(rerun_fr.x[rerun_fr.natoms - 1], x_init);
632                     }
633                     else
634                     {
635                         /* Determine the center of mass of the last molecule */
636                         clear_rvec(x_init);
637                         mass_tot = 0;
638                         for (i = 0; i < nat_cavity; i++)
639                         {
640                             for (d = 0; d < DIM; d++)
641                             {
642                                 x_init[d] += mass_cavity[i]
643                                              * rerun_fr.x[rerun_fr.natoms - nat_cavity + i][d];
644                             }
645                             mass_tot += mass_cavity[i];
646                         }
647                         for (d = 0; d < DIM; d++)
648                         {
649                             x_init[d] /= mass_tot;
650                         }
651                     }
652                 }
653             }
654
655             if (bNS)
656             {
657                 for (int a = a_tp0; a < a_tp1; a++)
658                 {
659                     x[a] = x_init;
660                 }
661
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);
665
666                 /* TODO: Avoid updating all atoms at every bNS step */
667                 fr->nbv->setAtomProperties(*mdatoms, fr->cginfo);
668
669                 fr->nbv->constructPairlist(InteractionLocality::Local, &top.excls, step, nrnb);
670
671                 bNS = FALSE;
672             }
673
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.
677              */
678             rvec x_tp;
679             if (bCavity || inputrec->nstlist > 1)
680             {
681                 /* Generate coordinates within |dx|=drmax of x_init */
682                 do
683                 {
684                     for (d = 0; d < DIM; d++)
685                     {
686                         dx[d] = (2 * dist(rng) - 1) * drmax;
687                     }
688                 } while (norm2(dx) > drmax * drmax);
689                 rvec_add(x_init, dx, x_tp);
690             }
691             else
692             {
693                 copy_rvec(x_init, x_tp);
694             }
695
696             if (a_tp1 - a_tp0 == 1)
697             {
698                 /* Insert a single atom, just copy the insertion location */
699                 copy_rvec(x_tp, x[a_tp0]);
700             }
701             else
702             {
703                 /* Copy the coordinates from the top file */
704                 for (i = a_tp0; i < a_tp1; i++)
705                 {
706                     copy_rvec(x_mol[i - a_tp0], x[i]);
707                 }
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,
713                             angleY, angleZ);
714                 /* Shift to the insertion location */
715                 for (i = a_tp0; i < a_tp1; i++)
716                 {
717                     rvec_inc(x[i], x_tp);
718                 }
719             }
720
721             /* Note: NonLocal refers to the inserted molecule */
722             fr->nbv->convertCoordinates(AtomLocality::NonLocal, false, x);
723
724             /* Clear some matrix variables  */
725             clear_mat(force_vir);
726             clear_mat(shake_vir);
727             clear_mat(vir);
728             clear_mat(pres);
729
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
735              * out of the box. */
736             /* Make do_force do a single node force calculation */
737             cr->nnodes = 1;
738
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);
753
754             cr->nnodes    = nnodes;
755             bStateChanged = FALSE;
756
757             if (fr->dispersionCorrection)
758             {
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;
767             }
768             else
769             {
770                 enerd->term[F_DISPCORR] = 0;
771             }
772             if (EEL_RF(fr->ic->eeltype))
773             {
774                 enerd->term[F_EPOT] += rfExclusionEnergy;
775             }
776
777             epot               = enerd->term[F_EPOT];
778             bEnergyOutOfBounds = FALSE;
779
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.
785              */
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.
793              */
794             if (epot != epot || epot > GMX_REAL_MAX)
795             {
796                 bEnergyOutOfBounds = TRUE;
797             }
798             if (bEnergyOutOfBounds)
799             {
800                 if (debug)
801                 {
802                     fprintf(debug,
803                             "\n  time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t,
804                             static_cast<int>(step), epot);
805                 }
806                 embU = 0;
807             }
808             else
809             {
810                 // Exponent argument is fine in SP range, but output can be in DP range
811                 embU = exp(static_cast<double>(-beta * epot));
812                 sum_embU += embU;
813                 /* Determine the weighted energy contributions of each energy group */
814                 e = 0;
815                 sum_UgembU[e++] += epot * embU;
816                 if (fr->bBHAM)
817                 {
818                     for (i = 0; i < ngid; i++)
819                     {
820                         sum_UgembU[e++] += enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)] * embU;
821                     }
822                 }
823                 else
824                 {
825                     for (i = 0; i < ngid; i++)
826                     {
827                         sum_UgembU[e++] += enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)] * embU;
828                     }
829                 }
830                 if (bDispCorr)
831                 {
832                     sum_UgembU[e++] += enerd->term[F_DISPCORR] * embU;
833                 }
834                 if (bCharge)
835                 {
836                     for (i = 0; i < ngid; i++)
837                     {
838                         sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
839                     }
840                     if (bRFExcl)
841                     {
842                         sum_UgembU[e++] += rfExclusionEnergy * embU;
843                     }
844                     if (EEL_FULL(fr->ic->eeltype))
845                     {
846                         sum_UgembU[e++] += enerd->term[F_COUL_RECIP] * embU;
847                     }
848                 }
849             }
850
851             if (embU == 0 || beta * epot > bU_bin_limit)
852             {
853                 bin[0]++;
854             }
855             else
856             {
857                 i = gmx::roundToInt((bU_logV_bin_limit - (beta * epot - logV + refvolshift)) * invbinw);
858                 if (i < 0)
859                 {
860                     i = 0;
861                 }
862                 if (i >= nbin)
863                 {
864                     realloc_bins(&bin, &nbin, i + 10);
865                 }
866                 bin[i]++;
867             }
868
869             if (debug)
870             {
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]);
873             }
874
875             if (dump_pdb && epot <= dump_ener)
876             {
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);
881             }
882
883             step++;
884             if ((step / stepblocksize) % cr->nnodes != cr->nodeid)
885             {
886                 /* Skip all steps assigned to the other MPI ranks */
887                 step += (cr->nnodes - 1) * stepblocksize;
888             }
889         }
890
891         if (PAR(cr))
892         {
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);
896         }
897
898         frame++;
899         V_all += V;
900         VembU_all += V * sum_embU / nsteps;
901
902         if (fp_tpi)
903         {
904             if (mdrunOptions.verbose || frame % 10 == 0 || frame < 10)
905             {
906                 fprintf(stderr, "mu %10.3e <mu> %10.3e\n", -log(sum_embU / nsteps) / beta,
907                         -log(VembU_all / V_all) / beta);
908             }
909
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++)
914             {
915                 fprintf(fp_tpi, " %12.5e", sum_UgembU[e] / nsteps);
916             }
917             fprintf(fp_tpi, "\n");
918             fflush(fp_tpi);
919         }
920
921         bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
922     } /* End of the loop  */
923     walltime_accounting_end_time(walltime_accounting);
924
925     close_trx(status);
926
927     if (fp_tpi != nullptr)
928     {
929         xvgrclose(fp_tpi);
930     }
931
932     if (fplog != nullptr)
933     {
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);
938
939         if (!std::isfinite(mu))
940         {
941             fprintf(fplog,
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");
944         }
945     }
946
947     /* Write the Boltzmann factor histogram */
948     if (PAR(cr))
949     {
950         /* When running in parallel sum the bins over the processes */
951         i = nbin;
952         global_max(cr, &i);
953         realloc_bins(&bin, &nbin, i);
954         gmx_sumd(nbin, bin, cr);
955     }
956     if (MASTER(cr))
957     {
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--)
964         {
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);
968         }
969         xvgrclose(fp_tpi);
970     }
971     sfree(bin);
972
973     sfree(sum_UgembU);
974
975     walltime_accounting_set_nsteps_done(walltime_accounting, frame * inputrec->nsteps);
976 }
977
978 } // namespace gmx