Introduce ObservablesHistory container
[gromacs.git] / src / gromacs / mdlib / tpi.cpp
blobb365d98270b0546a585927c5701c1ea3a1d5c8ff
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 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
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 gmx_edsam_t gmx_unused ed,
164 t_forcerec *fr,
165 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
166 gmx_membed_t gmx_unused *membed,
167 real gmx_unused cpt_period, real gmx_unused max_hours,
168 int gmx_unused imdport,
169 unsigned long gmx_unused Flags,
170 gmx_walltime_accounting_t walltime_accounting)
172 gmx_localtop_t *top;
173 gmx_groups_t *groups;
174 gmx_enerdata_t *enerd;
175 PaddedRVecVector f {};
176 real lambda, t, temp, beta, drmax, epot;
177 double embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
178 t_trxstatus *status;
179 t_trxframe rerun_fr;
180 gmx_bool bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
181 tensor force_vir, shake_vir, vir, pres;
182 int cg_tp, a_tp0, a_tp1, ngid, gid_tp, nener, e;
183 rvec *x_mol;
184 rvec mu_tot, x_init, dx, x_tp;
185 int nnodes, frame;
186 gmx_int64_t frame_step_prev, frame_step;
187 gmx_int64_t nsteps, stepblocksize = 0, step;
188 gmx_int64_t seed;
189 int i;
190 FILE *fp_tpi = nullptr;
191 char *ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
192 double dbl, dump_ener;
193 gmx_bool bCavity;
194 int nat_cavity = 0, d;
195 real *mass_cavity = nullptr, mass_tot;
196 int nbin;
197 double invbinw, *bin, refvolshift, logV, bUlogV;
198 real prescorr, enercorr, dvdlcorr;
199 gmx_bool bEnergyOutOfBounds;
200 const char *tpid_leg[2] = {"direct", "reweighted"};
202 GMX_UNUSED_VALUE(outputProvider);
204 /* Since there is no upper limit to the insertion energies,
205 * we need to set an upper limit for the distribution output.
207 real bU_bin_limit = 50;
208 real bU_logV_bin_limit = bU_bin_limit + 10;
210 if (inputrec->cutoff_scheme == ecutsVERLET)
212 gmx_fatal(FARGS, "TPI does not work (yet) with the Verlet cut-off scheme");
215 nnodes = cr->nnodes;
217 top = gmx_mtop_generate_local_top(top_global, inputrec->efep != efepNO);
219 groups = &top_global->groups;
221 bCavity = (inputrec->eI == eiTPIC);
222 if (bCavity)
224 ptr = getenv("GMX_TPIC_MASSES");
225 if (ptr == nullptr)
227 nat_cavity = 1;
229 else
231 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
232 * The center of mass of the last atoms is then used for TPIC.
234 nat_cavity = 0;
235 while (sscanf(ptr, "%20lf%n", &dbl, &i) > 0)
237 srenew(mass_cavity, nat_cavity+1);
238 mass_cavity[nat_cavity] = dbl;
239 fprintf(fplog, "mass[%d] = %f\n",
240 nat_cavity+1, mass_cavity[nat_cavity]);
241 nat_cavity++;
242 ptr += i;
244 if (nat_cavity == 0)
246 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
252 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
253 state_global->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
254 /* We never need full pbc for TPI */
255 fr->ePBC = epbcXYZ;
256 /* Determine the temperature for the Boltzmann weighting */
257 temp = inputrec->opts.ref_t[0];
258 if (fplog)
260 for (i = 1; (i < inputrec->opts.ngtc); i++)
262 if (inputrec->opts.ref_t[i] != temp)
264 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
265 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
268 fprintf(fplog,
269 "\n The temperature for test particle insertion is %.3f K\n\n",
270 temp);
272 beta = 1.0/(BOLTZ*temp);
274 /* Number of insertions per frame */
275 nsteps = inputrec->nsteps;
277 /* Use the same neighborlist with more insertions points
278 * in a sphere of radius drmax around the initial point
280 /* This should be a proper mdp parameter */
281 drmax = inputrec->rtpi;
283 /* An environment variable can be set to dump all configurations
284 * to pdb with an insertion energy <= this value.
286 dump_pdb = getenv("GMX_TPI_DUMP");
287 dump_ener = 0;
288 if (dump_pdb)
290 sscanf(dump_pdb, "%20lf", &dump_ener);
293 atoms2md(top_global, inputrec, -1, nullptr, top_global->natoms, mdatoms);
294 update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
296 snew(enerd, 1);
297 init_enerdata(groups->grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
298 /* We need to allocate one element extra, since we might use
299 * (unaligned) 4-wide SIMD loads to access rvec entries.
301 f.resize(top_global->natoms + 1);
303 /* Print to log file */
304 walltime_accounting_start(walltime_accounting);
305 wallcycle_start(wcycle, ewcRUN);
306 print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
308 /* The last charge group is the group to be inserted */
309 cg_tp = top->cgs.nr - 1;
310 a_tp0 = top->cgs.index[cg_tp];
311 a_tp1 = top->cgs.index[cg_tp+1];
312 if (debug)
314 fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
317 GMX_RELEASE_ASSERT(inputrec->rcoulomb <= inputrec->rlist && inputrec->rvdw <= inputrec->rlist, "Twin-range interactions are not supported with TPI");
319 snew(x_mol, a_tp1-a_tp0);
321 bDispCorr = (inputrec->eDispCorr != edispcNO);
322 bCharge = FALSE;
323 for (i = a_tp0; i < a_tp1; i++)
325 /* Copy the coordinates of the molecule to be insterted */
326 copy_rvec(state_global->x[i], x_mol[i-a_tp0]);
327 /* Check if we need to print electrostatic energies */
328 bCharge |= (mdatoms->chargeA[i] != 0 ||
329 (mdatoms->chargeB && mdatoms->chargeB[i] != 0));
331 bRFExcl = (bCharge && EEL_RF(fr->eeltype));
333 calc_cgcm(fplog, cg_tp, cg_tp+1, &(top->cgs), as_rvec_array(state_global->x.data()), fr->cg_cm);
334 if (bCavity)
336 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
338 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
339 fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
342 else
344 /* Center the molecule to be inserted at zero */
345 for (i = 0; i < a_tp1-a_tp0; i++)
347 rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
351 if (fplog)
353 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
354 a_tp1-a_tp0, bCharge ? "with" : "without");
356 fprintf(fplog, "\nWill insert %d times in each frame of %s\n",
357 (int)nsteps, opt2fn("-rerun", nfile, fnm));
360 if (!bCavity)
362 if (inputrec->nstlist > 1)
364 if (drmax == 0 && a_tp1-a_tp0 == 1)
366 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);
368 if (fplog)
370 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
374 else
376 if (fplog)
378 fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
382 ngid = groups->grps[egcENER].nr;
383 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
384 nener = 1 + ngid;
385 if (bDispCorr)
387 nener += 1;
389 if (bCharge)
391 nener += ngid;
392 if (bRFExcl)
394 nener += 1;
396 if (EEL_FULL(fr->eeltype))
398 nener += 1;
401 snew(sum_UgembU, nener);
403 /* Copy the random seed set by the user */
404 seed = inputrec->ld_seed;
406 gmx::ThreeFry2x64<16> rng(seed, gmx::RandomDomain::TestParticleInsertion); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
407 gmx::UniformRealDistribution<real> dist;
409 if (MASTER(cr))
411 fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
412 "TPI energies", "Time (ps)",
413 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
414 xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
415 snew(leg, 4+nener);
416 e = 0;
417 sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
418 leg[e++] = gmx_strdup(str);
419 sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
420 leg[e++] = gmx_strdup(str);
421 sprintf(str, "f. <e\\S-\\betaU\\N>");
422 leg[e++] = gmx_strdup(str);
423 sprintf(str, "f. V");
424 leg[e++] = gmx_strdup(str);
425 sprintf(str, "f. <Ue\\S-\\betaU\\N>");
426 leg[e++] = gmx_strdup(str);
427 for (i = 0; i < ngid; i++)
429 sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
430 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
431 leg[e++] = gmx_strdup(str);
433 if (bDispCorr)
435 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
436 leg[e++] = gmx_strdup(str);
438 if (bCharge)
440 for (i = 0; i < ngid; i++)
442 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
443 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
444 leg[e++] = gmx_strdup(str);
446 if (bRFExcl)
448 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
449 leg[e++] = gmx_strdup(str);
451 if (EEL_FULL(fr->eeltype))
453 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
454 leg[e++] = gmx_strdup(str);
457 xvgr_legend(fp_tpi, 4+nener, (const char**)leg, oenv);
458 for (i = 0; i < 4+nener; i++)
460 sfree(leg[i]);
462 sfree(leg);
464 clear_rvec(x_init);
465 V_all = 0;
466 VembU_all = 0;
468 invbinw = 10;
469 nbin = 10;
470 snew(bin, nbin);
472 /* Avoid frame step numbers <= -1 */
473 frame_step_prev = -1;
475 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
476 &rerun_fr, TRX_NEED_X);
477 frame = 0;
479 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
480 mdatoms->nr - (a_tp1 - a_tp0))
482 gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
483 "is not equal the number in the run input file (%d) "
484 "minus the number of atoms to insert (%d)\n",
485 rerun_fr.natoms, bCavity ? " minus one" : "",
486 mdatoms->nr, a_tp1-a_tp0);
489 refvolshift = log(det(rerun_fr.box));
491 switch (inputrec->eI)
493 case eiTPI:
494 stepblocksize = inputrec->nstlist;
495 break;
496 case eiTPIC:
497 stepblocksize = 1;
498 break;
499 default:
500 gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
503 while (bNotLastFrame)
505 frame_step = rerun_fr.step;
506 if (frame_step <= frame_step_prev)
508 /* We don't have step number in the trajectory file,
509 * or we have constant or decreasing step numbers.
510 * Ensure we have increasing step numbers, since we use
511 * the step numbers as a counter for random numbers.
513 frame_step = frame_step_prev + 1;
515 frame_step_prev = frame_step;
517 lambda = rerun_fr.lambda;
518 t = rerun_fr.time;
520 sum_embU = 0;
521 for (e = 0; e < nener; e++)
523 sum_UgembU[e] = 0;
526 /* Copy the coordinates from the input trajectory */
527 for (i = 0; i < rerun_fr.natoms; i++)
529 copy_rvec(rerun_fr.x[i], state_global->x[i]);
531 copy_mat(rerun_fr.box, state_global->box);
533 V = det(state_global->box);
534 logV = log(V);
536 bStateChanged = TRUE;
537 bNS = TRUE;
539 step = cr->nodeid*stepblocksize;
540 while (step < nsteps)
542 /* Restart random engine using the frame and insertion step
543 * as counters.
544 * Note that we need to draw several random values per iteration,
545 * but by using the internal subcounter functionality of ThreeFry2x64
546 * we can draw 131072 unique 64-bit values before exhausting
547 * the stream. This is a huge margin, and if something still goes
548 * wrong you will get an exception when the stream is exhausted.
550 rng.restart(frame_step, step);
551 dist.reset(); // erase any memory in the distribution
553 if (!bCavity)
555 /* Random insertion in the whole volume */
556 bNS = (step % inputrec->nstlist == 0);
557 if (bNS)
559 /* Generate a random position in the box */
560 for (d = 0; d < DIM; d++)
562 x_init[d] = dist(rng)*state_global->box[d][d];
566 if (inputrec->nstlist == 1)
568 copy_rvec(x_init, x_tp);
570 else
572 /* Generate coordinates within |dx|=drmax of x_init */
575 for (d = 0; d < DIM; d++)
577 dx[d] = (2*dist(rng) - 1)*drmax;
580 while (norm2(dx) > drmax*drmax);
581 rvec_add(x_init, dx, x_tp);
584 else
586 /* Random insertion around a cavity location
587 * given by the last coordinate of the trajectory.
589 if (step == 0)
591 if (nat_cavity == 1)
593 /* Copy the location of the cavity */
594 copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
596 else
598 /* Determine the center of mass of the last molecule */
599 clear_rvec(x_init);
600 mass_tot = 0;
601 for (i = 0; i < nat_cavity; i++)
603 for (d = 0; d < DIM; d++)
605 x_init[d] +=
606 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
608 mass_tot += mass_cavity[i];
610 for (d = 0; d < DIM; d++)
612 x_init[d] /= mass_tot;
616 /* Generate coordinates within |dx|=drmax of x_init */
619 for (d = 0; d < DIM; d++)
621 dx[d] = (2*dist(rng) - 1)*drmax;
624 while (norm2(dx) > drmax*drmax);
625 rvec_add(x_init, dx, x_tp);
628 if (a_tp1 - a_tp0 == 1)
630 /* Insert a single atom, just copy the insertion location */
631 copy_rvec(x_tp, state_global->x[a_tp0]);
633 else
635 /* Copy the coordinates from the top file */
636 for (i = a_tp0; i < a_tp1; i++)
638 copy_rvec(x_mol[i-a_tp0], state_global->x[i]);
640 /* Rotate the molecule randomly */
641 rotate_conf(a_tp1-a_tp0, as_rvec_array(state_global->x.data())+a_tp0, nullptr,
642 2*M_PI*dist(rng),
643 2*M_PI*dist(rng),
644 2*M_PI*dist(rng));
645 /* Shift to the insertion location */
646 for (i = a_tp0; i < a_tp1; i++)
648 rvec_inc(state_global->x[i], x_tp);
652 /* Clear some matrix variables */
653 clear_mat(force_vir);
654 clear_mat(shake_vir);
655 clear_mat(vir);
656 clear_mat(pres);
658 /* Set the charge group center of mass of the test particle */
659 copy_rvec(x_init, fr->cg_cm[top->cgs.nr-1]);
661 /* Calc energy (no forces) on new positions.
662 * Since we only need the intermolecular energy
663 * and the RF exclusion terms of the inserted molecule occur
664 * within a single charge group we can pass NULL for the graph.
665 * This also avoids shifts that would move charge groups
666 * out of the box. */
667 /* Make do_force do a single node force calculation */
668 cr->nnodes = 1;
669 do_force(fplog, cr, inputrec,
670 step, nrnb, wcycle, top, &top_global->groups,
671 state_global->box, &state_global->x, &state_global->hist,
672 &f, force_vir, mdatoms, enerd, fcd,
673 state_global->lambda,
674 nullptr, fr, nullptr, mu_tot, t, nullptr, FALSE,
675 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
676 (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS : 0) |
677 (bStateChanged ? GMX_FORCE_STATECHANGED : 0));
678 cr->nnodes = nnodes;
679 bStateChanged = FALSE;
680 bNS = FALSE;
682 /* Calculate long range corrections to pressure and energy */
683 calc_dispcorr(inputrec, fr, state_global->box,
684 lambda, pres, vir, &prescorr, &enercorr, &dvdlcorr);
685 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
686 enerd->term[F_DISPCORR] = enercorr;
687 enerd->term[F_EPOT] += enercorr;
688 enerd->term[F_PRES] += prescorr;
689 enerd->term[F_DVDL_VDW] += dvdlcorr;
691 epot = enerd->term[F_EPOT];
692 bEnergyOutOfBounds = FALSE;
694 /* If the compiler doesn't optimize this check away
695 * we catch the NAN energies.
696 * The epot>GMX_REAL_MAX check catches inf values,
697 * which should nicely result in embU=0 through the exp below,
698 * but it does not hurt to check anyhow.
700 /* Non-bonded Interaction usually diverge at r=0.
701 * With tabulated interaction functions the first few entries
702 * should be capped in a consistent fashion between
703 * repulsion, dispersion and Coulomb to avoid accidental
704 * negative values in the total energy.
705 * The table generation code in tables.c does this.
706 * With user tbales the user should take care of this.
708 if (epot != epot || epot > GMX_REAL_MAX)
710 bEnergyOutOfBounds = TRUE;
712 if (bEnergyOutOfBounds)
714 if (debug)
716 fprintf(debug, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, (int)step, epot);
718 embU = 0;
720 else
722 embU = exp(-beta*epot);
723 sum_embU += embU;
724 /* Determine the weighted energy contributions of each energy group */
725 e = 0;
726 sum_UgembU[e++] += epot*embU;
727 if (fr->bBHAM)
729 for (i = 0; i < ngid; i++)
731 sum_UgembU[e++] +=
732 enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)]*embU;
735 else
737 for (i = 0; i < ngid; i++)
739 sum_UgembU[e++] +=
740 enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)]*embU;
743 if (bDispCorr)
745 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
747 if (bCharge)
749 for (i = 0; i < ngid; i++)
751 sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
753 if (bRFExcl)
755 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
757 if (EEL_FULL(fr->eeltype))
759 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
764 if (embU == 0 || beta*epot > bU_bin_limit)
766 bin[0]++;
768 else
770 i = (int)((bU_logV_bin_limit
771 - (beta*epot - logV + refvolshift))*invbinw
772 + 0.5);
773 if (i < 0)
775 i = 0;
777 if (i >= nbin)
779 realloc_bins(&bin, &nbin, i+10);
781 bin[i]++;
784 if (debug)
786 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
787 (int)step, epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
790 if (dump_pdb && epot <= dump_ener)
792 sprintf(str, "t%g_step%d.pdb", t, (int)step);
793 sprintf(str2, "t: %f step %d ener: %f", t, (int)step, epot);
794 write_sto_conf_mtop(str, str2, top_global, as_rvec_array(state_global->x.data()), as_rvec_array(state_global->v.data()),
795 inputrec->ePBC, state_global->box);
798 step++;
799 if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
801 /* Skip all steps assigned to the other MPI ranks */
802 step += (cr->nnodes - 1)*stepblocksize;
806 if (PAR(cr))
808 /* When running in parallel sum the energies over the processes */
809 gmx_sumd(1, &sum_embU, cr);
810 gmx_sumd(nener, sum_UgembU, cr);
813 frame++;
814 V_all += V;
815 VembU_all += V*sum_embU/nsteps;
817 if (fp_tpi)
819 if (bVerbose || frame%10 == 0 || frame < 10)
821 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
822 -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
825 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
827 VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
828 sum_embU == 0 ? 20/beta : -log(sum_embU/nsteps)/beta,
829 sum_embU/nsteps, V);
830 for (e = 0; e < nener; e++)
832 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
834 fprintf(fp_tpi, "\n");
835 fflush(fp_tpi);
838 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
839 } /* End of the loop */
840 walltime_accounting_end(walltime_accounting);
842 close_trx(status);
844 if (fp_tpi != nullptr)
846 xvgrclose(fp_tpi);
849 if (fplog != nullptr)
851 fprintf(fplog, "\n");
852 fprintf(fplog, " <V> = %12.5e nm^3\n", V_all/frame);
853 fprintf(fplog, " <mu> = %12.5e kJ/mol\n", -log(VembU_all/V_all)/beta);
856 /* Write the Boltzmann factor histogram */
857 if (PAR(cr))
859 /* When running in parallel sum the bins over the processes */
860 i = nbin;
861 global_max(cr, &i);
862 realloc_bins(&bin, &nbin, i);
863 gmx_sumd(nbin, bin, cr);
865 if (MASTER(cr))
867 fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
868 "TPI energy distribution",
869 "\\betaU - log(V/<V>)", "count", oenv);
870 sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
871 xvgr_subtitle(fp_tpi, str, oenv);
872 xvgr_legend(fp_tpi, 2, (const char **)tpid_leg, oenv);
873 for (i = nbin-1; i > 0; i--)
875 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
876 fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
877 bUlogV,
878 (int)(bin[i]+0.5),
879 bin[i]*exp(-bUlogV)*V_all/VembU_all);
881 xvgrclose(fp_tpi);
883 sfree(bin);
885 sfree(sum_UgembU);
887 walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);
889 return 0;
892 } // namespace gmx