Fix membed with partial revert of 29943f
[gromacs/AngularHB.git] / src / gromacs / mdlib / tpi.cpp
blobaa1210a59fdc9641dc89ea5424e16311e0162ebb
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, 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/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 //! Global max algorithm
94 static void global_max(t_commrec *cr, int *n)
96 int *sum, i;
98 snew(sum, cr->nnodes);
99 sum[cr->nodeid] = *n;
100 gmx_sumi(cr->nnodes, sum, cr);
101 for (i = 0; i < cr->nnodes; i++)
103 *n = std::max(*n, sum[i]);
106 sfree(sum);
109 //! Reallocate arrays.
110 static void realloc_bins(double **bin, int *nbin, int nbin_new)
112 int i;
114 if (nbin_new != *nbin)
116 srenew(*bin, nbin_new);
117 for (i = *nbin; i < nbin_new; i++)
119 (*bin)[i] = 0;
121 *nbin = nbin_new;
125 namespace gmx
128 /*! \brief Do test particle insertion.
129 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
130 int nfile, const t_filenm fnm[],
131 const gmx_output_env_t *oenv, gmx_bool bVerbose,
132 int nstglobalcomm,
133 gmx_vsite_t *vsite, gmx_constr_t constr,
134 int stepout,
135 t_inputrec *inputrec,
136 gmx_mtop_t *top_global, t_fcdata *fcd,
137 t_state *state_global,
138 t_mdatoms *mdatoms,
139 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
140 gmx_edsam_t ed,
141 t_forcerec *fr,
142 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
143 real cpt_period, real max_hours,
144 int imdport,
145 unsigned long Flags,
146 gmx_walltime_accounting_t walltime_accounting)
148 double do_tpi(FILE *fplog, t_commrec *cr,
149 int nfile, const t_filenm fnm[],
150 const gmx_output_env_t *oenv, gmx_bool bVerbose,
151 int gmx_unused nstglobalcomm,
152 gmx_vsite_t gmx_unused *vsite, gmx_constr_t gmx_unused constr,
153 int gmx_unused stepout,
154 t_inputrec *inputrec,
155 gmx_mtop_t *top_global, t_fcdata *fcd,
156 t_state *state_global,
157 t_mdatoms *mdatoms,
158 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
159 gmx_edsam_t gmx_unused ed,
160 t_forcerec *fr,
161 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
162 gmx_membed_t gmx_unused *membed,
163 real gmx_unused cpt_period, real gmx_unused max_hours,
164 int gmx_unused imdport,
165 unsigned long gmx_unused Flags,
166 gmx_walltime_accounting_t walltime_accounting)
168 gmx_localtop_t *top;
169 gmx_groups_t *groups;
170 gmx_enerdata_t *enerd;
171 rvec *f;
172 real lambda, t, temp, beta, drmax, epot;
173 double embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
174 t_trxstatus *status;
175 t_trxframe rerun_fr;
176 gmx_bool bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
177 tensor force_vir, shake_vir, vir, pres;
178 int cg_tp, a_tp0, a_tp1, ngid, gid_tp, nener, e;
179 rvec *x_mol;
180 rvec mu_tot, x_init, dx, x_tp;
181 int nnodes, frame;
182 gmx_int64_t frame_step_prev, frame_step;
183 gmx_int64_t nsteps, stepblocksize = 0, step;
184 gmx_int64_t seed;
185 int i;
186 FILE *fp_tpi = NULL;
187 char *ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
188 double dbl, dump_ener;
189 gmx_bool bCavity;
190 int nat_cavity = 0, d;
191 real *mass_cavity = NULL, mass_tot;
192 int nbin;
193 double invbinw, *bin, refvolshift, logV, bUlogV;
194 real prescorr, enercorr, dvdlcorr;
195 gmx_bool bEnergyOutOfBounds;
196 const char *tpid_leg[2] = {"direct", "reweighted"};
198 /* Since there is no upper limit to the insertion energies,
199 * we need to set an upper limit for the distribution output.
201 real bU_bin_limit = 50;
202 real bU_logV_bin_limit = bU_bin_limit + 10;
204 if (inputrec->cutoff_scheme == ecutsVERLET)
206 gmx_fatal(FARGS, "TPI does not work (yet) with the Verlet cut-off scheme");
209 nnodes = cr->nnodes;
211 top = gmx_mtop_generate_local_top(top_global, inputrec->efep != efepNO);
213 groups = &top_global->groups;
215 bCavity = (inputrec->eI == eiTPIC);
216 if (bCavity)
218 ptr = getenv("GMX_TPIC_MASSES");
219 if (ptr == NULL)
221 nat_cavity = 1;
223 else
225 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
226 * The center of mass of the last atoms is then used for TPIC.
228 nat_cavity = 0;
229 while (sscanf(ptr, "%20lf%n", &dbl, &i) > 0)
231 srenew(mass_cavity, nat_cavity+1);
232 mass_cavity[nat_cavity] = dbl;
233 fprintf(fplog, "mass[%d] = %f\n",
234 nat_cavity+1, mass_cavity[nat_cavity]);
235 nat_cavity++;
236 ptr += i;
238 if (nat_cavity == 0)
240 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
246 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
247 state_global->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
248 /* We never need full pbc for TPI */
249 fr->ePBC = epbcXYZ;
250 /* Determine the temperature for the Boltzmann weighting */
251 temp = inputrec->opts.ref_t[0];
252 if (fplog)
254 for (i = 1; (i < inputrec->opts.ngtc); i++)
256 if (inputrec->opts.ref_t[i] != temp)
258 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
259 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
262 fprintf(fplog,
263 "\n The temperature for test particle insertion is %.3f K\n\n",
264 temp);
266 beta = 1.0/(BOLTZ*temp);
268 /* Number of insertions per frame */
269 nsteps = inputrec->nsteps;
271 /* Use the same neighborlist with more insertions points
272 * in a sphere of radius drmax around the initial point
274 /* This should be a proper mdp parameter */
275 drmax = inputrec->rtpi;
277 /* An environment variable can be set to dump all configurations
278 * to pdb with an insertion energy <= this value.
280 dump_pdb = getenv("GMX_TPI_DUMP");
281 dump_ener = 0;
282 if (dump_pdb)
284 sscanf(dump_pdb, "%20lf", &dump_ener);
287 atoms2md(top_global, inputrec, 0, NULL, top_global->natoms, mdatoms);
288 update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
290 snew(enerd, 1);
291 init_enerdata(groups->grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
292 snew(f, top_global->natoms);
294 /* Print to log file */
295 walltime_accounting_start(walltime_accounting);
296 wallcycle_start(wcycle, ewcRUN);
297 print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
299 /* The last charge group is the group to be inserted */
300 cg_tp = top->cgs.nr - 1;
301 a_tp0 = top->cgs.index[cg_tp];
302 a_tp1 = top->cgs.index[cg_tp+1];
303 if (debug)
305 fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
308 GMX_RELEASE_ASSERT(inputrec->rcoulomb <= inputrec->rlist && inputrec->rvdw <= inputrec->rlist, "Twin-range interactions are not supported with TPI");
310 snew(x_mol, a_tp1-a_tp0);
312 bDispCorr = (inputrec->eDispCorr != edispcNO);
313 bCharge = FALSE;
314 for (i = a_tp0; i < a_tp1; i++)
316 /* Copy the coordinates of the molecule to be insterted */
317 copy_rvec(state_global->x[i], x_mol[i-a_tp0]);
318 /* Check if we need to print electrostatic energies */
319 bCharge |= (mdatoms->chargeA[i] != 0 ||
320 (mdatoms->chargeB && mdatoms->chargeB[i] != 0));
322 bRFExcl = (bCharge && EEL_RF(fr->eeltype));
324 calc_cgcm(fplog, cg_tp, cg_tp+1, &(top->cgs), state_global->x, fr->cg_cm);
325 if (bCavity)
327 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
329 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
330 fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
333 else
335 /* Center the molecule to be inserted at zero */
336 for (i = 0; i < a_tp1-a_tp0; i++)
338 rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
342 if (fplog)
344 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
345 a_tp1-a_tp0, bCharge ? "with" : "without");
347 fprintf(fplog, "\nWill insert %d times in each frame of %s\n",
348 (int)nsteps, opt2fn("-rerun", nfile, fnm));
351 if (!bCavity)
353 if (inputrec->nstlist > 1)
355 if (drmax == 0 && a_tp1-a_tp0 == 1)
357 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);
359 if (fplog)
361 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
365 else
367 if (fplog)
369 fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
373 ngid = groups->grps[egcENER].nr;
374 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
375 nener = 1 + ngid;
376 if (bDispCorr)
378 nener += 1;
380 if (bCharge)
382 nener += ngid;
383 if (bRFExcl)
385 nener += 1;
387 if (EEL_FULL(fr->eeltype))
389 nener += 1;
392 snew(sum_UgembU, nener);
394 /* Copy the random seed set by the user */
395 seed = inputrec->ld_seed;
397 gmx::ThreeFry2x64<16> rng(seed, gmx::RandomDomain::TestParticleInsertion); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
398 gmx::UniformRealDistribution<real> dist;
400 if (MASTER(cr))
402 fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
403 "TPI energies", "Time (ps)",
404 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
405 xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
406 snew(leg, 4+nener);
407 e = 0;
408 sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
409 leg[e++] = gmx_strdup(str);
410 sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
411 leg[e++] = gmx_strdup(str);
412 sprintf(str, "f. <e\\S-\\betaU\\N>");
413 leg[e++] = gmx_strdup(str);
414 sprintf(str, "f. V");
415 leg[e++] = gmx_strdup(str);
416 sprintf(str, "f. <Ue\\S-\\betaU\\N>");
417 leg[e++] = gmx_strdup(str);
418 for (i = 0; i < ngid; i++)
420 sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
421 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
422 leg[e++] = gmx_strdup(str);
424 if (bDispCorr)
426 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
427 leg[e++] = gmx_strdup(str);
429 if (bCharge)
431 for (i = 0; i < ngid; i++)
433 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
434 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
435 leg[e++] = gmx_strdup(str);
437 if (bRFExcl)
439 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
440 leg[e++] = gmx_strdup(str);
442 if (EEL_FULL(fr->eeltype))
444 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
445 leg[e++] = gmx_strdup(str);
448 xvgr_legend(fp_tpi, 4+nener, (const char**)leg, oenv);
449 for (i = 0; i < 4+nener; i++)
451 sfree(leg[i]);
453 sfree(leg);
455 clear_rvec(x_init);
456 V_all = 0;
457 VembU_all = 0;
459 invbinw = 10;
460 nbin = 10;
461 snew(bin, nbin);
463 /* Avoid frame step numbers <= -1 */
464 frame_step_prev = -1;
466 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
467 &rerun_fr, TRX_NEED_X);
468 frame = 0;
470 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
471 mdatoms->nr - (a_tp1 - a_tp0))
473 gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
474 "is not equal the number in the run input file (%d) "
475 "minus the number of atoms to insert (%d)\n",
476 rerun_fr.natoms, bCavity ? " minus one" : "",
477 mdatoms->nr, a_tp1-a_tp0);
480 refvolshift = log(det(rerun_fr.box));
482 switch (inputrec->eI)
484 case eiTPI:
485 stepblocksize = inputrec->nstlist;
486 break;
487 case eiTPIC:
488 stepblocksize = 1;
489 break;
490 default:
491 gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
494 while (bNotLastFrame)
496 frame_step = rerun_fr.step;
497 if (frame_step <= frame_step_prev)
499 /* We don't have step number in the trajectory file,
500 * or we have constant or decreasing step numbers.
501 * Ensure we have increasing step numbers, since we use
502 * the step numbers as a counter for random numbers.
504 frame_step = frame_step_prev + 1;
506 frame_step_prev = frame_step;
508 lambda = rerun_fr.lambda;
509 t = rerun_fr.time;
511 sum_embU = 0;
512 for (e = 0; e < nener; e++)
514 sum_UgembU[e] = 0;
517 /* Copy the coordinates from the input trajectory */
518 for (i = 0; i < rerun_fr.natoms; i++)
520 copy_rvec(rerun_fr.x[i], state_global->x[i]);
522 copy_mat(rerun_fr.box, state_global->box);
524 V = det(state_global->box);
525 logV = log(V);
527 bStateChanged = TRUE;
528 bNS = TRUE;
530 step = cr->nodeid*stepblocksize;
531 while (step < nsteps)
533 /* Restart random engine using the frame and insertion step
534 * as counters.
535 * Note that we need to draw several random values per iteration,
536 * but by using the internal subcounter functionality of ThreeFry2x64
537 * we can draw 131072 unique 64-bit values before exhausting
538 * the stream. This is a huge margin, and if something still goes
539 * wrong you will get an exception when the stream is exhausted.
541 rng.restart(frame_step, step);
542 dist.reset(); // erase any memory in the distribution
544 if (!bCavity)
546 /* Random insertion in the whole volume */
547 bNS = (step % inputrec->nstlist == 0);
548 if (bNS)
550 /* Generate a random position in the box */
551 for (d = 0; d < DIM; d++)
553 x_init[d] = dist(rng)*state_global->box[d][d];
557 if (inputrec->nstlist == 1)
559 copy_rvec(x_init, x_tp);
561 else
563 /* Generate coordinates within |dx|=drmax of x_init */
566 for (d = 0; d < DIM; d++)
568 dx[d] = (2*dist(rng) - 1)*drmax;
571 while (norm2(dx) > drmax*drmax);
572 rvec_add(x_init, dx, x_tp);
575 else
577 /* Random insertion around a cavity location
578 * given by the last coordinate of the trajectory.
580 if (step == 0)
582 if (nat_cavity == 1)
584 /* Copy the location of the cavity */
585 copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
587 else
589 /* Determine the center of mass of the last molecule */
590 clear_rvec(x_init);
591 mass_tot = 0;
592 for (i = 0; i < nat_cavity; i++)
594 for (d = 0; d < DIM; d++)
596 x_init[d] +=
597 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
599 mass_tot += mass_cavity[i];
601 for (d = 0; d < DIM; d++)
603 x_init[d] /= mass_tot;
607 /* Generate coordinates within |dx|=drmax of x_init */
610 for (d = 0; d < DIM; d++)
612 dx[d] = (2*dist(rng) - 1)*drmax;
615 while (norm2(dx) > drmax*drmax);
616 rvec_add(x_init, dx, x_tp);
619 if (a_tp1 - a_tp0 == 1)
621 /* Insert a single atom, just copy the insertion location */
622 copy_rvec(x_tp, state_global->x[a_tp0]);
624 else
626 /* Copy the coordinates from the top file */
627 for (i = a_tp0; i < a_tp1; i++)
629 copy_rvec(x_mol[i-a_tp0], state_global->x[i]);
631 /* Rotate the molecule randomly */
632 rotate_conf(a_tp1-a_tp0, state_global->x+a_tp0, NULL,
633 2*M_PI*dist(rng),
634 2*M_PI*dist(rng),
635 2*M_PI*dist(rng));
636 /* Shift to the insertion location */
637 for (i = a_tp0; i < a_tp1; i++)
639 rvec_inc(state_global->x[i], x_tp);
643 /* Clear some matrix variables */
644 clear_mat(force_vir);
645 clear_mat(shake_vir);
646 clear_mat(vir);
647 clear_mat(pres);
649 /* Set the charge group center of mass of the test particle */
650 copy_rvec(x_init, fr->cg_cm[top->cgs.nr-1]);
652 /* Calc energy (no forces) on new positions.
653 * Since we only need the intermolecular energy
654 * and the RF exclusion terms of the inserted molecule occur
655 * within a single charge group we can pass NULL for the graph.
656 * This also avoids shifts that would move charge groups
657 * out of the box. */
658 /* Make do_force do a single node force calculation */
659 cr->nnodes = 1;
660 do_force(fplog, cr, inputrec,
661 step, nrnb, wcycle, top, &top_global->groups,
662 state_global->box, state_global->x, &state_global->hist,
663 f, force_vir, mdatoms, enerd, fcd,
664 state_global->lambda,
665 NULL, fr, NULL, mu_tot, t, NULL, NULL, FALSE,
666 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
667 (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS : 0) |
668 (bStateChanged ? GMX_FORCE_STATECHANGED : 0));
669 cr->nnodes = nnodes;
670 bStateChanged = FALSE;
671 bNS = FALSE;
673 /* Calculate long range corrections to pressure and energy */
674 calc_dispcorr(inputrec, fr, state_global->box,
675 lambda, pres, vir, &prescorr, &enercorr, &dvdlcorr);
676 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
677 enerd->term[F_DISPCORR] = enercorr;
678 enerd->term[F_EPOT] += enercorr;
679 enerd->term[F_PRES] += prescorr;
680 enerd->term[F_DVDL_VDW] += dvdlcorr;
682 epot = enerd->term[F_EPOT];
683 bEnergyOutOfBounds = FALSE;
685 /* If the compiler doesn't optimize this check away
686 * we catch the NAN energies.
687 * The epot>GMX_REAL_MAX check catches inf values,
688 * which should nicely result in embU=0 through the exp below,
689 * but it does not hurt to check anyhow.
691 /* Non-bonded Interaction usually diverge at r=0.
692 * With tabulated interaction functions the first few entries
693 * should be capped in a consistent fashion between
694 * repulsion, dispersion and Coulomb to avoid accidental
695 * negative values in the total energy.
696 * The table generation code in tables.c does this.
697 * With user tbales the user should take care of this.
699 if (epot != epot || epot > GMX_REAL_MAX)
701 bEnergyOutOfBounds = TRUE;
703 if (bEnergyOutOfBounds)
705 if (debug)
707 fprintf(debug, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, (int)step, epot);
709 embU = 0;
711 else
713 embU = exp(-beta*epot);
714 sum_embU += embU;
715 /* Determine the weighted energy contributions of each energy group */
716 e = 0;
717 sum_UgembU[e++] += epot*embU;
718 if (fr->bBHAM)
720 for (i = 0; i < ngid; i++)
722 sum_UgembU[e++] +=
723 enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)]*embU;
726 else
728 for (i = 0; i < ngid; i++)
730 sum_UgembU[e++] +=
731 enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)]*embU;
734 if (bDispCorr)
736 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
738 if (bCharge)
740 for (i = 0; i < ngid; i++)
742 sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
744 if (bRFExcl)
746 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
748 if (EEL_FULL(fr->eeltype))
750 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
755 if (embU == 0 || beta*epot > bU_bin_limit)
757 bin[0]++;
759 else
761 i = (int)((bU_logV_bin_limit
762 - (beta*epot - logV + refvolshift))*invbinw
763 + 0.5);
764 if (i < 0)
766 i = 0;
768 if (i >= nbin)
770 realloc_bins(&bin, &nbin, i+10);
772 bin[i]++;
775 if (debug)
777 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
778 (int)step, epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
781 if (dump_pdb && epot <= dump_ener)
783 sprintf(str, "t%g_step%d.pdb", t, (int)step);
784 sprintf(str2, "t: %f step %d ener: %f", t, (int)step, epot);
785 write_sto_conf_mtop(str, str2, top_global, state_global->x, state_global->v,
786 inputrec->ePBC, state_global->box);
789 step++;
790 if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
792 /* Skip all steps assigned to the other MPI ranks */
793 step += (cr->nnodes - 1)*stepblocksize;
797 if (PAR(cr))
799 /* When running in parallel sum the energies over the processes */
800 gmx_sumd(1, &sum_embU, cr);
801 gmx_sumd(nener, sum_UgembU, cr);
804 frame++;
805 V_all += V;
806 VembU_all += V*sum_embU/nsteps;
808 if (fp_tpi)
810 if (bVerbose || frame%10 == 0 || frame < 10)
812 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
813 -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
816 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
818 VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
819 sum_embU == 0 ? 20/beta : -log(sum_embU/nsteps)/beta,
820 sum_embU/nsteps, V);
821 for (e = 0; e < nener; e++)
823 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
825 fprintf(fp_tpi, "\n");
826 fflush(fp_tpi);
829 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
830 } /* End of the loop */
831 walltime_accounting_end(walltime_accounting);
833 close_trj(status);
835 if (fp_tpi != NULL)
837 xvgrclose(fp_tpi);
840 if (fplog != NULL)
842 fprintf(fplog, "\n");
843 fprintf(fplog, " <V> = %12.5e nm^3\n", V_all/frame);
844 fprintf(fplog, " <mu> = %12.5e kJ/mol\n", -log(VembU_all/V_all)/beta);
847 /* Write the Boltzmann factor histogram */
848 if (PAR(cr))
850 /* When running in parallel sum the bins over the processes */
851 i = nbin;
852 global_max(cr, &i);
853 realloc_bins(&bin, &nbin, i);
854 gmx_sumd(nbin, bin, cr);
856 if (MASTER(cr))
858 fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
859 "TPI energy distribution",
860 "\\betaU - log(V/<V>)", "count", oenv);
861 sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
862 xvgr_subtitle(fp_tpi, str, oenv);
863 xvgr_legend(fp_tpi, 2, (const char **)tpid_leg, oenv);
864 for (i = nbin-1; i > 0; i--)
866 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
867 fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
868 bUlogV,
869 (int)(bin[i]+0.5),
870 bin[i]*exp(-bUlogV)*V_all/VembU_all);
872 xvgrclose(fp_tpi);
874 sfree(bin);
876 sfree(sum_UgembU);
878 walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);
880 return 0;
883 } // namespace gmx