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