Restructure the load balancing timing
[gromacs.git] / src / gromacs / mdlib / tpi.cpp
blob0e99f3f59e7410de6ccc8a90e149b00f66116183
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, 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_mdlib
44 #include "gmxpre.h"
46 #include "tpi.h"
48 #include <cmath>
49 #include <cstdlib>
50 #include <cstring>
51 #include <ctime>
53 #include <algorithm>
55 #include "gromacs/commandline/filenm.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/ewald/pme.h"
58 #include "gromacs/fileio/confio.h"
59 #include "gromacs/fileio/trxio.h"
60 #include "gromacs/fileio/xvgr.h"
61 #include "gromacs/gmxlib/chargegroup.h"
62 #include "gromacs/gmxlib/conformation-utilities.h"
63 #include "gromacs/gmxlib/network.h"
64 #include "gromacs/gmxlib/nrnb.h"
65 #include "gromacs/math/units.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/mdlib/constr.h"
68 #include "gromacs/mdlib/force.h"
69 #include "gromacs/mdlib/mdatoms.h"
70 #include "gromacs/mdlib/mdebin.h"
71 #include "gromacs/mdlib/mdrun.h"
72 #include "gromacs/mdlib/ns.h"
73 #include "gromacs/mdlib/sim_util.h"
74 #include "gromacs/mdlib/tgroup.h"
75 #include "gromacs/mdlib/update.h"
76 #include "gromacs/mdlib/vsite.h"
77 #include "gromacs/mdtypes/commrec.h"
78 #include "gromacs/mdtypes/group.h"
79 #include "gromacs/mdtypes/inputrec.h"
80 #include "gromacs/mdtypes/md_enums.h"
81 #include "gromacs/mdtypes/state.h"
82 #include "gromacs/pbcutil/pbc.h"
83 #include "gromacs/random/threefry.h"
84 #include "gromacs/random/uniformrealdistribution.h"
85 #include "gromacs/timing/wallcycle.h"
86 #include "gromacs/timing/walltime_accounting.h"
87 #include "gromacs/topology/mtop_util.h"
88 #include "gromacs/trajectory/trajectoryframe.h"
89 #include "gromacs/utility/cstringutil.h"
90 #include "gromacs/utility/fatalerror.h"
91 #include "gromacs/utility/gmxassert.h"
92 #include "gromacs/utility/smalloc.h"
94 //! Global max algorithm
95 static void global_max(t_commrec *cr, int *n)
97 int *sum, i;
99 snew(sum, cr->nnodes);
100 sum[cr->nodeid] = *n;
101 gmx_sumi(cr->nnodes, sum, cr);
102 for (i = 0; i < cr->nnodes; i++)
104 *n = std::max(*n, sum[i]);
107 sfree(sum);
110 //! Reallocate arrays.
111 static void realloc_bins(double **bin, int *nbin, int nbin_new)
113 int i;
115 if (nbin_new != *nbin)
117 srenew(*bin, nbin_new);
118 for (i = *nbin; i < nbin_new; i++)
120 (*bin)[i] = 0;
122 *nbin = nbin_new;
126 namespace gmx
129 /*! \brief Do test particle insertion.
130 \copydoc integrator_t (FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
131 int nfile, const t_filenm fnm[],
132 const gmx_output_env_t *oenv, gmx_bool bVerbose,
133 int nstglobalcomm,
134 gmx_vsite_t *vsite, gmx_constr_t constr,
135 int stepout,
136 gmx::IMDOutputProvider *outputProvider,
137 t_inputrec *inputrec,
138 gmx_mtop_t *top_global, t_fcdata *fcd,
139 t_state *state_global,
140 t_mdatoms *mdatoms,
141 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
142 gmx_edsam_t ed,
143 t_forcerec *fr,
144 const ReplicaExchangeParameters &replExParams,
145 real cpt_period, real max_hours,
146 int imdport,
147 unsigned long Flags,
148 gmx_walltime_accounting_t walltime_accounting)
150 double do_tpi(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
151 int nfile, const t_filenm fnm[],
152 const gmx_output_env_t *oenv, gmx_bool bVerbose,
153 int gmx_unused nstglobalcomm,
154 gmx_vsite_t gmx_unused *vsite, gmx_constr_t gmx_unused constr,
155 int gmx_unused stepout,
156 gmx::IMDOutputProvider *outputProvider,
157 t_inputrec *inputrec,
158 gmx_mtop_t *top_global, t_fcdata *fcd,
159 t_state *state_global,
160 ObservablesHistory gmx_unused *observablesHistory,
161 t_mdatoms *mdatoms,
162 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
163 t_forcerec *fr,
164 const ReplicaExchangeParameters gmx_unused &replExParams,
165 gmx_membed_t gmx_unused *membed,
166 real gmx_unused cpt_period, real gmx_unused max_hours,
167 int gmx_unused imdport,
168 unsigned long gmx_unused Flags,
169 gmx_walltime_accounting_t walltime_accounting)
171 gmx_localtop_t *top;
172 gmx_groups_t *groups;
173 gmx_enerdata_t *enerd;
174 PaddedRVecVector f {};
175 real lambda, t, temp, beta, drmax, epot;
176 double embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
177 t_trxstatus *status;
178 t_trxframe rerun_fr;
179 gmx_bool bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
180 tensor force_vir, shake_vir, vir, pres;
181 int cg_tp, a_tp0, a_tp1, ngid, gid_tp, nener, e;
182 rvec *x_mol;
183 rvec mu_tot, x_init, dx, x_tp;
184 int nnodes, frame;
185 gmx_int64_t frame_step_prev, frame_step;
186 gmx_int64_t nsteps, stepblocksize = 0, step;
187 gmx_int64_t seed;
188 int i;
189 FILE *fp_tpi = nullptr;
190 char *ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
191 double dbl, dump_ener;
192 gmx_bool bCavity;
193 int nat_cavity = 0, d;
194 real *mass_cavity = nullptr, mass_tot;
195 int nbin;
196 double invbinw, *bin, refvolshift, logV, bUlogV;
197 real prescorr, enercorr, dvdlcorr;
198 gmx_bool bEnergyOutOfBounds;
199 const char *tpid_leg[2] = {"direct", "reweighted"};
201 GMX_UNUSED_VALUE(outputProvider);
203 /* Since there is no upper limit to the insertion energies,
204 * we need to set an upper limit for the distribution output.
206 real bU_bin_limit = 50;
207 real bU_logV_bin_limit = bU_bin_limit + 10;
209 if (inputrec->cutoff_scheme == ecutsVERLET)
211 gmx_fatal(FARGS, "TPI does not work (yet) with the Verlet cut-off scheme");
214 nnodes = cr->nnodes;
216 top = gmx_mtop_generate_local_top(top_global, inputrec->efep != efepNO);
218 groups = &top_global->groups;
220 bCavity = (inputrec->eI == eiTPIC);
221 if (bCavity)
223 ptr = getenv("GMX_TPIC_MASSES");
224 if (ptr == nullptr)
226 nat_cavity = 1;
228 else
230 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
231 * The center of mass of the last atoms is then used for TPIC.
233 nat_cavity = 0;
234 while (sscanf(ptr, "%20lf%n", &dbl, &i) > 0)
236 srenew(mass_cavity, nat_cavity+1);
237 mass_cavity[nat_cavity] = dbl;
238 fprintf(fplog, "mass[%d] = %f\n",
239 nat_cavity+1, mass_cavity[nat_cavity]);
240 nat_cavity++;
241 ptr += i;
243 if (nat_cavity == 0)
245 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
251 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
252 state_global->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
253 /* We never need full pbc for TPI */
254 fr->ePBC = epbcXYZ;
255 /* Determine the temperature for the Boltzmann weighting */
256 temp = inputrec->opts.ref_t[0];
257 if (fplog)
259 for (i = 1; (i < inputrec->opts.ngtc); i++)
261 if (inputrec->opts.ref_t[i] != temp)
263 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
264 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
267 fprintf(fplog,
268 "\n The temperature for test particle insertion is %.3f K\n\n",
269 temp);
271 beta = 1.0/(BOLTZ*temp);
273 /* Number of insertions per frame */
274 nsteps = inputrec->nsteps;
276 /* Use the same neighborlist with more insertions points
277 * in a sphere of radius drmax around the initial point
279 /* This should be a proper mdp parameter */
280 drmax = inputrec->rtpi;
282 /* An environment variable can be set to dump all configurations
283 * to pdb with an insertion energy <= this value.
285 dump_pdb = getenv("GMX_TPI_DUMP");
286 dump_ener = 0;
287 if (dump_pdb)
289 sscanf(dump_pdb, "%20lf", &dump_ener);
292 atoms2md(top_global, inputrec, -1, nullptr, top_global->natoms, mdatoms);
293 update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
295 snew(enerd, 1);
296 init_enerdata(groups->grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
297 /* We need to allocate one element extra, since we might use
298 * (unaligned) 4-wide SIMD loads to access rvec entries.
300 f.resize(top_global->natoms + 1);
302 /* Print to log file */
303 walltime_accounting_start(walltime_accounting);
304 wallcycle_start(wcycle, ewcRUN);
305 print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
307 /* The last charge group is the group to be inserted */
308 cg_tp = top->cgs.nr - 1;
309 a_tp0 = top->cgs.index[cg_tp];
310 a_tp1 = top->cgs.index[cg_tp+1];
311 if (debug)
313 fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
316 GMX_RELEASE_ASSERT(inputrec->rcoulomb <= inputrec->rlist && inputrec->rvdw <= inputrec->rlist, "Twin-range interactions are not supported with TPI");
318 snew(x_mol, a_tp1-a_tp0);
320 bDispCorr = (inputrec->eDispCorr != edispcNO);
321 bCharge = FALSE;
322 for (i = a_tp0; i < a_tp1; i++)
324 /* Copy the coordinates of the molecule to be insterted */
325 copy_rvec(state_global->x[i], x_mol[i-a_tp0]);
326 /* Check if we need to print electrostatic energies */
327 bCharge |= (mdatoms->chargeA[i] != 0 ||
328 (mdatoms->chargeB && mdatoms->chargeB[i] != 0));
330 bRFExcl = (bCharge && EEL_RF(fr->eeltype));
332 calc_cgcm(fplog, cg_tp, cg_tp+1, &(top->cgs), as_rvec_array(state_global->x.data()), fr->cg_cm);
333 if (bCavity)
335 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
337 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
338 fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
341 else
343 /* Center the molecule to be inserted at zero */
344 for (i = 0; i < a_tp1-a_tp0; i++)
346 rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
350 if (fplog)
352 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
353 a_tp1-a_tp0, bCharge ? "with" : "without");
355 fprintf(fplog, "\nWill insert %d times in each frame of %s\n",
356 (int)nsteps, opt2fn("-rerun", nfile, fnm));
359 if (!bCavity)
361 if (inputrec->nstlist > 1)
363 if (drmax == 0 && a_tp1-a_tp0 == 1)
365 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);
367 if (fplog)
369 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
373 else
375 if (fplog)
377 fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
381 ngid = groups->grps[egcENER].nr;
382 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
383 nener = 1 + ngid;
384 if (bDispCorr)
386 nener += 1;
388 if (bCharge)
390 nener += ngid;
391 if (bRFExcl)
393 nener += 1;
395 if (EEL_FULL(fr->eeltype))
397 nener += 1;
400 snew(sum_UgembU, nener);
402 /* Copy the random seed set by the user */
403 seed = inputrec->ld_seed;
405 gmx::ThreeFry2x64<16> rng(seed, gmx::RandomDomain::TestParticleInsertion); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
406 gmx::UniformRealDistribution<real> dist;
408 if (MASTER(cr))
410 fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
411 "TPI energies", "Time (ps)",
412 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
413 xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
414 snew(leg, 4+nener);
415 e = 0;
416 sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
417 leg[e++] = gmx_strdup(str);
418 sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
419 leg[e++] = gmx_strdup(str);
420 sprintf(str, "f. <e\\S-\\betaU\\N>");
421 leg[e++] = gmx_strdup(str);
422 sprintf(str, "f. V");
423 leg[e++] = gmx_strdup(str);
424 sprintf(str, "f. <Ue\\S-\\betaU\\N>");
425 leg[e++] = gmx_strdup(str);
426 for (i = 0; i < ngid; i++)
428 sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
429 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
430 leg[e++] = gmx_strdup(str);
432 if (bDispCorr)
434 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
435 leg[e++] = gmx_strdup(str);
437 if (bCharge)
439 for (i = 0; i < ngid; i++)
441 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
442 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
443 leg[e++] = gmx_strdup(str);
445 if (bRFExcl)
447 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
448 leg[e++] = gmx_strdup(str);
450 if (EEL_FULL(fr->eeltype))
452 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
453 leg[e++] = gmx_strdup(str);
456 xvgr_legend(fp_tpi, 4+nener, (const char**)leg, oenv);
457 for (i = 0; i < 4+nener; i++)
459 sfree(leg[i]);
461 sfree(leg);
463 clear_rvec(x_init);
464 V_all = 0;
465 VembU_all = 0;
467 invbinw = 10;
468 nbin = 10;
469 snew(bin, nbin);
471 /* Avoid frame step numbers <= -1 */
472 frame_step_prev = -1;
474 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
475 &rerun_fr, TRX_NEED_X);
476 frame = 0;
478 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
479 mdatoms->nr - (a_tp1 - a_tp0))
481 gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
482 "is not equal the number in the run input file (%d) "
483 "minus the number of atoms to insert (%d)\n",
484 rerun_fr.natoms, bCavity ? " minus one" : "",
485 mdatoms->nr, a_tp1-a_tp0);
488 refvolshift = log(det(rerun_fr.box));
490 switch (inputrec->eI)
492 case eiTPI:
493 stepblocksize = inputrec->nstlist;
494 break;
495 case eiTPIC:
496 stepblocksize = 1;
497 break;
498 default:
499 gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
502 while (bNotLastFrame)
504 frame_step = rerun_fr.step;
505 if (frame_step <= frame_step_prev)
507 /* We don't have step number in the trajectory file,
508 * or we have constant or decreasing step numbers.
509 * Ensure we have increasing step numbers, since we use
510 * the step numbers as a counter for random numbers.
512 frame_step = frame_step_prev + 1;
514 frame_step_prev = frame_step;
516 lambda = rerun_fr.lambda;
517 t = rerun_fr.time;
519 sum_embU = 0;
520 for (e = 0; e < nener; e++)
522 sum_UgembU[e] = 0;
525 /* Copy the coordinates from the input trajectory */
526 for (i = 0; i < rerun_fr.natoms; i++)
528 copy_rvec(rerun_fr.x[i], state_global->x[i]);
530 copy_mat(rerun_fr.box, state_global->box);
532 V = det(state_global->box);
533 logV = log(V);
535 bStateChanged = TRUE;
536 bNS = TRUE;
538 step = cr->nodeid*stepblocksize;
539 while (step < nsteps)
541 /* Restart random engine using the frame and insertion step
542 * as counters.
543 * Note that we need to draw several random values per iteration,
544 * but by using the internal subcounter functionality of ThreeFry2x64
545 * we can draw 131072 unique 64-bit values before exhausting
546 * the stream. This is a huge margin, and if something still goes
547 * wrong you will get an exception when the stream is exhausted.
549 rng.restart(frame_step, step);
550 dist.reset(); // erase any memory in the distribution
552 if (!bCavity)
554 /* Random insertion in the whole volume */
555 bNS = (step % inputrec->nstlist == 0);
556 if (bNS)
558 /* Generate a random position in the box */
559 for (d = 0; d < DIM; d++)
561 x_init[d] = dist(rng)*state_global->box[d][d];
565 if (inputrec->nstlist == 1)
567 copy_rvec(x_init, x_tp);
569 else
571 /* Generate coordinates within |dx|=drmax of x_init */
574 for (d = 0; d < DIM; d++)
576 dx[d] = (2*dist(rng) - 1)*drmax;
579 while (norm2(dx) > drmax*drmax);
580 rvec_add(x_init, dx, x_tp);
583 else
585 /* Random insertion around a cavity location
586 * given by the last coordinate of the trajectory.
588 if (step == 0)
590 if (nat_cavity == 1)
592 /* Copy the location of the cavity */
593 copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
595 else
597 /* Determine the center of mass of the last molecule */
598 clear_rvec(x_init);
599 mass_tot = 0;
600 for (i = 0; i < nat_cavity; i++)
602 for (d = 0; d < DIM; d++)
604 x_init[d] +=
605 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
607 mass_tot += mass_cavity[i];
609 for (d = 0; d < DIM; d++)
611 x_init[d] /= mass_tot;
615 /* Generate coordinates within |dx|=drmax of x_init */
618 for (d = 0; d < DIM; d++)
620 dx[d] = (2*dist(rng) - 1)*drmax;
623 while (norm2(dx) > drmax*drmax);
624 rvec_add(x_init, dx, x_tp);
627 if (a_tp1 - a_tp0 == 1)
629 /* Insert a single atom, just copy the insertion location */
630 copy_rvec(x_tp, state_global->x[a_tp0]);
632 else
634 /* Copy the coordinates from the top file */
635 for (i = a_tp0; i < a_tp1; i++)
637 copy_rvec(x_mol[i-a_tp0], state_global->x[i]);
639 /* Rotate the molecule randomly */
640 rotate_conf(a_tp1-a_tp0, as_rvec_array(state_global->x.data())+a_tp0, nullptr,
641 2*M_PI*dist(rng),
642 2*M_PI*dist(rng),
643 2*M_PI*dist(rng));
644 /* Shift to the insertion location */
645 for (i = a_tp0; i < a_tp1; i++)
647 rvec_inc(state_global->x[i], x_tp);
651 /* Clear some matrix variables */
652 clear_mat(force_vir);
653 clear_mat(shake_vir);
654 clear_mat(vir);
655 clear_mat(pres);
657 /* Set the charge group center of mass of the test particle */
658 copy_rvec(x_init, fr->cg_cm[top->cgs.nr-1]);
660 /* Calc energy (no forces) on new positions.
661 * Since we only need the intermolecular energy
662 * and the RF exclusion terms of the inserted molecule occur
663 * within a single charge group we can pass NULL for the graph.
664 * This also avoids shifts that would move charge groups
665 * out of the box. */
666 /* Make do_force do a single node force calculation */
667 cr->nnodes = 1;
668 do_force(fplog, cr, inputrec,
669 step, nrnb, wcycle, top, &top_global->groups,
670 state_global->box, &state_global->x, &state_global->hist,
671 &f, force_vir, mdatoms, enerd, fcd,
672 state_global->lambda,
673 nullptr, fr, nullptr, mu_tot, t, nullptr, FALSE,
674 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
675 (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS : 0) |
676 (bStateChanged ? GMX_FORCE_STATECHANGED : 0),
677 DdOpenBalanceRegionBeforeForceComputation::no,
678 DdCloseBalanceRegionAfterForceComputation::no);
679 cr->nnodes = nnodes;
680 bStateChanged = FALSE;
681 bNS = FALSE;
683 /* Calculate long range corrections to pressure and energy */
684 calc_dispcorr(inputrec, fr, state_global->box,
685 lambda, pres, vir, &prescorr, &enercorr, &dvdlcorr);
686 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
687 enerd->term[F_DISPCORR] = enercorr;
688 enerd->term[F_EPOT] += enercorr;
689 enerd->term[F_PRES] += prescorr;
690 enerd->term[F_DVDL_VDW] += dvdlcorr;
692 epot = enerd->term[F_EPOT];
693 bEnergyOutOfBounds = FALSE;
695 /* If the compiler doesn't optimize this check away
696 * we catch the NAN energies.
697 * The epot>GMX_REAL_MAX check catches inf values,
698 * which should nicely result in embU=0 through the exp below,
699 * but it does not hurt to check anyhow.
701 /* Non-bonded Interaction usually diverge at r=0.
702 * With tabulated interaction functions the first few entries
703 * should be capped in a consistent fashion between
704 * repulsion, dispersion and Coulomb to avoid accidental
705 * negative values in the total energy.
706 * The table generation code in tables.c does this.
707 * With user tbales the user should take care of this.
709 if (epot != epot || epot > GMX_REAL_MAX)
711 bEnergyOutOfBounds = TRUE;
713 if (bEnergyOutOfBounds)
715 if (debug)
717 fprintf(debug, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, (int)step, epot);
719 embU = 0;
721 else
723 embU = exp(-beta*epot);
724 sum_embU += embU;
725 /* Determine the weighted energy contributions of each energy group */
726 e = 0;
727 sum_UgembU[e++] += epot*embU;
728 if (fr->bBHAM)
730 for (i = 0; i < ngid; i++)
732 sum_UgembU[e++] +=
733 enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)]*embU;
736 else
738 for (i = 0; i < ngid; i++)
740 sum_UgembU[e++] +=
741 enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)]*embU;
744 if (bDispCorr)
746 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
748 if (bCharge)
750 for (i = 0; i < ngid; i++)
752 sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
754 if (bRFExcl)
756 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
758 if (EEL_FULL(fr->eeltype))
760 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
765 if (embU == 0 || beta*epot > bU_bin_limit)
767 bin[0]++;
769 else
771 i = (int)((bU_logV_bin_limit
772 - (beta*epot - logV + refvolshift))*invbinw
773 + 0.5);
774 if (i < 0)
776 i = 0;
778 if (i >= nbin)
780 realloc_bins(&bin, &nbin, i+10);
782 bin[i]++;
785 if (debug)
787 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
788 (int)step, epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
791 if (dump_pdb && epot <= dump_ener)
793 sprintf(str, "t%g_step%d.pdb", t, (int)step);
794 sprintf(str2, "t: %f step %d ener: %f", t, (int)step, epot);
795 write_sto_conf_mtop(str, str2, top_global, as_rvec_array(state_global->x.data()), as_rvec_array(state_global->v.data()),
796 inputrec->ePBC, state_global->box);
799 step++;
800 if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
802 /* Skip all steps assigned to the other MPI ranks */
803 step += (cr->nnodes - 1)*stepblocksize;
807 if (PAR(cr))
809 /* When running in parallel sum the energies over the processes */
810 gmx_sumd(1, &sum_embU, cr);
811 gmx_sumd(nener, sum_UgembU, cr);
814 frame++;
815 V_all += V;
816 VembU_all += V*sum_embU/nsteps;
818 if (fp_tpi)
820 if (bVerbose || frame%10 == 0 || frame < 10)
822 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
823 -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
826 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
828 VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
829 sum_embU == 0 ? 20/beta : -log(sum_embU/nsteps)/beta,
830 sum_embU/nsteps, V);
831 for (e = 0; e < nener; e++)
833 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
835 fprintf(fp_tpi, "\n");
836 fflush(fp_tpi);
839 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
840 } /* End of the loop */
841 walltime_accounting_end(walltime_accounting);
843 close_trx(status);
845 if (fp_tpi != nullptr)
847 xvgrclose(fp_tpi);
850 if (fplog != nullptr)
852 fprintf(fplog, "\n");
853 fprintf(fplog, " <V> = %12.5e nm^3\n", V_all/frame);
854 fprintf(fplog, " <mu> = %12.5e kJ/mol\n", -log(VembU_all/V_all)/beta);
857 /* Write the Boltzmann factor histogram */
858 if (PAR(cr))
860 /* When running in parallel sum the bins over the processes */
861 i = nbin;
862 global_max(cr, &i);
863 realloc_bins(&bin, &nbin, i);
864 gmx_sumd(nbin, bin, cr);
866 if (MASTER(cr))
868 fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
869 "TPI energy distribution",
870 "\\betaU - log(V/<V>)", "count", oenv);
871 sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
872 xvgr_subtitle(fp_tpi, str, oenv);
873 xvgr_legend(fp_tpi, 2, (const char **)tpid_leg, oenv);
874 for (i = nbin-1; i > 0; i--)
876 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
877 fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
878 bUlogV,
879 (int)(bin[i]+0.5),
880 bin[i]*exp(-bUlogV)*V_all/VembU_all);
882 xvgrclose(fp_tpi);
884 sfree(bin);
886 sfree(sum_UgembU);
888 walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);
890 return 0;
893 } // namespace gmx