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