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