Split txtdump.*, part 1
[gromacs.git] / src / gromacs / mdlib / tpi.cpp
blobb94c823d6e35c156917a1fe845439afd9d8b63e3
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, 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/trx.h"
60 #include "gromacs/fileio/trxio.h"
61 #include "gromacs/fileio/xvgr.h"
62 #include "gromacs/gmxlib/chargegroup.h"
63 #include "gromacs/gmxlib/conformation-utilities.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/gmxlib/nrnb.h"
66 #include "gromacs/math/units.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/mdlib/constr.h"
69 #include "gromacs/mdlib/force.h"
70 #include "gromacs/mdlib/mdatoms.h"
71 #include "gromacs/mdlib/mdebin.h"
72 #include "gromacs/mdlib/mdrun.h"
73 #include "gromacs/mdlib/ns.h"
74 #include "gromacs/mdlib/sim_util.h"
75 #include "gromacs/mdlib/tgroup.h"
76 #include "gromacs/mdlib/update.h"
77 #include "gromacs/mdlib/vsite.h"
78 #include "gromacs/mdtypes/commrec.h"
79 #include "gromacs/mdtypes/group.h"
80 #include "gromacs/mdtypes/md_enums.h"
81 #include "gromacs/pbcutil/pbc.h"
82 #include "gromacs/random/random.h"
83 #include "gromacs/timing/wallcycle.h"
84 #include "gromacs/timing/walltime_accounting.h"
85 #include "gromacs/topology/mtop_util.h"
86 #include "gromacs/utility/cstringutil.h"
87 #include "gromacs/utility/fatalerror.h"
88 #include "gromacs/utility/gmxassert.h"
89 #include "gromacs/utility/smalloc.h"
91 //! Global max algorithm
92 static void global_max(t_commrec *cr, int *n)
94 int *sum, i;
96 snew(sum, cr->nnodes);
97 sum[cr->nodeid] = *n;
98 gmx_sumi(cr->nnodes, sum, cr);
99 for (i = 0; i < cr->nnodes; i++)
101 *n = std::max(*n, sum[i]);
104 sfree(sum);
107 //! Reallocate arrays.
108 static void realloc_bins(double **bin, int *nbin, int nbin_new)
110 int i;
112 if (nbin_new != *nbin)
114 srenew(*bin, nbin_new);
115 for (i = *nbin; i < nbin_new; i++)
117 (*bin)[i] = 0;
119 *nbin = nbin_new;
123 //! Workaround to keep doxygen from generating warnings
124 static const gmx_uint64_t rnd_seed_tpi = RND_SEED_TPI;
126 namespace gmx
129 /*! \brief Do test particle insertion.
130 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
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 t_inputrec *inputrec,
137 gmx_mtop_t *top_global, t_fcdata *fcd,
138 t_state *state_global,
139 t_mdatoms *mdatoms,
140 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
141 gmx_edsam_t ed,
142 t_forcerec *fr,
143 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
144 gmx_membed_t *membed,
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,
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 t_inputrec *inputrec,
157 gmx_mtop_t *top_global, t_fcdata *fcd,
158 t_state *state_global,
159 t_mdatoms *mdatoms,
160 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
161 gmx_edsam_t gmx_unused ed,
162 t_forcerec *fr,
163 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
164 gmx_membed_t *membed,
165 real gmx_unused cpt_period, real gmx_unused max_hours,
166 int gmx_unused imdport,
167 unsigned long gmx_unused Flags,
168 gmx_walltime_accounting_t walltime_accounting)
170 gmx_localtop_t *top;
171 gmx_groups_t *groups;
172 gmx_enerdata_t *enerd;
173 rvec *f;
174 real lambda, t, temp, beta, drmax, epot;
175 double embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
176 t_trxstatus *status;
177 t_trxframe rerun_fr;
178 gmx_bool bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
179 tensor force_vir, shake_vir, vir, pres;
180 int cg_tp, a_tp0, a_tp1, ngid, gid_tp, nener, e;
181 rvec *x_mol;
182 rvec mu_tot, x_init, dx, x_tp;
183 int nnodes, frame;
184 gmx_int64_t frame_step_prev, frame_step;
185 gmx_int64_t nsteps, stepblocksize = 0, step;
186 gmx_int64_t rnd_count_stride, rnd_count;
187 gmx_int64_t seed;
188 double rnd[4];
189 int i;
190 FILE *fp_tpi = NULL;
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 = NULL, 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(membed);
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 == NULL)
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, 0, NULL, 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 snew(f, top_global->natoms);
299 /* Print to log file */
300 walltime_accounting_start(walltime_accounting);
301 wallcycle_start(wcycle, ewcRUN);
302 print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
304 /* The last charge group is the group to be inserted */
305 cg_tp = top->cgs.nr - 1;
306 a_tp0 = top->cgs.index[cg_tp];
307 a_tp1 = top->cgs.index[cg_tp+1];
308 if (debug)
310 fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
313 GMX_RELEASE_ASSERT(inputrec->rcoulomb <= inputrec->rlist && inputrec->rvdw <= inputrec->rlist, "Twin-range interactions are not supported with TPI");
315 snew(x_mol, a_tp1-a_tp0);
317 bDispCorr = (inputrec->eDispCorr != edispcNO);
318 bCharge = FALSE;
319 for (i = a_tp0; i < a_tp1; i++)
321 /* Copy the coordinates of the molecule to be insterted */
322 copy_rvec(state_global->x[i], x_mol[i-a_tp0]);
323 /* Check if we need to print electrostatic energies */
324 bCharge |= (mdatoms->chargeA[i] != 0 ||
325 (mdatoms->chargeB && mdatoms->chargeB[i] != 0));
327 bRFExcl = (bCharge && EEL_RF(fr->eeltype));
329 calc_cgcm(fplog, cg_tp, cg_tp+1, &(top->cgs), state_global->x, fr->cg_cm);
330 if (bCavity)
332 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
334 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
335 fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
338 else
340 /* Center the molecule to be inserted at zero */
341 for (i = 0; i < a_tp1-a_tp0; i++)
343 rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
347 if (fplog)
349 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
350 a_tp1-a_tp0, bCharge ? "with" : "without");
352 fprintf(fplog, "\nWill insert %d times in each frame of %s\n",
353 (int)nsteps, opt2fn("-rerun", nfile, fnm));
356 if (!bCavity)
358 if (inputrec->nstlist > 1)
360 if (drmax == 0 && a_tp1-a_tp0 == 1)
362 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);
364 if (fplog)
366 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
370 else
372 if (fplog)
374 fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
378 ngid = groups->grps[egcENER].nr;
379 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
380 nener = 1 + ngid;
381 if (bDispCorr)
383 nener += 1;
385 if (bCharge)
387 nener += ngid;
388 if (bRFExcl)
390 nener += 1;
392 if (EEL_FULL(fr->eeltype))
394 nener += 1;
397 snew(sum_UgembU, nener);
399 /* Copy the random seed set by the user */
400 seed = inputrec->ld_seed;
401 /* We use the frame step number as one random counter.
402 * The second counter use the insertion (step) count. But we
403 * need multiple random numbers per insertion. This number is
404 * not fixed, since we generate random locations in a sphere
405 * by putting locations in a cube and some of these fail.
406 * A count of 20 is already extremely unlikely, so 10000 is
407 * a safe margin for random numbers per insertion.
409 rnd_count_stride = 10000;
411 if (MASTER(cr))
413 fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
414 "TPI energies", "Time (ps)",
415 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
416 xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
417 snew(leg, 4+nener);
418 e = 0;
419 sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
420 leg[e++] = gmx_strdup(str);
421 sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
422 leg[e++] = gmx_strdup(str);
423 sprintf(str, "f. <e\\S-\\betaU\\N>");
424 leg[e++] = gmx_strdup(str);
425 sprintf(str, "f. V");
426 leg[e++] = gmx_strdup(str);
427 sprintf(str, "f. <Ue\\S-\\betaU\\N>");
428 leg[e++] = gmx_strdup(str);
429 for (i = 0; i < ngid; i++)
431 sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
432 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
433 leg[e++] = gmx_strdup(str);
435 if (bDispCorr)
437 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
438 leg[e++] = gmx_strdup(str);
440 if (bCharge)
442 for (i = 0; i < ngid; i++)
444 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
445 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
446 leg[e++] = gmx_strdup(str);
448 if (bRFExcl)
450 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
451 leg[e++] = gmx_strdup(str);
453 if (EEL_FULL(fr->eeltype))
455 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
456 leg[e++] = gmx_strdup(str);
459 xvgr_legend(fp_tpi, 4+nener, (const char**)leg, oenv);
460 for (i = 0; i < 4+nener; i++)
462 sfree(leg[i]);
464 sfree(leg);
466 clear_rvec(x_init);
467 V_all = 0;
468 VembU_all = 0;
470 invbinw = 10;
471 nbin = 10;
472 snew(bin, nbin);
474 /* Avoid frame step numbers <= -1 */
475 frame_step_prev = -1;
477 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
478 &rerun_fr, TRX_NEED_X);
479 frame = 0;
481 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
482 mdatoms->nr - (a_tp1 - a_tp0))
484 gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
485 "is not equal the number in the run input file (%d) "
486 "minus the number of atoms to insert (%d)\n",
487 rerun_fr.natoms, bCavity ? " minus one" : "",
488 mdatoms->nr, a_tp1-a_tp0);
491 refvolshift = log(det(rerun_fr.box));
493 switch (inputrec->eI)
495 case eiTPI:
496 stepblocksize = inputrec->nstlist;
497 break;
498 case eiTPIC:
499 stepblocksize = 1;
500 break;
501 default:
502 gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
505 while (bNotLastFrame)
507 frame_step = rerun_fr.step;
508 if (frame_step <= frame_step_prev)
510 /* We don't have step number in the trajectory file,
511 * or we have constant or decreasing step numbers.
512 * Ensure we have increasing step numbers, since we use
513 * the step numbers as a counter for random numbers.
515 frame_step = frame_step_prev + 1;
517 frame_step_prev = frame_step;
519 lambda = rerun_fr.lambda;
520 t = rerun_fr.time;
522 sum_embU = 0;
523 for (e = 0; e < nener; e++)
525 sum_UgembU[e] = 0;
528 /* Copy the coordinates from the input trajectory */
529 for (i = 0; i < rerun_fr.natoms; i++)
531 copy_rvec(rerun_fr.x[i], state_global->x[i]);
533 copy_mat(rerun_fr.box, state_global->box);
535 V = det(state_global->box);
536 logV = log(V);
538 bStateChanged = TRUE;
539 bNS = TRUE;
541 step = cr->nodeid*stepblocksize;
542 while (step < nsteps)
544 /* Initialize the second counter for random numbers using
545 * the insertion step index. This ensures that we get
546 * the same random numbers independently of how many
547 * MPI ranks we use. Also for the same seed, we get
548 * the same initial random sequence for different nsteps.
550 rnd_count = step*rnd_count_stride;
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 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, rnd_seed_tpi, rnd);
560 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, rnd_seed_tpi, rnd+2);
561 for (d = 0; d < DIM; d++)
563 x_init[d] = rnd[d]*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 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, rnd_seed_tpi, rnd);
576 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, rnd_seed_tpi, rnd+2);
577 for (d = 0; d < DIM; d++)
579 dx[d] = (2*rnd[d] - 1)*drmax;
582 while (norm2(dx) > drmax*drmax);
583 rvec_add(x_init, dx, x_tp);
586 else
588 /* Random insertion around a cavity location
589 * given by the last coordinate of the trajectory.
591 if (step == 0)
593 if (nat_cavity == 1)
595 /* Copy the location of the cavity */
596 copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
598 else
600 /* Determine the center of mass of the last molecule */
601 clear_rvec(x_init);
602 mass_tot = 0;
603 for (i = 0; i < nat_cavity; i++)
605 for (d = 0; d < DIM; d++)
607 x_init[d] +=
608 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
610 mass_tot += mass_cavity[i];
612 for (d = 0; d < DIM; d++)
614 x_init[d] /= mass_tot;
618 /* Generate coordinates within |dx|=drmax of x_init */
621 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, rnd_seed_tpi, rnd);
622 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, rnd_seed_tpi, rnd+2);
623 for (d = 0; d < DIM; d++)
625 dx[d] = (2*rnd[d] - 1)*drmax;
628 while (norm2(dx) > drmax*drmax);
629 rvec_add(x_init, dx, x_tp);
632 if (a_tp1 - a_tp0 == 1)
634 /* Insert a single atom, just copy the insertion location */
635 copy_rvec(x_tp, state_global->x[a_tp0]);
637 else
639 /* Copy the coordinates from the top file */
640 for (i = a_tp0; i < a_tp1; i++)
642 copy_rvec(x_mol[i-a_tp0], state_global->x[i]);
644 /* Rotate the molecule randomly */
645 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, rnd_seed_tpi, rnd);
646 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, rnd_seed_tpi, rnd+2);
647 rotate_conf(a_tp1-a_tp0, state_global->x+a_tp0, NULL,
648 2*M_PI*rnd[0],
649 2*M_PI*rnd[1],
650 2*M_PI*rnd[2]);
651 /* Shift to the insertion location */
652 for (i = a_tp0; i < a_tp1; i++)
654 rvec_inc(state_global->x[i], x_tp);
658 /* Clear some matrix variables */
659 clear_mat(force_vir);
660 clear_mat(shake_vir);
661 clear_mat(vir);
662 clear_mat(pres);
664 /* Set the charge group center of mass of the test particle */
665 copy_rvec(x_init, fr->cg_cm[top->cgs.nr-1]);
667 /* Calc energy (no forces) on new positions.
668 * Since we only need the intermolecular energy
669 * and the RF exclusion terms of the inserted molecule occur
670 * within a single charge group we can pass NULL for the graph.
671 * This also avoids shifts that would move charge groups
672 * out of the box. */
673 /* Make do_force do a single node force calculation */
674 cr->nnodes = 1;
675 do_force(fplog, cr, inputrec,
676 step, nrnb, wcycle, top, &top_global->groups,
677 state_global->box, state_global->x, &state_global->hist,
678 f, force_vir, mdatoms, enerd, fcd,
679 state_global->lambda,
680 NULL, fr, NULL, mu_tot, t, NULL, NULL, FALSE,
681 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
682 (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS : 0) |
683 (bStateChanged ? GMX_FORCE_STATECHANGED : 0));
684 cr->nnodes = nnodes;
685 bStateChanged = FALSE;
686 bNS = FALSE;
688 /* Calculate long range corrections to pressure and energy */
689 calc_dispcorr(inputrec, fr, top_global->natoms, state_global->box,
690 lambda, pres, vir, &prescorr, &enercorr, &dvdlcorr);
691 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
692 enerd->term[F_DISPCORR] = enercorr;
693 enerd->term[F_EPOT] += enercorr;
694 enerd->term[F_PRES] += prescorr;
695 enerd->term[F_DVDL_VDW] += dvdlcorr;
697 epot = enerd->term[F_EPOT];
698 bEnergyOutOfBounds = FALSE;
700 /* If the compiler doesn't optimize this check away
701 * we catch the NAN energies.
702 * The epot>GMX_REAL_MAX check catches inf values,
703 * which should nicely result in embU=0 through the exp below,
704 * but it does not hurt to check anyhow.
706 /* Non-bonded Interaction usually diverge at r=0.
707 * With tabulated interaction functions the first few entries
708 * should be capped in a consistent fashion between
709 * repulsion, dispersion and Coulomb to avoid accidental
710 * negative values in the total energy.
711 * The table generation code in tables.c does this.
712 * With user tbales the user should take care of this.
714 if (epot != epot || epot > GMX_REAL_MAX)
716 bEnergyOutOfBounds = TRUE;
718 if (bEnergyOutOfBounds)
720 if (debug)
722 fprintf(debug, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, (int)step, epot);
724 embU = 0;
726 else
728 embU = exp(-beta*epot);
729 sum_embU += embU;
730 /* Determine the weighted energy contributions of each energy group */
731 e = 0;
732 sum_UgembU[e++] += epot*embU;
733 if (fr->bBHAM)
735 for (i = 0; i < ngid; i++)
737 sum_UgembU[e++] +=
738 enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)]*embU;
741 else
743 for (i = 0; i < ngid; i++)
745 sum_UgembU[e++] +=
746 enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)]*embU;
749 if (bDispCorr)
751 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
753 if (bCharge)
755 for (i = 0; i < ngid; i++)
757 sum_UgembU[e++] += enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] * embU;
759 if (bRFExcl)
761 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
763 if (EEL_FULL(fr->eeltype))
765 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
770 if (embU == 0 || beta*epot > bU_bin_limit)
772 bin[0]++;
774 else
776 i = (int)((bU_logV_bin_limit
777 - (beta*epot - logV + refvolshift))*invbinw
778 + 0.5);
779 if (i < 0)
781 i = 0;
783 if (i >= nbin)
785 realloc_bins(&bin, &nbin, i+10);
787 bin[i]++;
790 if (debug)
792 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
793 (int)step, epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
796 if (dump_pdb && epot <= dump_ener)
798 sprintf(str, "t%g_step%d.pdb", t, (int)step);
799 sprintf(str2, "t: %f step %d ener: %f", t, (int)step, epot);
800 write_sto_conf_mtop(str, str2, top_global, state_global->x, state_global->v,
801 inputrec->ePBC, state_global->box);
804 step++;
805 if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
807 /* Skip all steps assigned to the other MPI ranks */
808 step += (cr->nnodes - 1)*stepblocksize;
812 if (PAR(cr))
814 /* When running in parallel sum the energies over the processes */
815 gmx_sumd(1, &sum_embU, cr);
816 gmx_sumd(nener, sum_UgembU, cr);
819 frame++;
820 V_all += V;
821 VembU_all += V*sum_embU/nsteps;
823 if (fp_tpi)
825 if (bVerbose || frame%10 == 0 || frame < 10)
827 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
828 -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
831 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
833 VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
834 sum_embU == 0 ? 20/beta : -log(sum_embU/nsteps)/beta,
835 sum_embU/nsteps, V);
836 for (e = 0; e < nener; e++)
838 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
840 fprintf(fp_tpi, "\n");
841 fflush(fp_tpi);
844 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
845 } /* End of the loop */
846 walltime_accounting_end(walltime_accounting);
848 close_trj(status);
850 if (fp_tpi != NULL)
852 xvgrclose(fp_tpi);
855 if (fplog != NULL)
857 fprintf(fplog, "\n");
858 fprintf(fplog, " <V> = %12.5e nm^3\n", V_all/frame);
859 fprintf(fplog, " <mu> = %12.5e kJ/mol\n", -log(VembU_all/V_all)/beta);
862 /* Write the Boltzmann factor histogram */
863 if (PAR(cr))
865 /* When running in parallel sum the bins over the processes */
866 i = nbin;
867 global_max(cr, &i);
868 realloc_bins(&bin, &nbin, i);
869 gmx_sumd(nbin, bin, cr);
871 if (MASTER(cr))
873 fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
874 "TPI energy distribution",
875 "\\betaU - log(V/<V>)", "count", oenv);
876 sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
877 xvgr_subtitle(fp_tpi, str, oenv);
878 xvgr_legend(fp_tpi, 2, (const char **)tpid_leg, oenv);
879 for (i = nbin-1; i > 0; i--)
881 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
882 fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
883 bUlogV,
884 (int)(bin[i]+0.5),
885 bin[i]*exp(-bUlogV)*V_all/VembU_all);
887 xvgrclose(fp_tpi);
889 sfree(bin);
891 sfree(sum_UgembU);
893 walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);
895 return 0;
898 } // namespace gmx