4c808700fe5724721552dd7d96de899176ed21ff
[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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /*! \internal \file
39  *
40  * \brief This file defines the integrator for test particle insertion
41  *
42  * \author Berk Hess <hess@kth.se>
43  * \ingroup module_mdrun
44  */
45 #include "gmxpre.h"
46
47 #include <cmath>
48 #include <cstdlib>
49 #include <cstring>
50 #include <ctime>
51
52 #include <algorithm>
53
54 #include <cfenv>
55
56 #include "gromacs/commandline/filenm.h"
57 #include "gromacs/domdec/dlbtiming.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/ewald/pme.h"
60 #include "gromacs/fileio/confio.h"
61 #include "gromacs/fileio/trxio.h"
62 #include "gromacs/fileio/xvgr.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/gmx_omp_nthreads.h"
74 #include "gromacs/mdlib/mdatoms.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/nbnxm/nbnxm.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"
99
100 #include "legacysimulator.h"
101
102 //! Global max algorithm
103 static void global_max(t_commrec* cr, int* n)
104 {
105     int *sum, i;
106
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++)
111     {
112         *n = std::max(*n, sum[i]);
113     }
114
115     sfree(sum);
116 }
117
118 //! Reallocate arrays.
119 static void realloc_bins(double** bin, int* nbin, int nbin_new)
120 {
121     int i;
122
123     if (nbin_new != *nbin)
124     {
125         srenew(*bin, nbin_new);
126         for (i = *nbin; i < nbin_new; i++)
127         {
128             (*bin)[i] = 0;
129         }
130         *nbin = nbin_new;
131     }
132 }
133
134 //! Computes and returns the RF exclusion energy for the last molecule starting at \p beginAtom
135 static real reactionFieldExclusionCorrection(gmx::ArrayRef<const gmx::RVec> x,
136                                              const t_mdatoms&               mdatoms,
137                                              const interaction_const_t&     ic,
138                                              const int                      beginAtom)
139 {
140     real energy = 0;
141
142     for (int i = beginAtom; i < mdatoms.homenr; i++)
143     {
144         const real qi = mdatoms.chargeA[i];
145         energy -= 0.5 * qi * qi * ic.c_rf;
146
147         for (int j = i + 1; j < mdatoms.homenr; j++)
148         {
149             const real qj  = mdatoms.chargeA[j];
150             const real rsq = distance2(x[i], x[j]);
151             energy += qi * qj * (ic.k_rf * rsq - ic.c_rf);
152         }
153     }
154
155     return ic.epsfac * energy;
156 }
157
158 namespace gmx
159 {
160
161 // TODO: Convert to use the nbnxm kernels by putting the system and the teset molecule on two separate search grids
162 void LegacySimulator::do_tpi()
163 {
164     GMX_RELEASE_ASSERT(gmx_omp_nthreads_get(emntDefault) == 1, "TPI does not support OpenMP");
165
166     gmx_localtop_t              top;
167     PaddedHostVector<gmx::RVec> f{};
168     real                        lambda, t, temp, beta, drmax, epot;
169     double                      embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
170     t_trxstatus*                status;
171     t_trxframe                  rerun_fr;
172     gmx_bool                    bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
173     tensor                      force_vir, shake_vir, vir, pres;
174     int                         a_tp0, a_tp1, ngid, gid_tp, nener, e;
175     rvec*                       x_mol;
176     rvec                        mu_tot, x_init, dx;
177     int                         nnodes, frame;
178     int64_t                     frame_step_prev, frame_step;
179     int64_t                     nsteps, stepblocksize = 0, step;
180     int64_t                     seed;
181     int                         i;
182     FILE*                       fp_tpi = nullptr;
183     char *                      ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
184     double                      dbl, dump_ener;
185     gmx_bool                    bCavity;
186     int                         nat_cavity  = 0, d;
187     real *                      mass_cavity = nullptr, mass_tot;
188     int                         nbin;
189     double                      invbinw, *bin, refvolshift, logV, bUlogV;
190     gmx_bool                    bEnergyOutOfBounds;
191     const char*                 tpid_leg[2] = { "direct", "reweighted" };
192     auto                        mdatoms     = mdAtoms->mdatoms();
193
194     GMX_UNUSED_VALUE(outputProvider);
195
196     GMX_LOG(mdlog.info)
197             .asParagraph()
198             .appendText(
199                     "Note that it is planned to change the command gmx mdrun -tpi "
200                     "(and -tpic) to make the functionality available in a different "
201                     "form in a future version of GROMACS, e.g. gmx test-particle-insertion.");
202
203     /* Since there is no upper limit to the insertion energies,
204      * we need to set an upper limit for the distribution output.
205      */
206     real bU_bin_limit      = 50;
207     real bU_logV_bin_limit = bU_bin_limit + 10;
208
209     nnodes = cr->nnodes;
210
211     gmx_mtop_generate_local_top(*top_global, &top, inputrec->efep != efepNO);
212
213     SimulationGroups* groups = &top_global->groups;
214
215     bCavity = (inputrec->eI == eiTPIC);
216     if (bCavity)
217     {
218         ptr = getenv("GMX_TPIC_MASSES");
219         if (ptr == nullptr)
220         {
221             nat_cavity = 1;
222         }
223         else
224         {
225             /* Read (multiple) masses from env var GMX_TPIC_MASSES,
226              * The center of mass of the last atoms is then used for TPIC.
227              */
228             nat_cavity = 0;
229             while (sscanf(ptr, "%20lf%n", &dbl, &i) > 0)
230             {
231                 srenew(mass_cavity, nat_cavity + 1);
232                 mass_cavity[nat_cavity] = dbl;
233                 fprintf(fplog, "mass[%d] = %f\n", nat_cavity + 1, mass_cavity[nat_cavity]);
234                 nat_cavity++;
235                 ptr += i;
236             }
237             if (nat_cavity == 0)
238             {
239                 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
240             }
241         }
242     }
243
244     /*
245        init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
246        state_global->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
247     /* We never need full pbc for TPI */
248     fr->ePBC = epbcXYZ;
249     /* Determine the temperature for the Boltzmann weighting */
250     temp = inputrec->opts.ref_t[0];
251     if (fplog)
252     {
253         for (i = 1; (i < inputrec->opts.ngtc); i++)
254         {
255             if (inputrec->opts.ref_t[i] != temp)
256             {
257                 fprintf(fplog,
258                         "\nWARNING: The temperatures of the different temperature coupling groups "
259                         "are not identical\n\n");
260                 fprintf(stderr,
261                         "\nWARNING: The temperatures of the different temperature coupling groups "
262                         "are not identical\n\n");
263             }
264         }
265         fprintf(fplog, "\n  The temperature for test particle insertion is %.3f K\n\n", temp);
266     }
267     beta = 1.0 / (BOLTZ * temp);
268
269     /* Number of insertions per frame */
270     nsteps = inputrec->nsteps;
271
272     /* Use the same neighborlist with more insertions points
273      * in a sphere of radius drmax around the initial point
274      */
275     /* This should be a proper mdp parameter */
276     drmax = inputrec->rtpi;
277
278     /* An environment variable can be set to dump all configurations
279      * to pdb with an insertion energy <= this value.
280      */
281     dump_pdb  = getenv("GMX_TPI_DUMP");
282     dump_ener = 0;
283     if (dump_pdb)
284     {
285         sscanf(dump_pdb, "%20lf", &dump_ener);
286     }
287
288     atoms2md(top_global, inputrec, -1, nullptr, top_global->natoms, mdAtoms);
289     update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
290
291     f.resizeWithPadding(top_global->natoms);
292
293     /* Print to log file  */
294     walltime_accounting_start_time(walltime_accounting);
295     wallcycle_start(wcycle, ewcRUN);
296     print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
297
298     /* The last charge group is the group to be inserted */
299     const t_atoms& atomsToInsert = top_global->moltype[top_global->molblock.back().type].atoms;
300     a_tp0                        = top_global->natoms - atomsToInsert.nr;
301     a_tp1                        = top_global->natoms;
302     if (debug)
303     {
304         fprintf(debug, "TPI atoms %d-%d\n", a_tp0, a_tp1);
305     }
306
307     auto x = makeArrayRef(state_global->x);
308
309     if (EEL_PME(fr->ic->eeltype))
310     {
311         gmx_pme_reinit_atoms(fr->pmedata, a_tp0, nullptr);
312     }
313
314     /* With reacion-field we have distance dependent potentials
315      * between excluded atoms, we need to add these separately
316      * for the inserted molecule.
317      */
318     real rfExclusionEnergy = 0;
319     if (EEL_RF(fr->ic->eeltype))
320     {
321         rfExclusionEnergy = reactionFieldExclusionCorrection(x, *mdatoms, *fr->ic, a_tp0);
322         if (debug)
323         {
324             fprintf(debug, "RF exclusion correction for inserted molecule: %f kJ/mol\n", rfExclusionEnergy);
325         }
326     }
327
328     snew(x_mol, a_tp1 - a_tp0);
329
330     bDispCorr = (inputrec->eDispCorr != edispcNO);
331     bCharge   = FALSE;
332     for (i = a_tp0; i < a_tp1; i++)
333     {
334         /* Copy the coordinates of the molecule to be insterted */
335         copy_rvec(x[i], x_mol[i - a_tp0]);
336         /* Check if we need to print electrostatic energies */
337         bCharge |= (mdatoms->chargeA[i] != 0
338                     || ((mdatoms->chargeB != nullptr) && mdatoms->chargeB[i] != 0));
339     }
340     bRFExcl = (bCharge && EEL_RF(fr->ic->eeltype));
341
342     // Calculate the center of geometry of the molecule to insert
343     rvec cog = { 0, 0, 0 };
344     for (int a = a_tp0; a < a_tp1; a++)
345     {
346         rvec_inc(cog, x[a]);
347     }
348     svmul(1.0_real / (a_tp1 - a_tp0), cog, cog);
349     real molRadius = 0;
350     for (int a = a_tp0; a < a_tp1; a++)
351     {
352         molRadius = std::max(molRadius, distance2(x[a], cog));
353     }
354     molRadius = std::sqrt(molRadius);
355
356     const real maxCutoff = std::max(inputrec->rvdw, inputrec->rcoulomb);
357     if (bCavity)
358     {
359         if (norm(cog) > 0.5 * maxCutoff && fplog)
360         {
361             fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
362             fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
363         }
364     }
365     else
366     {
367         /* Center the molecule to be inserted at zero */
368         for (i = 0; i < a_tp1 - a_tp0; i++)
369         {
370             rvec_dec(x_mol[i], cog);
371         }
372     }
373
374     if (fplog)
375     {
376         fprintf(fplog, "\nWill insert %d atoms %s partial charges\n", a_tp1 - a_tp0,
377                 bCharge ? "with" : "without");
378
379         fprintf(fplog, "\nWill insert %" PRId64 " times in each frame of %s\n", nsteps,
380                 opt2fn("-rerun", nfile, fnm));
381     }
382
383     if (!bCavity)
384     {
385         if (inputrec->nstlist > 1)
386         {
387
388             /* With the same pair list we insert in a sphere of radius rtpi  in different orientations */
389             if (drmax == 0 && a_tp1 - a_tp0 == 1)
390             {
391                 gmx_fatal(FARGS,
392                           "Re-using the neighborlist %d times for insertions of a single atom in a "
393                           "sphere of radius %f does not make sense",
394                           inputrec->nstlist, drmax);
395             }
396             if (fplog)
397             {
398                 fprintf(fplog,
399                         "Will use the same neighborlist for %d insertions in a sphere of radius "
400                         "%f\n",
401                         inputrec->nstlist, drmax);
402             }
403         }
404     }
405     else
406     {
407         if (fplog)
408         {
409             fprintf(fplog,
410                     "Will insert randomly in a sphere of radius %f around the center of the "
411                     "cavity\n",
412                     drmax);
413         }
414     }
415
416     /* With the same pair list we insert in a sphere of radius rtpi
417      * in different orientations. We generate the pairlist with all
418      * inserted atoms located in the center of the sphere, so we need
419      * a buffer of size of the sphere and molecule radius.
420      */
421     inputrec->rlist = maxCutoff + 2 * inputrec->rtpi + 2 * molRadius;
422     fr->rlist       = inputrec->rlist;
423     fr->nbv->changePairlistRadii(inputrec->rlist, inputrec->rlist);
424
425     ngid   = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
426     gid_tp = GET_CGINFO_GID(fr->cginfo[a_tp0]);
427     for (int a = a_tp0 + 1; a < a_tp1; a++)
428     {
429         if (GET_CGINFO_GID(fr->cginfo[a]) != gid_tp)
430         {
431             fprintf(fplog,
432                     "NOTE: Atoms in the molecule to insert belong to different energy groups.\n"
433                     "      Only contributions to the group of the first atom will be reported.\n");
434             break;
435         }
436     }
437     nener = 1 + ngid;
438     if (bDispCorr)
439     {
440         nener += 1;
441     }
442     if (bCharge)
443     {
444         nener += ngid;
445         if (bRFExcl)
446         {
447             nener += 1;
448         }
449         if (EEL_FULL(fr->ic->eeltype))
450         {
451             nener += 1;
452         }
453     }
454     snew(sum_UgembU, nener);
455
456     /* Copy the random seed set by the user */
457     seed = inputrec->ld_seed;
458
459     gmx::ThreeFry2x64<16> rng(
460             seed, gmx::RandomDomain::TestParticleInsertion); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
461     gmx::UniformRealDistribution<real> dist;
462
463     if (MASTER(cr))
464     {
465         fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm), "TPI energies", "Time (ps)",
466                           "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
467         xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
468         snew(leg, 4 + nener);
469         e = 0;
470         sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
471         leg[e++] = gmx_strdup(str);
472         sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
473         leg[e++] = gmx_strdup(str);
474         sprintf(str, "f. <e\\S-\\betaU\\N>");
475         leg[e++] = gmx_strdup(str);
476         sprintf(str, "f. V");
477         leg[e++] = gmx_strdup(str);
478         sprintf(str, "f. <Ue\\S-\\betaU\\N>");
479         leg[e++] = gmx_strdup(str);
480         for (i = 0; i < ngid; i++)
481         {
482             sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
483                     *(groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput][i]]));
484             leg[e++] = gmx_strdup(str);
485         }
486         if (bDispCorr)
487         {
488             sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
489             leg[e++] = gmx_strdup(str);
490         }
491         if (bCharge)
492         {
493             for (i = 0; i < ngid; i++)
494             {
495                 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
496                         *(groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput][i]]));
497                 leg[e++] = gmx_strdup(str);
498             }
499             if (bRFExcl)
500             {
501                 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
502                 leg[e++] = gmx_strdup(str);
503             }
504             if (EEL_FULL(fr->ic->eeltype))
505             {
506                 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
507                 leg[e++] = gmx_strdup(str);
508             }
509         }
510         xvgr_legend(fp_tpi, 4 + nener, leg, oenv);
511         for (i = 0; i < 4 + nener; i++)
512         {
513             sfree(leg[i]);
514         }
515         sfree(leg);
516     }
517     clear_rvec(x_init);
518     V_all     = 0;
519     VembU_all = 0;
520
521     invbinw = 10;
522     nbin    = 10;
523     snew(bin, nbin);
524
525     /* Avoid frame step numbers <= -1 */
526     frame_step_prev = -1;
527
528     bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
529     frame         = 0;
530
531     if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) != mdatoms->nr - (a_tp1 - a_tp0))
532     {
533         gmx_fatal(FARGS,
534                   "Number of atoms in trajectory (%d)%s "
535                   "is not equal the number in the run input file (%d) "
536                   "minus the number of atoms to insert (%d)\n",
537                   rerun_fr.natoms, bCavity ? " minus one" : "", mdatoms->nr, a_tp1 - a_tp0);
538     }
539
540     refvolshift = log(det(rerun_fr.box));
541
542     switch (inputrec->eI)
543     {
544         case eiTPI: stepblocksize = inputrec->nstlist; break;
545         case eiTPIC: stepblocksize = 1; break;
546         default: gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
547     }
548
549     while (bNotLastFrame)
550     {
551         frame_step = rerun_fr.step;
552         if (frame_step <= frame_step_prev)
553         {
554             /* We don't have step number in the trajectory file,
555              * or we have constant or decreasing step numbers.
556              * Ensure we have increasing step numbers, since we use
557              * the step numbers as a counter for random numbers.
558              */
559             frame_step = frame_step_prev + 1;
560         }
561         frame_step_prev = frame_step;
562
563         lambda = rerun_fr.lambda;
564         t      = rerun_fr.time;
565
566         sum_embU = 0;
567         for (e = 0; e < nener; e++)
568         {
569             sum_UgembU[e] = 0;
570         }
571
572         /* Copy the coordinates from the input trajectory */
573         auto x = makeArrayRef(state_global->x);
574         for (i = 0; i < rerun_fr.natoms; i++)
575         {
576             copy_rvec(rerun_fr.x[i], x[i]);
577         }
578         copy_mat(rerun_fr.box, state_global->box);
579         const matrix& box = state_global->box;
580
581         V    = det(box);
582         logV = log(V);
583
584         bStateChanged = TRUE;
585         bNS           = TRUE;
586
587         put_atoms_in_box(fr->ePBC, box, x);
588
589         /* Put all atoms except for the inserted ones on the grid */
590         rvec vzero       = { 0, 0, 0 };
591         rvec boxDiagonal = { box[XX][XX], box[YY][YY], box[ZZ][ZZ] };
592         nbnxn_put_on_grid(fr->nbv.get(), box, 0, vzero, boxDiagonal, nullptr, { 0, a_tp0 }, -1,
593                           fr->cginfo, x, 0, nullptr);
594
595         step = cr->nodeid * stepblocksize;
596         while (step < nsteps)
597         {
598             /* Restart random engine using the frame and insertion step
599              * as counters.
600              * Note that we need to draw several random values per iteration,
601              * but by using the internal subcounter functionality of ThreeFry2x64
602              * we can draw 131072 unique 64-bit values before exhausting
603              * the stream. This is a huge margin, and if something still goes
604              * wrong you will get an exception when the stream is exhausted.
605              */
606             rng.restart(frame_step, step);
607             dist.reset(); // erase any memory in the distribution
608
609             if (!bCavity)
610             {
611                 /* Random insertion in the whole volume */
612                 bNS = (step % inputrec->nstlist == 0);
613                 if (bNS)
614                 {
615                     /* Generate a random position in the box */
616                     for (d = 0; d < DIM; d++)
617                     {
618                         x_init[d] = dist(rng) * state_global->box[d][d];
619                     }
620                 }
621             }
622             else
623             {
624                 /* Random insertion around a cavity location
625                  * given by the last coordinate of the trajectory.
626                  */
627                 if (step == 0)
628                 {
629                     if (nat_cavity == 1)
630                     {
631                         /* Copy the location of the cavity */
632                         copy_rvec(rerun_fr.x[rerun_fr.natoms - 1], x_init);
633                     }
634                     else
635                     {
636                         /* Determine the center of mass of the last molecule */
637                         clear_rvec(x_init);
638                         mass_tot = 0;
639                         for (i = 0; i < nat_cavity; i++)
640                         {
641                             for (d = 0; d < DIM; d++)
642                             {
643                                 x_init[d] += mass_cavity[i]
644                                              * rerun_fr.x[rerun_fr.natoms - nat_cavity + i][d];
645                             }
646                             mass_tot += mass_cavity[i];
647                         }
648                         for (d = 0; d < DIM; d++)
649                         {
650                             x_init[d] /= mass_tot;
651                         }
652                     }
653                 }
654             }
655
656             if (bNS)
657             {
658                 for (int a = a_tp0; a < a_tp1; a++)
659                 {
660                     x[a] = x_init;
661                 }
662
663                 /* Put the inserted molecule on it's own search grid */
664                 nbnxn_put_on_grid(fr->nbv.get(), box, 1, x_init, x_init, nullptr, { a_tp0, a_tp1 },
665                                   -1, fr->cginfo, x, 0, nullptr);
666
667                 /* TODO: Avoid updating all atoms at every bNS step */
668                 fr->nbv->setAtomProperties(*mdatoms, fr->cginfo);
669
670                 fr->nbv->constructPairlist(InteractionLocality::Local, top.excls, step, nrnb);
671
672                 bNS = FALSE;
673             }
674
675             /* Add random displacement uniformly distributed in a sphere
676              * of radius rtpi. We don't need to do this is we generate
677              * a new center location every step.
678              */
679             rvec x_tp;
680             if (bCavity || inputrec->nstlist > 1)
681             {
682                 /* Generate coordinates within |dx|=drmax of x_init */
683                 do
684                 {
685                     for (d = 0; d < DIM; d++)
686                     {
687                         dx[d] = (2 * dist(rng) - 1) * drmax;
688                     }
689                 } while (norm2(dx) > drmax * drmax);
690                 rvec_add(x_init, dx, x_tp);
691             }
692             else
693             {
694                 copy_rvec(x_init, x_tp);
695             }
696
697             if (a_tp1 - a_tp0 == 1)
698             {
699                 /* Insert a single atom, just copy the insertion location */
700                 copy_rvec(x_tp, x[a_tp0]);
701             }
702             else
703             {
704                 /* Copy the coordinates from the top file */
705                 for (i = a_tp0; i < a_tp1; i++)
706                 {
707                     copy_rvec(x_mol[i - a_tp0], x[i]);
708                 }
709                 /* Rotate the molecule randomly */
710                 real angleX = 2 * M_PI * dist(rng);
711                 real angleY = 2 * M_PI * dist(rng);
712                 real angleZ = 2 * M_PI * dist(rng);
713                 rotate_conf(a_tp1 - a_tp0, state_global->x.rvec_array() + a_tp0, nullptr, angleX,
714                             angleY, angleZ);
715                 /* Shift to the insertion location */
716                 for (i = a_tp0; i < a_tp1; i++)
717                 {
718                     rvec_inc(x[i], x_tp);
719                 }
720             }
721
722             /* Note: NonLocal refers to the inserted molecule */
723             fr->nbv->convertCoordinates(AtomLocality::NonLocal, false, x);
724
725             /* Clear some matrix variables  */
726             clear_mat(force_vir);
727             clear_mat(shake_vir);
728             clear_mat(vir);
729             clear_mat(pres);
730
731             /* Calc energy (no forces) on new positions.
732              * Since we only need the intermolecular energy
733              * and the RF exclusion terms of the inserted molecule occur
734              * within a single charge group we can pass NULL for the graph.
735              * This also avoids shifts that would move charge groups
736              * out of the box. */
737             /* Make do_force do a single node force calculation */
738             cr->nnodes = 1;
739
740             // TPI might place a particle so close that the potential
741             // is infinite. Since this is intended to happen, we
742             // temporarily suppress any exceptions that the processor
743             // might raise, then restore the old behaviour.
744             std::fenv_t floatingPointEnvironment;
745             std::feholdexcept(&floatingPointEnvironment);
746             do_force(fplog, cr, ms, inputrec, nullptr, nullptr, imdSession, pull_work, step, nrnb,
747                      wcycle, &top, state_global->box, state_global->x.arrayRefWithPadding(),
748                      &state_global->hist, f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
749                      state_global->lambda, nullptr, fr, runScheduleWork, nullptr, mu_tot, t, nullptr,
750                      GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY | (bStateChanged ? GMX_FORCE_STATECHANGED : 0),
751                      DDBalanceRegionHandler(nullptr));
752             std::feclearexcept(FE_DIVBYZERO | FE_INVALID | FE_OVERFLOW);
753             std::feupdateenv(&floatingPointEnvironment);
754
755             cr->nnodes    = nnodes;
756             bStateChanged = FALSE;
757
758             if (fr->dispersionCorrection)
759             {
760                 /* Calculate long range corrections to pressure and energy */
761                 const DispersionCorrection::Correction correction =
762                         fr->dispersionCorrection->calculate(state_global->box, lambda);
763                 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
764                 enerd->term[F_DISPCORR] = correction.energy;
765                 enerd->term[F_EPOT] += correction.energy;
766                 enerd->term[F_PRES] += correction.pressure;
767                 enerd->term[F_DVDL] += correction.dvdl;
768             }
769             else
770             {
771                 enerd->term[F_DISPCORR] = 0;
772             }
773             if (EEL_RF(fr->ic->eeltype))
774             {
775                 enerd->term[F_EPOT] += rfExclusionEnergy;
776             }
777
778             epot               = enerd->term[F_EPOT];
779             bEnergyOutOfBounds = FALSE;
780
781             /* If the compiler doesn't optimize this check away
782              * we catch the NAN energies.
783              * The epot>GMX_REAL_MAX check catches inf values,
784              * which should nicely result in embU=0 through the exp below,
785              * but it does not hurt to check anyhow.
786              */
787             /* Non-bonded Interaction usually diverge at r=0.
788              * With tabulated interaction functions the first few entries
789              * should be capped in a consistent fashion between
790              * repulsion, dispersion and Coulomb to avoid accidental
791              * negative values in the total energy.
792              * The table generation code in tables.c does this.
793              * With user tbales the user should take care of this.
794              */
795             if (epot != epot || epot > GMX_REAL_MAX)
796             {
797                 bEnergyOutOfBounds = TRUE;
798             }
799             if (bEnergyOutOfBounds)
800             {
801                 if (debug)
802                 {
803                     fprintf(debug,
804                             "\n  time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t,
805                             static_cast<int>(step), epot);
806                 }
807                 embU = 0;
808             }
809             else
810             {
811                 // Exponent argument is fine in SP range, but output can be in DP range
812                 embU = exp(static_cast<double>(-beta * epot));
813                 sum_embU += embU;
814                 /* Determine the weighted energy contributions of each energy group */
815                 e = 0;
816                 sum_UgembU[e++] += epot * embU;
817                 if (fr->bBHAM)
818                 {
819                     for (i = 0; i < ngid; i++)
820                     {
821                         sum_UgembU[e++] += enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)] * embU;
822                     }
823                 }
824                 else
825                 {
826                     for (i = 0; i < ngid; i++)
827                     {
828                         sum_UgembU[e++] += enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)] * embU;
829                     }
830                 }
831                 if (bDispCorr)
832                 {
833                     sum_UgembU[e++] += enerd->term[F_DISPCORR] * embU;
834                 }
835                 if (bCharge)
836                 {
837                     for (i = 0; i < ngid; i++)
838                     {
839                         sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
840                     }
841                     if (bRFExcl)
842                     {
843                         sum_UgembU[e++] += rfExclusionEnergy * embU;
844                     }
845                     if (EEL_FULL(fr->ic->eeltype))
846                     {
847                         sum_UgembU[e++] += enerd->term[F_COUL_RECIP] * embU;
848                     }
849                 }
850             }
851
852             if (embU == 0 || beta * epot > bU_bin_limit)
853             {
854                 bin[0]++;
855             }
856             else
857             {
858                 i = gmx::roundToInt((bU_logV_bin_limit - (beta * epot - logV + refvolshift)) * invbinw);
859                 if (i < 0)
860                 {
861                     i = 0;
862                 }
863                 if (i >= nbin)
864                 {
865                     realloc_bins(&bin, &nbin, i + 10);
866                 }
867                 bin[i]++;
868             }
869
870             if (debug)
871             {
872                 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n", static_cast<int>(step),
873                         epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
874             }
875
876             if (dump_pdb && epot <= dump_ener)
877             {
878                 sprintf(str, "t%g_step%d.pdb", t, static_cast<int>(step));
879                 sprintf(str2, "t: %f step %d ener: %f", t, static_cast<int>(step), epot);
880                 write_sto_conf_mtop(str, str2, top_global, state_global->x.rvec_array(),
881                                     state_global->v.rvec_array(), inputrec->ePBC, state_global->box);
882             }
883
884             step++;
885             if ((step / stepblocksize) % cr->nnodes != cr->nodeid)
886             {
887                 /* Skip all steps assigned to the other MPI ranks */
888                 step += (cr->nnodes - 1) * stepblocksize;
889             }
890         }
891
892         if (PAR(cr))
893         {
894             /* When running in parallel sum the energies over the processes */
895             gmx_sumd(1, &sum_embU, cr);
896             gmx_sumd(nener, sum_UgembU, cr);
897         }
898
899         frame++;
900         V_all += V;
901         VembU_all += V * sum_embU / nsteps;
902
903         if (fp_tpi)
904         {
905             if (mdrunOptions.verbose || frame % 10 == 0 || frame < 10)
906             {
907                 fprintf(stderr, "mu %10.3e <mu> %10.3e\n", -log(sum_embU / nsteps) / beta,
908                         -log(VembU_all / V_all) / beta);
909             }
910
911             fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e", t,
912                     VembU_all == 0 ? 20 / beta : -log(VembU_all / V_all) / beta,
913                     sum_embU == 0 ? 20 / beta : -log(sum_embU / nsteps) / beta, sum_embU / nsteps, V);
914             for (e = 0; e < nener; e++)
915             {
916                 fprintf(fp_tpi, " %12.5e", sum_UgembU[e] / nsteps);
917             }
918             fprintf(fp_tpi, "\n");
919             fflush(fp_tpi);
920         }
921
922         bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
923     } /* End of the loop  */
924     walltime_accounting_end_time(walltime_accounting);
925
926     close_trx(status);
927
928     if (fp_tpi != nullptr)
929     {
930         xvgrclose(fp_tpi);
931     }
932
933     if (fplog != nullptr)
934     {
935         fprintf(fplog, "\n");
936         fprintf(fplog, "  <V>  = %12.5e nm^3\n", V_all / frame);
937         const double mu = -log(VembU_all / V_all) / beta;
938         fprintf(fplog, "  <mu> = %12.5e kJ/mol\n", mu);
939
940         if (!std::isfinite(mu))
941         {
942             fprintf(fplog,
943                     "\nThe computed chemical potential is not finite - consider increasing the "
944                     "number of steps and/or the number of frames to insert into.\n");
945         }
946     }
947
948     /* Write the Boltzmann factor histogram */
949     if (PAR(cr))
950     {
951         /* When running in parallel sum the bins over the processes */
952         i = nbin;
953         global_max(cr, &i);
954         realloc_bins(&bin, &nbin, i);
955         gmx_sumd(nbin, bin, cr);
956     }
957     if (MASTER(cr))
958     {
959         fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm), "TPI energy distribution",
960                           "\\betaU - log(V/<V>)", "count", oenv);
961         sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
962         xvgr_subtitle(fp_tpi, str, oenv);
963         xvgr_legend(fp_tpi, 2, tpid_leg, oenv);
964         for (i = nbin - 1; i > 0; i--)
965         {
966             bUlogV = -i / invbinw + bU_logV_bin_limit - refvolshift + log(V_all / frame);
967             fprintf(fp_tpi, "%6.2f %10d %12.5e\n", bUlogV, roundToInt(bin[i]),
968                     bin[i] * exp(-bUlogV) * V_all / VembU_all);
969         }
970         xvgrclose(fp_tpi);
971     }
972     sfree(bin);
973
974     sfree(sum_UgembU);
975
976     walltime_accounting_set_nsteps_done(walltime_accounting, frame * inputrec->nsteps);
977 }
978
979 } // namespace gmx