2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief This file defines the integrator for test particle insertion
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_mdlib
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
)
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
]);
110 //! Reallocate arrays.
111 static void realloc_bins(double **bin
, int *nbin
, int nbin_new
)
115 if (nbin_new
!= *nbin
)
117 srenew(*bin
, nbin_new
);
118 for (i
= *nbin
; i
< nbin_new
; i
++)
129 /*! \brief Do test particle insertion.
130 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
131 const gmx_multi_sim_t *,
132 const gmx::MDLogger &mdlog,
133 int nfile, const t_filenm fnm[],
134 const gmx_output_env_t *oenv,
135 const MdrunOptions &mdrunOptions,
136 gmx_vsite_t *vsite, gmx_constr_t constr,
137 gmx::IMDOutputProvider *outputProvider,
138 t_inputrec *inputrec,
139 gmx_mtop_t *top_global, t_fcdata *fcd,
140 t_state *state_global,
141 gmx::MDAtoms *mdAtoms,
142 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
145 const ReplicaExchangeParameters &replExParams,
146 gmx_membed_t gmx_unused *membed,
147 gmx_walltime_accounting_t walltime_accounting)
149 double do_tpi(FILE *fplog
, t_commrec
*cr
,
150 const gmx_multisim_t
*ms
,
151 const gmx::MDLogger gmx_unused
&mdlog
,
152 int nfile
, const t_filenm fnm
[],
153 const gmx_output_env_t
*oenv
,
154 const MdrunOptions
&mdrunOptions
,
155 gmx_vsite_t gmx_unused
*vsite
, gmx_constr_t gmx_unused constr
,
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 gmx::MDAtoms
*mdAtoms
,
162 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
,
164 const ReplicaExchangeParameters gmx_unused
&replExParams
,
165 gmx_membed_t gmx_unused
*membed
,
166 gmx_walltime_accounting_t walltime_accounting
)
169 gmx_groups_t
*groups
;
170 gmx_enerdata_t
*enerd
;
171 PaddedRVecVector f
{};
172 real lambda
, t
, temp
, beta
, drmax
, epot
;
173 double embU
, sum_embU
, *sum_UgembU
, V
, V_all
, VembU_all
;
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
;
180 rvec mu_tot
, x_init
, dx
, x_tp
;
182 gmx_int64_t frame_step_prev
, frame_step
;
183 gmx_int64_t nsteps
, stepblocksize
= 0, step
;
186 FILE *fp_tpi
= nullptr;
187 char *ptr
, *dump_pdb
, **leg
, str
[STRLEN
], str2
[STRLEN
];
188 double dbl
, dump_ener
;
190 int nat_cavity
= 0, d
;
191 real
*mass_cavity
= nullptr, mass_tot
;
193 double invbinw
, *bin
, refvolshift
, logV
, bUlogV
;
194 real prescorr
, enercorr
, dvdlcorr
;
195 gmx_bool bEnergyOutOfBounds
;
196 const char *tpid_leg
[2] = {"direct", "reweighted"};
197 auto mdatoms
= mdAtoms
->mdatoms();
199 GMX_UNUSED_VALUE(outputProvider
);
201 /* Since there is no upper limit to the insertion energies,
202 * we need to set an upper limit for the distribution output.
204 real bU_bin_limit
= 50;
205 real bU_logV_bin_limit
= bU_bin_limit
+ 10;
207 if (inputrec
->cutoff_scheme
== ecutsVERLET
)
209 gmx_fatal(FARGS
, "TPI does not work (yet) with the Verlet cut-off scheme");
214 top
= gmx_mtop_generate_local_top(top_global
, inputrec
->efep
!= efepNO
);
216 groups
= &top_global
->groups
;
218 bCavity
= (inputrec
->eI
== eiTPIC
);
221 ptr
= getenv("GMX_TPIC_MASSES");
228 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
229 * The center of mass of the last atoms is then used for TPIC.
232 while (sscanf(ptr
, "%20lf%n", &dbl
, &i
) > 0)
234 srenew(mass_cavity
, nat_cavity
+1);
235 mass_cavity
[nat_cavity
] = dbl
;
236 fprintf(fplog
, "mass[%d] = %f\n",
237 nat_cavity
+1, mass_cavity
[nat_cavity
]);
243 gmx_fatal(FARGS
, "Found %d masses in GMX_TPIC_MASSES", nat_cavity
);
249 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
250 state_global->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
251 /* We never need full pbc for TPI */
253 /* Determine the temperature for the Boltzmann weighting */
254 temp
= inputrec
->opts
.ref_t
[0];
257 for (i
= 1; (i
< inputrec
->opts
.ngtc
); i
++)
259 if (inputrec
->opts
.ref_t
[i
] != temp
)
261 fprintf(fplog
, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
262 fprintf(stderr
, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
266 "\n The temperature for test particle insertion is %.3f K\n\n",
269 beta
= 1.0/(BOLTZ
*temp
);
271 /* Number of insertions per frame */
272 nsteps
= inputrec
->nsteps
;
274 /* Use the same neighborlist with more insertions points
275 * in a sphere of radius drmax around the initial point
277 /* This should be a proper mdp parameter */
278 drmax
= inputrec
->rtpi
;
280 /* An environment variable can be set to dump all configurations
281 * to pdb with an insertion energy <= this value.
283 dump_pdb
= getenv("GMX_TPI_DUMP");
287 sscanf(dump_pdb
, "%20lf", &dump_ener
);
290 atoms2md(top_global
, inputrec
, -1, nullptr, top_global
->natoms
, mdAtoms
);
291 update_mdatoms(mdatoms
, inputrec
->fepvals
->init_lambda
);
294 init_enerdata(groups
->grps
[egcENER
].nr
, inputrec
->fepvals
->n_lambda
, enerd
);
295 /* We need to allocate one element extra, since we might use
296 * (unaligned) 4-wide SIMD loads to access rvec entries.
298 f
.resize(gmx::paddedRVecVectorSize(top_global
->natoms
));
300 /* Print to log file */
301 walltime_accounting_start(walltime_accounting
);
302 wallcycle_start(wcycle
, ewcRUN
);
303 print_start(fplog
, cr
, walltime_accounting
, "Test Particle Insertion");
305 /* The last charge group is the group to be inserted */
306 cg_tp
= top
->cgs
.nr
- 1;
307 a_tp0
= top
->cgs
.index
[cg_tp
];
308 a_tp1
= top
->cgs
.index
[cg_tp
+1];
311 fprintf(debug
, "TPI cg %d, atoms %d-%d\n", cg_tp
, a_tp0
, a_tp1
);
314 GMX_RELEASE_ASSERT(inputrec
->rcoulomb
<= inputrec
->rlist
&& inputrec
->rvdw
<= inputrec
->rlist
, "Twin-range interactions are not supported with TPI");
316 snew(x_mol
, a_tp1
-a_tp0
);
318 bDispCorr
= (inputrec
->eDispCorr
!= edispcNO
);
320 for (i
= a_tp0
; i
< a_tp1
; i
++)
322 /* Copy the coordinates of the molecule to be insterted */
323 copy_rvec(state_global
->x
[i
], x_mol
[i
-a_tp0
]);
324 /* Check if we need to print electrostatic energies */
325 bCharge
|= (mdatoms
->chargeA
[i
] != 0 ||
326 (mdatoms
->chargeB
&& mdatoms
->chargeB
[i
] != 0));
328 bRFExcl
= (bCharge
&& EEL_RF(fr
->ic
->eeltype
));
330 calc_cgcm(fplog
, cg_tp
, cg_tp
+1, &(top
->cgs
), as_rvec_array(state_global
->x
.data()), fr
->cg_cm
);
333 if (norm(fr
->cg_cm
[cg_tp
]) > 0.5*inputrec
->rlist
&& fplog
)
335 fprintf(fplog
, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
336 fprintf(stderr
, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
341 /* Center the molecule to be inserted at zero */
342 for (i
= 0; i
< a_tp1
-a_tp0
; i
++)
344 rvec_dec(x_mol
[i
], fr
->cg_cm
[cg_tp
]);
350 fprintf(fplog
, "\nWill insert %d atoms %s partial charges\n",
351 a_tp1
-a_tp0
, bCharge
? "with" : "without");
353 fprintf(fplog
, "\nWill insert %d times in each frame of %s\n",
354 (int)nsteps
, opt2fn("-rerun", nfile
, fnm
));
359 if (inputrec
->nstlist
> 1)
361 if (drmax
== 0 && a_tp1
-a_tp0
== 1)
363 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 fprintf(fplog
, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec
->nstlist
, drmax
);
375 fprintf(fplog
, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax
);
379 ngid
= groups
->grps
[egcENER
].nr
;
380 gid_tp
= GET_CGINFO_GID(fr
->cginfo
[cg_tp
]);
393 if (EEL_FULL(fr
->ic
->eeltype
))
398 snew(sum_UgembU
, nener
);
400 /* Copy the random seed set by the user */
401 seed
= inputrec
->ld_seed
;
403 gmx::ThreeFry2x64
<16> rng(seed
, gmx::RandomDomain::TestParticleInsertion
); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
404 gmx::UniformRealDistribution
<real
> dist
;
408 fp_tpi
= xvgropen(opt2fn("-tpi", nfile
, fnm
),
409 "TPI energies", "Time (ps)",
410 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv
);
411 xvgr_subtitle(fp_tpi
, "f. are averages over one frame", oenv
);
414 sprintf(str
, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
415 leg
[e
++] = gmx_strdup(str
);
416 sprintf(str
, "f. -kT log<e\\S-\\betaU\\N>");
417 leg
[e
++] = gmx_strdup(str
);
418 sprintf(str
, "f. <e\\S-\\betaU\\N>");
419 leg
[e
++] = gmx_strdup(str
);
420 sprintf(str
, "f. V");
421 leg
[e
++] = gmx_strdup(str
);
422 sprintf(str
, "f. <Ue\\S-\\betaU\\N>");
423 leg
[e
++] = gmx_strdup(str
);
424 for (i
= 0; i
< ngid
; i
++)
426 sprintf(str
, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
427 *(groups
->grpname
[groups
->grps
[egcENER
].nm_ind
[i
]]));
428 leg
[e
++] = gmx_strdup(str
);
432 sprintf(str
, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
433 leg
[e
++] = gmx_strdup(str
);
437 for (i
= 0; i
< ngid
; i
++)
439 sprintf(str
, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
440 *(groups
->grpname
[groups
->grps
[egcENER
].nm_ind
[i
]]));
441 leg
[e
++] = gmx_strdup(str
);
445 sprintf(str
, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
446 leg
[e
++] = gmx_strdup(str
);
448 if (EEL_FULL(fr
->ic
->eeltype
))
450 sprintf(str
, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
451 leg
[e
++] = gmx_strdup(str
);
454 xvgr_legend(fp_tpi
, 4+nener
, (const char**)leg
, oenv
);
455 for (i
= 0; i
< 4+nener
; i
++)
469 /* Avoid frame step numbers <= -1 */
470 frame_step_prev
= -1;
472 bNotLastFrame
= read_first_frame(oenv
, &status
, opt2fn("-rerun", nfile
, fnm
),
473 &rerun_fr
, TRX_NEED_X
);
476 if (rerun_fr
.natoms
- (bCavity
? nat_cavity
: 0) !=
477 mdatoms
->nr
- (a_tp1
- a_tp0
))
479 gmx_fatal(FARGS
, "Number of atoms in trajectory (%d)%s "
480 "is not equal the number in the run input file (%d) "
481 "minus the number of atoms to insert (%d)\n",
482 rerun_fr
.natoms
, bCavity
? " minus one" : "",
483 mdatoms
->nr
, a_tp1
-a_tp0
);
486 refvolshift
= log(det(rerun_fr
.box
));
488 switch (inputrec
->eI
)
491 stepblocksize
= inputrec
->nstlist
;
497 gmx_fatal(FARGS
, "Unknown integrator %s", ei_names
[inputrec
->eI
]);
500 while (bNotLastFrame
)
502 frame_step
= rerun_fr
.step
;
503 if (frame_step
<= frame_step_prev
)
505 /* We don't have step number in the trajectory file,
506 * or we have constant or decreasing step numbers.
507 * Ensure we have increasing step numbers, since we use
508 * the step numbers as a counter for random numbers.
510 frame_step
= frame_step_prev
+ 1;
512 frame_step_prev
= frame_step
;
514 lambda
= rerun_fr
.lambda
;
518 for (e
= 0; e
< nener
; e
++)
523 /* Copy the coordinates from the input trajectory */
524 for (i
= 0; i
< rerun_fr
.natoms
; i
++)
526 copy_rvec(rerun_fr
.x
[i
], state_global
->x
[i
]);
528 copy_mat(rerun_fr
.box
, state_global
->box
);
530 V
= det(state_global
->box
);
533 bStateChanged
= TRUE
;
536 step
= cr
->nodeid
*stepblocksize
;
537 while (step
< nsteps
)
539 /* Restart random engine using the frame and insertion step
541 * Note that we need to draw several random values per iteration,
542 * but by using the internal subcounter functionality of ThreeFry2x64
543 * we can draw 131072 unique 64-bit values before exhausting
544 * the stream. This is a huge margin, and if something still goes
545 * wrong you will get an exception when the stream is exhausted.
547 rng
.restart(frame_step
, step
);
548 dist
.reset(); // erase any memory in the distribution
552 /* Random insertion in the whole volume */
553 bNS
= (step
% inputrec
->nstlist
== 0);
556 /* Generate a random position in the box */
557 for (d
= 0; d
< DIM
; d
++)
559 x_init
[d
] = dist(rng
)*state_global
->box
[d
][d
];
563 if (inputrec
->nstlist
== 1)
565 copy_rvec(x_init
, x_tp
);
569 /* Generate coordinates within |dx|=drmax of x_init */
572 for (d
= 0; d
< DIM
; d
++)
574 dx
[d
] = (2*dist(rng
) - 1)*drmax
;
577 while (norm2(dx
) > drmax
*drmax
);
578 rvec_add(x_init
, dx
, x_tp
);
583 /* Random insertion around a cavity location
584 * given by the last coordinate of the trajectory.
590 /* Copy the location of the cavity */
591 copy_rvec(rerun_fr
.x
[rerun_fr
.natoms
-1], x_init
);
595 /* Determine the center of mass of the last molecule */
598 for (i
= 0; i
< nat_cavity
; i
++)
600 for (d
= 0; d
< DIM
; d
++)
603 mass_cavity
[i
]*rerun_fr
.x
[rerun_fr
.natoms
-nat_cavity
+i
][d
];
605 mass_tot
+= mass_cavity
[i
];
607 for (d
= 0; d
< DIM
; d
++)
609 x_init
[d
] /= mass_tot
;
613 /* Generate coordinates within |dx|=drmax of x_init */
616 for (d
= 0; d
< DIM
; d
++)
618 dx
[d
] = (2*dist(rng
) - 1)*drmax
;
621 while (norm2(dx
) > drmax
*drmax
);
622 rvec_add(x_init
, dx
, x_tp
);
625 if (a_tp1
- a_tp0
== 1)
627 /* Insert a single atom, just copy the insertion location */
628 copy_rvec(x_tp
, state_global
->x
[a_tp0
]);
632 /* Copy the coordinates from the top file */
633 for (i
= a_tp0
; i
< a_tp1
; i
++)
635 copy_rvec(x_mol
[i
-a_tp0
], state_global
->x
[i
]);
637 /* Rotate the molecule randomly */
638 rotate_conf(a_tp1
-a_tp0
, as_rvec_array(state_global
->x
.data())+a_tp0
, nullptr,
642 /* Shift to the insertion location */
643 for (i
= a_tp0
; i
< a_tp1
; i
++)
645 rvec_inc(state_global
->x
[i
], x_tp
);
649 /* Clear some matrix variables */
650 clear_mat(force_vir
);
651 clear_mat(shake_vir
);
655 /* Set the charge group center of mass of the test particle */
656 copy_rvec(x_init
, fr
->cg_cm
[top
->cgs
.nr
-1]);
658 /* Calc energy (no forces) on new positions.
659 * Since we only need the intermolecular energy
660 * and the RF exclusion terms of the inserted molecule occur
661 * within a single charge group we can pass NULL for the graph.
662 * This also avoids shifts that would move charge groups
664 /* Make do_force do a single node force calculation */
666 do_force(fplog
, cr
, ms
, inputrec
,
667 step
, nrnb
, wcycle
, top
, &top_global
->groups
,
668 state_global
->box
, state_global
->x
, &state_global
->hist
,
669 f
, force_vir
, mdatoms
, enerd
, fcd
,
670 state_global
->lambda
,
671 nullptr, fr
, nullptr, mu_tot
, t
, nullptr,
672 GMX_FORCE_NONBONDED
| GMX_FORCE_ENERGY
|
673 (bNS
? GMX_FORCE_DYNAMICBOX
| GMX_FORCE_NS
: 0) |
674 (bStateChanged
? GMX_FORCE_STATECHANGED
: 0),
675 DdOpenBalanceRegionBeforeForceComputation::no
,
676 DdCloseBalanceRegionAfterForceComputation::no
);
678 bStateChanged
= FALSE
;
681 /* Calculate long range corrections to pressure and energy */
682 calc_dispcorr(inputrec
, fr
, state_global
->box
,
683 lambda
, pres
, vir
, &prescorr
, &enercorr
, &dvdlcorr
);
684 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
685 enerd
->term
[F_DISPCORR
] = enercorr
;
686 enerd
->term
[F_EPOT
] += enercorr
;
687 enerd
->term
[F_PRES
] += prescorr
;
688 enerd
->term
[F_DVDL_VDW
] += dvdlcorr
;
690 epot
= enerd
->term
[F_EPOT
];
691 bEnergyOutOfBounds
= FALSE
;
693 /* If the compiler doesn't optimize this check away
694 * we catch the NAN energies.
695 * The epot>GMX_REAL_MAX check catches inf values,
696 * which should nicely result in embU=0 through the exp below,
697 * but it does not hurt to check anyhow.
699 /* Non-bonded Interaction usually diverge at r=0.
700 * With tabulated interaction functions the first few entries
701 * should be capped in a consistent fashion between
702 * repulsion, dispersion and Coulomb to avoid accidental
703 * negative values in the total energy.
704 * The table generation code in tables.c does this.
705 * With user tbales the user should take care of this.
707 if (epot
!= epot
|| epot
> GMX_REAL_MAX
)
709 bEnergyOutOfBounds
= TRUE
;
711 if (bEnergyOutOfBounds
)
715 fprintf(debug
, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t
, (int)step
, epot
);
721 embU
= exp(-beta
*epot
);
723 /* Determine the weighted energy contributions of each energy group */
725 sum_UgembU
[e
++] += epot
*embU
;
728 for (i
= 0; i
< ngid
; i
++)
731 enerd
->grpp
.ener
[egBHAMSR
][GID(i
, gid_tp
, ngid
)]*embU
;
736 for (i
= 0; i
< ngid
; i
++)
739 enerd
->grpp
.ener
[egLJSR
][GID(i
, gid_tp
, ngid
)]*embU
;
744 sum_UgembU
[e
++] += enerd
->term
[F_DISPCORR
]*embU
;
748 for (i
= 0; i
< ngid
; i
++)
750 sum_UgembU
[e
++] += enerd
->grpp
.ener
[egCOULSR
][GID(i
, gid_tp
, ngid
)] * embU
;
754 sum_UgembU
[e
++] += enerd
->term
[F_RF_EXCL
]*embU
;
756 if (EEL_FULL(fr
->ic
->eeltype
))
758 sum_UgembU
[e
++] += enerd
->term
[F_COUL_RECIP
]*embU
;
763 if (embU
== 0 || beta
*epot
> bU_bin_limit
)
769 i
= (int)((bU_logV_bin_limit
770 - (beta
*epot
- logV
+ refvolshift
))*invbinw
778 realloc_bins(&bin
, &nbin
, i
+10);
785 fprintf(debug
, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
786 (int)step
, epot
, x_tp
[XX
], x_tp
[YY
], x_tp
[ZZ
]);
789 if (dump_pdb
&& epot
<= dump_ener
)
791 sprintf(str
, "t%g_step%d.pdb", t
, (int)step
);
792 sprintf(str2
, "t: %f step %d ener: %f", t
, (int)step
, epot
);
793 write_sto_conf_mtop(str
, str2
, top_global
, as_rvec_array(state_global
->x
.data()), as_rvec_array(state_global
->v
.data()),
794 inputrec
->ePBC
, state_global
->box
);
798 if ((step
/stepblocksize
) % cr
->nnodes
!= cr
->nodeid
)
800 /* Skip all steps assigned to the other MPI ranks */
801 step
+= (cr
->nnodes
- 1)*stepblocksize
;
807 /* When running in parallel sum the energies over the processes */
808 gmx_sumd(1, &sum_embU
, cr
);
809 gmx_sumd(nener
, sum_UgembU
, cr
);
814 VembU_all
+= V
*sum_embU
/nsteps
;
818 if (mdrunOptions
.verbose
|| frame
%10 == 0 || frame
< 10)
820 fprintf(stderr
, "mu %10.3e <mu> %10.3e\n",
821 -log(sum_embU
/nsteps
)/beta
, -log(VembU_all
/V_all
)/beta
);
824 fprintf(fp_tpi
, "%10.3f %12.5e %12.5e %12.5e %12.5e",
826 VembU_all
== 0 ? 20/beta
: -log(VembU_all
/V_all
)/beta
,
827 sum_embU
== 0 ? 20/beta
: -log(sum_embU
/nsteps
)/beta
,
829 for (e
= 0; e
< nener
; e
++)
831 fprintf(fp_tpi
, " %12.5e", sum_UgembU
[e
]/nsteps
);
833 fprintf(fp_tpi
, "\n");
837 bNotLastFrame
= read_next_frame(oenv
, status
, &rerun_fr
);
838 } /* End of the loop */
839 walltime_accounting_end(walltime_accounting
);
843 if (fp_tpi
!= nullptr)
848 if (fplog
!= nullptr)
850 fprintf(fplog
, "\n");
851 fprintf(fplog
, " <V> = %12.5e nm^3\n", V_all
/frame
);
852 fprintf(fplog
, " <mu> = %12.5e kJ/mol\n", -log(VembU_all
/V_all
)/beta
);
855 /* Write the Boltzmann factor histogram */
858 /* When running in parallel sum the bins over the processes */
861 realloc_bins(&bin
, &nbin
, i
);
862 gmx_sumd(nbin
, bin
, cr
);
866 fp_tpi
= xvgropen(opt2fn("-tpid", nfile
, fnm
),
867 "TPI energy distribution",
868 "\\betaU - log(V/<V>)", "count", oenv
);
869 sprintf(str
, "number \\betaU > %g: %9.3e", bU_bin_limit
, bin
[0]);
870 xvgr_subtitle(fp_tpi
, str
, oenv
);
871 xvgr_legend(fp_tpi
, 2, (const char **)tpid_leg
, oenv
);
872 for (i
= nbin
-1; i
> 0; i
--)
874 bUlogV
= -i
/invbinw
+ bU_logV_bin_limit
- refvolshift
+ log(V_all
/frame
);
875 fprintf(fp_tpi
, "%6.2f %10d %12.5e\n",
878 bin
[i
]*exp(-bUlogV
)*V_all
/VembU_all
);
886 walltime_accounting_set_nsteps_done(walltime_accounting
, frame
*inputrec
->nsteps
);