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,2019, 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_mdrun
55 #include "gromacs/commandline/filenm.h"
56 #include "gromacs/domdec/dlbtiming.h"
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/ewald/pme.h"
59 #include "gromacs/fileio/confio.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/dispersioncorrection.h"
70 #include "gromacs/mdlib/energyoutput.h"
71 #include "gromacs/mdlib/force.h"
72 #include "gromacs/mdlib/force_flags.h"
73 #include "gromacs/mdlib/mdatoms.h"
74 #include "gromacs/mdlib/tgroup.h"
75 #include "gromacs/mdlib/update.h"
76 #include "gromacs/mdlib/vsite.h"
77 #include "gromacs/mdrunutility/printtime.h"
78 #include "gromacs/mdtypes/commrec.h"
79 #include "gromacs/mdtypes/forcerec.h"
80 #include "gromacs/mdtypes/group.h"
81 #include "gromacs/mdtypes/inputrec.h"
82 #include "gromacs/mdtypes/md_enums.h"
83 #include "gromacs/mdtypes/mdrunoptions.h"
84 #include "gromacs/mdtypes/state.h"
85 #include "gromacs/pbcutil/pbc.h"
86 #include "gromacs/random/threefry.h"
87 #include "gromacs/random/uniformrealdistribution.h"
88 #include "gromacs/timing/wallcycle.h"
89 #include "gromacs/timing/walltime_accounting.h"
90 #include "gromacs/topology/mtop_util.h"
91 #include "gromacs/trajectory/trajectoryframe.h"
92 #include "gromacs/utility/cstringutil.h"
93 #include "gromacs/utility/fatalerror.h"
94 #include "gromacs/utility/gmxassert.h"
95 #include "gromacs/utility/logger.h"
96 #include "gromacs/utility/smalloc.h"
98 #include "integrator.h"
100 //! Global max algorithm
101 static void global_max(t_commrec
*cr
, int *n
)
105 snew(sum
, cr
->nnodes
);
106 sum
[cr
->nodeid
] = *n
;
107 gmx_sumi(cr
->nnodes
, sum
, cr
);
108 for (i
= 0; i
< cr
->nnodes
; i
++)
110 *n
= std::max(*n
, sum
[i
]);
116 //! Reallocate arrays.
117 static void realloc_bins(double **bin
, int *nbin
, int nbin_new
)
121 if (nbin_new
!= *nbin
)
123 srenew(*bin
, nbin_new
);
124 for (i
= *nbin
; i
< nbin_new
; i
++)
135 // TODO: Convert to use the nbnxm kernels by putting the system and the teset molecule on two separate search grids
140 PaddedVector
<gmx::RVec
> f
{};
141 real lambda
, t
, temp
, beta
, drmax
, epot
;
142 double embU
, sum_embU
, *sum_UgembU
, V
, V_all
, VembU_all
;
145 gmx_bool bDispCorr
, bCharge
, bRFExcl
, bNotLastFrame
, bStateChanged
, bNS
;
146 tensor force_vir
, shake_vir
, vir
, pres
;
147 int a_tp0
, a_tp1
, ngid
, gid_tp
, nener
, e
;
149 rvec mu_tot
, x_init
, dx
, x_tp
;
151 int64_t frame_step_prev
, frame_step
;
152 int64_t nsteps
, stepblocksize
= 0, step
;
155 FILE *fp_tpi
= nullptr;
156 char *ptr
, *dump_pdb
, **leg
, str
[STRLEN
], str2
[STRLEN
];
157 double dbl
, dump_ener
;
159 int nat_cavity
= 0, d
;
160 real
*mass_cavity
= nullptr, mass_tot
;
162 double invbinw
, *bin
, refvolshift
, logV
, bUlogV
;
163 real prescorr
, enercorr
, dvdlcorr
;
164 gmx_bool bEnergyOutOfBounds
;
165 const char *tpid_leg
[2] = {"direct", "reweighted"};
166 auto mdatoms
= mdAtoms
->mdatoms();
168 GMX_UNUSED_VALUE(outputProvider
);
170 GMX_LOG(mdlog
.info
).asParagraph().
171 appendText("Note that it is planned to change the command gmx mdrun -tpi "
172 "(and -tpic) to make the functionality available in a different "
173 "form in a future version of GROMACS, e.g. gmx test-particle-insertion.");
175 /* Since there is no upper limit to the insertion energies,
176 * we need to set an upper limit for the distribution output.
178 real bU_bin_limit
= 50;
179 real bU_logV_bin_limit
= bU_bin_limit
+ 10;
181 if (inputrec
->cutoff_scheme
== ecutsVERLET
)
183 gmx_fatal(FARGS
, "TPI does not work (yet) with the Verlet cut-off scheme");
188 gmx_mtop_generate_local_top(*top_global
, &top
, inputrec
->efep
!= efepNO
);
190 SimulationGroups
*groups
= &top_global
->groups
;
192 bCavity
= (inputrec
->eI
== eiTPIC
);
195 ptr
= getenv("GMX_TPIC_MASSES");
202 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
203 * The center of mass of the last atoms is then used for TPIC.
206 while (sscanf(ptr
, "%20lf%n", &dbl
, &i
) > 0)
208 srenew(mass_cavity
, nat_cavity
+1);
209 mass_cavity
[nat_cavity
] = dbl
;
210 fprintf(fplog
, "mass[%d] = %f\n",
211 nat_cavity
+1, mass_cavity
[nat_cavity
]);
217 gmx_fatal(FARGS
, "Found %d masses in GMX_TPIC_MASSES", nat_cavity
);
223 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
224 state_global->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
225 /* We never need full pbc for TPI */
227 /* Determine the temperature for the Boltzmann weighting */
228 temp
= inputrec
->opts
.ref_t
[0];
231 for (i
= 1; (i
< inputrec
->opts
.ngtc
); i
++)
233 if (inputrec
->opts
.ref_t
[i
] != temp
)
235 fprintf(fplog
, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
236 fprintf(stderr
, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
240 "\n The temperature for test particle insertion is %.3f K\n\n",
243 beta
= 1.0/(BOLTZ
*temp
);
245 /* Number of insertions per frame */
246 nsteps
= inputrec
->nsteps
;
248 /* Use the same neighborlist with more insertions points
249 * in a sphere of radius drmax around the initial point
251 /* This should be a proper mdp parameter */
252 drmax
= inputrec
->rtpi
;
254 /* An environment variable can be set to dump all configurations
255 * to pdb with an insertion energy <= this value.
257 dump_pdb
= getenv("GMX_TPI_DUMP");
261 sscanf(dump_pdb
, "%20lf", &dump_ener
);
264 atoms2md(top_global
, inputrec
, -1, nullptr, top_global
->natoms
, mdAtoms
);
265 update_mdatoms(mdatoms
, inputrec
->fepvals
->init_lambda
);
267 f
.resizeWithPadding(top_global
->natoms
);
269 /* Print to log file */
270 walltime_accounting_start_time(walltime_accounting
);
271 wallcycle_start(wcycle
, ewcRUN
);
272 print_start(fplog
, cr
, walltime_accounting
, "Test Particle Insertion");
274 /* The last charge group is the group to be inserted */
275 const t_atoms
&atomsToInsert
= top_global
->moltype
[top_global
->molblock
.back().type
].atoms
;
276 a_tp0
= top_global
->natoms
- atomsToInsert
.nr
;
277 a_tp1
= top_global
->natoms
;
280 fprintf(debug
, "TPI atoms %d-%d\n", a_tp0
, a_tp1
);
283 GMX_RELEASE_ASSERT(inputrec
->rcoulomb
<= inputrec
->rlist
&& inputrec
->rvdw
<= inputrec
->rlist
, "Twin-range interactions are not supported with TPI");
285 snew(x_mol
, a_tp1
-a_tp0
);
287 bDispCorr
= (inputrec
->eDispCorr
!= edispcNO
);
289 auto x
= makeArrayRef(state_global
->x
);
290 for (i
= a_tp0
; i
< a_tp1
; i
++)
292 /* Copy the coordinates of the molecule to be insterted */
293 copy_rvec(x
[i
], x_mol
[i
-a_tp0
]);
294 /* Check if we need to print electrostatic energies */
295 bCharge
|= (mdatoms
->chargeA
[i
] != 0 ||
296 ((mdatoms
->chargeB
!= nullptr) && mdatoms
->chargeB
[i
] != 0));
298 bRFExcl
= (bCharge
&& EEL_RF(fr
->ic
->eeltype
));
300 // TODO: Calculate the center of geometry of the molecule to insert
302 calc_cgcm(fplog
, cg_tp
, cg_tp
+1, &(top
.cgs
), state_global
->x
.rvec_array(), fr
->cg_cm
);
305 if (norm(fr
->cg_cm
[cg_tp
]) > 0.5*inputrec
->rlist
&& fplog
)
307 fprintf(fplog
, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
308 fprintf(stderr
, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
313 /* Center the molecule to be inserted at zero */
314 for (i
= 0; i
< a_tp1
-a_tp0
; i
++)
316 rvec_dec(x_mol
[i
], fr
->cg_cm
[cg_tp
]);
323 fprintf(fplog
, "\nWill insert %d atoms %s partial charges\n",
324 a_tp1
-a_tp0
, bCharge
? "with" : "without");
326 fprintf(fplog
, "\nWill insert %" PRId64
" times in each frame of %s\n",
327 nsteps
, opt2fn("-rerun", nfile
, fnm
));
332 if (inputrec
->nstlist
> 1)
334 if (drmax
== 0 && a_tp1
-a_tp0
== 1)
336 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
);
340 fprintf(fplog
, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec
->nstlist
, drmax
);
348 fprintf(fplog
, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax
);
352 ngid
= groups
->groups
[SimulationAtomGroupType::EnergyOutput
].nr
;
353 // TODO: Figure out which energy group to use
355 gid_tp
= GET_CGINFO_GID(fr
->cginfo
[cg_tp
]);
370 if (EEL_FULL(fr
->ic
->eeltype
))
375 snew(sum_UgembU
, nener
);
377 /* Copy the random seed set by the user */
378 seed
= inputrec
->ld_seed
;
380 gmx::ThreeFry2x64
<16> rng(seed
, gmx::RandomDomain::TestParticleInsertion
); // 16 bits internal counter => 2^16 * 2 = 131072 values per stream
381 gmx::UniformRealDistribution
<real
> dist
;
385 fp_tpi
= xvgropen(opt2fn("-tpi", nfile
, fnm
),
386 "TPI energies", "Time (ps)",
387 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv
);
388 xvgr_subtitle(fp_tpi
, "f. are averages over one frame", oenv
);
391 sprintf(str
, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
392 leg
[e
++] = gmx_strdup(str
);
393 sprintf(str
, "f. -kT log<e\\S-\\betaU\\N>");
394 leg
[e
++] = gmx_strdup(str
);
395 sprintf(str
, "f. <e\\S-\\betaU\\N>");
396 leg
[e
++] = gmx_strdup(str
);
397 sprintf(str
, "f. V");
398 leg
[e
++] = gmx_strdup(str
);
399 sprintf(str
, "f. <Ue\\S-\\betaU\\N>");
400 leg
[e
++] = gmx_strdup(str
);
401 for (i
= 0; i
< ngid
; i
++)
403 sprintf(str
, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
404 *(groups
->groupNames
[groups
->groups
[SimulationAtomGroupType::EnergyOutput
].nm_ind
[i
]]));
405 leg
[e
++] = gmx_strdup(str
);
409 sprintf(str
, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
410 leg
[e
++] = gmx_strdup(str
);
414 for (i
= 0; i
< ngid
; i
++)
416 sprintf(str
, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
417 *(groups
->groupNames
[groups
->groups
[SimulationAtomGroupType::EnergyOutput
].nm_ind
[i
]]));
418 leg
[e
++] = gmx_strdup(str
);
422 sprintf(str
, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
423 leg
[e
++] = gmx_strdup(str
);
425 if (EEL_FULL(fr
->ic
->eeltype
))
427 sprintf(str
, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
428 leg
[e
++] = gmx_strdup(str
);
431 xvgr_legend(fp_tpi
, 4+nener
, leg
, oenv
);
432 for (i
= 0; i
< 4+nener
; i
++)
446 /* Avoid frame step numbers <= -1 */
447 frame_step_prev
= -1;
449 bNotLastFrame
= read_first_frame(oenv
, &status
, opt2fn("-rerun", nfile
, fnm
),
450 &rerun_fr
, TRX_NEED_X
);
453 if (rerun_fr
.natoms
- (bCavity
? nat_cavity
: 0) !=
454 mdatoms
->nr
- (a_tp1
- a_tp0
))
456 gmx_fatal(FARGS
, "Number of atoms in trajectory (%d)%s "
457 "is not equal the number in the run input file (%d) "
458 "minus the number of atoms to insert (%d)\n",
459 rerun_fr
.natoms
, bCavity
? " minus one" : "",
460 mdatoms
->nr
, a_tp1
-a_tp0
);
463 refvolshift
= log(det(rerun_fr
.box
));
465 switch (inputrec
->eI
)
468 stepblocksize
= inputrec
->nstlist
;
474 gmx_fatal(FARGS
, "Unknown integrator %s", ei_names
[inputrec
->eI
]);
477 while (bNotLastFrame
)
479 frame_step
= rerun_fr
.step
;
480 if (frame_step
<= frame_step_prev
)
482 /* We don't have step number in the trajectory file,
483 * or we have constant or decreasing step numbers.
484 * Ensure we have increasing step numbers, since we use
485 * the step numbers as a counter for random numbers.
487 frame_step
= frame_step_prev
+ 1;
489 frame_step_prev
= frame_step
;
491 lambda
= rerun_fr
.lambda
;
495 for (e
= 0; e
< nener
; e
++)
500 /* Copy the coordinates from the input trajectory */
501 auto x
= makeArrayRef(state_global
->x
);
502 for (i
= 0; i
< rerun_fr
.natoms
; i
++)
504 copy_rvec(rerun_fr
.x
[i
], x
[i
]);
506 copy_mat(rerun_fr
.box
, state_global
->box
);
508 V
= det(state_global
->box
);
511 bStateChanged
= TRUE
;
514 step
= cr
->nodeid
*stepblocksize
;
515 while (step
< nsteps
)
517 /* Restart random engine using the frame and insertion step
519 * Note that we need to draw several random values per iteration,
520 * but by using the internal subcounter functionality of ThreeFry2x64
521 * we can draw 131072 unique 64-bit values before exhausting
522 * the stream. This is a huge margin, and if something still goes
523 * wrong you will get an exception when the stream is exhausted.
525 rng
.restart(frame_step
, step
);
526 dist
.reset(); // erase any memory in the distribution
530 /* Random insertion in the whole volume */
531 bNS
= (step
% inputrec
->nstlist
== 0);
534 /* Generate a random position in the box */
535 for (d
= 0; d
< DIM
; d
++)
537 x_init
[d
] = dist(rng
)*state_global
->box
[d
][d
];
541 if (inputrec
->nstlist
== 1)
543 copy_rvec(x_init
, x_tp
);
547 /* Generate coordinates within |dx|=drmax of x_init */
550 for (d
= 0; d
< DIM
; d
++)
552 dx
[d
] = (2*dist(rng
) - 1)*drmax
;
555 while (norm2(dx
) > drmax
*drmax
);
556 rvec_add(x_init
, dx
, x_tp
);
561 /* Random insertion around a cavity location
562 * given by the last coordinate of the trajectory.
568 /* Copy the location of the cavity */
569 copy_rvec(rerun_fr
.x
[rerun_fr
.natoms
-1], x_init
);
573 /* Determine the center of mass of the last molecule */
576 for (i
= 0; i
< nat_cavity
; i
++)
578 for (d
= 0; d
< DIM
; d
++)
581 mass_cavity
[i
]*rerun_fr
.x
[rerun_fr
.natoms
-nat_cavity
+i
][d
];
583 mass_tot
+= mass_cavity
[i
];
585 for (d
= 0; d
< DIM
; d
++)
587 x_init
[d
] /= mass_tot
;
591 /* Generate coordinates within |dx|=drmax of x_init */
594 for (d
= 0; d
< DIM
; d
++)
596 dx
[d
] = (2*dist(rng
) - 1)*drmax
;
599 while (norm2(dx
) > drmax
*drmax
);
600 rvec_add(x_init
, dx
, x_tp
);
603 if (a_tp1
- a_tp0
== 1)
605 /* Insert a single atom, just copy the insertion location */
606 copy_rvec(x_tp
, x
[a_tp0
]);
610 /* Copy the coordinates from the top file */
611 for (i
= a_tp0
; i
< a_tp1
; i
++)
613 copy_rvec(x_mol
[i
-a_tp0
], x
[i
]);
615 /* Rotate the molecule randomly */
616 real angleX
= 2*M_PI
*dist(rng
);
617 real angleY
= 2*M_PI
*dist(rng
);
618 real angleZ
= 2*M_PI
*dist(rng
);
619 rotate_conf(a_tp1
-a_tp0
, state_global
->x
.rvec_array()+a_tp0
, nullptr,
620 angleX
, angleY
, angleZ
);
621 /* Shift to the insertion location */
622 for (i
= a_tp0
; i
< a_tp1
; i
++)
624 rvec_inc(x
[i
], x_tp
);
628 /* Clear some matrix variables */
629 clear_mat(force_vir
);
630 clear_mat(shake_vir
);
634 /* Set the center of geometry mass of the test molecule */
635 // TODO: Compute and set the COG
637 copy_rvec(x_init
, fr
->cg_cm
[top
.cgs
.nr
-1]);
640 /* Calc energy (no forces) on new positions.
641 * Since we only need the intermolecular energy
642 * and the RF exclusion terms of the inserted molecule occur
643 * within a single charge group we can pass NULL for the graph.
644 * This also avoids shifts that would move charge groups
646 /* Make do_force do a single node force calculation */
649 // TPI might place a particle so close that the potential
650 // is infinite. Since this is intended to happen, we
651 // temporarily suppress any exceptions that the processor
652 // might raise, then restore the old behaviour.
653 std::fenv_t floatingPointEnvironment
;
654 std::feholdexcept(&floatingPointEnvironment
);
655 do_force(fplog
, cr
, ms
, inputrec
, nullptr, nullptr, imdSession
,
657 step
, nrnb
, wcycle
, &top
,
658 state_global
->box
, state_global
->x
.arrayRefWithPadding(), &state_global
->hist
,
659 f
.arrayRefWithPadding(), force_vir
, mdatoms
, enerd
, fcd
,
660 state_global
->lambda
,
661 nullptr, fr
, ppForceWorkload
, nullptr, mu_tot
, t
, nullptr,
662 GMX_FORCE_NONBONDED
| GMX_FORCE_ENERGY
|
663 (bNS
? GMX_FORCE_DYNAMICBOX
| GMX_FORCE_NS
: 0) |
664 (bStateChanged
? GMX_FORCE_STATECHANGED
: 0),
665 DDBalanceRegionHandler(nullptr));
666 std::feclearexcept(FE_DIVBYZERO
| FE_INVALID
| FE_OVERFLOW
);
667 std::feupdateenv(&floatingPointEnvironment
);
670 bStateChanged
= 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
)
707 fprintf(debug
, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t
, static_cast<int>(step
), epot
);
713 // Exponent argument is fine in SP range, but output can be in DP range
714 embU
= exp(static_cast<double>(-beta
*epot
));
716 /* Determine the weighted energy contributions of each energy group */
718 sum_UgembU
[e
++] += epot
*embU
;
721 for (i
= 0; i
< ngid
; i
++)
724 enerd
->grpp
.ener
[egBHAMSR
][GID(i
, gid_tp
, ngid
)]*embU
;
729 for (i
= 0; i
< ngid
; i
++)
732 enerd
->grpp
.ener
[egLJSR
][GID(i
, gid_tp
, ngid
)]*embU
;
737 sum_UgembU
[e
++] += enerd
->term
[F_DISPCORR
]*embU
;
741 for (i
= 0; i
< ngid
; i
++)
743 sum_UgembU
[e
++] += enerd
->grpp
.ener
[egCOULSR
][GID(i
, gid_tp
, ngid
)] * embU
;
747 sum_UgembU
[e
++] += enerd
->term
[F_RF_EXCL
]*embU
;
749 if (EEL_FULL(fr
->ic
->eeltype
))
751 sum_UgembU
[e
++] += enerd
->term
[F_COUL_RECIP
]*embU
;
756 if (embU
== 0 || beta
*epot
> bU_bin_limit
)
762 i
= gmx::roundToInt((bU_logV_bin_limit
763 - (beta
*epot
- logV
+ refvolshift
))*invbinw
);
770 realloc_bins(&bin
, &nbin
, i
+10);
777 fprintf(debug
, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
778 static_cast<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
, static_cast<int>(step
));
784 sprintf(str2
, "t: %f step %d ener: %f", t
, static_cast<int>(step
), epot
);
785 write_sto_conf_mtop(str
, str2
, top_global
, state_global
->x
.rvec_array(), state_global
->v
.rvec_array(),
786 inputrec
->ePBC
, state_global
->box
);
790 if ((step
/stepblocksize
) % cr
->nnodes
!= cr
->nodeid
)
792 /* Skip all steps assigned to the other MPI ranks */
793 step
+= (cr
->nnodes
- 1)*stepblocksize
;
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
);
806 VembU_all
+= V
*sum_embU
/nsteps
;
810 if (mdrunOptions
.verbose
|| 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
,
821 for (e
= 0; e
< nener
; e
++)
823 fprintf(fp_tpi
, " %12.5e", sum_UgembU
[e
]/nsteps
);
825 fprintf(fp_tpi
, "\n");
829 bNotLastFrame
= read_next_frame(oenv
, status
, &rerun_fr
);
830 } /* End of the loop */
831 walltime_accounting_end_time(walltime_accounting
);
835 if (fp_tpi
!= nullptr)
840 if (fplog
!= nullptr)
842 fprintf(fplog
, "\n");
843 fprintf(fplog
, " <V> = %12.5e nm^3\n", V_all
/frame
);
844 const double mu
= -log(VembU_all
/V_all
)/beta
;
845 fprintf(fplog
, " <mu> = %12.5e kJ/mol\n", mu
);
847 if (!std::isfinite(mu
))
849 fprintf(fplog
, "\nThe computed chemical potential is not finite - consider increasing the number of steps and/or the number of frames to insert into.\n");
853 /* Write the Boltzmann factor histogram */
856 /* When running in parallel sum the bins over the processes */
859 realloc_bins(&bin
, &nbin
, i
);
860 gmx_sumd(nbin
, bin
, cr
);
864 fp_tpi
= xvgropen(opt2fn("-tpid", nfile
, fnm
),
865 "TPI energy distribution",
866 "\\betaU - log(V/<V>)", "count", oenv
);
867 sprintf(str
, "number \\betaU > %g: %9.3e", bU_bin_limit
, bin
[0]);
868 xvgr_subtitle(fp_tpi
, str
, oenv
);
869 xvgr_legend(fp_tpi
, 2, tpid_leg
, oenv
);
870 for (i
= nbin
-1; i
> 0; i
--)
872 bUlogV
= -i
/invbinw
+ bU_logV_bin_limit
- refvolshift
+ log(V_all
/frame
);
873 fprintf(fp_tpi
, "%6.2f %10d %12.5e\n",
876 bin
[i
]*exp(-bUlogV
)*V_all
/VembU_all
);
884 walltime_accounting_set_nsteps_done(walltime_accounting
, frame
*inputrec
->nsteps
);