Use gmx::Range in Nbnxm gridding functions
[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
135 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
163 LegacySimulator::do_tpi()
164 {
165     GMX_RELEASE_ASSERT(gmx_omp_nthreads_get(emntDefault) == 1, "TPI does not support OpenMP");
166
167     gmx_localtop_t              top;
168     PaddedHostVector<gmx::RVec> f {};
169     real                        lambda, t, temp, beta, drmax, epot;
170     double                      embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
171     t_trxstatus                *status;
172     t_trxframe                  rerun_fr;
173     gmx_bool                    bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
174     tensor                      force_vir, shake_vir, vir, pres;
175     int                         a_tp0, a_tp1, ngid, gid_tp, nener, e;
176     rvec                       *x_mol;
177     rvec                        mu_tot, x_init, dx;
178     int                         nnodes, frame;
179     int64_t                     frame_step_prev, frame_step;
180     int64_t                     nsteps, stepblocksize = 0, step;
181     int64_t                     seed;
182     int                         i;
183     FILE                       *fp_tpi = nullptr;
184     char                       *ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
185     double                      dbl, dump_ener;
186     gmx_bool                    bCavity;
187     int                         nat_cavity  = 0, d;
188     real                       *mass_cavity = nullptr, mass_tot;
189     int                         nbin;
190     double                      invbinw, *bin, refvolshift, logV, bUlogV;
191     gmx_bool                    bEnergyOutOfBounds;
192     const char                 *tpid_leg[2] = {"direct", "reweighted"};
193     auto                        mdatoms     = mdAtoms->mdatoms();
194
195     GMX_UNUSED_VALUE(outputProvider);
196
197     GMX_LOG(mdlog.info).asParagraph().
198         appendText("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",
233                         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, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
258                 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
259             }
260         }
261         fprintf(fplog,
262                 "\n  The temperature for test particle insertion is %.3f K\n\n",
263                 temp);
264     }
265     beta = 1.0/(BOLTZ*temp);
266
267     /* Number of insertions per frame */
268     nsteps = inputrec->nsteps;
269
270     /* Use the same neighborlist with more insertions points
271      * in a sphere of radius drmax around the initial point
272      */
273     /* This should be a proper mdp parameter */
274     drmax = inputrec->rtpi;
275
276     /* An environment variable can be set to dump all configurations
277      * to pdb with an insertion energy <= this value.
278      */
279     dump_pdb  = getenv("GMX_TPI_DUMP");
280     dump_ener = 0;
281     if (dump_pdb)
282     {
283         sscanf(dump_pdb, "%20lf", &dump_ener);
284     }
285
286     atoms2md(top_global, inputrec, -1, nullptr, top_global->natoms, mdAtoms);
287     update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
288
289     f.resizeWithPadding(top_global->natoms);
290
291     /* Print to log file  */
292     walltime_accounting_start_time(walltime_accounting);
293     wallcycle_start(wcycle, ewcRUN);
294     print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
295
296     /* The last charge group is the group to be inserted */
297     const t_atoms &atomsToInsert = top_global->moltype[top_global->molblock.back().type].atoms;
298     a_tp0 = top_global->natoms - atomsToInsert.nr;
299     a_tp1 = top_global->natoms;
300     if (debug)
301     {
302         fprintf(debug, "TPI atoms %d-%d\n", a_tp0, a_tp1);
303     }
304
305     auto x = makeArrayRef(state_global->x);
306
307     if (EEL_PME(fr->ic->eeltype))
308     {
309         gmx_pme_reinit_atoms(fr->pmedata, a_tp0, nullptr);
310     }
311
312     /* With reacion-field we have distance dependent potentials
313      * between excluded atoms, we need to add these separately
314      * for the inserted molecule.
315      */
316     real rfExclusionEnergy = 0;
317     if (EEL_RF(fr->ic->eeltype))
318     {
319         rfExclusionEnergy = reactionFieldExclusionCorrection(x, *mdatoms, *fr->ic, a_tp0);
320         if (debug)
321         {
322             fprintf(debug, "RF exclusion correction for inserted molecule: %f kJ/mol\n", rfExclusionEnergy);
323         }
324     }
325
326     snew(x_mol, a_tp1-a_tp0);
327
328     bDispCorr = (inputrec->eDispCorr != edispcNO);
329     bCharge   = FALSE;
330     for (i = a_tp0; i < a_tp1; i++)
331     {
332         /* Copy the coordinates of the molecule to be insterted */
333         copy_rvec(x[i], x_mol[i-a_tp0]);
334         /* Check if we need to print electrostatic energies */
335         bCharge |= (mdatoms->chargeA[i] != 0 ||
336                     ((mdatoms->chargeB != nullptr) && mdatoms->chargeB[i] != 0));
337     }
338     bRFExcl = (bCharge && EEL_RF(fr->ic->eeltype));
339
340     // Calculate the center of geometry of the molecule to insert
341     rvec cog = { 0, 0, 0 };
342     for (int a = a_tp0; a < a_tp1; a++)
343     {
344         rvec_inc(cog, x[a]);
345     }
346     svmul(1.0_real/(a_tp1 - a_tp0), cog, cog);
347     real molRadius = 0;
348     for (int a = a_tp0; a < a_tp1; a++)
349     {
350         molRadius = std::max(molRadius, distance2(x[a], cog));
351     }
352     molRadius = std::sqrt(molRadius);
353
354     const real maxCutoff =  std::max(inputrec->rvdw, inputrec->rcoulomb);
355     if (bCavity)
356     {
357         if (norm(cog) > 0.5*maxCutoff && fplog)
358         {
359             fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
360             fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
361         }
362     }
363     else
364     {
365         /* Center the molecule to be inserted at zero */
366         for (i = 0; i < a_tp1-a_tp0; i++)
367         {
368             rvec_dec(x_mol[i], cog);
369         }
370     }
371
372     if (fplog)
373     {
374         fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
375                 a_tp1-a_tp0, bCharge ? "with" : "without");
376
377         fprintf(fplog, "\nWill insert %" PRId64 " times in each frame of %s\n",
378                 nsteps, opt2fn("-rerun", nfile, fnm));
379     }
380
381     if (!bCavity)
382     {
383         if (inputrec->nstlist > 1)
384         {
385
386             /* With the same pair list we insert in a sphere of radius rtpi  in different orientations */
387             if (drmax == 0 && a_tp1-a_tp0 == 1)
388             {
389                 gmx_fatal(FARGS, "Re-using the neighborlist %d times for insertions of a single atom in a sphere of radius %f does not make sense", inputrec->nstlist, drmax);
390             }
391             if (fplog)
392             {
393                 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
394             }
395         }
396     }
397     else
398     {
399         if (fplog)
400         {
401             fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
402         }
403     }
404
405     /* With the same pair list we insert in a sphere of radius rtpi
406      * in different orientations. We generate the pairlist with all
407      * inserted atoms located in the center of the sphere, so we need
408      * a buffer of size of the sphere and molecule radius.
409      */
410     inputrec->rlist = maxCutoff + 2*inputrec->rtpi + 2*molRadius;
411     fr->rlist       = inputrec->rlist;
412     fr->nbv->changePairlistRadii(inputrec->rlist, inputrec->rlist);
413
414     ngid   = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
415     gid_tp = GET_CGINFO_GID(fr->cginfo[a_tp0]);
416     for (int a = a_tp0 + 1; a < a_tp1; a++)
417     {
418         if (GET_CGINFO_GID(fr->cginfo[a]) != gid_tp)
419         {
420             fprintf(fplog,
421                     "NOTE: Atoms in the molecule to insert belong to different energy groups.\n"
422                     "      Only contributions to the group of the first atom will be reported.\n");
423             break;
424         }
425     }
426     nener  = 1 + ngid;
427     if (bDispCorr)
428     {
429         nener += 1;
430     }
431     if (bCharge)
432     {
433         nener += ngid;
434         if (bRFExcl)
435         {
436             nener += 1;
437         }
438         if (EEL_FULL(fr->ic->eeltype))
439         {
440             nener += 1;
441         }
442     }
443     snew(sum_UgembU, nener);
444
445     /* Copy the random seed set by the user */
446     seed = inputrec->ld_seed;
447
448     gmx::ThreeFry2x64<16>                rng(seed, gmx::RandomDomain::TestParticleInsertion); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
449     gmx::UniformRealDistribution<real>   dist;
450
451     if (MASTER(cr))
452     {
453         fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
454                           "TPI energies", "Time (ps)",
455                           "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
456         xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
457         snew(leg, 4+nener);
458         e = 0;
459         sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
460         leg[e++] = gmx_strdup(str);
461         sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
462         leg[e++] = gmx_strdup(str);
463         sprintf(str, "f. <e\\S-\\betaU\\N>");
464         leg[e++] = gmx_strdup(str);
465         sprintf(str, "f. V");
466         leg[e++] = gmx_strdup(str);
467         sprintf(str, "f. <Ue\\S-\\betaU\\N>");
468         leg[e++] = gmx_strdup(str);
469         for (i = 0; i < ngid; i++)
470         {
471             sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
472                     *(groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput][i]]));
473             leg[e++] = gmx_strdup(str);
474         }
475         if (bDispCorr)
476         {
477             sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
478             leg[e++] = gmx_strdup(str);
479         }
480         if (bCharge)
481         {
482             for (i = 0; i < ngid; i++)
483             {
484                 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
485                         *(groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput][i]]));
486                 leg[e++] = gmx_strdup(str);
487             }
488             if (bRFExcl)
489             {
490                 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
491                 leg[e++] = gmx_strdup(str);
492             }
493             if (EEL_FULL(fr->ic->eeltype))
494             {
495                 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
496                 leg[e++] = gmx_strdup(str);
497             }
498         }
499         xvgr_legend(fp_tpi, 4+nener, leg, oenv);
500         for (i = 0; i < 4+nener; i++)
501         {
502             sfree(leg[i]);
503         }
504         sfree(leg);
505     }
506     clear_rvec(x_init);
507     V_all     = 0;
508     VembU_all = 0;
509
510     invbinw = 10;
511     nbin    = 10;
512     snew(bin, nbin);
513
514     /* Avoid frame step numbers <= -1 */
515     frame_step_prev = -1;
516
517     bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
518                                      &rerun_fr, TRX_NEED_X);
519     frame = 0;
520
521     if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
522         mdatoms->nr - (a_tp1 - a_tp0))
523     {
524         gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
525                   "is not equal the number in the run input file (%d) "
526                   "minus the number of atoms to insert (%d)\n",
527                   rerun_fr.natoms, bCavity ? " minus one" : "",
528                   mdatoms->nr, a_tp1-a_tp0);
529     }
530
531     refvolshift = log(det(rerun_fr.box));
532
533     switch (inputrec->eI)
534     {
535         case eiTPI:
536             stepblocksize = inputrec->nstlist;
537             break;
538         case eiTPIC:
539             stepblocksize = 1;
540             break;
541         default:
542             gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
543     }
544
545     while (bNotLastFrame)
546     {
547         frame_step      = rerun_fr.step;
548         if (frame_step <= frame_step_prev)
549         {
550             /* We don't have step number in the trajectory file,
551              * or we have constant or decreasing step numbers.
552              * Ensure we have increasing step numbers, since we use
553              * the step numbers as a counter for random numbers.
554              */
555             frame_step  = frame_step_prev + 1;
556         }
557         frame_step_prev = frame_step;
558
559         lambda = rerun_fr.lambda;
560         t      = rerun_fr.time;
561
562         sum_embU = 0;
563         for (e = 0; e < nener; e++)
564         {
565             sum_UgembU[e] = 0;
566         }
567
568         /* Copy the coordinates from the input trajectory */
569         auto x = makeArrayRef(state_global->x);
570         for (i = 0; i < rerun_fr.natoms; i++)
571         {
572             copy_rvec(rerun_fr.x[i], x[i]);
573         }
574         copy_mat(rerun_fr.box, state_global->box);
575         const matrix &box = state_global->box;
576
577         V    = det(box);
578         logV = log(V);
579
580         bStateChanged = TRUE;
581         bNS           = TRUE;
582
583         put_atoms_in_box(fr->ePBC, box, x);
584
585         /* Put all atoms except for the inserted ones on the grid */
586         rvec vzero       = { 0, 0, 0 };
587         rvec boxDiagonal = { box[XX][XX], box[YY][YY], box[ZZ][ZZ] };
588         nbnxn_put_on_grid(fr->nbv.get(), box,
589                           0, vzero, boxDiagonal,
590                           nullptr, { 0, a_tp0 }, -1,
591                           fr->cginfo, x,
592                           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] +=
643                                     mass_cavity[i]*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,
664                                   1, x_init, x_init,
665                                   nullptr, { a_tp0, a_tp1 }, -1,
666                                   fr->cginfo, x,
667                                   0, nullptr);
668
669                 /* TODO: Avoid updating all atoms at every bNS step */
670                 fr->nbv->setAtomProperties(*mdatoms, fr->cginfo);
671
672                 fr->nbv->constructPairlist(Nbnxm::InteractionLocality::Local,
673                                            &top.excls, step, nrnb);
674
675                 bNS = FALSE;
676             }
677
678             /* Add random displacement uniformly distributed in a sphere
679              * of radius rtpi. We don't need to do this is we generate
680              * a new center location every step.
681              */
682             rvec x_tp;
683             if (bCavity || inputrec->nstlist > 1)
684             {
685                 /* Generate coordinates within |dx|=drmax of x_init */
686                 do
687                 {
688                     for (d = 0; d < DIM; d++)
689                     {
690                         dx[d] = (2*dist(rng) - 1)*drmax;
691                     }
692                 }
693                 while (norm2(dx) > drmax*drmax);
694                 rvec_add(x_init, dx, x_tp);
695             }
696             else
697             {
698                 copy_rvec(x_init, x_tp);
699             }
700
701             if (a_tp1 - a_tp0 == 1)
702             {
703                 /* Insert a single atom, just copy the insertion location */
704                 copy_rvec(x_tp, x[a_tp0]);
705             }
706             else
707             {
708                 /* Copy the coordinates from the top file */
709                 for (i = a_tp0; i < a_tp1; i++)
710                 {
711                     copy_rvec(x_mol[i-a_tp0], x[i]);
712                 }
713                 /* Rotate the molecule randomly */
714                 real angleX = 2*M_PI*dist(rng);
715                 real angleY = 2*M_PI*dist(rng);
716                 real angleZ = 2*M_PI*dist(rng);
717                 rotate_conf(a_tp1-a_tp0, state_global->x.rvec_array()+a_tp0, nullptr,
718                             angleX, angleY, angleZ);
719                 /* Shift to the insertion location */
720                 for (i = a_tp0; i < a_tp1; i++)
721                 {
722                     rvec_inc(x[i], x_tp);
723                 }
724             }
725
726             /* Note: NonLocal refers to the inserted molecule */
727             fr->nbv->convertCoordinates(Nbnxm::AtomLocality::NonLocal, false, x);
728
729             /* Clear some matrix variables  */
730             clear_mat(force_vir);
731             clear_mat(shake_vir);
732             clear_mat(vir);
733             clear_mat(pres);
734
735             /* Calc energy (no forces) on new positions.
736              * Since we only need the intermolecular energy
737              * and the RF exclusion terms of the inserted molecule occur
738              * within a single charge group we can pass NULL for the graph.
739              * This also avoids shifts that would move charge groups
740              * out of the box. */
741             /* Make do_force do a single node force calculation */
742             cr->nnodes = 1;
743
744             // TPI might place a particle so close that the potential
745             // is infinite. Since this is intended to happen, we
746             // temporarily suppress any exceptions that the processor
747             // might raise, then restore the old behaviour.
748             std::fenv_t floatingPointEnvironment;
749             std::feholdexcept(&floatingPointEnvironment);
750             do_force(fplog, cr, ms, inputrec, nullptr, nullptr, imdSession,
751                      pull_work,
752                      step, nrnb, wcycle, &top,
753                      state_global->box, state_global->x.arrayRefWithPadding(), &state_global->hist,
754                      f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
755                      state_global->lambda,
756                      nullptr, fr, runScheduleWork, nullptr, mu_tot, t, nullptr,
757                      GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
758                      (bStateChanged ? GMX_FORCE_STATECHANGED : 0),
759                      DDBalanceRegionHandler(nullptr));
760             std::feclearexcept(FE_DIVBYZERO | FE_INVALID | FE_OVERFLOW);
761             std::feupdateenv(&floatingPointEnvironment);
762
763             cr->nnodes    = nnodes;
764             bStateChanged = FALSE;
765
766             if (fr->dispersionCorrection)
767             {
768                 /* Calculate long range corrections to pressure and energy */
769                 const DispersionCorrection::Correction correction =
770                     fr->dispersionCorrection->calculate(state_global->box, lambda);
771                 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
772                 enerd->term[F_DISPCORR] = correction.energy;
773                 enerd->term[F_EPOT]    += correction.energy;
774                 enerd->term[F_PRES]    += correction.pressure;
775                 enerd->term[F_DVDL]    += correction.dvdl;
776             }
777             else
778             {
779                 enerd->term[F_DISPCORR]  = 0;
780             }
781             if (EEL_RF(fr->ic->eeltype))
782             {
783                 enerd->term[F_EPOT]    += rfExclusionEnergy;
784             }
785
786             epot               = enerd->term[F_EPOT];
787             bEnergyOutOfBounds = FALSE;
788
789             /* If the compiler doesn't optimize this check away
790              * we catch the NAN energies.
791              * The epot>GMX_REAL_MAX check catches inf values,
792              * which should nicely result in embU=0 through the exp below,
793              * but it does not hurt to check anyhow.
794              */
795             /* Non-bonded Interaction usually diverge at r=0.
796              * With tabulated interaction functions the first few entries
797              * should be capped in a consistent fashion between
798              * repulsion, dispersion and Coulomb to avoid accidental
799              * negative values in the total energy.
800              * The table generation code in tables.c does this.
801              * With user tbales the user should take care of this.
802              */
803             if (epot != epot || epot > GMX_REAL_MAX)
804             {
805                 bEnergyOutOfBounds = TRUE;
806             }
807             if (bEnergyOutOfBounds)
808             {
809                 if (debug)
810                 {
811                     fprintf(debug, "\n  time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, static_cast<int>(step), epot);
812                 }
813                 embU = 0;
814             }
815             else
816             {
817                 // Exponent argument is fine in SP range, but output can be in DP range
818                 embU      = exp(static_cast<double>(-beta*epot));
819                 sum_embU += embU;
820                 /* Determine the weighted energy contributions of each energy group */
821                 e                = 0;
822                 sum_UgembU[e++] += epot*embU;
823                 if (fr->bBHAM)
824                 {
825                     for (i = 0; i < ngid; i++)
826                     {
827                         sum_UgembU[e++] +=
828                             enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)]*embU;
829                     }
830                 }
831                 else
832                 {
833                     for (i = 0; i < ngid; i++)
834                     {
835                         sum_UgembU[e++] +=
836                             enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)]*embU;
837                     }
838                 }
839                 if (bDispCorr)
840                 {
841                     sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
842                 }
843                 if (bCharge)
844                 {
845                     for (i = 0; i < ngid; i++)
846                     {
847                         sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
848                     }
849                     if (bRFExcl)
850                     {
851                         sum_UgembU[e++] += rfExclusionEnergy*embU;
852                     }
853                     if (EEL_FULL(fr->ic->eeltype))
854                     {
855                         sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
856                     }
857                 }
858             }
859
860             if (embU == 0 || beta*epot > bU_bin_limit)
861             {
862                 bin[0]++;
863             }
864             else
865             {
866                 i = gmx::roundToInt((bU_logV_bin_limit
867                                      - (beta*epot - logV + refvolshift))*invbinw);
868                 if (i < 0)
869                 {
870                     i = 0;
871                 }
872                 if (i >= nbin)
873                 {
874                     realloc_bins(&bin, &nbin, i+10);
875                 }
876                 bin[i]++;
877             }
878
879             if (debug)
880             {
881                 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
882                         static_cast<int>(step), epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
883             }
884
885             if (dump_pdb && epot <= dump_ener)
886             {
887                 sprintf(str, "t%g_step%d.pdb", t, static_cast<int>(step));
888                 sprintf(str2, "t: %f step %d ener: %f", t, static_cast<int>(step), epot);
889                 write_sto_conf_mtop(str, str2, top_global, state_global->x.rvec_array(), state_global->v.rvec_array(),
890                                     inputrec->ePBC, state_global->box);
891             }
892
893             step++;
894             if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
895             {
896                 /* Skip all steps assigned to the other MPI ranks */
897                 step += (cr->nnodes - 1)*stepblocksize;
898             }
899         }
900
901         if (PAR(cr))
902         {
903             /* When running in parallel sum the energies over the processes */
904             gmx_sumd(1,    &sum_embU, cr);
905             gmx_sumd(nener, sum_UgembU, cr);
906         }
907
908         frame++;
909         V_all     += V;
910         VembU_all += V*sum_embU/nsteps;
911
912         if (fp_tpi)
913         {
914             if (mdrunOptions.verbose || frame%10 == 0 || frame < 10)
915             {
916                 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
917                         -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
918             }
919
920             fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
921                     t,
922                     VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
923                     sum_embU == 0  ? 20/beta : -log(sum_embU/nsteps)/beta,
924                     sum_embU/nsteps, V);
925             for (e = 0; e < nener; e++)
926             {
927                 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
928             }
929             fprintf(fp_tpi, "\n");
930             fflush(fp_tpi);
931         }
932
933         bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
934     }   /* End of the loop  */
935     walltime_accounting_end_time(walltime_accounting);
936
937     close_trx(status);
938
939     if (fp_tpi != nullptr)
940     {
941         xvgrclose(fp_tpi);
942     }
943
944     if (fplog != nullptr)
945     {
946         fprintf(fplog, "\n");
947         fprintf(fplog, "  <V>  = %12.5e nm^3\n", V_all/frame);
948         const double mu = -log(VembU_all/V_all)/beta;
949         fprintf(fplog, "  <mu> = %12.5e kJ/mol\n", mu);
950
951         if (!std::isfinite(mu))
952         {
953             fprintf(fplog, "\nThe computed chemical potential is not finite - consider increasing the number of steps and/or the number of frames to insert into.\n");
954         }
955     }
956
957     /* Write the Boltzmann factor histogram */
958     if (PAR(cr))
959     {
960         /* When running in parallel sum the bins over the processes */
961         i = nbin;
962         global_max(cr, &i);
963         realloc_bins(&bin, &nbin, i);
964         gmx_sumd(nbin, bin, cr);
965     }
966     if (MASTER(cr))
967     {
968         fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
969                           "TPI energy distribution",
970                           "\\betaU - log(V/<V>)", "count", oenv);
971         sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
972         xvgr_subtitle(fp_tpi, str, oenv);
973         xvgr_legend(fp_tpi, 2, tpid_leg, oenv);
974         for (i = nbin-1; i > 0; i--)
975         {
976             bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
977             fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
978                     bUlogV,
979                     roundToInt(bin[i]),
980                     bin[i]*exp(-bUlogV)*V_all/VembU_all);
981         }
982         xvgrclose(fp_tpi);
983     }
984     sfree(bin);
985
986     sfree(sum_UgembU);
987
988     walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);
989 }
990
991 } // namespace gmx