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