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