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