2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
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"
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"
79 static void global_max(t_commrec *cr, int *n)
83 snew(sum, cr->nnodes);
85 gmx_sumi(cr->nnodes, sum, cr);
86 for (i = 0; i < cr->nnodes; i++)
94 static void realloc_bins(double **bin, int *nbin, int nbin_new)
98 if (nbin_new != *nbin)
100 srenew(*bin, nbin_new);
101 for (i = *nbin; i < nbin_new; i++)
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,
119 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
120 gmx_edsam_t gmx_unused ed,
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)
130 const char *TPI = "Test Particle Insertion";
132 gmx_groups_t *groups;
133 gmx_enerdata_t *enerd;
135 real lambda, t, temp, beta, drmax, epot;
136 double embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
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;
143 rvec mu_tot, x_init, dx, x_tp;
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;
152 char *ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
153 double dbl, dump_ener;
155 int nat_cavity = 0, d;
156 real *mass_cavity = NULL, mass_tot;
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"};
163 /* Since there is no upper limit to the insertion energies,
164 * we need to set an upper limit for the distribution output.
166 real bU_bin_limit = 50;
167 real bU_logV_bin_limit = bU_bin_limit + 10;
169 if (inputrec->cutoff_scheme == ecutsVERLET)
171 gmx_fatal(FARGS, "TPI does not work (yet) with the Verlet cut-off scheme");
176 top = gmx_mtop_generate_local_top(top_global, inputrec);
178 groups = &top_global->groups;
180 bCavity = (inputrec->eI == eiTPIC);
183 ptr = getenv("GMX_TPIC_MASSES");
190 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
191 * The center of mass of the last atoms is then used for TPIC.
194 while (sscanf(ptr, "%lf%n", &dbl, &i) > 0)
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]);
205 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
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 */
215 /* Determine the temperature for the Boltzmann weighting */
216 temp = inputrec->opts.ref_t[0];
219 for (i = 1; (i < inputrec->opts.ngtc); i++)
221 if (inputrec->opts.ref_t[i] != temp)
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");
228 "\n The temperature for test particle insertion is %.3f K\n\n",
231 beta = 1.0/(BOLTZ*temp);
233 /* Number of insertions per frame */
234 nsteps = inputrec->nsteps;
236 /* Use the same neighborlist with more insertions points
237 * in a sphere of radius drmax around the initial point
239 /* This should be a proper mdp parameter */
240 drmax = inputrec->rtpi;
242 /* An environment variable can be set to dump all configurations
243 * to pdb with an insertion energy <= this value.
245 dump_pdb = getenv("GMX_TPI_DUMP");
249 sscanf(dump_pdb, "%lf", &dump_ener);
252 atoms2md(top_global, inputrec, 0, NULL, top_global->natoms, mdatoms);
253 update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
256 init_enerdata(groups->grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
257 snew(f, top_global->natoms);
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");
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];
270 fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
272 if (a_tp1 - a_tp0 > 1 &&
273 (inputrec->rlist < inputrec->rcoulomb ||
274 inputrec->rlist < inputrec->rvdw))
276 gmx_fatal(FARGS, "Can not do TPI for multi-atom molecule with a twin-range cut-off");
278 snew(x_mol, a_tp1-a_tp0);
280 bDispCorr = (inputrec->eDispCorr != edispcNO);
282 for (i = a_tp0; i < a_tp1; i++)
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));
290 bRFExcl = (bCharge && EEL_RF(fr->eeltype) && fr->eeltype != eelRF_NEC);
292 calc_cgcm(fplog, cg_tp, cg_tp+1, &(top->cgs), state->x, fr->cg_cm);
295 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
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");
303 /* Center the molecule to be inserted at zero */
304 for (i = 0; i < a_tp1-a_tp0; i++)
306 rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
312 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
313 a_tp1-a_tp0, bCharge ? "with" : "without");
315 fprintf(fplog, "\nWill insert %d times in each frame of %s\n",
316 (int)nsteps, opt2fn("-rerun", nfile, fnm));
321 if (inputrec->nstlist > 1)
323 if (drmax == 0 && a_tp1-a_tp0 == 1)
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);
329 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
337 fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
341 ngid = groups->grps[egcENER].nr;
342 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
355 if (EEL_FULL(fr->eeltype))
360 snew(sum_UgembU, nener);
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.
372 rnd_count_stride = 10000;
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);
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++)
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);
400 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
401 leg[e++] = gmx_strdup(str);
405 for (i = 0; i < ngid; i++)
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);
413 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
414 leg[e++] = gmx_strdup(str);
416 if (EEL_FULL(fr->eeltype))
418 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
419 leg[e++] = gmx_strdup(str);
422 xvgr_legend(fp_tpi, 4+nener, (const char**)leg, oenv);
423 for (i = 0; i < 4+nener; i++)
437 /* Avoid frame step numbers <= -1 */
438 frame_step_prev = -1;
440 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
441 &rerun_fr, TRX_NEED_X);
444 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
445 mdatoms->nr - (a_tp1 - a_tp0))
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);
454 refvolshift = log(det(rerun_fr.box));
456 switch (inputrec->eI)
459 stepblocksize = inputrec->nstlist;
465 gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
469 /* Make sure we don't detect SIMD overflow generated before this point */
470 gmx_simd_check_and_reset_overflow();
473 while (bNotLastFrame)
475 frame_step = rerun_fr.step;
476 if (frame_step <= frame_step_prev)
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.
483 frame_step = frame_step_prev + 1;
485 frame_step_prev = frame_step;
487 lambda = rerun_fr.lambda;
491 for (e = 0; e < nener; e++)
496 /* Copy the coordinates from the input trajectory */
497 for (i = 0; i < rerun_fr.natoms; i++)
499 copy_rvec(rerun_fr.x[i], state->x[i]);
501 copy_mat(rerun_fr.box, state->box);
506 bStateChanged = TRUE;
509 step = cr->nodeid*stepblocksize;
510 while (step < nsteps)
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.
518 rnd_count = step*rnd_count_stride;
522 /* Random insertion in the whole volume */
523 bNS = (step % inputrec->nstlist == 0);
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++)
531 x_init[d] = rnd[d]*state->box[d][d];
534 if (inputrec->nstlist == 1)
536 copy_rvec(x_init, x_tp);
540 /* Generate coordinates within |dx|=drmax of x_init */
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++)
547 dx[d] = (2*rnd[d] - 1)*drmax;
550 while (norm2(dx) > drmax*drmax);
551 rvec_add(x_init, dx, x_tp);
556 /* Random insertion around a cavity location
557 * given by the last coordinate of the trajectory.
563 /* Copy the location of the cavity */
564 copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
568 /* Determine the center of mass of the last molecule */
571 for (i = 0; i < nat_cavity; i++)
573 for (d = 0; d < DIM; d++)
576 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
578 mass_tot += mass_cavity[i];
580 for (d = 0; d < DIM; d++)
582 x_init[d] /= mass_tot;
586 /* Generate coordinates within |dx|=drmax of x_init */
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++)
593 dx[d] = (2*rnd[d] - 1)*drmax;
596 while (norm2(dx) > drmax*drmax);
597 rvec_add(x_init, dx, x_tp);
600 if (a_tp1 - a_tp0 == 1)
602 /* Insert a single atom, just copy the insertion location */
603 copy_rvec(x_tp, state->x[a_tp0]);
607 /* Copy the coordinates from the top file */
608 for (i = a_tp0; i < a_tp1; i++)
610 copy_rvec(x_mol[i-a_tp0], state->x[i]);
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,
619 /* Shift to the insertion location */
620 for (i = a_tp0; i < a_tp1; i++)
622 rvec_inc(state->x[i], x_tp);
626 /* Clear some matrix variables */
627 clear_mat(force_vir);
628 clear_mat(shake_vir);
632 /* Set the charge group center of mass of the test particle */
633 copy_rvec(x_init, fr->cg_cm[top->cgs.nr-1]);
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
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.
646 /* Make do_force do a single node force calculation */
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,
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));
658 bStateChanged = FALSE;
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;
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())
678 fprintf(debug, "Found an SSE overflow, assuming the energy is out of bounds\n");
680 bEnergyOutOfBounds = TRUE;
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.
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.
697 if (epot != epot || epot > GMX_REAL_MAX)
699 bEnergyOutOfBounds = TRUE;
701 if (bEnergyOutOfBounds)
705 fprintf(debug, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, (int)step, epot);
711 embU = exp(-beta*epot);
713 /* Determine the weighted energy contributions of each energy group */
715 sum_UgembU[e++] += epot*embU;
718 for (i = 0; i < ngid; i++)
721 (enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)] +
722 enerd->grpp.ener[egBHAMLR][GID(i, gid_tp, ngid)])*embU;
727 for (i = 0; i < ngid; i++)
730 (enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)] +
731 enerd->grpp.ener[egLJLR][GID(i, gid_tp, ngid)])*embU;
736 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
740 for (i = 0; i < ngid; i++)
743 (enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] +
744 enerd->grpp.ener[egCOULLR][GID(i, gid_tp, ngid)])*embU;
748 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
750 if (EEL_FULL(fr->eeltype))
752 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
757 if (embU == 0 || beta*epot > bU_bin_limit)
763 i = (int)((bU_logV_bin_limit
764 - (beta*epot - logV + refvolshift))*invbinw
772 realloc_bins(&bin, &nbin, i+10);
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]);
783 if (dump_pdb && epot <= dump_ener)
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);
792 if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
794 /* Skip all steps assigned to the other MPI ranks */
795 step += (cr->nnodes - 1)*stepblocksize;
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);
808 VembU_all += V*sum_embU/nsteps;
812 if (bVerbose || frame%10 == 0 || frame < 10)
814 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
815 -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
818 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
820 VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
821 sum_embU == 0 ? 20/beta : -log(sum_embU/nsteps)/beta,
823 for (e = 0; e < nener; e++)
825 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
827 fprintf(fp_tpi, "\n");
831 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
832 } /* End of the loop */
833 walltime_accounting_end(walltime_accounting);
839 gmx_fio_fclose(fp_tpi);
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);
849 /* Write the Boltzmann factor histogram */
852 /* When running in parallel sum the bins over the processes */
855 realloc_bins(&bin, &nbin, i);
856 gmx_sumd(nbin, bin, cr);
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--)
868 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
869 fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
872 bin[i]*exp(-bUlogV)*V_all/VembU_all);
874 gmx_fio_fclose(fp_tpi);
880 walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);