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