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