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