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