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