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