Update Awh initialization and lifetime management
[gromacs.git] / src / gromacs / mdrun / tpi.cpp
blobb3d28a63daf1a0b2195795e61132005928ed4d91
1 /*
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,2015,2016,2017,2018, 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.
37 /*! \internal \file
39 * \brief This file defines the integrator for test particle insertion
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_mdrun
44 #include "gmxpre.h"
46 #include <cmath>
47 #include <cstdlib>
48 #include <cstring>
49 #include <ctime>
51 #include <algorithm>
53 #include "gromacs/commandline/filenm.h"
54 #include "gromacs/domdec/domdec.h"
55 #include "gromacs/ewald/pme.h"
56 #include "gromacs/fileio/confio.h"
57 #include "gromacs/fileio/trxio.h"
58 #include "gromacs/fileio/xvgr.h"
59 #include "gromacs/gmxlib/chargegroup.h"
60 #include "gromacs/gmxlib/conformation-utilities.h"
61 #include "gromacs/gmxlib/network.h"
62 #include "gromacs/gmxlib/nrnb.h"
63 #include "gromacs/math/units.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/mdlib/constr.h"
66 #include "gromacs/mdlib/force.h"
67 #include "gromacs/mdlib/force_flags.h"
68 #include "gromacs/mdlib/mdatoms.h"
69 #include "gromacs/mdlib/mdebin.h"
70 #include "gromacs/mdlib/mdrun.h"
71 #include "gromacs/mdlib/ns.h"
72 #include "gromacs/mdlib/sim_util.h"
73 #include "gromacs/mdlib/tgroup.h"
74 #include "gromacs/mdlib/update.h"
75 #include "gromacs/mdlib/vsite.h"
76 #include "gromacs/mdtypes/commrec.h"
77 #include "gromacs/mdtypes/group.h"
78 #include "gromacs/mdtypes/inputrec.h"
79 #include "gromacs/mdtypes/md_enums.h"
80 #include "gromacs/mdtypes/state.h"
81 #include "gromacs/pbcutil/pbc.h"
82 #include "gromacs/random/threefry.h"
83 #include "gromacs/random/uniformrealdistribution.h"
84 #include "gromacs/timing/wallcycle.h"
85 #include "gromacs/timing/walltime_accounting.h"
86 #include "gromacs/topology/mtop_util.h"
87 #include "gromacs/trajectory/trajectoryframe.h"
88 #include "gromacs/utility/cstringutil.h"
89 #include "gromacs/utility/fatalerror.h"
90 #include "gromacs/utility/gmxassert.h"
91 #include "gromacs/utility/smalloc.h"
93 #include "integrator.h"
95 //! Global max algorithm
96 static void global_max(t_commrec *cr, int *n)
98 int *sum, i;
100 snew(sum, cr->nnodes);
101 sum[cr->nodeid] = *n;
102 gmx_sumi(cr->nnodes, sum, cr);
103 for (i = 0; i < cr->nnodes; i++)
105 *n = std::max(*n, sum[i]);
108 sfree(sum);
111 //! Reallocate arrays.
112 static void realloc_bins(double **bin, int *nbin, int nbin_new)
114 int i;
116 if (nbin_new != *nbin)
118 srenew(*bin, nbin_new);
119 for (i = *nbin; i < nbin_new; i++)
121 (*bin)[i] = 0;
123 *nbin = nbin_new;
127 namespace gmx
130 void
131 Integrator::do_tpi()
133 gmx_localtop_t *top;
134 gmx_groups_t *groups;
135 gmx_enerdata_t *enerd;
136 PaddedRVecVector 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 seed;
150 int i;
151 FILE *fp_tpi = nullptr;
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 = nullptr, mass_tot;
157 int nbin;
158 double invbinw, *bin, refvolshift, logV, bUlogV;
159 real prescorr, enercorr, dvdlcorr;
160 gmx_bool bEnergyOutOfBounds;
161 const char *tpid_leg[2] = {"direct", "reweighted"};
162 auto mdatoms = mdAtoms->mdatoms();
164 GMX_UNUSED_VALUE(outputProvider);
166 /* Since there is no upper limit to the insertion energies,
167 * we need to set an upper limit for the distribution output.
169 real bU_bin_limit = 50;
170 real bU_logV_bin_limit = bU_bin_limit + 10;
172 if (inputrec->cutoff_scheme == ecutsVERLET)
174 gmx_fatal(FARGS, "TPI does not work (yet) with the Verlet cut-off scheme");
177 nnodes = cr->nnodes;
179 top = gmx_mtop_generate_local_top(top_global, inputrec->efep != efepNO);
181 groups = &top_global->groups;
183 bCavity = (inputrec->eI == eiTPIC);
184 if (bCavity)
186 ptr = getenv("GMX_TPIC_MASSES");
187 if (ptr == nullptr)
189 nat_cavity = 1;
191 else
193 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
194 * The center of mass of the last atoms is then used for TPIC.
196 nat_cavity = 0;
197 while (sscanf(ptr, "%20lf%n", &dbl, &i) > 0)
199 srenew(mass_cavity, nat_cavity+1);
200 mass_cavity[nat_cavity] = dbl;
201 fprintf(fplog, "mass[%d] = %f\n",
202 nat_cavity+1, mass_cavity[nat_cavity]);
203 nat_cavity++;
204 ptr += i;
206 if (nat_cavity == 0)
208 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
214 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
215 state_global->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
216 /* We never need full pbc for TPI */
217 fr->ePBC = epbcXYZ;
218 /* Determine the temperature for the Boltzmann weighting */
219 temp = inputrec->opts.ref_t[0];
220 if (fplog)
222 for (i = 1; (i < inputrec->opts.ngtc); i++)
224 if (inputrec->opts.ref_t[i] != temp)
226 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
227 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
230 fprintf(fplog,
231 "\n The temperature for test particle insertion is %.3f K\n\n",
232 temp);
234 beta = 1.0/(BOLTZ*temp);
236 /* Number of insertions per frame */
237 nsteps = inputrec->nsteps;
239 /* Use the same neighborlist with more insertions points
240 * in a sphere of radius drmax around the initial point
242 /* This should be a proper mdp parameter */
243 drmax = inputrec->rtpi;
245 /* An environment variable can be set to dump all configurations
246 * to pdb with an insertion energy <= this value.
248 dump_pdb = getenv("GMX_TPI_DUMP");
249 dump_ener = 0;
250 if (dump_pdb)
252 sscanf(dump_pdb, "%20lf", &dump_ener);
255 atoms2md(top_global, inputrec, -1, nullptr, top_global->natoms, mdAtoms);
256 update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
258 snew(enerd, 1);
259 init_enerdata(groups->grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
260 /* We need to allocate one element extra, since we might use
261 * (unaligned) 4-wide SIMD loads to access rvec entries.
263 f.resize(gmx::paddedRVecVectorSize(top_global->natoms));
265 /* Print to log file */
266 walltime_accounting_start(walltime_accounting);
267 wallcycle_start(wcycle, ewcRUN);
268 print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
270 /* The last charge group is the group to be inserted */
271 cg_tp = top->cgs.nr - 1;
272 a_tp0 = top->cgs.index[cg_tp];
273 a_tp1 = top->cgs.index[cg_tp+1];
274 if (debug)
276 fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
279 GMX_RELEASE_ASSERT(inputrec->rcoulomb <= inputrec->rlist && inputrec->rvdw <= inputrec->rlist, "Twin-range interactions are not supported with TPI");
281 snew(x_mol, a_tp1-a_tp0);
283 bDispCorr = (inputrec->eDispCorr != edispcNO);
284 bCharge = FALSE;
285 for (i = a_tp0; i < a_tp1; i++)
287 /* Copy the coordinates of the molecule to be insterted */
288 copy_rvec(state_global->x[i], x_mol[i-a_tp0]);
289 /* Check if we need to print electrostatic energies */
290 bCharge |= (mdatoms->chargeA[i] != 0 ||
291 (mdatoms->chargeB && mdatoms->chargeB[i] != 0));
293 bRFExcl = (bCharge && EEL_RF(fr->ic->eeltype));
295 calc_cgcm(fplog, cg_tp, cg_tp+1, &(top->cgs), as_rvec_array(state_global->x.data()), fr->cg_cm);
296 if (bCavity)
298 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
300 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
301 fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
304 else
306 /* Center the molecule to be inserted at zero */
307 for (i = 0; i < a_tp1-a_tp0; i++)
309 rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
313 if (fplog)
315 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
316 a_tp1-a_tp0, bCharge ? "with" : "without");
318 fprintf(fplog, "\nWill insert %d times in each frame of %s\n",
319 (int)nsteps, opt2fn("-rerun", nfile, fnm));
322 if (!bCavity)
324 if (inputrec->nstlist > 1)
326 if (drmax == 0 && a_tp1-a_tp0 == 1)
328 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);
330 if (fplog)
332 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
336 else
338 if (fplog)
340 fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
344 ngid = groups->grps[egcENER].nr;
345 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
346 nener = 1 + ngid;
347 if (bDispCorr)
349 nener += 1;
351 if (bCharge)
353 nener += ngid;
354 if (bRFExcl)
356 nener += 1;
358 if (EEL_FULL(fr->ic->eeltype))
360 nener += 1;
363 snew(sum_UgembU, nener);
365 /* Copy the random seed set by the user */
366 seed = inputrec->ld_seed;
368 gmx::ThreeFry2x64<16> rng(seed, gmx::RandomDomain::TestParticleInsertion); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
369 gmx::UniformRealDistribution<real> dist;
371 if (MASTER(cr))
373 fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
374 "TPI energies", "Time (ps)",
375 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
376 xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
377 snew(leg, 4+nener);
378 e = 0;
379 sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
380 leg[e++] = gmx_strdup(str);
381 sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
382 leg[e++] = gmx_strdup(str);
383 sprintf(str, "f. <e\\S-\\betaU\\N>");
384 leg[e++] = gmx_strdup(str);
385 sprintf(str, "f. V");
386 leg[e++] = gmx_strdup(str);
387 sprintf(str, "f. <Ue\\S-\\betaU\\N>");
388 leg[e++] = gmx_strdup(str);
389 for (i = 0; i < ngid; i++)
391 sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
392 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
393 leg[e++] = gmx_strdup(str);
395 if (bDispCorr)
397 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
398 leg[e++] = gmx_strdup(str);
400 if (bCharge)
402 for (i = 0; i < ngid; i++)
404 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
405 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
406 leg[e++] = gmx_strdup(str);
408 if (bRFExcl)
410 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
411 leg[e++] = gmx_strdup(str);
413 if (EEL_FULL(fr->ic->eeltype))
415 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
416 leg[e++] = gmx_strdup(str);
419 xvgr_legend(fp_tpi, 4+nener, (const char**)leg, oenv);
420 for (i = 0; i < 4+nener; i++)
422 sfree(leg[i]);
424 sfree(leg);
426 clear_rvec(x_init);
427 V_all = 0;
428 VembU_all = 0;
430 invbinw = 10;
431 nbin = 10;
432 snew(bin, nbin);
434 /* Avoid frame step numbers <= -1 */
435 frame_step_prev = -1;
437 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
438 &rerun_fr, TRX_NEED_X);
439 frame = 0;
441 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
442 mdatoms->nr - (a_tp1 - a_tp0))
444 gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
445 "is not equal the number in the run input file (%d) "
446 "minus the number of atoms to insert (%d)\n",
447 rerun_fr.natoms, bCavity ? " minus one" : "",
448 mdatoms->nr, a_tp1-a_tp0);
451 refvolshift = log(det(rerun_fr.box));
453 switch (inputrec->eI)
455 case eiTPI:
456 stepblocksize = inputrec->nstlist;
457 break;
458 case eiTPIC:
459 stepblocksize = 1;
460 break;
461 default:
462 gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
465 while (bNotLastFrame)
467 frame_step = rerun_fr.step;
468 if (frame_step <= frame_step_prev)
470 /* We don't have step number in the trajectory file,
471 * or we have constant or decreasing step numbers.
472 * Ensure we have increasing step numbers, since we use
473 * the step numbers as a counter for random numbers.
475 frame_step = frame_step_prev + 1;
477 frame_step_prev = frame_step;
479 lambda = rerun_fr.lambda;
480 t = rerun_fr.time;
482 sum_embU = 0;
483 for (e = 0; e < nener; e++)
485 sum_UgembU[e] = 0;
488 /* Copy the coordinates from the input trajectory */
489 for (i = 0; i < rerun_fr.natoms; i++)
491 copy_rvec(rerun_fr.x[i], state_global->x[i]);
493 copy_mat(rerun_fr.box, state_global->box);
495 V = det(state_global->box);
496 logV = log(V);
498 bStateChanged = TRUE;
499 bNS = TRUE;
501 step = cr->nodeid*stepblocksize;
502 while (step < nsteps)
504 /* Restart random engine using the frame and insertion step
505 * as counters.
506 * Note that we need to draw several random values per iteration,
507 * but by using the internal subcounter functionality of ThreeFry2x64
508 * we can draw 131072 unique 64-bit values before exhausting
509 * the stream. This is a huge margin, and if something still goes
510 * wrong you will get an exception when the stream is exhausted.
512 rng.restart(frame_step, step);
513 dist.reset(); // erase any memory in the distribution
515 if (!bCavity)
517 /* Random insertion in the whole volume */
518 bNS = (step % inputrec->nstlist == 0);
519 if (bNS)
521 /* Generate a random position in the box */
522 for (d = 0; d < DIM; d++)
524 x_init[d] = dist(rng)*state_global->box[d][d];
528 if (inputrec->nstlist == 1)
530 copy_rvec(x_init, x_tp);
532 else
534 /* Generate coordinates within |dx|=drmax of x_init */
537 for (d = 0; d < DIM; d++)
539 dx[d] = (2*dist(rng) - 1)*drmax;
542 while (norm2(dx) > drmax*drmax);
543 rvec_add(x_init, dx, x_tp);
546 else
548 /* Random insertion around a cavity location
549 * given by the last coordinate of the trajectory.
551 if (step == 0)
553 if (nat_cavity == 1)
555 /* Copy the location of the cavity */
556 copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
558 else
560 /* Determine the center of mass of the last molecule */
561 clear_rvec(x_init);
562 mass_tot = 0;
563 for (i = 0; i < nat_cavity; i++)
565 for (d = 0; d < DIM; d++)
567 x_init[d] +=
568 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
570 mass_tot += mass_cavity[i];
572 for (d = 0; d < DIM; d++)
574 x_init[d] /= mass_tot;
578 /* Generate coordinates within |dx|=drmax of x_init */
581 for (d = 0; d < DIM; d++)
583 dx[d] = (2*dist(rng) - 1)*drmax;
586 while (norm2(dx) > drmax*drmax);
587 rvec_add(x_init, dx, x_tp);
590 if (a_tp1 - a_tp0 == 1)
592 /* Insert a single atom, just copy the insertion location */
593 copy_rvec(x_tp, state_global->x[a_tp0]);
595 else
597 /* Copy the coordinates from the top file */
598 for (i = a_tp0; i < a_tp1; i++)
600 copy_rvec(x_mol[i-a_tp0], state_global->x[i]);
602 /* Rotate the molecule randomly */
603 real angleX = 2*M_PI*dist(rng);
604 real angleY = 2*M_PI*dist(rng);
605 real angleZ = 2*M_PI*dist(rng);
606 rotate_conf(a_tp1-a_tp0, as_rvec_array(state_global->x.data())+a_tp0, nullptr,
607 angleX, angleY, angleZ);
608 /* Shift to the insertion location */
609 for (i = a_tp0; i < a_tp1; i++)
611 rvec_inc(state_global->x[i], x_tp);
615 /* Clear some matrix variables */
616 clear_mat(force_vir);
617 clear_mat(shake_vir);
618 clear_mat(vir);
619 clear_mat(pres);
621 /* Set the charge group center of mass of the test particle */
622 copy_rvec(x_init, fr->cg_cm[top->cgs.nr-1]);
624 /* Calc energy (no forces) on new positions.
625 * Since we only need the intermolecular energy
626 * and the RF exclusion terms of the inserted molecule occur
627 * within a single charge group we can pass NULL for the graph.
628 * This also avoids shifts that would move charge groups
629 * out of the box. */
630 /* Make do_force do a single node force calculation */
631 cr->nnodes = 1;
632 do_force(fplog, cr, ms, inputrec, nullptr,
633 step, nrnb, wcycle, top, &top_global->groups,
634 state_global->box, state_global->x, &state_global->hist,
635 f, force_vir, mdatoms, enerd, fcd,
636 state_global->lambda,
637 nullptr, fr, nullptr, mu_tot, t, nullptr,
638 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
639 (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS : 0) |
640 (bStateChanged ? GMX_FORCE_STATECHANGED : 0),
641 DdOpenBalanceRegionBeforeForceComputation::no,
642 DdCloseBalanceRegionAfterForceComputation::no);
643 cr->nnodes = nnodes;
644 bStateChanged = FALSE;
645 bNS = FALSE;
647 /* Calculate long range corrections to pressure and energy */
648 calc_dispcorr(inputrec, fr, state_global->box,
649 lambda, pres, vir, &prescorr, &enercorr, &dvdlcorr);
650 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
651 enerd->term[F_DISPCORR] = enercorr;
652 enerd->term[F_EPOT] += enercorr;
653 enerd->term[F_PRES] += prescorr;
654 enerd->term[F_DVDL_VDW] += dvdlcorr;
656 epot = enerd->term[F_EPOT];
657 bEnergyOutOfBounds = FALSE;
659 /* If the compiler doesn't optimize this check away
660 * we catch the NAN energies.
661 * The epot>GMX_REAL_MAX check catches inf values,
662 * which should nicely result in embU=0 through the exp below,
663 * but it does not hurt to check anyhow.
665 /* Non-bonded Interaction usually diverge at r=0.
666 * With tabulated interaction functions the first few entries
667 * should be capped in a consistent fashion between
668 * repulsion, dispersion and Coulomb to avoid accidental
669 * negative values in the total energy.
670 * The table generation code in tables.c does this.
671 * With user tbales the user should take care of this.
673 if (epot != epot || epot > GMX_REAL_MAX)
675 bEnergyOutOfBounds = TRUE;
677 if (bEnergyOutOfBounds)
679 if (debug)
681 fprintf(debug, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, (int)step, epot);
683 embU = 0;
685 else
687 // Exponent argument is fine in SP range, but output can be in DP range
688 embU = exp(static_cast<double>(-beta*epot));
689 sum_embU += embU;
690 /* Determine the weighted energy contributions of each energy group */
691 e = 0;
692 sum_UgembU[e++] += epot*embU;
693 if (fr->bBHAM)
695 for (i = 0; i < ngid; i++)
697 sum_UgembU[e++] +=
698 enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)]*embU;
701 else
703 for (i = 0; i < ngid; i++)
705 sum_UgembU[e++] +=
706 enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)]*embU;
709 if (bDispCorr)
711 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
713 if (bCharge)
715 for (i = 0; i < ngid; i++)
717 sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
719 if (bRFExcl)
721 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
723 if (EEL_FULL(fr->ic->eeltype))
725 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
730 if (embU == 0 || beta*epot > bU_bin_limit)
732 bin[0]++;
734 else
736 i = (int)((bU_logV_bin_limit
737 - (beta*epot - logV + refvolshift))*invbinw
738 + 0.5);
739 if (i < 0)
741 i = 0;
743 if (i >= nbin)
745 realloc_bins(&bin, &nbin, i+10);
747 bin[i]++;
750 if (debug)
752 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
753 (int)step, epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
756 if (dump_pdb && epot <= dump_ener)
758 sprintf(str, "t%g_step%d.pdb", t, (int)step);
759 sprintf(str2, "t: %f step %d ener: %f", t, (int)step, epot);
760 write_sto_conf_mtop(str, str2, top_global, as_rvec_array(state_global->x.data()), as_rvec_array(state_global->v.data()),
761 inputrec->ePBC, state_global->box);
764 step++;
765 if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
767 /* Skip all steps assigned to the other MPI ranks */
768 step += (cr->nnodes - 1)*stepblocksize;
772 if (PAR(cr))
774 /* When running in parallel sum the energies over the processes */
775 gmx_sumd(1, &sum_embU, cr);
776 gmx_sumd(nener, sum_UgembU, cr);
779 frame++;
780 V_all += V;
781 VembU_all += V*sum_embU/nsteps;
783 if (fp_tpi)
785 if (mdrunOptions.verbose || frame%10 == 0 || frame < 10)
787 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
788 -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
791 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
793 VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
794 sum_embU == 0 ? 20/beta : -log(sum_embU/nsteps)/beta,
795 sum_embU/nsteps, V);
796 for (e = 0; e < nener; e++)
798 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
800 fprintf(fp_tpi, "\n");
801 fflush(fp_tpi);
804 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
805 } /* End of the loop */
806 walltime_accounting_end(walltime_accounting);
808 close_trx(status);
810 if (fp_tpi != nullptr)
812 xvgrclose(fp_tpi);
815 if (fplog != nullptr)
817 fprintf(fplog, "\n");
818 fprintf(fplog, " <V> = %12.5e nm^3\n", V_all/frame);
819 const double mu = -log(VembU_all/V_all)/beta;
820 fprintf(fplog, " <mu> = %12.5e kJ/mol\n", mu);
822 if (!std::isfinite(mu))
824 fprintf(fplog, "\nThe computed chemical potential is not finite - consider increasing the number of steps and/or the number of frames to insert into.\n");
828 /* Write the Boltzmann factor histogram */
829 if (PAR(cr))
831 /* When running in parallel sum the bins over the processes */
832 i = nbin;
833 global_max(cr, &i);
834 realloc_bins(&bin, &nbin, i);
835 gmx_sumd(nbin, bin, cr);
837 if (MASTER(cr))
839 fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
840 "TPI energy distribution",
841 "\\betaU - log(V/<V>)", "count", oenv);
842 sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
843 xvgr_subtitle(fp_tpi, str, oenv);
844 xvgr_legend(fp_tpi, 2, (const char **)tpid_leg, oenv);
845 for (i = nbin-1; i > 0; i--)
847 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
848 fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
849 bUlogV,
850 (int)(bin[i]+0.5),
851 bin[i]*exp(-bUlogV)*V_all/VembU_all);
853 xvgrclose(fp_tpi);
855 sfree(bin);
857 sfree(sum_UgembU);
859 walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);
862 } // namespace gmx